fvFieldDecomposerTemplates.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 "fvFieldDecomposer.H"
27 #include "processorFvPatchField.H"
28 #include "processorFvsPatchField.H"
31 #include "emptyFvPatchFields.H"
32 #include "stringOps.H"
33 
34 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35 
36 template<class Type>
37 Foam::tmp<Foam::Field<Type>> Foam::fvFieldDecomposer::mapCellToFace
38 (
39  const labelUList& owner,
40  const labelUList& neighbour,
41  const Field<Type>& field,
42  const labelUList& addressing
43 )
44 {
45  tmp<Field<Type>> tfld(new Field<Type>(addressing.size()));
46  Field<Type>& fld = tfld.ref();
47 
48  forAll(addressing, i)
49  {
50  fld[i] =
51  field
52  [
53  addressing[i] > 0
54  ? neighbour[addressing[i] - 1]
55  : owner[- addressing[i] - 1]
56  ];
57  }
58 
59  return tfld;
60 }
61 
62 
63 template<class Type>
64 Foam::tmp<Foam::Field<Type>> Foam::fvFieldDecomposer::mapFaceToFace
65 (
66  const Field<Type>& field,
67  const labelUList& addressing,
68  const bool isFlux
69 )
70 {
71  tmp<Field<Type>> tfld(new Field<Type>(addressing.size()));
72  Field<Type>& fld = tfld.ref();
73 
74  forAll(addressing, i)
75  {
76  fld[i] =
77  (isFlux && addressing[i] < 0 ? -1 : +1)
78  *field[mag(addressing[i]) - 1];
79  }
80 
81  return tfld;
82 }
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
87 template<class Type>
90 (
91  const IOobject& fieldIoObject
92 ) const
93 {
94  // Read the field
95  const typename VolField<Type>::Internal field
96  (
97  IOobject
98  (
99  fieldIoObject.name(),
100  completeMesh_.time().name(),
101  completeMesh_,
104  false
105  ),
106  completeMesh_
107  );
108 
109  // Construct the processor fields
110  PtrList<typename VolField<Type>::Internal> procFields(procMeshes_.size());
111  forAll(procMeshes_, proci)
112  {
113  // Create the processor field with the dummy patch fields
114  procFields.set
115  (
116  proci,
117  new typename VolField<Type>::Internal
118  (
119  IOobject
120  (
121  fieldIoObject.name(),
122  procMeshes_[proci].time().name(),
123  procMeshes_[proci],
126  false
127  ),
128  procMeshes_[proci],
129  field.dimensions(),
130  Field<Type>(field.primitiveField(), cellProcAddressing_[proci])
131  )
132  );
133  }
134 
135  return procFields;
136 }
137 
138 
139 template<class Type>
142 (
143  const IOobject& fieldIoObject
144 ) const
145 {
146  // Read the field
147  const VolField<Type> field
148  (
149  IOobject
150  (
151  fieldIoObject.name(),
152  completeMesh_.time().name(),
153  completeMesh_,
156  false
157  ),
158  completeMesh_
159  );
160 
161  // Construct the processor fields
162  PtrList<VolField<Type>> procFields(procMeshes_.size());
163  forAll(procMeshes_, proci)
164  {
165  // Create dummy patch fields
166  PtrList<fvPatchField<Type>> patchFields
167  (
168  procMeshes_[proci].boundary().size()
169  );
170  forAll(procMeshes_[proci].boundary(), procPatchi)
171  {
172  patchFields.set
173  (
174  procPatchi,
176  (
178  procMeshes_[proci].boundary()[procPatchi],
180  )
181  );
182  }
183 
184  // Create the processor field with the dummy patch fields
185  procFields.set
186  (
187  proci,
188  new VolField<Type>
189  (
190  IOobject
191  (
192  fieldIoObject.name(),
193  procMeshes_[proci].time().name(),
194  procMeshes_[proci],
197  false
198  ),
199  procMeshes_[proci],
200  field.dimensions(),
201  Field<Type>(field.primitiveField(), cellProcAddressing_[proci]),
202  patchFields,
203  field.sources().table()
204  )
205  );
206 
207  // Alias the created proc field
208  VolField<Type>& vf = procFields[proci];
209 
210  // Change the patch fields to the correct type using a mapper
211  // constructor (with reference to the now correct internal field)
212  typename VolField<Type>::Boundary& bf = vf.boundaryFieldRef();
213  forAll(bf, procPatchi)
214  {
215  const fvPatch& procPatch =
216  procMeshes_[proci].boundary()[procPatchi];
217 
218  const label completePatchi = completePatchID(proci, procPatchi);
219 
220  if (completePatchi == procPatchi)
221  {
222  bf.set
223  (
224  procPatchi,
226  (
227  field.boundaryField()[completePatchi],
228  procPatch,
229  vf(),
230  patchFieldDecomposers_[proci][procPatchi]
231  )
232  );
233  }
234  else if (isA<processorCyclicFvPatch>(procPatch))
235  {
236  if (field.boundaryField()[completePatchi].overridesConstraint())
237  {
238  OStringStream str;
239  str << "\nThe field \"" << field.name()
240  << "\" on cyclic patch \""
241  << field.boundaryField()[completePatchi].patch().name()
242  << "\" cannot be decomposed as it is not a cyclic "
243  << "patch field. A \"patchType cyclic;\" setting has "
244  << "been used to override the cyclic patch type.\n\n"
245  << "Cyclic patches like this with non-cyclic boundary "
246  << "conditions should be confined to a single "
247  << "processor using decomposition constraints.";
249  << stringOps::breakIntoIndentedLines(str.str()).c_str()
250  << exit(FatalError);
251  }
252 
253  const label nbrCompletePatchi =
254  refCast<const processorCyclicFvPatch>(procPatch)
255  .referPatch().nbrPatchIndex();
256 
257  // Use `fvPatchField<Type>::New` rather than
258  // `new processorCyclicFvPatchField<Type>` so that derivations
259  // (such as non-conformal processor cyclics) are constructed
260  bf.set
261  (
262  procPatchi,
264  (
265  procPatch.type(),
266  procPatch,
267  vf()
268  )
269  );
270 
271  bf[procPatchi] =
272  mapCellToFace
273  (
274  labelUList(),
275  completeMesh_.lduAddr().patchAddr(nbrCompletePatchi),
276  field.primitiveField(),
277  faceProcAddressingBf_[proci][procPatchi]
278  );
279  }
280  else if (isA<processorFvPatch>(procPatch))
281  {
282  bf.set
283  (
284  procPatchi,
286  (
287  procPatch.type(),
288  procPatch,
289  vf()
290  )
291  );
292 
293  bf[procPatchi] =
294  mapCellToFace
295  (
296  completeMesh_.owner(),
297  completeMesh_.neighbour(),
298  field.primitiveField(),
299  faceProcAddressingBf_[proci][procPatchi]
300  );
301  }
302  else
303  {
305  << "Unknown type." << abort(FatalError);
306  }
307  }
308  }
309 
310  return procFields;
311 }
312 
313 
314 template<class Type>
317 (
318  const IOobject& fieldIoObject
319 ) const
320 {
321  // Read the field
322  const SurfaceField<Type> field
323  (
324  IOobject
325  (
326  fieldIoObject.name(),
327  completeMesh_.time().name(),
328  completeMesh_,
331  false
332  ),
333  completeMesh_
334  );
335 
336  // Construct the processor fields
337  PtrList<SurfaceField<Type>> procFields(procMeshes_.size());
338  forAll(procMeshes_, proci)
339  {
340  const SubList<label> faceAddressingIf
341  (
342  faceProcAddressing_[proci],
343  procMeshes_[proci].nInternalFaces()
344  );
345 
346  // Create dummy patch fields
347  PtrList<fvsPatchField<Type>> patchFields
348  (
349  procMeshes_[proci].boundary().size()
350  );
351  forAll(procMeshes_[proci].boundary(), procPatchi)
352  {
353  patchFields.set
354  (
355  procPatchi,
357  (
359  procMeshes_[proci].boundary()[procPatchi],
361  )
362  );
363  }
364 
365  // Create the processor field with the dummy patch fields
366  procFields.set
367  (
368  proci,
370  (
371  IOobject
372  (
373  field.name(),
374  procMeshes_[proci].time().name(),
375  procMeshes_[proci],
378  false
379  ),
380  procMeshes_[proci],
381  field.dimensions(),
382  mapFaceToFace
383  (
384  field,
385  faceAddressingIf,
386  isFlux(field)
387  ),
388  patchFields
389  )
390  );
391 
392  // Alias the created proc field
393  SurfaceField<Type>& sf = procFields[proci];
394 
395  // Change the patch fields to the correct type using a mapper
396  // constructor (with reference to the now correct internal field)
398  forAll(procMeshes_[proci].boundary(), procPatchi)
399  {
400  const fvPatch& procPatch =
401  procMeshes_[proci].boundary()[procPatchi];
402 
403  const label completePatchi = completePatchID(proci, procPatchi);
404 
405  if (completePatchi == procPatchi)
406  {
407  bf.set
408  (
409  procPatchi,
411  (
412  field.boundaryField()[procPatchi],
413  procPatch,
414  sf(),
415  patchFieldDecomposers_[proci][procPatchi]
416  )
417  );
418  }
419  else if (isA<processorCyclicFvPatch>(procPatch))
420  {
421  bf.set
422  (
423  procPatchi,
425  (
426  procPatch,
427  sf(),
428  mapFaceToFace
429  (
430  field.boundaryField()[completePatchi],
431  faceProcAddressingBf_[proci][procPatchi],
432  isFlux(field)
433  )
434  )
435  );
436  }
437  else if (isA<processorFvPatch>(procPatch))
438  {
439  bf.set
440  (
441  procPatchi,
443  (
444  procPatch,
445  sf(),
446  mapFaceToFace
447  (
448  field.primitiveField(),
449  faceProcAddressingBf_[proci][procPatchi],
450  isFlux(field)
451  )
452  )
453  );
454  }
455  else
456  {
458  << "Unknown type." << abort(FatalError);
459  }
460  }
461  }
462 
463  return procFields;
464 }
465 
466 
467 template<class Type>
469 (
470  const IOobjectList& objects
471 )
472 {
473  const word& fieldClassName = VolField<Type>::Internal::typeName;
474 
475  IOobjectList fields = objects.lookupClass(fieldClassName);
476 
477  if (fields.size())
478  {
479  Info<< nl << " Decomposing " << fieldClassName << "s" << nl << endl;
480 
481  forAllConstIter(IOobjectList, fields, fieldIter)
482  {
483  Info<< " " << fieldIter()->name() << endl;
484 
486  decomposeVolInternalField<Type>(*fieldIter());
487 
488  forAll(procFields, proci)
489  {
490  procFields[proci].write();
491  }
492  }
493  }
494 }
495 
496 
497 template<class Type>
499 (
500  const IOobjectList& objects
501 )
502 {
503  const word& fieldClassName = VolField<Type>::typeName;
504 
505  IOobjectList fields = objects.lookupClass(fieldClassName);
506 
507  if (fields.size())
508  {
509  Info<< nl << " Decomposing " << fieldClassName << "s" << nl << endl;
510 
511  forAllConstIter(IOobjectList, fields, fieldIter)
512  {
513  Info<< " " << fieldIter()->name() << endl;
514 
515  const PtrList<VolField<Type>> procFields =
516  decomposeVolField<Type>(*fieldIter());
517 
518  forAll(procFields, proci)
519  {
520  procFields[proci].write();
521  }
522  }
523  }
524 }
525 
526 
527 template<class Type>
529 (
530  const IOobjectList& objects
531 )
532 {
533  const word& fieldClassName = SurfaceField<Type>::typeName;
534 
535  IOobjectList fields = objects.lookupClass(fieldClassName);
536 
537  if (fields.size())
538  {
539  Info<< nl << " Decomposing " << fieldClassName << "s" << nl << endl;
540 
541  forAllConstIter(IOobjectList, fields, fieldIter)
542  {
543  Info<< " " << fieldIter()->name() << endl;
544 
545  const PtrList<SurfaceField<Type>> procFields =
546  decomposeFvSurfaceField<Type>(*fieldIter());
547 
548  forAll(procFields, proci)
549  {
550  procFields[proci].write();
551  }
552  }
553  }
554 }
555 
556 
557 // ************************************************************************* //
#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 dimensionSet & dimensions() const
Return dimensions.
const PrimitiveField< Type > & primitiveField() const
Return a const-reference to the primitive field.
Pre-declare SubField and related Field type.
Definition: Field.H:83
Generic GeometricBoundaryField class.
const HashPtrTable< Source > & table() const
Access the underlying field table.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const Sources & sources() const
Return const-reference to the sources.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const Internal::FieldType & primitiveField() const
Return a const-reference to the primitive field.
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
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
Foam::calculatedFvsPatchField.
PtrList< SurfaceField< Type > > decomposeFvSurfaceField(const IOobject &fieldIoObject) const
Decompose a surface field.
void decomposeVolFields(const IOobjectList &objects)
Read, decompose and write all volume fields.
PtrList< VolField< Type > > decomposeVolField(const IOobject &fieldIoObject) const
Decompose a volume field.
void decomposeFvSurfaceFields(const IOobjectList &objects)
Read, decompose and write all surface fields.
void decomposeVolInternalFields(const IOobjectList &objects)
Read, decompose and write all volume internal fields.
PtrList< typename VolField< Type >::Internal > decomposeVolInternalField(const IOobject &fieldIoObject) const
Decompose a volume internal field.
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
Foam::processorCyclicFvsPatchField.
Foam::processorFvsPatchField.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
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
volScalarField sf(fieldObject, mesh)
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
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
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
errorManip< error > abort(error &err)
Definition: errorManip.H:131
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)
error FatalError
UList< label > labelUList
Definition: UList.H:65
static const char nl
Definition: Ostream.H:267
faceListList boundary(nPatches)
objects