fvFieldReconstructorTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fvFieldReconstructor.H"
27 #include "Time.H"
28 #include "PtrList.H"
29 #include "fvPatchFields.H"
30 #include "emptyFvPatch.H"
31 #include "emptyFvPatchField.H"
32 #include "emptyFvsPatchField.H"
33 #include "processorCyclicFvPatch.H"
34 #include "reverseFieldMapper.H"
35 #include "setSizeFieldMapper.H"
36 #include "stringOps.H"
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
40 template<class FieldType>
41 bool Foam::fvFieldReconstructor::reconstructs
42 (
43  const IOobjectList& objects,
44  const HashSet<word>& selectedFields
45 )
46 {
47  IOobjectList fields = objects.lookupClass(FieldType::typeName);
48 
49  if (fields.size() && selectedFields.empty())
50  {
51  return true;
52  }
53 
54  forAllConstIter(IOobjectList, fields, fieldIter)
55  {
56  if (selectedFields.found(fieldIter()->name()))
57  {
58  return true;
59  }
60  }
61 
62  return false;
63 }
64 
65 
66 template<class Type>
67 void Foam::fvFieldReconstructor::rmapFaceToFace
68 (
69  Field<Type>& toField,
70  const Field<Type>& fromField,
71  const labelUList& addressing,
72  const bool isFlux
73 )
74 {
75  forAll(addressing, i)
76  {
77  toField[mag(addressing[i]) - 1] =
78  (isFlux && addressing[i] < 0 ? -1 : +1)*fromField[i];
79  }
80 }
81 
82 
83 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
84 
85 template<class Type>
88 (
89  const IOobject& fieldIoObject
90 ) const
91 {
92  // Read the field for all the processors
93  PtrList<DimensionedField<Type, volMesh>> procFields(procMeshes_.size());
94  forAll(procMeshes_, proci)
95  {
96  procFields.set
97  (
98  proci,
100  (
101  IOobject
102  (
103  fieldIoObject.name(),
104  procMeshes_[proci].time().name(),
105  procMeshes_[proci],
108  false
109  ),
110  procMeshes_[proci]
111  )
112  );
113  }
114 
115  // Create the internalField
116  Field<Type> internalField(completeMesh_.nCells());
117 
118  forAll(procMeshes_, proci)
119  {
120  const DimensionedField<Type, volMesh>& procField = procFields[proci];
121 
122  // Set the cell values in the reconstructed field
123  internalField.rmap
124  (
125  procField.primitiveField(),
126  cellProcAddressing_[proci]
127  );
128  }
129 
131  (
133  (
134  IOobject
135  (
136  fieldIoObject.name(),
137  completeMesh_.time().name(),
138  completeMesh_,
141  false
142  ),
143  completeMesh_,
144  procFields[0].dimensions(),
145  internalField
146  )
147  );
148 }
149 
150 
151 template<class Type>
154 (
155  const IOobject& fieldIoObject
156 ) const
157 {
158  // Read the field for all the processors
159  PtrList<VolField<Type>> procFields(procMeshes_.size());
160  forAll(procMeshes_, proci)
161  {
162  procFields.set
163  (
164  proci,
165  new VolField<Type>
166  (
167  IOobject
168  (
169  fieldIoObject.name(),
170  procMeshes_[proci].time().name(),
171  procMeshes_[proci],
174  false
175  ),
176  procMeshes_[proci]
177  )
178  );
179  }
180 
181  // Create the internalField
182  Field<Type> internalField(completeMesh_.nCells());
183 
184  // Create the patch fields
185  PtrList<fvPatchField<Type>> patchFields(completeMesh_.boundary().size());
186 
187  forAll(procFields, proci)
188  {
189  const VolField<Type>& procField =
190  procFields[proci];
191 
192  // Set the cell values in the reconstructed field
193  internalField.rmap
194  (
195  procField.primitiveField(),
196  cellProcAddressing_[proci]
197  );
198 
199  // Set the boundary patch values in the reconstructed field
200  forAll(procField.boundaryField(), procPatchi)
201  {
202  const fvPatch& procPatch =
203  procMeshes_[proci].boundary()[procPatchi];
204 
205  const label completePatchi = completePatchID(proci, procPatchi);
206 
207  if (completePatchi == procPatchi)
208  {
209  if (!patchFields(completePatchi))
210  {
211  patchFields.set
212  (
213  completePatchi,
215  (
216  procField.boundaryField()[procPatchi],
217  completeMesh_.boundary()[completePatchi],
220  (
221  completeMesh_.boundary()[completePatchi].size()
222  )
223  )
224  );
225  }
226 
227  patchFields[completePatchi].map
228  (
229  procField.boundaryField()[procPatchi],
231  (
232  faceProcAddressingBf_[proci][procPatchi] - 1
233  )
234  );
235  }
236  else if (isA<processorCyclicFvPatch>(procPatch))
237  {
238  if (!patchFields(completePatchi))
239  {
240  patchFields.set
241  (
242  completePatchi,
244  (
245  completeMesh_.boundary()[completePatchi].type(),
246  completeMesh_.boundary()[completePatchi],
248  )
249  );
250  }
251 
252  if (patchFields[completePatchi].overridesConstraint())
253  {
254  OStringStream str;
255  str << "\nThe field \"" << procFields[0].name()
256  << "\" on cyclic patch \""
257  << patchFields[completePatchi].patch().name()
258  << "\" cannot be reconstructed as it is not a cyclic "
259  << "patch field. A \"patchType cyclic;\" setting has "
260  << "been used to override the cyclic patch type.\n\n"
261  << "Cyclic patches like this with non-cyclic boundary "
262  << "conditions should be confined to a single "
263  << "processor using decomposition constraints.";
265  << stringOps::breakIntoIndentedLines(str.str()).c_str()
266  << exit(FatalError);
267  }
268 
269  patchFields[completePatchi].map
270  (
271  procField.boundaryField()[procPatchi],
273  (
274  faceProcAddressingBf_[proci][procPatchi] - 1
275  )
276  );
277  }
278  }
279  }
280 
281  // Construct and return the field
282  return tmp<VolField<Type>>
283  (
284  new VolField<Type>
285  (
286  IOobject
287  (
288  fieldIoObject.name(),
289  completeMesh_.time().name(),
290  completeMesh_,
293  false
294  ),
295  completeMesh_,
296  procFields[0].dimensions(),
297  internalField,
298  patchFields,
299  procFields[0].sources().table()
300  )
301  );
302 }
303 
304 
305 template<class Type>
308 (
309  const IOobject& fieldIoObject
310 ) const
311 {
312  // Read the field for all the processors
313  PtrList<SurfaceField<Type>> procFields(procMeshes_.size());
314  forAll(procMeshes_, proci)
315  {
316  procFields.set
317  (
318  proci,
320  (
321  IOobject
322  (
323  fieldIoObject.name(),
324  procMeshes_[proci].time().name(),
325  procMeshes_[proci],
328  false
329  ),
330  procMeshes_[proci]
331  )
332  );
333  }
334 
335  // Create the internalField
336  Field<Type> internalField(completeMesh_.nInternalFaces());
337 
338  // Create the patch fields
339  PtrList<fvsPatchField<Type>> patchFields(completeMesh_.boundary().size());
340 
341  forAll(procMeshes_, proci)
342  {
343  const SurfaceField<Type>& procField =
344  procFields[proci];
345 
346  // Set the internal face values in the reconstructed field
347  rmapFaceToFace
348  (
349  internalField,
350  procField.primitiveField(),
352  (
353  faceProcAddressing_[proci],
354  procMeshes_[proci].nInternalFaces()
355  ),
356  isFlux(procFields[proci])
357  );
358 
359  // Set the boundary patch values in the reconstructed field
360  forAll(procField.boundaryField(), procPatchi)
361  {
362  const fvPatch& procPatch =
363  procMeshes_[proci].boundary()[procPatchi];
364 
365  const label completePatchi = completePatchID(proci, procPatchi);
366 
367  if (completePatchi == procPatchi)
368  {
369  if (!patchFields(completePatchi))
370  {
371  patchFields.set
372  (
373  completePatchi,
375  (
376  procField.boundaryField()[procPatchi],
377  completeMesh_.boundary()[completePatchi],
380  (
381  completeMesh_.boundary()[completePatchi].size()
382  )
383  )
384  );
385  }
386 
387  patchFields[completePatchi].map
388  (
389  procField.boundaryField()[procPatchi],
391  (
392  faceProcAddressingBf_[proci][procPatchi] - 1
393  )
394  );
395  }
396  else if (isA<processorCyclicFvPatch>(procPatch))
397  {
398  if (!patchFields(completePatchi))
399  {
400  patchFields.set
401  (
402  completePatchi,
404  (
405  completeMesh_.boundary()[completePatchi].type(),
406  completeMesh_.boundary()[completePatchi],
408  )
409  );
410  }
411 
412  patchFields[completePatchi].map
413  (
414  procField.boundaryField()[procPatchi],
416  (
417  faceProcAddressingBf_[proci][procPatchi] - 1
418  )
419  );
420  }
421  else if (isA<processorFvPatch>(procPatch))
422  {
423  rmapFaceToFace
424  (
425  internalField,
426  procField.boundaryField()[procPatchi],
427  faceProcAddressingBf_[proci][procPatchi],
428  isFlux(procFields[proci])
429  );
430  }
431  }
432  }
433 
434  // Construct and return the field
435  return tmp<SurfaceField<Type>>
436  (
438  (
439  IOobject
440  (
441  fieldIoObject.name(),
442  completeMesh_.time().name(),
443  completeMesh_,
446  false
447  ),
448  completeMesh_,
449  procFields[0].dimensions(),
450  internalField,
451  patchFields
452  )
453  );
454 }
455 
456 
457 template<class Type>
459 (
460  const IOobjectList& objects,
461  const HashSet<word>& selectedFields
462 )
463 {
464  const word& fieldClassName = DimensionedField<Type, volMesh>::typeName;
465 
466  IOobjectList fields = objects.lookupClass(fieldClassName);
467 
468  if (fields.size())
469  {
470  Info<< nl << " Reconstructing " << fieldClassName << "s"
471  << nl << endl;
472 
473  forAllConstIter(IOobjectList, fields, fieldIter)
474  {
475  if
476  (
477  selectedFields.empty()
478  || selectedFields.found(fieldIter()->name())
479  )
480  {
481  Info<< " " << fieldIter()->name() << endl;
482 
483  reconstructVolInternalField<Type>(*fieldIter())().write();
484  }
485  }
486  }
487 }
488 
489 
490 template<class Type>
492 (
493  const IOobjectList& objects,
494  const HashSet<word>& selectedFields
495 )
496 {
497  const word& fieldClassName =
499 
500  IOobjectList fields = objects.lookupClass(fieldClassName);
501 
502  if (fields.size())
503  {
504  Info<< nl << " Reconstructing " << fieldClassName << "s"
505  << nl << endl;
506 
507  forAllConstIter(IOobjectList, fields, fieldIter)
508  {
509  if
510  (
511  selectedFields.empty()
512  || selectedFields.found(fieldIter()->name())
513  )
514  {
515  Info<< " " << fieldIter()->name() << endl;
516 
517  reconstructVolField<Type>(*fieldIter())().write();
518  }
519  }
520  }
521 }
522 
523 
524 template<class Type>
526 (
527  const IOobjectList& objects,
528  const HashSet<word>& selectedFields
529 )
530 {
531  const word& fieldClassName =
533 
534  IOobjectList fields = objects.lookupClass(fieldClassName);
535 
536  if (fields.size())
537  {
538  Info<< nl << " Reconstructing " << fieldClassName << "s"
539  << nl << endl;
540 
541  forAllConstIter(IOobjectList, fields, fieldIter)
542  {
543  if
544  (
545  selectedFields.empty()
546  || selectedFields.found(fieldIter()->name())
547  )
548  {
549  Info<< " " << fieldIter()->name() << endl;
550 
551  reconstructFvSurfaceField<Type>(*fieldIter())().write();
552  }
553  }
554  }
555 }
556 
557 
558 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const PrimitiveField< Type > & primitiveField() const
Return a const-reference to the primitive field.
Pre-declare SubField and related Field type.
Definition: Field.H:83
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:386
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const Internal::FieldType & primitiveField() const
Return a const-reference to the primitive field.
A HashTable with keys but without contents.
Definition: HashSet.H:62
bool empty() const
Return true if the hash table is empty.
Definition: HashTableI.H:72
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:307
Output to memory buffer stream.
Definition: OStringStream.H:52
string str() const
Return the string.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
A List obtained as a section of another List.
Definition: SubList.H:56
tmp< VolField< Type > > reconstructVolField(const IOobject &fieldIoObject) const
Read and reconstruct a volume field.
tmp< DimensionedField< Type, volMesh > > reconstructVolInternalField(const IOobject &fieldIoObject) const
Read and reconstruct a volume internal field.
void reconstructVolInternalFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Read, reconstruct and write all/selected volume internal fields.
void reconstructVolFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Read, reconstruct and write all/selected volume fields.
tmp< SurfaceField< Type > > reconstructFvSurfaceField(const IOobject &fieldIoObject) const
Read and reconstruct a surface field.
void reconstructFvSurfaceFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Read, reconstruct and write all/selected surface fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:91
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:86
Reverse field mapper.
Mapper which sets the field size. It does not actually map values.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:234
string breakIntoIndentedLines(const string &str, const string::size_type nLength=80, const string::size_type nIndent=0)
Break a string up into indented lines.
Definition: stringOps.C:963
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
bool isFlux(const DimensionedField< Type, surfaceMesh > &df)
Check if surfaceField is a flux.
Definition: surfaceFields.H:55
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
UList< label > labelUList
Definition: UList.H:65
static const char nl
Definition: Ostream.H:267
objects