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-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
33template<class ZoneType, class ZoneMesh>
34void 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
80template<class GeoField>
82{
83 typedef GeometricField
84 <
85 typename GeoField::value_type,
88 > excludeType;
89
91 (
92 mesh.objectRegistry::lookupClass<GeoField>()
93 );
94
95 forAllConstIters(flds, iter)
96 {
97 const GeoField& fld = *iter();
98 if (!isA<excludeType>(fld))
99 {
100 Pout<< "Field:" << iter.key() << " internalsize:" << fld.size()
101 //<< " value:" << fld
102 << endl;
103 }
104 }
105}
106
107
108template<class GeoField>
110{
112 (
113 mesh.objectRegistry::lookupClass<GeoField>()
114 );
115
116 forAllConstIters(flds, iter)
117 {
118 const GeoField& fld = *iter();
119
120 Pout<< "Field:" << iter.key() << " internalsize:" << fld.size()
121 //<< " value:" << fld
122 << endl;
123
124 for (const auto& patchFld : fld.boundaryField())
125 {
126 Pout<< " " << patchFld.patch().index()
127 << ' ' << patchFld.patch().name()
128 << ' ' << patchFld.type()
129 << ' ' << patchFld.size()
130 << nl;
131 }
132 }
133}
134
135
136template<class T, class Mesh>
137void Foam::fvMeshDistribute::saveBoundaryFields
138(
140) const
141{
142 // Save whole boundary field
143
145
147 (
148 mesh_.objectRegistry::lookupClass<const fldType>()
149 );
150
151 bflds.resize(flds.size());
152
153 label i = 0;
154 forAllConstIters(flds, iter)
155 {
156 const fldType& fld = *iter();
157
158 bflds.set(i, fld.boundaryField().clone());
159
160 ++i;
161 }
162}
163
164
165template<class T, class Mesh>
166void Foam::fvMeshDistribute::mapBoundaryFields
167(
168 const mapPolyMesh& map,
169 const PtrList<FieldField<fvsPatchField, T>>& oldBflds
170)
171{
172 // Map boundary field
173
174 const labelList& oldPatchStarts = map.oldPatchStarts();
175 const labelList& faceMap = map.faceMap();
176
177 typedef GeometricField<T, fvsPatchField, Mesh> fldType;
178
179 HashTable<fldType*> flds
180 (
181 mesh_.objectRegistry::lookupClass<fldType>()
182 );
183
184 if (flds.size() != oldBflds.size())
185 {
187 << abort(FatalError);
188 }
189
190 label fieldi = 0;
191
192 forAllIters(flds, iter)
193 {
194 fldType& fld = *iter();
195 auto& bfld = fld.boundaryFieldRef();
196
197 const FieldField<fvsPatchField, T>& oldBfld = oldBflds[fieldi++];
198
199 // Pull from old boundary field into bfld.
200
201 forAll(bfld, patchi)
202 {
203 fvsPatchField<T>& patchFld = bfld[patchi];
204 label facei = patchFld.patch().start();
205
206 forAll(patchFld, i)
207 {
208 label oldFacei = faceMap[facei++];
209
210 // Find patch and local patch face oldFacei was in.
211 forAll(oldPatchStarts, oldPatchi)
212 {
213 label oldLocalI = oldFacei - oldPatchStarts[oldPatchi];
214
215 if (oldLocalI >= 0 && oldLocalI < oldBfld[oldPatchi].size())
216 {
217 patchFld[i] = oldBfld[oldPatchi][oldLocalI];
218 }
219 }
220 }
221 }
222 }
223}
224
225
226template<class T>
227void Foam::fvMeshDistribute::saveInternalFields
228(
229 PtrList<Field<T>>& iflds
230) const
231{
232 typedef GeometricField<T, fvsPatchField, surfaceMesh> fldType;
233
234 HashTable<const fldType*> flds
235 (
236 mesh_.objectRegistry::lookupClass<const fldType>()
237 );
238
239 iflds.resize(flds.size());
240
241 label i = 0;
242
243 forAllConstIters(flds, iter)
244 {
245 const fldType& fld = *iter();
246
247 iflds.set(i, fld.primitiveField().clone());
248
249 ++i;
250 }
251}
252
253
254template<class T>
255void Foam::fvMeshDistribute::mapExposedFaces
256(
257 const mapPolyMesh& map,
258 const PtrList<Field<T>>& oldFlds
259)
260{
261 // Set boundary values of exposed internal faces
262
263 const labelList& faceMap = map.faceMap();
264
265 typedef GeometricField<T, fvsPatchField, surfaceMesh> fldType;
266
267 HashTable<fldType*> flds
268 (
269 mesh_.objectRegistry::lookupClass<fldType>()
270 );
271
272 if (flds.size() != oldFlds.size())
273 {
275 << "problem"
276 << abort(FatalError);
277 }
278
279
280 label fieldI = 0;
281
282 forAllIters(flds, iter)
283 {
284 fldType& fld = *iter();
285 const bool oriented = fld.oriented()();
286
287 typename fldType::Boundary& bfld = fld.boundaryFieldRef();
288
289 const Field<T>& oldInternal = oldFlds[fieldI++];
290
291 // Pull from old internal field into bfld.
292
293 forAll(bfld, patchi)
294 {
295 fvsPatchField<T>& patchFld = bfld[patchi];
296
297 forAll(patchFld, i)
298 {
299 const label faceI = patchFld.patch().start()+i;
300
301 label oldFaceI = faceMap[faceI];
302
303 if (oldFaceI < oldInternal.size())
304 {
305 patchFld[i] = oldInternal[oldFaceI];
306
307 if (oriented && map.flipFaceFlux().found(faceI))
308 {
309 patchFld[i] = flipOp()(patchFld[i]);
310 }
311 }
312 }
313 }
314 }
315}
316
317
318template<class GeoField, class PatchFieldType>
319void Foam::fvMeshDistribute::initPatchFields
320(
321 const typename GeoField::value_type& initVal
322)
323{
324 // Init patch fields of certain type
325
326 HashTable<GeoField*> flds
327 (
328 mesh_.objectRegistry::lookupClass<GeoField>()
329 );
330
331 forAllIters(flds, iter)
332 {
333 GeoField& fld = *iter();
334
335 auto& bfld = fld.boundaryFieldRef();
336
337 forAll(bfld, patchi)
338 {
339 if (isA<PatchFieldType>(bfld[patchi]))
340 {
341 bfld[patchi] == initVal;
342 }
343 }
344 }
345}
346
347
348//template<class GeoField>
349//void Foam::fvMeshDistribute::correctBoundaryConditions()
350//{
351// // CorrectBoundaryConditions patch fields of certain type
352//
353// HashTable<GeoField*> flds
354// (
355// mesh_.objectRegistry::lookupClass<GeoField>()
356// );
357//
358// forAllIters(flds, iter)
359// {
360// GeoField& fld = *iter();
361// fld.correctBoundaryConditions();
362// }
363//}
364
365
366template<class GeoField>
367void Foam::fvMeshDistribute::getFieldNames
368(
369 const fvMesh& mesh,
370 HashTable<wordList>& allFieldNames,
371 const word& excludeType,
372 const bool syncPar
373)
374{
375 wordList& list = allFieldNames(GeoField::typeName);
376 list = mesh.sortedNames<GeoField>();
377
378 if (!excludeType.empty())
379 {
380 const wordList& excludeList = allFieldNames(excludeType);
381
382 DynamicList<word> newList(list.size());
383 for(const auto& name : list)
384 {
385 if (!excludeList.found(name))
386 {
387 newList.append(name);
388 }
389 }
390 if (newList.size() < list.size())
391 {
392 list = std::move(newList);
393 }
394 }
395
396
397 // Check all procs have same names
398 if (syncPar && Pstream::parRun())
399 {
400 // Check and report error(s) on master
401
402 const globalIndex procAddr
403 (
404 // Don't need to collect master itself
405 (Pstream::master() ? 0 : list.size()),
406 globalIndex::gatherOnly{}
407 );
408
409 const wordList allNames(procAddr.gather(list));
410
411 // Automatically restricted to master
412 for (const int proci : procAddr.subProcs())
413 {
414 const auto procNames(allNames.slice(procAddr.range(proci)));
415
416 if (procNames != list)
417 {
419 << "When checking for equal " << GeoField::typeName
420 << " :" << nl
421 << "processor0 has:" << list << nl
422 << "processor" << proci << " has:" << procNames << nl
423 << GeoField::typeName
424 << " need to be synchronised on all processors."
425 << exit(FatalError);
426 break;
427 }
428 }
429 }
430}
431
432
433template<class GeoField>
434void Foam::fvMeshDistribute::sendFields
435(
436 const label domain,
437 const HashTable<wordList>& allFieldNames,
438 const fvMeshSubset& subsetter,
439 Ostream& toNbr
440)
441{
442 // Send fields. Note order supplied so we can receive in exactly the same
443 // order.
444 // Note that field gets written as entry in dictionary so we
445 // can construct from subdictionary.
446 // (since otherwise the reading as-a-dictionary mixes up entries from
447 // consecutive fields)
448 // The dictionary constructed is:
449 // volScalarField
450 // {
451 // p {internalField ..; boundaryField ..;}
452 // k {internalField ..; boundaryField ..;}
453 // }
454 // volVectorField
455 // {
456 // U {internalField ... }
457 // }
458
459 // volVectorField {U {internalField ..; boundaryField ..;}}
460
461 const wordList& fieldNames =
462 allFieldNames.lookup(GeoField::typeName, wordList::null());
463
464 toNbr << GeoField::typeName << token::NL << token::BEGIN_BLOCK << token::NL;
465
466 for (const word& fieldName : fieldNames)
467 {
468 if (debug)
469 {
470 Pout<< "Subsetting " << GeoField::typeName
471 << " field " << fieldName
472 << " for domain:" << domain << endl;
473 }
474
475 // Send all fieldNames. This has to be exactly the same set as is
476 // being received!
477 const GeoField& fld =
478 subsetter.baseMesh().lookupObject<GeoField>(fieldName);
479
480 // Note: use subsetter to get sub field. Override default behaviour
481 // to warn for unset fields since they will be reset later on
482 tmp<GeoField> tsubfld = subsetter.interpolate(fld, true);
483
484 toNbr
485 << fieldName << token::NL << token::BEGIN_BLOCK
486 << tsubfld
488 }
489 toNbr << token::END_BLOCK << token::NL;
490}
491
492
493template<class GeoField>
494void Foam::fvMeshDistribute::receiveFields
495(
496 const label domain,
497 const HashTable<wordList>& allFieldNames,
498 fvMesh& mesh,
499 PtrList<GeoField>& fields,
500 const dictionary& allFieldsDict
501)
502{
503 // Opposite of sendFields
504
505 const wordList& fieldNames =
506 allFieldNames.lookup(GeoField::typeName, wordList::null());
507
508 const dictionary& fieldDicts =
509 allFieldsDict.subDict(GeoField::typeName);
510
511
512 if (debug)
513 {
514 Pout<< "Receiving:" << GeoField::typeName
515 << " fields:" << fieldNames
516 << " from domain:" << domain << endl;
517 }
518
519 fields.resize(fieldNames.size());
520
521 label fieldi = 0;
522 for (const word& fieldName : fieldNames)
523 {
524 if (debug)
525 {
526 Pout<< "Constructing type:" << GeoField::typeName
527 << " field:" << fieldName
528 << " from domain:" << domain << endl;
529 }
530
531 fields.set
532 (
533 fieldi++,
534 new GeoField
535 (
536 IOobject
537 (
538 fieldName,
539 mesh.time().timeName(),
540 mesh,
543 ),
544 mesh,
545 fieldDicts.subDict(fieldName)
546 )
547 );
548 }
549}
550
551
552// ************************************************************************* //
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
globalIndex procAddr(aMesh.nFaces())
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:80
Generic GeometricField class.
A HashTable similar to std::unordered_map.
Definition: HashTable.H:123
static const List< word > & null()
Return a null List.
Definition: ListI.H:109
virtual const fileName & name() const
Get the name of the stream.
Definition: OSstream.H:107
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
Type type(bool followLink=true, bool checkGzip=false) const
Definition: fileName.C:360
static void printIntFieldInfo(const fvMesh &)
Print some field info.
static void printFieldInfo(const fvMesh &)
Print some field info.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
labelRange range() const
Return start/size range of local processor data.
Definition: globalIndexI.H:232
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]).
labelRange subProcs() const noexcept
Range of process indices for addressed sub-offsets (processes)
Definition: globalIndexI.H:159
wordList sortedNames() const
The sorted names of all objects.
splitCell * master() const
Definition: splitCell.H:113
@ BEGIN_BLOCK
Begin block [isseparator].
Definition: token.H:159
@ END_BLOCK
End block [isseparator].
Definition: token.H:160
@ NL
Newline [isspace].
Definition: token.H:124
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:54
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
List< word > wordList
A List of words.
Definition: fileName.H:63
List< label > labelList
A List of labels.
Definition: List.H:66
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:260
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278