fvFieldReconstructorReconstructFields.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-2022 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 
35 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
36 
37 template<class Type>
38 void Foam::fvFieldReconstructor::rmapFaceToFace
39 (
40  Field<Type>& toField,
41  const Field<Type>& fromField,
42  const labelUList& addressing,
43  const bool isFlux
44 )
45 {
46  forAll(addressing, i)
47  {
48  toField[mag(addressing[i]) - 1] =
49  (isFlux && addressing[i] < 0 ? -1 : +1)*fromField[i];
50  }
51 }
52 
53 
54 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
55 
56 template<class Type>
59 (
60  const IOobject& fieldIoObject,
61  const PtrList<DimensionedField<Type, volMesh>>& procFields
62 ) const
63 {
64  // Create the internalField
65  Field<Type> internalField(completeMesh_.nCells());
66 
67  forAll(procMeshes_, proci)
68  {
69  const DimensionedField<Type, volMesh>& procField = procFields[proci];
70 
71  // Set the cell values in the reconstructed field
72  internalField.rmap
73  (
74  procField.field(),
75  cellProcAddressing_[proci]
76  );
77  }
78 
79  return tmp<DimensionedField<Type, volMesh>>
80  (
81  new DimensionedField<Type, volMesh>
82  (
83  fieldIoObject,
84  completeMesh_,
85  procFields[0].dimensions(),
86  internalField
87  )
88  );
89 }
90 
91 
92 template<class Type>
95 (
96  const IOobject& fieldIoObject
97 ) const
98 {
99  PtrList<DimensionedField<Type, volMesh>>
100  procFields(procMeshes_.size());
101 
102  forAll(procMeshes_, proci)
103  {
104  procFields.set
105  (
106  proci,
107  new DimensionedField<Type, volMesh>
108  (
109  IOobject
110  (
111  fieldIoObject.name(),
112  procMeshes_[proci].time().timeName(),
113  procMeshes_[proci],
116  false
117  ),
118  procMeshes_[proci]
119  )
120  );
121  }
122 
124  (
125  IOobject
126  (
127  fieldIoObject.name(),
128  completeMesh_.time().timeName(),
129  completeMesh_,
132  false
133  ),
134  procFields
135  );
136 }
137 
138 
139 template<class Type>
142 (
143  const IOobject& fieldIoObject,
144  const PtrList<GeometricField<Type, fvPatchField, volMesh>>& procFields
145 ) const
146 {
147  // Create the internalField
148  Field<Type> internalField(completeMesh_.nCells());
149 
150  // Create the patch fields
151  PtrList<fvPatchField<Type>> patchFields(completeMesh_.boundary().size());
152 
153  forAll(procFields, proci)
154  {
155  const GeometricField<Type, fvPatchField, volMesh>& procField =
156  procFields[proci];
157 
158  // Set the cell values in the reconstructed field
159  internalField.rmap
160  (
161  procField.primitiveField(),
162  cellProcAddressing_[proci]
163  );
164 
165  // Set the boundary patch values in the reconstructed field
166  forAll(procField.boundaryField(), procPatchi)
167  {
168  const fvPatch& procPatch =
169  procMeshes_[proci].boundary()[procPatchi];
170 
171  const label completePatchi = completePatchID(proci, procPatchi);
172 
173  if (completePatchi == procPatchi)
174  {
175  if (!patchFields(completePatchi))
176  {
177  patchFields.set
178  (
179  completePatchi,
181  (
182  procField.boundaryField()[procPatchi],
183  completeMesh_.boundary()[completePatchi],
185  fvPatchFieldReconstructor
186  (
187  completeMesh_.boundary()[completePatchi].size()
188  )
189  )
190  );
191  }
192 
193  patchFields[completePatchi].rmap
194  (
195  procField.boundaryField()[procPatchi],
196  faceProcAddressingBf_[proci][procPatchi] - 1
197  );
198  }
199  else if (isA<processorCyclicFvPatch>(procPatch))
200  {
201  if (!patchFields(completePatchi))
202  {
203  patchFields.set
204  (
205  completePatchi,
207  (
208  completeMesh_.boundary()[completePatchi].type(),
209  completeMesh_.boundary()[completePatchi],
211  )
212  );
213  }
214 
215  if (patchFields[completePatchi].overridesConstraint())
216  {
218  << "\nThe field \"" << procFields[0].name()
219  << "\" on cyclic patch \""
220  << patchFields[completePatchi].patch().name()
221  << "\" cannot be reconstructed as it is not a cyclic "
222  << "patch field. A \"patchType cyclic;\" setting has "
223  << "been used to override the cyclic patch type.\n\n"
224  << "Cyclic patches like this with non-cyclic boundary "
225  << "conditions should be confined to a single "
226  << "processor using decomposition constraints."
227  << exit(FatalError);
228  }
229 
230  patchFields[completePatchi].rmap
231  (
232  procField.boundaryField()[procPatchi],
233  faceProcAddressingBf_[proci][procPatchi] - 1
234  );
235  }
236  }
237  }
238 
239  // Construct and return the field
240  return tmp<GeometricField<Type, fvPatchField, volMesh>>
241  (
242  new GeometricField<Type, fvPatchField, volMesh>
243  (
244  fieldIoObject,
245  completeMesh_,
246  procFields[0].dimensions(),
247  internalField,
248  patchFields
249  )
250  );
251 }
252 
253 
254 template<class Type>
257 (
258  const IOobject& fieldIoObject
259 ) const
260 {
261  PtrList<GeometricField<Type, fvPatchField, volMesh>>
262  procFields(procMeshes_.size());
263 
264  forAll(procMeshes_, proci)
265  {
266  procFields.set
267  (
268  proci,
269  new GeometricField<Type, fvPatchField, volMesh>
270  (
271  IOobject
272  (
273  fieldIoObject.name(),
274  procMeshes_[proci].time().timeName(),
275  procMeshes_[proci],
278  false
279  ),
280  procMeshes_[proci]
281  )
282  );
283  }
284 
286  (
287  IOobject
288  (
289  fieldIoObject.name(),
290  completeMesh_.time().timeName(),
291  completeMesh_,
294  false
295  ),
296  procFields
297  );
298 }
299 
300 
301 template<class Type>
304 (
305  const IOobject& fieldIoObject,
306  const PtrList<GeometricField<Type, fvsPatchField, surfaceMesh>>& procFields
307 ) const
308 {
309  // Create the internalField
310  Field<Type> internalField(completeMesh_.nInternalFaces());
311 
312  // Create the patch fields
313  PtrList<fvsPatchField<Type>> patchFields(completeMesh_.boundary().size());
314 
315  forAll(procMeshes_, proci)
316  {
317  const GeometricField<Type, fvsPatchField, surfaceMesh>& procField =
318  procFields[proci];
319 
320  // Set the internal face values in the reconstructed field
321  rmapFaceToFace
322  (
323  internalField,
324  procField.primitiveField(),
325  SubList<label>
326  (
327  faceProcAddressing_[proci],
328  procMeshes_[proci].nInternalFaces()
329  ),
330  isFlux(procFields[proci])
331  );
332 
333  // Set the boundary patch values in the reconstructed field
334  forAll(procField.boundaryField(), procPatchi)
335  {
336  const fvPatch& procPatch =
337  procMeshes_[proci].boundary()[procPatchi];
338 
339  const label completePatchi = completePatchID(proci, procPatchi);
340 
341  if (completePatchi == procPatchi)
342  {
343  if (!patchFields(completePatchi))
344  {
345  patchFields.set
346  (
347  completePatchi,
349  (
350  procField.boundaryField()[procPatchi],
351  completeMesh_.boundary()[completePatchi],
353  fvPatchFieldReconstructor
354  (
355  completeMesh_.boundary()[completePatchi].size()
356  )
357  )
358  );
359  }
360 
361  patchFields[completePatchi].rmap
362  (
363  procField.boundaryField()[procPatchi],
364  faceProcAddressingBf_[proci][procPatchi] - 1
365  );
366  }
367  else if (isA<processorCyclicFvPatch>(procPatch))
368  {
369  if (!patchFields(completePatchi))
370  {
371  patchFields.set
372  (
373  completePatchi,
375  (
376  completeMesh_.boundary()[completePatchi].type(),
377  completeMesh_.boundary()[completePatchi],
379  )
380  );
381  }
382 
383  patchFields[completePatchi].rmap
384  (
385  procField.boundaryField()[procPatchi],
386  faceProcAddressingBf_[proci][procPatchi] - 1
387  );
388  }
389  else if (isA<processorFvPatch>(procPatch))
390  {
391  rmapFaceToFace
392  (
393  internalField,
394  procField.boundaryField()[procPatchi],
395  faceProcAddressingBf_[proci][procPatchi],
396  isFlux(procFields[proci])
397  );
398  }
399  }
400  }
401 
402  // Construct and return the field
403  return tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
404  (
405  new GeometricField<Type, fvsPatchField, surfaceMesh>
406  (
407  fieldIoObject,
408  completeMesh_,
409  procFields[0].dimensions(),
410  internalField,
411  patchFields
412  )
413  );
414 }
415 
416 
417 template<class Type>
420 (
421  const IOobject& fieldIoObject
422 ) const
423 {
424  PtrList<GeometricField<Type, fvsPatchField, surfaceMesh>>
425  procFields(procMeshes_.size());
426 
427  forAll(procMeshes_, proci)
428  {
429  procFields.set
430  (
431  proci,
432  new GeometricField<Type, fvsPatchField, surfaceMesh>
433  (
434  IOobject
435  (
436  fieldIoObject.name(),
437  procMeshes_[proci].time().timeName(),
438  procMeshes_[proci],
441  false
442  ),
443  procMeshes_[proci]
444  )
445  );
446  }
447 
449  (
450  IOobject
451  (
452  fieldIoObject.name(),
453  completeMesh_.time().timeName(),
454  completeMesh_,
457  false
458  ),
459  procFields
460  );
461 }
462 
463 
464 template<class Type>
466 (
467  const IOobjectList& objects,
468  const HashSet<word>& selectedFields
469 )
470 {
471  const word& fieldClassName = DimensionedField<Type, volMesh>::typeName;
472 
473  IOobjectList fields = objects.lookupClass(fieldClassName);
474 
475  if (fields.size())
476  {
477  Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
478 
479  forAllConstIter(IOobjectList, fields, fieldIter)
480  {
481  if
482  (
483  selectedFields.empty()
484  || selectedFields.found(fieldIter()->name())
485  )
486  {
487  Info<< " " << fieldIter()->name() << endl;
488 
489  reconstructFvVolumeInternalField<Type>(*fieldIter())().write();
490 
491  nReconstructed_++;
492  }
493  }
494  Info<< endl;
495  }
496 }
497 
498 
499 template<class Type>
501 (
502  const IOobjectList& objects,
503  const HashSet<word>& selectedFields
504 )
505 {
506  const word& fieldClassName =
508 
509  IOobjectList fields = objects.lookupClass(fieldClassName);
510 
511  if (fields.size())
512  {
513  Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
514 
515  forAllConstIter(IOobjectList, fields, fieldIter)
516  {
517  if
518  (
519  selectedFields.empty()
520  || selectedFields.found(fieldIter()->name())
521  )
522  {
523  Info<< " " << fieldIter()->name() << endl;
524 
525  reconstructFvVolumeField<Type>(*fieldIter())().write();
526 
527  nReconstructed_++;
528  }
529  }
530  Info<< endl;
531  }
532 }
533 
534 
535 template<class Type>
537 (
538  const IOobjectList& objects,
539  const HashSet<word>& selectedFields
540 )
541 {
542  const word& fieldClassName =
544 
545  IOobjectList fields = objects.lookupClass(fieldClassName);
546 
547  if (fields.size())
548  {
549  Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
550 
551  forAllConstIter(IOobjectList, fields, fieldIter)
552  {
553  if
554  (
555  selectedFields.empty()
556  || selectedFields.found(fieldIter()->name())
557  )
558  {
559  Info<< " " << fieldIter()->name() << endl;
560 
561  reconstructFvSurfaceField<Type>(*fieldIter())().write();
562 
563  nReconstructed_++;
564  }
565  }
566  Info<< endl;
567  }
568 }
569 
570 
571 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual Ostream & write(const char)=0
Write character.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
tmp< GeometricField< Type, fvPatchField, volMesh > > reconstructFvVolumeField(const IOobject &fieldIoObject, const PtrList< GeometricField< Type, fvPatchField, volMesh >> &) const
Reconstruct volume field.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
label nInternalFaces() const
static const char *const typeName
Definition: Field.H:105
void reconstructFvSurfaceFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Read, reconstruct and write all/selected surface fields.
label nCells() const
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
UList< label > labelUList
Definition: UList.H:65
static tmp< fvsPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, surfaceMesh > &)
Return a pointer to a new patchField created on freestore given.
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
static tmp< fvPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
bool isFlux(const DimensionedField< Type, surfaceMesh > &df)
Check if surfaceField is a flux.
Definition: surfaceFields.H:55
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< DimensionedField< Type, volMesh > > reconstructFvVolumeInternalField(const IOobject &fieldIoObject, const PtrList< DimensionedField< Type, volMesh >> &procFields) const
Reconstruct volume internal field.
void reconstructFvVolumeFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Read, reconstruct and write all/selected volume fields.
A class for managing temporary objects.
Definition: PtrList.H:53
void reconstructFvVolumeInternalFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Read, reconstruct and write all/selected volume internal fields.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:800
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > reconstructFvSurfaceField(const IOobject &fieldIoObject, const PtrList< GeometricField< Type, fvsPatchField, surfaceMesh >> &) const
Reconstruct surface field.