findRefCell.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2018 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "findRefCell.H"
30 
31 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
32 
34 (
35  const volScalarField& field,
36  const volScalarField& fieldRef,
37  const dictionary& dict,
38  label& refCelli,
39  scalar& refValue,
40  const bool forceReference
41 )
42 {
43  if (fieldRef.needReference() || forceReference)
44  {
45  const word refCellName = field.name() + "RefCell";
46  const word refPointName = field.name() + "RefPoint";
47  const word refValueName = field.name() + "RefValue";
48 
49  if (dict.found(refCellName))
50  {
51  if (Pstream::master())
52  {
53  dict.readEntry(refCellName, refCelli);
54 
55  if (refCelli < 0 || refCelli >= field.mesh().nCells())
56  {
58  << "Illegal master cellID " << refCelli
59  << ". Should be 0.." << field.mesh().nCells()
60  << exit(FatalIOError);
61  }
62  }
63  else
64  {
65  refCelli = -1;
66  }
67  }
68  else if (dict.found(refPointName))
69  {
70  point refPointi(dict.get<point>(refPointName));
71 
72  // Try fast approximate search avoiding octree construction
73  refCelli = field.mesh().findCell(refPointi, polyMesh::FACE_PLANES);
74 
75  label hasRef = (refCelli >= 0 ? 1 : 0);
76  label sumHasRef = returnReduce<label>(hasRef, sumOp<label>());
77 
78  // If reference cell not found, use octree search
79  // with cell tet-decompositoin
80  if (sumHasRef != 1)
81  {
82  refCelli = field.mesh().findCell(refPointi);
83 
84  hasRef = (refCelli >= 0 ? 1 : 0);
85  sumHasRef = returnReduce<label>(hasRef, sumOp<label>());
86  }
87 
88  if (sumHasRef != 1)
89  {
91  << "Unable to set reference cell for field " << field.name()
92  << nl << " Reference point " << refPointName
93  << " " << refPointi
94  << " found on " << sumHasRef << " domains (should be one)"
95  << nl << exit(FatalIOError);
96  }
97  }
98  else
99  {
101  << "Unable to set reference cell for field " << field.name()
102  << nl
103  << " Please supply either " << refCellName
104  << " or " << refPointName << nl << exit(FatalIOError);
105  }
106 
107  dict.readEntry(refValueName, refValue);
108 
109  return true;
110  }
111  else
112  {
113  refCelli = -1;
114  }
115 
116  return false;
117 }
118 
119 
120 bool Foam::setRefCell
121 (
122  const volScalarField& field,
123  const dictionary& dict,
124  label& refCelli,
125  scalar& refValue,
126  const bool forceReference
127 )
128 {
129  return setRefCell(field, field, dict, refCelli, refValue, forceReference);
130 }
131 
132 
133 Foam::scalar Foam::getRefCellValue
134 (
135  const volScalarField& field,
136  const label refCelli
137 )
138 {
139  #ifdef FULLDEBUG
140  if (refCelli >= field.mesh().nCells())
141  {
143  << "Illegal reference cellID " << refCelli
144  << ". Mesh has " << field.mesh().nCells() << ". cells."
145  << exit(FatalError);
146  }
147  #endif
148  scalar refCellValue = (refCelli >= 0 ? field[refCelli] : 0.0);
149 
150  // Currently distributing the value to all processors. This is generally
151  // not needed since only the processor holding the reference cell needs
152  // it. Tdb.
153  return returnReduce(refCellValue, sumOp<scalar>());
154 }
155 
156 
157 // ************************************************************************* //
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
findRefCell.H
Find the reference cell nearest (in index) to the given cell but which is not on a cyclic,...
Foam::FatalIOError
IOerror FatalIOError
Foam::GeometricField::needReference
bool needReference() const
Does the field need a reference level for solution.
Definition: GeometricField.C:949
Foam::sumOp
Definition: ops.H:213
Foam::setRefCell
bool setRefCell(const volScalarField &field, const volScalarField &fieldRef, const dictionary &dict, label &refCelli, scalar &refValue, const bool forceReference=false)
If the field fieldRef needs referencing find the reference cell nearest.
Definition: findRefCell.C:34
field
rDeltaTY field()
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::getRefCellValue
scalar getRefCellValue(const volScalarField &field, const label refCelli)
Return the current value of field in the reference cell.
Definition: findRefCell.C:134
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::Vector< scalar >
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::GeometricField< scalar, fvPatchField, volMesh >