fvMeshDistributeTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "mapPolyMesh.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class ZoneType, class ZoneMesh>
34 void Foam::fvMeshDistribute::reorderZones
35 (
36  const wordList& zoneNames,
37  ZoneMesh& zones
38 )
39 {
40  zones.clearAddressing();
41 
42  // Shift old ones to new position
43  UPtrList<ZoneType> newZonePtrs(zoneNames.size());
44  forAll(zones, zonei)
45  {
46  auto* zonePtr = zones.get(zonei);
47  if (!zonePtr)
48  {
49  FatalErrorInFunction << "Problem with zones " << zones.names()
50  << exit(FatalError);
51  }
52  const label newIndex = zoneNames.find(zonePtr->name());
53  zonePtr->index() = newIndex;
54  newZonePtrs.set(newIndex, zonePtr);
55  }
56 
57  // Add empty zones for unknown ones
58  forAll(newZonePtrs, i)
59  {
60  if (!newZonePtrs.get(i))
61  {
62  newZonePtrs.set
63  (
64  i,
65  new ZoneType
66  (
67  zoneNames[i],
68  i,
69  zones
70  )
71  );
72  }
73  }
74 
75  // Transfer
76  zones.swap(newZonePtrs);
77 }
78 
79 
80 template<class GeoField>
82 {
83  typedef GeometricField
84  <
85  typename GeoField::value_type,
87  volMesh
88  > excludeType;
89 
90  for (const GeoField& field : mesh.objectRegistry::csorted<GeoField>())
91  {
92  if (!isA<excludeType>(field))
93  {
94  Pout<< "Field:" << field.name() << " size:" << field.size()
95  //<< " value:" << field
96  << endl;
97  }
98  }
99 }
100 
101 
102 template<class GeoField>
104 {
105  for (const GeoField& field : mesh.objectRegistry::csorted<GeoField>())
106  {
107  Pout<< "Field:" << field.name() << " size:" << field.size()
108  //<< " value:" << field
109  << endl;
110 
111  for (const auto& patchFld : field.boundaryField())
112  {
113  Pout<< " " << patchFld.patch().index()
114  << ' ' << patchFld.patch().name()
115  << ' ' << patchFld.type()
116  << ' ' << patchFld.size()
117  << nl;
118  }
119  }
120 }
121 
122 
123 template<class T, class Mesh>
124 void Foam::fvMeshDistribute::saveBoundaryFields
125 (
126  PtrList<FieldField<fvsPatchField, T>>& bflds
127 ) const
128 {
129  // Save whole boundary field
130 
131  typedef GeometricField<T, fvsPatchField, Mesh> fldType;
132 
133  const UPtrList<const fldType> flds
134  (
135  mesh_.objectRegistry::csorted<fldType>()
136  );
137 
138  bflds.resize_null(flds.size());
139 
140  label fieldi = 0;
141  for (const fldType& fld : flds)
142  {
143  bflds.set(fieldi, fld.boundaryField().clone());
144 
145  ++fieldi;
146  }
147 }
148 
149 
150 template<class T, class Mesh>
151 void Foam::fvMeshDistribute::mapBoundaryFields
152 (
153  const mapPolyMesh& map,
154  const PtrList<FieldField<fvsPatchField, T>>& oldBflds
155 )
156 {
157  // Map boundary field
158 
159  const labelList& oldPatchStarts = map.oldPatchStarts();
160  const labelList& faceMap = map.faceMap();
161 
162  typedef GeometricField<T, fvsPatchField, Mesh> fldType;
163 
164  UPtrList<fldType> flds
165  (
166  mesh_.objectRegistry::sorted<fldType>()
167  );
168 
169  if (flds.size() != oldBflds.size())
170  {
172  << abort(FatalError);
173  }
174 
175  forAll(flds, fieldi)
176  {
177  auto& bfld = flds[fieldi].boundaryFieldRef();
178  const auto& oldBfld = oldBflds[fieldi];
179 
180  // Pull from old boundary field into bfld.
181 
182  forAll(bfld, patchi)
183  {
184  fvsPatchField<T>& patchFld = bfld[patchi];
185  label facei = patchFld.patch().start();
186 
187  forAll(patchFld, i)
188  {
189  label oldFacei = faceMap[facei++];
190 
191  // Find patch and local patch face oldFacei was in.
192  forAll(oldPatchStarts, oldPatchi)
193  {
194  label oldLocalI = oldFacei - oldPatchStarts[oldPatchi];
195 
196  if (oldLocalI >= 0 && oldLocalI < oldBfld[oldPatchi].size())
197  {
198  patchFld[i] = oldBfld[oldPatchi][oldLocalI];
199  }
200  }
201  }
202  }
203  }
204 }
205 
206 
207 template<class T>
208 void Foam::fvMeshDistribute::saveInternalFields
209 (
210  PtrList<Field<T>>& iflds
211 ) const
212 {
213  typedef GeometricField<T, fvsPatchField, surfaceMesh> fldType;
214 
215  const UPtrList<const fldType> fields
216  (
217  mesh_.objectRegistry::csorted<fldType>()
218  );
219 
220  iflds.resize_null(fields.size());
221 
222  forAll(fields, fieldi)
223  {
224  iflds.set(fieldi, fields[fieldi].primitiveField().clone());
225  }
226 }
227 
228 
229 template<class T>
230 void Foam::fvMeshDistribute::mapExposedFaces
231 (
232  const mapPolyMesh& map,
233  const PtrList<Field<T>>& oldFlds
234 )
235 {
236  // Set boundary values of exposed internal faces
237 
238  const labelList& faceMap = map.faceMap();
239 
240  typedef GeometricField<T, fvsPatchField, surfaceMesh> fldType;
241 
242  UPtrList<fldType> flds
243  (
244  mesh_.objectRegistry::sorted<fldType>()
245  );
246 
247  if (flds.size() != oldFlds.size())
248  {
250  << "problem"
251  << abort(FatalError);
252  }
253 
254 
255  forAll(flds, fieldi)
256  {
257  auto& fld = flds[fieldi];
258  const auto& oldInternal = oldFlds[fieldi];
259 
260  const bool oriented = fld.is_oriented();
261 
262  auto& bfld = fld.boundaryFieldRef();
263 
264 
265  // Pull from old internal field into bfld.
266 
267  forAll(bfld, patchi)
268  {
269  fvsPatchField<T>& patchFld = bfld[patchi];
270 
271  forAll(patchFld, i)
272  {
273  const label faceI = patchFld.patch().start()+i;
274 
275  label oldFaceI = faceMap[faceI];
276 
277  if (oldFaceI < oldInternal.size())
278  {
279  patchFld[i] = oldInternal[oldFaceI];
280 
281  if (oriented && map.flipFaceFlux().found(faceI))
282  {
283  patchFld[i] = flipOp()(patchFld[i]);
284  }
285  }
286  }
287  }
288  }
289 }
290 
291 
292 template<class GeoField, class PatchFieldType>
293 void Foam::fvMeshDistribute::initPatchFields
294 (
295  const typename GeoField::value_type& initVal
296 )
297 {
298  // Init patch fields of certain type
299  // - field order is irrelevant
300 
301  for (GeoField& fld : mesh_.objectRegistry::objects<GeoField>())
302  {
303  auto& bfld = fld.boundaryFieldRef();
304 
305  forAll(bfld, patchi)
306  {
307  if (isA<PatchFieldType>(bfld[patchi]))
308  {
309  bfld[patchi] == initVal;
310  }
311  }
312  }
313 }
314 
315 
316 //template<class GeoField>
317 //void Foam::fvMeshDistribute::correctBoundaryConditions()
318 //{
319 // // CorrectBoundaryConditions patch fields of certain type
320 //
321 // for (GeoField& fld : mesh_.objectRegistry::sorted<GeoField>())
322 // {
323 // fld.correctBoundaryConditions();
324 // }
325 //}
326 
327 
328 template<class GeoField>
329 void Foam::fvMeshDistribute::getFieldNames
330 (
331  const fvMesh& mesh,
332  HashTable<wordList>& allFieldNames,
333  const word& excludeType,
334  const bool syncPar
335 )
336 {
337  wordList& list = allFieldNames(GeoField::typeName);
338  list = mesh.sortedNames<GeoField>();
339 
340  if (!excludeType.empty())
341  {
342  const wordList& excludeList =
343  allFieldNames.lookup(excludeType, wordList::null());
344 
345  if (!excludeList.empty())
346  {
347  DynamicList<word> newList(list.size());
348  for (const auto& name : list)
349  {
350  if (!excludeList.contains(name))
351  {
352  newList.push_back(name);
353  }
354  }
355  if (newList.size() < list.size())
356  {
357  list = std::move(newList);
358  }
359  }
360  }
361 
362 
363  // Check all procs have same names
364  if (syncPar && Pstream::parRun())
365  {
366  // Check and report error(s) on master
367  // - don't need indexing on master itself
368 
369  const globalIndex procAddr(globalIndex::gatherNonLocal{}, list.size());
370 
371  const wordList allNames(procAddr.gather(list));
372 
373  // Automatically restricted to master
374  for (const int proci : procAddr.subProcs())
375  {
376  const auto procNames(allNames.slice(procAddr.range(proci)));
377 
378  if (procNames != list)
379  {
381  << "When checking for equal " << GeoField::typeName
382  << " :" << nl
383  << "processor0 has:" << list << nl
384  << "processor" << proci << " has:" << procNames << nl
385  << GeoField::typeName
386  << " need to be synchronised on all processors."
387  << exit(FatalError);
388  break;
389  }
390  }
391  }
392 }
393 
394 
395 template<class GeoField>
396 void Foam::fvMeshDistribute::sendFields
397 (
398  const label domain,
399  const HashTable<wordList>& allFieldNames,
400  const fvMeshSubset& subsetter,
401  Ostream& toNbr
402 )
403 {
404  // Send fields. Note order supplied so we can receive in exactly the same
405  // order.
406  // Note that field gets written as entry in dictionary so we
407  // can construct from subdictionary.
408  // (since otherwise the reading as-a-dictionary mixes up entries from
409  // consecutive fields)
410  // The dictionary constructed is:
411  // volScalarField
412  // {
413  // p {internalField ..; boundaryField ..;}
414  // k {internalField ..; boundaryField ..;}
415  // }
416  // volVectorField
417  // {
418  // U {internalField ... }
419  // }
420 
421  // volVectorField {U {internalField ..; boundaryField ..;}}
422 
423  const wordList& fieldNames =
424  allFieldNames.lookup(GeoField::typeName, wordList::null());
425 
426  toNbr << GeoField::typeName << token::NL << token::BEGIN_BLOCK << token::NL;
427 
428  for (const word& fieldName : fieldNames)
429  {
430  if (debug)
431  {
432  Pout<< "Subsetting " << GeoField::typeName
433  << " field " << fieldName
434  << " for domain:" << domain << endl;
435  }
436 
437  // Send all fieldNames. This has to be exactly the same set as is
438  // being received!
439  const GeoField& fld =
440  subsetter.baseMesh().lookupObject<GeoField>(fieldName);
441 
442  // Note: use subsetter to get sub field. Override default behaviour
443  // to warn for unset fields since they will be reset later on
444  tmp<GeoField> tsubfld = subsetter.interpolate(fld, true);
445 
446  toNbr
447  << fieldName << token::NL << token::BEGIN_BLOCK
448  << tsubfld
450  }
451  toNbr << token::END_BLOCK << token::NL;
452 }
453 
454 
455 template<class GeoField>
456 void Foam::fvMeshDistribute::receiveFields
457 (
458  const label domain,
459  const HashTable<wordList>& allFieldNames,
460  fvMesh& mesh,
461  PtrList<GeoField>& fields,
462  const dictionary& allFieldsDict
463 )
464 {
465  // Opposite of sendFields
466 
467  const wordList& fieldNames =
468  allFieldNames.lookup(GeoField::typeName, wordList::null());
469 
470  const dictionary& fieldDicts =
471  allFieldsDict.subDict(GeoField::typeName);
472 
473 
474  if (debug)
475  {
476  Pout<< "Receiving:" << GeoField::typeName
477  << " fields:" << fieldNames
478  << " from domain:" << domain << endl;
479  }
480 
481  fields.resize(fieldNames.size());
482 
483  label fieldi = 0;
484  for (const word& fieldName : fieldNames)
485  {
486  if (debug)
487  {
488  Pout<< "Constructing type:" << GeoField::typeName
489  << " field:" << fieldName
490  << " from domain:" << domain << endl;
491  }
492 
493  fields.set
494  (
495  fieldi++,
496  new GeoField
497  (
498  IOobject
499  (
500  fieldName,
501  mesh.time().timeName(),
502  mesh,
505  ),
506  mesh,
507  fieldDicts.subDict(fieldName)
508  )
509  );
510  }
511 }
512 
513 
514 // ************************************************************************* //
Begin block [isseparator].
Definition: token.H:165
labelRange range(const label proci) const
Return start/size range of proci data.
Definition: globalIndexI.H:282
rDeltaTY field()
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Type type(bool followLink=true, bool checkGzip=false) const
Return the directory entry type: UNDEFINED, FILE, DIRECTORY (or SYMLINK).
Definition: fileName.C:353
Newline [isspace].
Definition: token.H:130
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
static const List< word > & null()
Return a null List.
Definition: ListI.H:130
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:128
void push_back(const T &val)
Append an element at the end of the list.
Definition: ListI.H:227
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:45
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
wordList sortedNames() const
The sorted names of all objects.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
static void printFieldInfo(const fvMesh &)
Print some field info.
int debug
Static debugging option.
labelRange subProcs() const noexcept
Range of process indices for addressed sub-offsets (processes)
Definition: globalIndexI.H:197
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
List< word > wordList
List of word.
Definition: fileName.H:59
globalIndex procAddr(aMesh.nFaces())
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
Automatically write from objectRegistry::writeObject()
static void gather(const labelUList &offsets, const label comm, const ProcIDsContainer &procIDs, const UList< Type > &fld, List< Type > &allFld, const int tag=UPstream::msgType(), const UPstream::commsTypes=UPstream::commsTypes::nonBlocking)
Collect data in processor order on master (== procIDs[0]).
List< label > labelList
A List of labels.
Definition: List.H:62
End block [isseparator].
Definition: token.H:166
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
static void printIntFieldInfo(const fvMesh &)
Print some field info.