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-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 "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 template<class Type>
87 Foam::fvFieldDecomposer::decomposeVolInternalField
88 (
89  const IOobject& fieldIoObject
90 ) const
91 {
92  // Read the field
93  const typename VolField<Type>::Internal field
94  (
95  IOobject
96  (
97  fieldIoObject.name(),
98  completeMesh_.time().name(),
99  completeMesh_,
102  false
103  ),
104  completeMesh_
105  );
106 
107  // Construct the processor fields
108  PtrList<typename VolField<Type>::Internal> procFields(procMeshes_.size());
109  forAll(procMeshes_, proci)
110  {
111  // Create the processor field with the dummy patch fields
112  procFields.set
113  (
114  proci,
115  new typename VolField<Type>::Internal
116  (
117  IOobject
118  (
119  fieldIoObject.name(),
120  procMeshes_[proci].time().name(),
121  procMeshes_[proci],
124  false
125  ),
126  procMeshes_[proci],
127  field.dimensions(),
128  Field<Type>(field.primitiveField(), cellProcAddressing_[proci])
129  )
130  );
131  }
132 
133  return procFields;
134 }
135 
136 
137 template<class Type>
139 Foam::fvFieldDecomposer::decomposeVolField
140 (
141  const IOobject& fieldIoObject
142 ) const
143 {
144  // Read the field
145  const VolField<Type> field
146  (
147  IOobject
148  (
149  fieldIoObject.name(),
150  completeMesh_.time().name(),
151  completeMesh_,
154  false
155  ),
156  completeMesh_
157  );
158 
159  // Construct the processor fields
160  PtrList<VolField<Type>> procFields(procMeshes_.size());
161  forAll(procMeshes_, proci)
162  {
163  // Create dummy patch fields
164  PtrList<fvPatchField<Type>> patchFields
165  (
166  procMeshes_[proci].boundary().size()
167  );
168  forAll(procMeshes_[proci].boundary(), procPatchi)
169  {
170  patchFields.set
171  (
172  procPatchi,
174  (
176  procMeshes_[proci].boundary()[procPatchi],
178  )
179  );
180  }
181 
182  // Create the processor field with the dummy patch fields
183  procFields.set
184  (
185  proci,
186  new VolField<Type>
187  (
188  IOobject
189  (
190  fieldIoObject.name(),
191  procMeshes_[proci].time().name(),
192  procMeshes_[proci],
195  false
196  ),
197  procMeshes_[proci],
198  field.dimensions(),
199  Field<Type>(field.primitiveField(), cellProcAddressing_[proci]),
200  patchFields,
201  field.sources().table()
202  )
203  );
204 
205  // Alias the created proc field
206  VolField<Type>& vf = procFields[proci];
207 
208  // Change the patch fields to the correct type using a mapper
209  // constructor (with reference to the now correct internal field)
210  typename VolField<Type>::Boundary& bf = vf.boundaryFieldRef();
211  forAll(bf, procPatchi)
212  {
213  const fvPatch& procPatch =
214  procMeshes_[proci].boundary()[procPatchi];
215 
216  const label completePatchi = completePatchID(proci, procPatchi);
217 
218  if (completePatchi == procPatchi)
219  {
220  bf.set
221  (
222  procPatchi,
224  (
225  field.boundaryField()[completePatchi],
226  procPatch,
227  vf(),
228  patchFieldDecomposers_[proci][procPatchi]
229  )
230  );
231  }
232  else if (isA<processorCyclicFvPatch>(procPatch))
233  {
234  if (field.boundaryField()[completePatchi].overridesConstraint())
235  {
236  OStringStream str;
237  str << "\nThe field \"" << field.name()
238  << "\" on cyclic patch \""
239  << field.boundaryField()[completePatchi].patch().name()
240  << "\" cannot be decomposed as it is not a cyclic "
241  << "patch field. A \"patchType cyclic;\" setting has "
242  << "been used to override the cyclic patch type.\n\n"
243  << "Cyclic patches like this with non-cyclic boundary "
244  << "conditions should be confined to a single "
245  << "processor using decomposition constraints.";
247  << stringOps::breakIntoIndentedLines(str.str()).c_str()
248  << exit(FatalError);
249  }
250 
251  const label nbrCompletePatchi =
252  refCast<const processorCyclicFvPatch>(procPatch)
253  .referPatch().nbrPatchIndex();
254 
255  // Use `fvPatchField<Type>::New` rather than
256  // `new processorCyclicFvPatchField<Type>` so that derivations
257  // (such as non-conformal processor cyclics) are constructed
258  bf.set
259  (
260  procPatchi,
262  (
263  procPatch.type(),
264  procPatch,
265  vf()
266  )
267  );
268 
269  bf[procPatchi] =
270  mapCellToFace
271  (
272  labelUList(),
273  completeMesh_.lduAddr().patchAddr(nbrCompletePatchi),
274  field.primitiveField(),
275  faceProcAddressingBf_[proci][procPatchi]
276  );
277  }
278  else if (isA<processorFvPatch>(procPatch))
279  {
280  bf.set
281  (
282  procPatchi,
284  (
285  procPatch.type(),
286  procPatch,
287  vf()
288  )
289  );
290 
291  bf[procPatchi] =
292  mapCellToFace
293  (
294  completeMesh_.owner(),
295  completeMesh_.neighbour(),
296  field.primitiveField(),
297  faceProcAddressingBf_[proci][procPatchi]
298  );
299  }
300  else
301  {
303  << "Unknown type." << abort(FatalError);
304  }
305  }
306  }
307 
308  return procFields;
309 }
310 
311 
312 template<class Type>
314 Foam::fvFieldDecomposer::decomposeFvSurfaceField
315 (
316  const IOobject& fieldIoObject
317 ) const
318 {
319  // Read the field
320  const SurfaceField<Type> field
321  (
322  IOobject
323  (
324  fieldIoObject.name(),
325  completeMesh_.time().name(),
326  completeMesh_,
329  false
330  ),
331  completeMesh_
332  );
333 
334  // Construct the processor fields
335  PtrList<SurfaceField<Type>> procFields(procMeshes_.size());
336  forAll(procMeshes_, proci)
337  {
338  const SubList<label> faceAddressingIf
339  (
340  faceProcAddressing_[proci],
341  procMeshes_[proci].nInternalFaces()
342  );
343 
344  // Create dummy patch fields
345  PtrList<fvsPatchField<Type>> patchFields
346  (
347  procMeshes_[proci].boundary().size()
348  );
349  forAll(procMeshes_[proci].boundary(), procPatchi)
350  {
351  patchFields.set
352  (
353  procPatchi,
355  (
357  procMeshes_[proci].boundary()[procPatchi],
359  )
360  );
361  }
362 
363  // Create the processor field with the dummy patch fields
364  procFields.set
365  (
366  proci,
367  new SurfaceField<Type>
368  (
369  IOobject
370  (
371  field.name(),
372  procMeshes_[proci].time().name(),
373  procMeshes_[proci],
376  false
377  ),
378  procMeshes_[proci],
379  field.dimensions(),
380  mapFaceToFace
381  (
382  field,
383  faceAddressingIf,
384  isFlux(field)
385  ),
386  patchFields
387  )
388  );
389 
390  // Alias the created proc field
391  SurfaceField<Type>& sf = procFields[proci];
392 
393  // Change the patch fields to the correct type using a mapper
394  // constructor (with reference to the now correct internal field)
396  forAll(procMeshes_[proci].boundary(), procPatchi)
397  {
398  const fvPatch& procPatch =
399  procMeshes_[proci].boundary()[procPatchi];
400 
401  const label completePatchi = completePatchID(proci, procPatchi);
402 
403  if (completePatchi == procPatchi)
404  {
405  bf.set
406  (
407  procPatchi,
409  (
410  field.boundaryField()[procPatchi],
411  procPatch,
412  sf(),
413  patchFieldDecomposers_[proci][procPatchi]
414  )
415  );
416  }
417  else if (isA<processorCyclicFvPatch>(procPatch))
418  {
419  bf.set
420  (
421  procPatchi,
422  new processorCyclicFvsPatchField<Type>
423  (
424  procPatch,
425  sf(),
426  mapFaceToFace
427  (
428  field.boundaryField()[completePatchi],
429  faceProcAddressingBf_[proci][procPatchi],
430  isFlux(field)
431  )
432  )
433  );
434  }
435  else if (isA<processorFvPatch>(procPatch))
436  {
437  bf.set
438  (
439  procPatchi,
440  new processorFvsPatchField<Type>
441  (
442  procPatch,
443  sf(),
444  mapFaceToFace
445  (
446  field.primitiveField(),
447  faceProcAddressingBf_[proci][procPatchi],
448  isFlux(field)
449  )
450  )
451  );
452  }
453  else
454  {
456  << "Unknown type." << abort(FatalError);
457  }
458  }
459  }
460 
461  return procFields;
462 }
463 
464 
465 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
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 
485  PtrList<typename VolField<Type>::Internal> procFields =
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  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  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: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
DimensionedField< Type, GeoMesh > Internal
Type of the internal field from which this GeometricField is derived.
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
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
void decomposeVolFields(const IOobjectList &objects)
Read, decompose and write all volume fields.
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.
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
volScalarField sf(fieldObject, mesh)
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(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(), 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: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
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
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
dimensioned< scalar > mag(const dimensioned< Type > &)
error FatalError
UList< label > labelUList
Definition: UList.H:65
static const char nl
Definition: Ostream.H:266
faceListList boundary(nPatches)
objects