mapFields.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) 2016-2023 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "mapFields.H"
29 #include "meshToMesh.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace functionObjects
37 {
38  defineTypeNameAndDebug(mapFields, 0);
39  addToRunTimeSelectionTable(functionObject, mapFields, dictionary);
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 void Foam::functionObjects::mapFields::createInterpolation
47 (
48  const dictionary& dict
49 )
50 {
51  const fvMesh& meshTarget = mesh_;
52  const word mapRegionName(dict.get<word>("mapRegion"));
53 
54  Info<< name() << ':' << nl
55  << " Reading mesh " << mapRegionName << endl;
56 
57  mapRegionPtr_.reset
58  (
59  new fvMesh
60  (
61  IOobject
62  (
63  mapRegionName,
64  meshTarget.time().constant(),
65  meshTarget.time(),
67  )
68  )
69  );
70 
71  const fvMesh& mapRegion = mapRegionPtr_();
72  const word mapMethodName(dict.get<word>("mapMethod"));
73  if (!meshToMesh::interpolationMethodNames_.found(mapMethodName))
74  {
76  << type() << " " << name() << ": unknown map method "
77  << mapMethodName << nl
78  << "Available methods include: "
80  << exit(FatalError);
81  }
82 
84  (
86  );
87 
88  // Lookup corresponding AMI method
89  word patchMapMethodName = meshToMesh::interpolationMethodAMI(mapMethod);
90 
91  // Optionally override
92  if (dict.readIfPresent("patchMapMethod", patchMapMethodName))
93  {
94  Info<< " Patch mapping method: " << patchMapMethodName << endl;
95  }
96 
97  Info<< " Creating mesh to mesh interpolation" << endl;
98 
99  if (dict.get<bool>("consistent"))
100  {
101  interpPtr_.reset
102  (
103  new meshToMesh
104  (
105  mapRegion,
106  meshTarget,
107  mapMethodName,
108  patchMapMethodName
109  )
110  );
111  }
112  else
113  {
114  HashTable<word> patchMap;
115 
116  bool createPatchMap = dict.getOrDefault<bool>("createPatchMap", false);
117 
118  const entry* ePtr = dict.findEntry("patchMap");
119 
120  if (createPatchMap && ePtr)
121  {
123  << "Requested 'createPatchMap' but 'patchMap' entry provided. "
124  << "Please remove one of these entries"
125  << exit(FatalIOError);
126  }
127 
128  if (!createPatchMap && !ePtr)
129  {
131  << "Either the 'createPatchMap' or 'patchMap' entry must be "
132  << "provided when not using the 'consistent' mapping option"
133  << exit(FatalIOError);
134  }
135 
136  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
137  const polyBoundaryMesh& pbmT = mapRegion.boundaryMesh();
138  DynamicList<label> cuttingIndices;
139 
140  if (createPatchMap)
141  {
142  Info<< " Creating patchMap and cuttingPatches" << endl;
143 
144  if (dict.found("cuttingPatches"))
145  {
147  << "Ignoring user supplied cuttingPatches when "
148  << "createPatchMap option is active"
149  << endl;
150  }
151 
152  forAll(pbmT, patchiT)
153  {
154  const polyPatch& ppT = pbmT[patchiT];
155 
156  if (!polyPatch::constraintType(ppT.type()))
157  {
158  const word& patchNameT = ppT.name();
159  const label patchi = pbm.findPatchID(patchNameT);
160 
161  if (patchi == -1)
162  {
163  Info<< " Adding to cuttingPatches: "
164  << ppT.name() << endl;
165 
166  cuttingIndices.push_back(ppT.index());
167  }
168  else
169  {
170  if (returnReduce(ppT.size(), sumOp<label>()) > 0)
171  {
172  Info<< " Adding to patchMap: " << ppT.name()
173  << endl;
174 
175  patchMap.set(patchNameT, patchNameT);
176  }
177  }
178  }
179  }
180  }
181  else
182  {
183  // Read patch map
184  patchMap = HashTable<word>(ePtr->stream());
185 
186  // Read cutting patches using wordRe's
187  const wordRes cuttingPatchRes = dict.get<wordRes>("cuttingPatches");
188  cuttingIndices.push_back(pbmT.indices(cuttingPatchRes));
189  }
190 
191  const wordList cuttingPatches(pbmT.names(), cuttingIndices);
192 
193  // Checks
194  {
195  // Find patch names that exist on target mesh that are not included
196  // in the patchMap
197  wordHashSet unknownPatchNames;
198  for (const auto& ppT : pbmT)
199  {
200  if
201  (
202  !polyPatch::constraintType(ppT.type())
203  && !patchMap.found(ppT.name())
204  && returnReduce(ppT.size(), sumOp<label>()) > 0
205  )
206  {
207  unknownPatchNames.insert(ppT.name());
208  }
209  }
210 
211  for (const label patchiT : cuttingIndices)
212  {
213  const word& patchName = pbmT[patchiT].name();
214 
215  unknownPatchNames.erase(patchName);
216 
217  if (patchMap.found(patchName))
218  {
219  Info<< " Removing cutting patch from patchMap: "
220  << patchName << endl;
221 
222  patchMap.erase(patchName);
223  }
224  }
225 
226  if (unknownPatchNames.size())
227  {
229  << "Patches not present in source mesh. "
230  << "Add to cuttingPatches? Patches: "
231  << unknownPatchNames
232  << exit(FatalError);
233  }
234  }
235 
236  interpPtr_.reset
237  (
238  new meshToMesh
239  (
240  mapRegion,
241  meshTarget,
242  mapMethodName,
243  patchMapMethodName,
244  patchMap,
245  cuttingPatches
246  )
247  );
248  }
249 }
250 
251 
252 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
253 
255 (
256  const word& name,
257  const Time& runTime,
258  const dictionary& dict
259 )
260 :
261  fvMeshFunctionObject(name, runTime, dict),
262  mapRegionPtr_(),
263  interpPtr_(),
264  fieldNames_()
265 {
266  read(dict);
267 }
268 
269 
270 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
271 
273 {
275  {
276  dict.readEntry("fields", fieldNames_);
277  createInterpolation(dict);
278 
279  return true;
280  }
281 
282  return false;
283 }
284 
285 
287 {
288  Log << type() << " " << name() << " execute:" << nl;
289 
290  bool ok = false;
291 
292  ok = mapFieldType<scalar>() || ok;
293  ok = mapFieldType<vector>() || ok;
294  ok = mapFieldType<sphericalTensor>() || ok;
295  ok = mapFieldType<symmTensor>() || ok;
296  ok = mapFieldType<tensor>() || ok;
297 
298  if (log)
299  {
300  if (!ok)
301  {
302  Info<< " none" << nl;
303  }
305  Info<< endl;
306  }
307  return true;
308 }
309 
310 
312 {
313  Log << type() << " " << name() << " write:" << nl;
314 
315  bool ok = false;
316 
317  ok = writeFieldType<scalar>() || ok;
318  ok = writeFieldType<vector>() || ok;
319  ok = writeFieldType<sphericalTensor>() || ok;
320  ok = writeFieldType<symmTensor>() || ok;
321  ok = writeFieldType<tensor>() || ok;
322 
323  if (log)
324  {
325  if (!ok)
326  {
327  Info<< " none" << nl;
328  }
329 
330  Info<< endl;
331  }
332 
333  return true;
334 }
335 
336 
337 // ************************************************************************* //
const polyBoundaryMesh & pbm
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
static const Enum< interpolationMethod > interpolationMethodNames_
Definition: meshToMesh.H:77
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
mapFields(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: mapFields.C:248
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:232
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
const word & name() const noexcept
Return the name of this functionObject.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
virtual bool read(const dictionary &)
Read the mapFields data.
Definition: mapFields.C:265
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual const word & type() const =0
Runtime type information.
virtual bool write()
Calculate the mapFields and write.
Definition: mapFields.C:304
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition: HashSet.H:73
static word interpolationMethodAMI(const interpolationMethod method)
Conversion between mesh and patch interpolation methods.
Definition: meshToMesh.C:639
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
List< word > wordList
List of word.
Definition: fileName.H:59
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition: polyPatch.C:269
#define WarningInFunction
Report a warning using Foam::Warning.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
#define Log
Definition: PDRblock.C:28
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual bool read(const dictionary &dict)
Read optional controls.
virtual bool execute()
Execute, currently does nothing.
Definition: mapFields.C:279
interpolationMethod
Enumeration specifying interpolation method.
Definition: meshToMesh.H:69
bool found
const fvMesh & mesh_
Reference to the fvMesh.
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...