fvMeshDistributeTemplates.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 "polyTopoChangeMap.H"
27 #include "processorFvPatchField.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 template<class GeoField>
33 {
35  (
36  mesh.fields<GeoField>(false, fvMesh::curGeometryFields)
37  );
38 
39  forAll(fields, i)
40  {
41  const GeoField& field = fields[i];
42 
43  Pout<< "Field:" << field.name() << " internal size:" << field.size()
44  << endl;
45 
46  forAll(field.boundaryField(), patchi)
47  {
48  Pout<< " " << patchi
49  << ' ' << field.boundaryField()[patchi].patch().name()
50  << ' ' << field.boundaryField()[patchi].type()
51  << ' ' << field.boundaryField()[patchi].size()
52  << endl;
53  }
54  }
55 }
56 
57 
58 template<class Type, class Mesh>
59 void Foam::fvMeshDistribute::saveBoundaryFields
60 (
62 ) const
63 {
64  // Save whole boundary field
65 
67  (
68  mesh_.fields<SurfaceField<Type>>(false, fvMesh::curGeometryFields)
69  );
70 
71  bfields.setSize(fields.size());
72 
73  forAll(fields, i)
74  {
75  const SurfaceField<Type>& field = fields[i];
76 
77  bfields.set(i, field.boundaryField().clone().ptr());
78  }
79 }
80 
81 
82 template<class Type, class Mesh>
83 void Foam::fvMeshDistribute::mapBoundaryFields
84 (
85  const polyTopoChangeMap& map,
86  const PtrList<FieldField<surfaceMesh::PatchField, Type>>& oldBfields
87 )
88 {
89  // Map boundary field
90 
91  const labelList& oldPatchStarts = map.oldPatchStarts();
92  const labelList& faceMap = map.faceMap();
93 
94  UPtrList<SurfaceField<Type>> fields
95  (
96  mesh_.fields<SurfaceField<Type>>(false, fvMesh::curGeometryFields)
97  );
98 
99  forAll(fields, i)
100  {
101  SurfaceField<Type>& field = fields[i];
102 
103  typename SurfaceField<Type>::Boundary& bfield =
104  field.boundaryFieldRef();
105 
106  const FieldField<surfaceMesh::PatchField, Type>& oldBfield =
107  oldBfields[i];
108 
109  // Pull from old boundary field into bfield.
110 
111  forAll(bfield, patchi)
112  {
113  fvsPatchField<Type>& patchField = bfield[patchi];
114  label facei = patchField.patch().start();
115 
116  forAll(patchField, i)
117  {
118  label oldFacei = faceMap[facei++];
119 
120  // Find patch and local patch face oldFacei was in.
121  forAll(oldPatchStarts, oldPatchi)
122  {
123  label oldLocalI = oldFacei - oldPatchStarts[oldPatchi];
124 
125  if
126  (
127  oldLocalI >= 0
128  && oldLocalI < oldBfield[oldPatchi].size()
129  )
130  {
131  patchField[i] = oldBfield[oldPatchi][oldLocalI];
132  }
133  }
134  }
135  }
136  }
137 }
138 
139 
140 template<class Type>
141 void Foam::fvMeshDistribute::initMapExposedFaces
142 (
143  PtrList<Field<Type>>& ifields
144 ) const
145 {
146  const UPtrList<SurfaceField<Type>> fields
147  (
148  mesh_.fields<SurfaceField<Type>>(false, fvMesh::curGeometryFields)
149  );
150 
151  ifields.setSize(fields.size());
152 
153  forAll(fields, i)
154  {
155  ifields.set(i, Field<Type>(mesh_.nFaces()));
156 
157  const SurfaceField<Type>& field = fields[i];
158 
159  SubList<Type>(ifields[i], field.primitiveField().size()) =
160  field.primitiveField();
161 
162  forAll(field.boundaryField(), patchi)
163  {
164  const fvsPatchField<Type>& pfield = field.boundaryField()[patchi];
165 
166  SubList<Type>(ifields[i], pfield.size(), pfield.patch().start()) =
167  pfield;
168  }
169  }
170 }
171 
172 
173 template<class Type>
174 void Foam::fvMeshDistribute::mapExposedFaces
175 (
176  const polyTopoChangeMap& map,
177  const PtrList<Field<Type>>& oldFields
178 )
179 {
180  UPtrList<SurfaceField<Type>> fields
181  (
182  mesh_.fields<SurfaceField<Type>>(false, fvMesh::curGeometryFields)
183  );
184 
185  forAll(fields, i)
186  {
187  SurfaceField<Type>& field = fields[i];
188 
189  const Field<Type>& oldField = oldFields[i];
190 
191  const bool negateIfFlipped = isFlux(field);
192 
193  forAll(field.boundaryField(), patchi)
194  {
195  fvsPatchField<Type>& patchField = field.boundaryFieldRef()[patchi];
196 
197  forAll(patchField, i)
198  {
199  const label facei = patchField.patch().start()+i;
200  const label oldFacei = map.faceMap()[facei];
201 
202  if (oldFacei < map.nOldInternalFaces())
203  {
204  if (negateIfFlipped && map.flipFaceFlux().found(facei))
205  {
206  patchField[i] = flipOp()(oldField[oldFacei]);
207  }
208  else
209  {
210  patchField[i] = oldField[oldFacei];
211  }
212  }
213  else
214  {
215  patchField[i] = oldField[oldFacei];
216  }
217  }
218  }
219  }
220 }
221 
222 
223 template<class GeoField>
224 void Foam::fvMeshDistribute::correctCoupledPatchFields()
225 {
226  UPtrList<GeoField> fields
227  (
228  mesh_.fields<GeoField>(false, fvMesh::curGeometryFields)
229  );
230 
231  // Ensure the deltaCoeffs are available for constraint patch evaluation
232  mesh_.deltaCoeffs();
233 
234  forAll(fields, i)
235  {
236  GeoField& field = fields[i];
237 
238  typename GeoField::Boundary& bfield = field.boundaryFieldRef();
239 
240  if
241  (
244  )
245  {
246  label nReq = Pstream::nRequests();
247 
248  forAll(bfield, patchi)
249  {
250  if (bfield[patchi].coupled())
251  {
252  bfield[patchi].initEvaluate(Pstream::defaultCommsType);
253  }
254  }
255 
256  // Block for any outstanding requests
257  if
258  (
261  )
262  {
263  Pstream::waitRequests(nReq);
264  }
265 
266  forAll(bfield, patchi)
267  {
268  if (bfield[patchi].coupled())
269  {
270  bfield[patchi].evaluate(Pstream::defaultCommsType);
271  }
272  }
273  }
275  {
276  const lduSchedule& patchSchedule =
277  mesh_.globalData().patchSchedule();
278 
279  forAll(patchSchedule, patchEvali)
280  {
281  if (bfield[patchEvali].coupled())
282  {
283  if (patchSchedule[patchEvali].init)
284  {
285  bfield[patchSchedule[patchEvali].patch]
286  .initEvaluate(Pstream::commsTypes::scheduled);
287  }
288  else
289  {
290  bfield[patchSchedule[patchEvali].patch]
292  }
293  }
294  }
295  }
296  }
297 }
298 
299 
300 template<class GeoField>
301 void Foam::fvMeshDistribute::sendFields
302 (
303  const label domain,
304  const wordList& fieldNames,
305  const fvMeshSubset& subsetter,
306  Ostream& toNbr
307 )
308 {
309  // Send fields. Note order supplied so we can receive in exactly the same
310  // order.
311  // Note that field gets written as entry in dictionary so we
312  // can construct from subdictionary.
313  // (since otherwise the reading as-a-dictionary mixes up entries from
314  // consecutive fields)
315  // The dictionary constructed is:
316  // volScalarField
317  // {
318  // p {internalField ..; boundaryField ..;}
319  // k {internalField ..; boundaryField ..;}
320  // }
321  // volVectorField
322  // {
323  // U {internalField ... }
324  // }
325 
326  // volVectorField {U {internalField ..; boundaryField ..;}}
327 
328  toNbr << GeoField::typeName << token::NL << token::BEGIN_BLOCK << token::NL;
329  forAll(fieldNames, i)
330  {
331  if (debug)
332  {
333  Pout<< "Subsetting field " << fieldNames[i]
334  << " for domain:" << domain << endl;
335  }
336 
337  // Send all fieldNames. This has to be exactly the same set as is
338  // being received!
339  const GeoField& field =
340  subsetter.baseMesh().lookupObject<GeoField>(fieldNames[i]);
341 
342  tmp<GeoField> tsubfield = subsetter.interpolate(field);
343 
344  toNbr
346  << tsubfield
348  }
349  toNbr << token::END_BLOCK << token::NL;
350 }
351 
352 
353 template<class GeoField>
354 void Foam::fvMeshDistribute::receiveFields
355 (
356  const label domain,
357  const wordList& fieldNames,
358  typename GeoField::Mesh& mesh,
359  PtrList<GeoField>& fields,
360  const dictionary& fieldDicts
361 )
362 {
363  if (debug)
364  {
365  Pout<< "Receiving fields " << fieldNames
366  << " from domain:" << domain << endl;
367  }
368 
369  fields.setSize(fieldNames.size());
370 
371  forAll(fieldNames, i)
372  {
373  if (debug)
374  {
375  Pout<< "Constructing field " << fieldNames[i]
376  << " from domain:" << domain << endl;
377  }
378 
379  fields.set
380  (
381  i,
382  new GeoField
383  (
384  IOobject
385  (
386  fieldNames[i],
387  mesh.thisDb().time().name(),
388  mesh.thisDb(),
391  ),
392  mesh,
393  fieldDicts.subDict(fieldNames[i])
394  )
395  );
396  }
397 }
398 
399 
400 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Generic field type.
Definition: FieldField.H:77
tmp< FieldField< Field, Type > > clone() const
Clone.
Definition: FieldField.C:188
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
GeometricBoundaryField< Type, GeoMesh, PrimitiveField > Boundary
Type of the boundary field.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:85
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
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:137
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:147
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:272
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:66
const word & name() const
Return const reference to name.
static void printFieldInfo(const fvMesh &)
Print some field info.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
static const HashSet< word > curGeometryFields
Set of names of registered current-time geometric fields.
Definition: fvMesh.H:326
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:426
UPtrList< GeoField > fields(bool strict=false, const HashSet< word > &geometryFields=fvMesh::geometryFields) const
Return the list of fields of type GeoField.
const Time & time() const
Return time.
@ BEGIN_BLOCK
Definition: token.H:113
@ END_BLOCK
Definition: token.H:114
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
label patchi
static List< word > fieldNames
Definition: globalFoam.H:46
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
List< word > wordList
A List of words.
Definition: fileName.H:54
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
bool isFlux(const DimensionedField< Type, surfaceMesh > &df)
Check if surfaceField is a flux.
Definition: surfaceFields.H:55
List< lduScheduleEntry > lduSchedule
Definition: lduSchedule.H:74
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53