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-2024 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 template<class Type>
85 Foam::fvFieldReconstructor::reconstructVolInternalField
86 (
87  const IOobject& fieldIoObject
88 ) const
89 {
90  // Read the field for all the processors
91  PtrList<DimensionedField<Type, volMesh>> procFields(procMeshes_.size());
92  forAll(procMeshes_, proci)
93  {
94  procFields.set
95  (
96  proci,
97  new DimensionedField<Type, volMesh>
98  (
99  IOobject
100  (
101  fieldIoObject.name(),
102  procMeshes_[proci].time().name(),
103  procMeshes_[proci],
106  false
107  ),
108  procMeshes_[proci]
109  )
110  );
111  }
112 
113  // Create the internalField
114  Field<Type> internalField(completeMesh_.nCells());
115 
116  forAll(procMeshes_, proci)
117  {
118  const DimensionedField<Type, volMesh>& procField = procFields[proci];
119 
120  // Set the cell values in the reconstructed field
121  internalField.rmap
122  (
123  procField.primitiveField(),
124  cellProcAddressing_[proci]
125  );
126  }
127 
128  return tmp<DimensionedField<Type, volMesh>>
129  (
130  new DimensionedField<Type, volMesh>
131  (
132  IOobject
133  (
134  fieldIoObject.name(),
135  completeMesh_.time().name(),
136  completeMesh_,
139  false
140  ),
141  completeMesh_,
142  procFields[0].dimensions(),
143  internalField
144  )
145  );
146 }
147 
148 
149 template<class Type>
151 Foam::fvFieldReconstructor::reconstructVolField
152 (
153  const IOobject& fieldIoObject
154 ) const
155 {
156  // Read the field for all the processors
157  PtrList<VolField<Type>> procFields(procMeshes_.size());
158  forAll(procMeshes_, proci)
159  {
160  procFields.set
161  (
162  proci,
163  new VolField<Type>
164  (
165  IOobject
166  (
167  fieldIoObject.name(),
168  procMeshes_[proci].time().name(),
169  procMeshes_[proci],
172  false
173  ),
174  procMeshes_[proci]
175  )
176  );
177  }
178 
179  // Create the internalField
180  Field<Type> internalField(completeMesh_.nCells());
181 
182  // Create the patch fields
183  PtrList<fvPatchField<Type>> patchFields(completeMesh_.boundary().size());
184 
185  forAll(procFields, proci)
186  {
187  const VolField<Type>& procField =
188  procFields[proci];
189 
190  // Set the cell values in the reconstructed field
191  internalField.rmap
192  (
193  procField.primitiveField(),
194  cellProcAddressing_[proci]
195  );
196 
197  // Set the boundary patch values in the reconstructed field
198  forAll(procField.boundaryField(), procPatchi)
199  {
200  const fvPatch& procPatch =
201  procMeshes_[proci].boundary()[procPatchi];
202 
203  const label completePatchi = completePatchID(proci, procPatchi);
204 
205  if (completePatchi == procPatchi)
206  {
207  if (!patchFields(completePatchi))
208  {
209  patchFields.set
210  (
211  completePatchi,
213  (
214  procField.boundaryField()[procPatchi],
215  completeMesh_.boundary()[completePatchi],
217  setSizeFieldMapper
218  (
219  completeMesh_.boundary()[completePatchi].size()
220  )
221  )
222  );
223  }
224 
225  patchFields[completePatchi].map
226  (
227  procField.boundaryField()[procPatchi],
228  reverseFieldMapper
229  (
230  faceProcAddressingBf_[proci][procPatchi] - 1
231  )
232  );
233  }
234  else if (isA<processorCyclicFvPatch>(procPatch))
235  {
236  if (!patchFields(completePatchi))
237  {
238  patchFields.set
239  (
240  completePatchi,
242  (
243  completeMesh_.boundary()[completePatchi].type(),
244  completeMesh_.boundary()[completePatchi],
246  )
247  );
248  }
249 
250  if (patchFields[completePatchi].overridesConstraint())
251  {
252  OStringStream str;
253  str << "\nThe field \"" << procFields[0].name()
254  << "\" on cyclic patch \""
255  << patchFields[completePatchi].patch().name()
256  << "\" cannot be reconstructed as it is not a cyclic "
257  << "patch field. A \"patchType cyclic;\" setting has "
258  << "been used to override the cyclic patch type.\n\n"
259  << "Cyclic patches like this with non-cyclic boundary "
260  << "conditions should be confined to a single "
261  << "processor using decomposition constraints.";
263  << stringOps::breakIntoIndentedLines(str.str()).c_str()
264  << exit(FatalError);
265  }
266 
267  patchFields[completePatchi].map
268  (
269  procField.boundaryField()[procPatchi],
270  reverseFieldMapper
271  (
272  faceProcAddressingBf_[proci][procPatchi] - 1
273  )
274  );
275  }
276  }
277  }
278 
279  // Construct and return the field
280  return tmp<VolField<Type>>
281  (
282  new VolField<Type>
283  (
284  IOobject
285  (
286  fieldIoObject.name(),
287  completeMesh_.time().name(),
288  completeMesh_,
291  false
292  ),
293  completeMesh_,
294  procFields[0].dimensions(),
295  internalField,
296  patchFields,
297  procFields[0].sources().table()
298  )
299  );
300 }
301 
302 
303 template<class Type>
305 Foam::fvFieldReconstructor::reconstructFvSurfaceField
306 (
307  const IOobject& fieldIoObject
308 ) const
309 {
310  // Read the field for all the processors
311  PtrList<SurfaceField<Type>> procFields(procMeshes_.size());
312  forAll(procMeshes_, proci)
313  {
314  procFields.set
315  (
316  proci,
317  new SurfaceField<Type>
318  (
319  IOobject
320  (
321  fieldIoObject.name(),
322  procMeshes_[proci].time().name(),
323  procMeshes_[proci],
326  false
327  ),
328  procMeshes_[proci]
329  )
330  );
331  }
332 
333  // Create the internalField
334  Field<Type> internalField(completeMesh_.nInternalFaces());
335 
336  // Create the patch fields
337  PtrList<fvsPatchField<Type>> patchFields(completeMesh_.boundary().size());
338 
339  forAll(procMeshes_, proci)
340  {
341  const SurfaceField<Type>& procField =
342  procFields[proci];
343 
344  // Set the internal face values in the reconstructed field
345  rmapFaceToFace
346  (
347  internalField,
348  procField.primitiveField(),
349  SubList<label>
350  (
351  faceProcAddressing_[proci],
352  procMeshes_[proci].nInternalFaces()
353  ),
354  isFlux(procFields[proci])
355  );
356 
357  // Set the boundary patch values in the reconstructed field
358  forAll(procField.boundaryField(), procPatchi)
359  {
360  const fvPatch& procPatch =
361  procMeshes_[proci].boundary()[procPatchi];
362 
363  const label completePatchi = completePatchID(proci, procPatchi);
364 
365  if (completePatchi == procPatchi)
366  {
367  if (!patchFields(completePatchi))
368  {
369  patchFields.set
370  (
371  completePatchi,
373  (
374  procField.boundaryField()[procPatchi],
375  completeMesh_.boundary()[completePatchi],
377  setSizeFieldMapper
378  (
379  completeMesh_.boundary()[completePatchi].size()
380  )
381  )
382  );
383  }
384 
385  patchFields[completePatchi].map
386  (
387  procField.boundaryField()[procPatchi],
388  reverseFieldMapper
389  (
390  faceProcAddressingBf_[proci][procPatchi] - 1
391  )
392  );
393  }
394  else if (isA<processorCyclicFvPatch>(procPatch))
395  {
396  if (!patchFields(completePatchi))
397  {
398  patchFields.set
399  (
400  completePatchi,
402  (
403  completeMesh_.boundary()[completePatchi].type(),
404  completeMesh_.boundary()[completePatchi],
406  )
407  );
408  }
409 
410  patchFields[completePatchi].map
411  (
412  procField.boundaryField()[procPatchi],
413  reverseFieldMapper
414  (
415  faceProcAddressingBf_[proci][procPatchi] - 1
416  )
417  );
418  }
419  else if (isA<processorFvPatch>(procPatch))
420  {
421  rmapFaceToFace
422  (
423  internalField,
424  procField.boundaryField()[procPatchi],
425  faceProcAddressingBf_[proci][procPatchi],
426  isFlux(procFields[proci])
427  );
428  }
429  }
430  }
431 
432  // Construct and return the field
433  return tmp<SurfaceField<Type>>
434  (
435  new SurfaceField<Type>
436  (
437  IOobject
438  (
439  fieldIoObject.name(),
440  completeMesh_.time().name(),
441  completeMesh_,
444  false
445  ),
446  completeMesh_,
447  procFields[0].dimensions(),
448  internalField,
449  patchFields
450  )
451  );
452 }
453 
454 
455 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
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:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
static const char *const typeName
Definition: Field.H:106
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.
void reconstructFvSurfaceFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Read, reconstruct and write all/selected surface fields.
static tmp< fvPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
static tmp< fvsPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, surfaceMesh > &)
Return a pointer to a new patchField created on freestore given.
A class for managing temporary objects.
Definition: tmp.H:55
#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:228
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:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
messageStream Info
bool isFlux(const DimensionedField< Type, surfaceMesh > &df)
Check if surfaceField is a flux.
Definition: surfaceFields.H:55
dimensioned< scalar > mag(const dimensioned< Type > &)
error FatalError
UList< label > labelUList
Definition: UList.H:65
static const char nl
Definition: Ostream.H:266
objects