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-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 "polyTopoChangeMap.H"
27 #include "processorFvPatchField.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 template<class GeoField>
33 {
35  (
36  mesh.objectRegistry::lookupClass<GeoField>()
37  );
38 
39  forAllConstIter(typename HashTable<const GeoField*>, flds, iter)
40  {
41  const GeoField& fld = *iter();
42 
43  Pout<< "Field:" << iter.key() << " internal size:" << fld.size()
44  << endl;
45 
46  forAll(fld.boundaryField(), patchi)
47  {
48  Pout<< " " << patchi
49  << ' ' << fld.boundaryField()[patchi].patch().name()
50  << ' ' << fld.boundaryField()[patchi].type()
51  << ' ' << fld.boundaryField()[patchi].size()
52  << endl;
53  }
54  }
55 }
56 
57 
58 template<class T, class Mesh>
59 void Foam::fvMeshDistribute::saveBoundaryFields
60 (
62 ) const
63 {
64  // Save whole boundary field
65 
67 
69  (
70  static_cast<const fvMesh&>(mesh_).objectRegistry::lookupClass<fldType>()
71  );
72 
73  bflds.setSize(flds.size());
74 
75  label i = 0;
76 
77  forAllConstIter(typename HashTable<const fldType*>, flds, iter)
78  {
79  const fldType& fld = *iter();
80 
81  bflds.set(i, fld.boundaryField().clone().ptr());
82 
83  i++;
84  }
85 }
86 
87 
88 template<class T, class Mesh>
89 void Foam::fvMeshDistribute::mapBoundaryFields
90 (
91  const polyTopoChangeMap& map,
92  const PtrList<FieldField<fvsPatchField, T>>& oldBflds
93 )
94 {
95  // Map boundary field
96 
97  const labelList& oldPatchStarts = map.oldPatchStarts();
98  const labelList& faceMap = map.faceMap();
99 
101 
103  (
104  mesh_.objectRegistry::lookupClass<fldType>()
105  );
106 
107  if (flds.size() != oldBflds.size())
108  {
110  << abort(FatalError);
111  }
112 
113  label fieldi = 0;
114 
115  forAllIter(typename HashTable<fldType*>, flds, iter)
116  {
117  fldType& fld = *iter();
118  typename fldType::Boundary& bfld =
119  fld.boundaryFieldRef();
120 
121  const FieldField<fvsPatchField, T>& oldBfld = oldBflds[fieldi++];
122 
123  // Pull from old boundary field into bfld.
124 
125  forAll(bfld, patchi)
126  {
127  fvsPatchField<T>& patchFld = bfld[patchi];
128  label facei = patchFld.patch().start();
129 
130  forAll(patchFld, i)
131  {
132  label oldFacei = faceMap[facei++];
133 
134  // Find patch and local patch face oldFacei was in.
135  forAll(oldPatchStarts, oldPatchi)
136  {
137  label oldLocalI = oldFacei - oldPatchStarts[oldPatchi];
138 
139  if (oldLocalI >= 0 && oldLocalI < oldBfld[oldPatchi].size())
140  {
141  patchFld[i] = oldBfld[oldPatchi][oldLocalI];
142  }
143  }
144  }
145  }
146  }
147 }
148 
149 
150 template<class T>
151 void Foam::fvMeshDistribute::saveInternalFields
152 (
153  PtrList<Field<T>>& iflds
154 ) const
155 {
157 
159  (
160  static_cast<const fvMesh&>(mesh_).objectRegistry::lookupClass<fldType>()
161  );
162 
163  iflds.setSize(flds.size());
164 
165  label i = 0;
166 
167  forAllConstIter(typename HashTable<const fldType*>, flds, iter)
168  {
169  const fldType& fld = *iter();
170 
171  iflds.set(i, fld.primitiveField().clone());
172 
173  i++;
174  }
175 }
176 
177 
178 template<class T>
179 void Foam::fvMeshDistribute::mapExposedFaces
180 (
181  const polyTopoChangeMap& map,
182  const PtrList<Field<T>>& oldFlds
183 )
184 {
185  // Set boundary values of exposed internal faces
186 
187  const labelList& faceMap = map.faceMap();
188 
190 
192  (
193  mesh_.objectRegistry::lookupClass<fldType>()
194  );
195 
196  if (flds.size() != oldFlds.size())
197  {
198  FatalErrorIn("fvMeshDistribute::mapExposedFaces(..)") << "problem"
199  << abort(FatalError);
200  }
201 
202 
203  label fieldI = 0;
204 
205  forAllIter(typename HashTable<fldType*>, flds, iter)
206  {
207  fldType& fld = *iter();
208  typename fldType::Boundary& bfld = fld.boundaryFieldRef();
209  const bool negateIfFlipped = isFlux(fld);
210 
211  const Field<T>& oldInternal = oldFlds[fieldI++];
212 
213  // Pull from old internal field into bfld.
214 
215  forAll(bfld, patchi)
216  {
217  fvsPatchField<T>& patchFld = bfld[patchi];
218 
219  forAll(patchFld, i)
220  {
221  const label facei = patchFld.patch().start()+i;
222  const label oldFacei = faceMap[facei];
223 
224  if (oldFacei < oldInternal.size())
225  {
226  patchFld[i] = oldInternal[oldFacei];
227 
228  if (negateIfFlipped && map.flipFaceFlux().found(facei))
229  {
230  patchFld[i] = flipOp()(patchFld[i]);
231  }
232  }
233  }
234  }
235  }
236 }
237 
238 
239 template<class GeoField>
240 void Foam::fvMeshDistribute::correctProcessorPatchFields()
241 {
243  processorPatchFieldType;
244 
246  (
247  mesh_.objectRegistry::lookupClass<GeoField>()
248  );
249 
250  forAllIter(typename HashTable<GeoField*>, flds, iter)
251  {
252  GeoField& fld = *iter();
253 
254  typename GeoField::Boundary& bfld = fld.boundaryFieldRef();
255 
256  if
257  (
260  )
261  {
262  label nReq = Pstream::nRequests();
263 
264  forAll(bfld, patchi)
265  {
266  if (isA<processorPatchFieldType>(bfld[patchi]))
267  {
268  bfld[patchi].initEvaluate(Pstream::defaultCommsType);
269  }
270  }
271 
272  // Block for any outstanding requests
273  if
274  (
277  )
278  {
279  Pstream::waitRequests(nReq);
280  }
281 
282  forAll(bfld, patchi)
283  {
284  if (isA<processorPatchFieldType>(bfld[patchi]))
285  {
286  bfld[patchi].evaluate(Pstream::defaultCommsType);
287  }
288  }
289  }
291  {
292  const lduSchedule& patchSchedule =
293  mesh_.globalData().patchSchedule();
294 
295  forAll(patchSchedule, patchEvali)
296  {
297  if (isA<processorPatchFieldType>(bfld[patchEvali]))
298  {
299  if (patchSchedule[patchEvali].init)
300  {
301  bfld[patchSchedule[patchEvali].patch]
302  .initEvaluate(Pstream::commsTypes::scheduled);
303  }
304  else
305  {
306  bfld[patchSchedule[patchEvali].patch]
308  }
309  }
310  }
311  }
312  }
313 }
314 
315 
316 template<class GeoField>
317 void Foam::fvMeshDistribute::sendFields
318 (
319  const label domain,
320  const wordList& fieldNames,
321  const fvMeshSubset& subsetter,
322  Ostream& toNbr
323 )
324 {
325  // Send fields. Note order supplied so we can receive in exactly the same
326  // order.
327  // Note that field gets written as entry in dictionary so we
328  // can construct from subdictionary.
329  // (since otherwise the reading as-a-dictionary mixes up entries from
330  // consecutive fields)
331  // The dictionary constructed is:
332  // volScalarField
333  // {
334  // p {internalField ..; boundaryField ..;}
335  // k {internalField ..; boundaryField ..;}
336  // }
337  // volVectorField
338  // {
339  // U {internalField ... }
340  // }
341 
342  // volVectorField {U {internalField ..; boundaryField ..;}}
343 
344  toNbr << GeoField::typeName << token::NL << token::BEGIN_BLOCK << token::NL;
345  forAll(fieldNames, i)
346  {
347  if (debug)
348  {
349  Pout<< "Subsetting field " << fieldNames[i]
350  << " for domain:" << domain << endl;
351  }
352 
353  // Send all fieldNames. This has to be exactly the same set as is
354  // being received!
355  const GeoField& fld =
356  subsetter.baseMesh().lookupObject<GeoField>(fieldNames[i]);
357 
358  tmp<GeoField> tsubfld = subsetter.interpolate(fld);
359 
360  toNbr
361  << fieldNames[i] << token::NL << token::BEGIN_BLOCK
362  << tsubfld
363  << token::NL << token::END_BLOCK << token::NL;
364  }
365  toNbr << token::END_BLOCK << token::NL;
366 }
367 
368 
369 template<class GeoField>
370 void Foam::fvMeshDistribute::receiveFields
371 (
372  const label domain,
373  const wordList& fieldNames,
374  typename GeoField::Mesh& mesh,
375  PtrList<GeoField>& fields,
376  const dictionary& fieldDicts
377 )
378 {
379  if (debug)
380  {
381  Pout<< "Receiving fields " << fieldNames
382  << " from domain:" << domain << endl;
383  }
384 
385  fields.setSize(fieldNames.size());
386 
387  forAll(fieldNames, i)
388  {
389  if (debug)
390  {
391  Pout<< "Constructing field " << fieldNames[i]
392  << " from domain:" << domain << endl;
393  }
394 
395  fields.set
396  (
397  i,
398  new GeoField
399  (
400  IOobject
401  (
402  fieldNames[i],
403  mesh.thisDb().time().timeName(),
404  mesh.thisDb(),
407  ),
408  mesh,
409  fieldDicts.subDict(fieldNames[i])
410  )
411  );
412  }
413 }
414 
415 
416 // ************************************************************************* //
const labelHashSet & flipFaceFlux() const
Map of flipped face flux faces.
Class containing functor to negate primitives. Dummy for all other types.
Definition: flipOp.H:50
const labelList & oldPatchStarts() const
Return list of the old patch start labels.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
const labelList & faceMap() const
Old face map.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#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
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const lduSchedule & patchSchedule() const
Order in which the patches should be initialised/evaluated.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:137
Generic GeometricField class.
fvMesh & mesh
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Generic field type.
Definition: FieldField.H:51
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1435
An STL-conforming hash table.
Definition: HashTable.H:61
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static void printFieldInfo(const fvMesh &)
Print some field info.
Post-processing mesh subset tool. Given the original mesh and the list of selected cells...
Definition: fvMeshSubset.H:73
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:147
bool isFlux(const DimensionedField< Type, surfaceMesh > &df)
Check if surfaceField is a flux.
Definition: surfaceFields.H:55
const fvPatch & patch() const
Return patch.
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:272
static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &cellMap, const labelList &faceMap)
Map volume field.
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
label patchi
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
This boundary condition enables processor communication across patches.
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:301
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
const fvMesh & baseMesh() const
Original mesh.
Definition: fvMeshSubset.H:222