Loading...
Searching...
No Matches
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-2025 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 ZoneMesh>
34void Foam::fvMeshDistribute::reorderZones
35(
36 const wordList& zoneNames,
37 ZoneMesh& zones
38)
39{
40 zones.clearAddressing();
41
42 typedef typename ZoneMesh::zone_type zone_type;
43
44 // Shift old ones to new position
45 PtrList<zone_type>& oldPtrs = zones;
46 PtrList<zone_type> newPtrs(zoneNames.size());
47
48 forAll(zones, zonei)
49 {
50 if (!zones.get(zonei))
51 {
53 << "Problem with zone[" << zonei << "] " << zones.names()
54 << exit(FatalError);
55 }
56 }
57
58 forAll(zones, zonei)
59 {
60 // Take ownership
61 auto ptr = oldPtrs.release(zonei);
62 if (!ptr)
63 {
65 << "Problem with zone[" << zonei << "] == nullptr" << nl
66 << exit(FatalError);
67 }
68
69 const label newIndex = zoneNames.find(ptr->name());
70 // Check if not found?
71 ptr->index() = newIndex;
72 newPtrs.set(newIndex, std::move(ptr));
73 }
74
75 // Add empty zones for unknown ones
76 forAll(newPtrs, i)
77 {
78 if (!newPtrs.get(i))
79 {
80 newPtrs.emplace_set(i, zoneNames[i], i, zones);
81 }
82 }
84 // Transfer
85 oldPtrs = std::move(newPtrs);
86}
87
88
89template<class GeoField>
91{
92 typedef GeometricField
93 <
94 typename GeoField::value_type,
97 > excludeType;
98
99 for (const GeoField& field : mesh.objectRegistry::csorted<GeoField>())
100 {
102 {
103 Pout<< "Field:" << field.name() << " size:" << field.size()
104 //<< " value:" << field
106 }
107 }
108}
109
110
111template<class GeoField>
113{
114 for (const GeoField& field : mesh.objectRegistry::csorted<GeoField>())
115 {
116 Pout<< "Field:" << field.name() << " size:" << field.size()
117 //<< " value:" << field
118 << endl;
119
120 for (const auto& patchFld : field.boundaryField())
121 {
122 Pout<< " " << patchFld.patch().index()
123 << ' ' << patchFld.patch().name()
124 << ' ' << patchFld.type()
125 << ' ' << patchFld.size()
126 << nl;
127 }
128 }
129}
130
131
132template<class T, class Mesh>
133void Foam::fvMeshDistribute::saveBoundaryFields
134(
136) const
137{
138 // Save whole boundary field
139
141
142 const UPtrList<const fldType> flds
143 (
144 mesh_.objectRegistry::csorted<fldType>()
145 );
146
147 bflds.resize_null(flds.size());
148
149 label fieldi = 0;
150 for (const fldType& fld : flds)
151 {
152 bflds.set(fieldi, fld.boundaryField().clone());
153
154 ++fieldi;
155 }
156}
157
158
159template<class T, class Mesh>
160void Foam::fvMeshDistribute::mapBoundaryFields
161(
162 const mapPolyMesh& map,
164)
165{
166 // Map boundary field
167
168 const labelList& oldPatchStarts = map.oldPatchStarts();
169 const labelList& faceMap = map.faceMap();
170
172
174 (
175 mesh_.objectRegistry::sorted<fldType>()
176 );
177
178 if (flds.size() != oldBflds.size())
179 {
181 << abort(FatalError);
182 }
183
184 forAll(flds, fieldi)
185 {
186 auto& bfld = flds[fieldi].boundaryFieldRef();
187 const auto& oldBfld = oldBflds[fieldi];
188
189 // Pull from old boundary field into bfld.
190
191 forAll(bfld, patchi)
192 {
193 fvsPatchField<T>& patchFld = bfld[patchi];
194 label facei = patchFld.patch().start();
195
196 forAll(patchFld, i)
197 {
198 label oldFacei = faceMap[facei++];
199
200 // Find patch and local patch face oldFacei was in.
201 forAll(oldPatchStarts, oldPatchi)
202 {
203 label oldLocalI = oldFacei - oldPatchStarts[oldPatchi];
204
205 if (oldLocalI >= 0 && oldLocalI < oldBfld[oldPatchi].size())
206 {
207 patchFld[i] = oldBfld[oldPatchi][oldLocalI];
208 }
209 }
210 }
211 }
212 }
213}
214
215
216template<class T>
217void Foam::fvMeshDistribute::saveInternalFields
218(
219 PtrList<Field<T>>& iflds
220) const
221{
223
225 (
226 mesh_.objectRegistry::csorted<fldType>()
227 );
228
229 iflds.resize_null(fields.size());
230
231 forAll(fields, fieldi)
232 {
233 // Bypass DynamicField level, save as Field
234 iflds.set
235 (
236 fieldi,
237 static_cast<const Field<T>&>(fields[fieldi]).clone()
238 );
239 }
240}
241
242
243template<class T>
244void Foam::fvMeshDistribute::mapExposedFaces
245(
246 const mapPolyMesh& map,
247 const PtrList<Field<T>>& oldFlds
248)
249{
250 // Set boundary values of exposed internal faces
251
252 const labelList& faceMap = map.faceMap();
253
255
257 (
258 mesh_.objectRegistry::sorted<fldType>()
259 );
260
261 if (flds.size() != oldFlds.size())
262 {
264 << "problem"
265 << abort(FatalError);
266 }
267
268
269 forAll(flds, fieldi)
270 {
271 auto& fld = flds[fieldi];
272 const auto& oldInternal = oldFlds[fieldi];
273
274 const bool oriented = fld.is_oriented();
275
276 auto& bfld = fld.boundaryFieldRef();
277
278
279 // Pull from old internal field into bfld.
280
281 forAll(bfld, patchi)
282 {
283 fvsPatchField<T>& patchFld = bfld[patchi];
284
285 forAll(patchFld, i)
286 {
287 const label faceI = patchFld.patch().start()+i;
288
289 label oldFaceI = faceMap[faceI];
290
291 if (oldFaceI < oldInternal.size())
292 {
293 patchFld[i] = oldInternal[oldFaceI];
294
295 if (oriented && map.flipFaceFlux().found(faceI))
296 {
297 patchFld[i] = flipOp()(patchFld[i]);
298 }
299 }
300 }
301 }
302 }
303}
304
305
306template<class GeoField, class PatchFieldType>
307void Foam::fvMeshDistribute::initPatchFields
308(
309 const typename GeoField::value_type& initVal
310)
311{
312 // Init patch fields of certain type
313 // - field order is irrelevant
314
315 for (GeoField& fld : mesh_.objectRegistry::objects<GeoField>())
316 {
317 auto& bfld = fld.boundaryFieldRef();
318
319 forAll(bfld, patchi)
320 {
321 if (isA<PatchFieldType>(bfld[patchi]))
322 {
323 bfld[patchi] == initVal;
324 }
325 }
326 }
327}
328
329
330//template<class GeoField>
331//void Foam::fvMeshDistribute::correctBoundaryConditions()
332//{
333// // CorrectBoundaryConditions patch fields of certain type
334//
335// for (GeoField& fld : mesh_.objectRegistry::sorted<GeoField>())
336// {
337// fld.correctBoundaryConditions();
338// }
339//}
340
341
342template<class GeoField>
343void Foam::fvMeshDistribute::getFieldNames
344(
345 const fvMesh& mesh,
346 HashTable<wordList>& allFieldNames,
347 const word& excludeType,
348 const bool syncPar
349)
350{
351 wordList& list = allFieldNames(GeoField::typeName);
352 list = mesh.sortedNames<GeoField>();
353
354 if (!excludeType.empty())
355 {
356 const wordList& excludeList =
357 allFieldNames.lookup(excludeType, wordList::null());
358
359 if (!excludeList.empty())
360 {
361 DynamicList<word> newList(list.size());
362 for (const auto& name : list)
363 {
364 if (!excludeList.contains(name))
365 {
366 newList.push_back(name);
367 }
368 }
369 if (newList.size() < list.size())
370 {
371 list = std::move(newList);
372 }
373 }
374 }
375
376
377 // Check all procs have same names
378 if (syncPar && Pstream::parRun())
379 {
380 // Check and report error(s) on master
381 // - don't need indexing on master itself
382
383 const globalIndex procAddr(globalIndex::gatherNonLocal{}, list.size());
384
385 const wordList allNames(procAddr.gather(list));
386
387 // Automatically restricted to master
388 for (const int proci : procAddr.subProcs())
389 {
390 const auto procNames(allNames.slice(procAddr.range(proci)));
391
392 if (procNames != list)
393 {
395 << "When checking for equal " << GeoField::typeName
396 << " :" << nl
397 << "processor0 has:" << list << nl
398 << "processor" << proci << " has:" << procNames << nl
399 << GeoField::typeName
400 << " need to be synchronised on all processors."
401 << exit(FatalError);
402 break;
403 }
404 }
405 }
406}
407
408
409template<class GeoField>
410void Foam::fvMeshDistribute::sendFields
411(
412 const label domain,
413 const HashTable<wordList>& allFieldNames,
414 const fvMeshSubset& subsetter,
415 Ostream& toNbr
416)
417{
418 // Send fields. Note order supplied so we can receive in exactly the same
419 // order.
420 // Note that field gets written as entry in dictionary so we
421 // can construct from subdictionary.
422 // (since otherwise the reading as-a-dictionary mixes up entries from
423 // consecutive fields)
424 // The dictionary constructed is:
425 // volScalarField
426 // {
427 // p {internalField ..; boundaryField ..;}
428 // k {internalField ..; boundaryField ..;}
429 // }
430 // volVectorField
431 // {
432 // U {internalField ... }
433 // }
434
435 // volVectorField {U {internalField ..; boundaryField ..;}}
436
437 const wordList& fieldNames =
438 allFieldNames.lookup(GeoField::typeName, wordList::null());
439
440 toNbr << GeoField::typeName << token::NL << token::BEGIN_BLOCK << token::NL;
441
442 for (const word& fieldName : fieldNames)
443 {
444 if (debug)
445 {
446 Pout<< "Subsetting " << GeoField::typeName
447 << " field " << fieldName
448 << " for domain:" << domain << endl;
449 }
450
451 // Send all fieldNames. This has to be exactly the same set as is
452 // being received!
453 const GeoField& fld =
454 subsetter.baseMesh().lookupObject<GeoField>(fieldName);
455
456 // Note: use subsetter to get sub field. Override default behaviour
457 // to warn for unset fields since they will be reset later on
458 tmp<GeoField> tsubfld = subsetter.interpolate(fld, true);
459
460 toNbr
461 << fieldName << token::NL << token::BEGIN_BLOCK
462 << tsubfld
464 }
465 toNbr << token::END_BLOCK << token::NL;
466}
467
468
469template<class GeoField>
470void Foam::fvMeshDistribute::receiveFields
471(
472 const label domain,
473 const HashTable<wordList>& allFieldNames,
474 fvMesh& mesh,
476 const dictionary& allFieldsDict
477)
478{
479 // Opposite of sendFields
480
481 const wordList& fieldNames =
482 allFieldNames.lookup(GeoField::typeName, wordList::null());
483
484 const dictionary& fieldDicts =
485 allFieldsDict.subDict(GeoField::typeName);
486
487
488 if (debug)
489 {
490 Pout<< "Receiving:" << GeoField::typeName
491 << " fields:" << fieldNames
492 << " from domain:" << domain << endl;
493 }
494
495 fields.resize(fieldNames.size());
496
497 label fieldi = 0;
498 for (const word& fieldName : fieldNames)
499 {
500 if (debug)
501 {
502 Pout<< "Constructing type:" << GeoField::typeName
503 << " field:" << fieldName
504 << " from domain:" << domain << endl;
505 }
506
507 fields.set
508 (
509 fieldi++,
510 new GeoField
511 (
513 (
514 fieldName,
515 mesh.time().timeName(),
516 mesh,
519 ),
520 mesh,
521 fieldDicts.subDict(fieldName)
522 )
523 );
524 }
525}
526
527
528// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Generic GeometricField class.
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
@ NO_READ
Nothing to be read.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static const List< word > & null() noexcept
Definition List.H:138
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
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 bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
ZoneType zone_type
The zone type. Same as PtrList<ZoneType>::value_type.
Definition ZoneMesh.H:142
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
static void printIntFieldInfo(const fvMesh &)
Print some field info.
static void printFieldInfo(const fvMesh &)
Print some field info.
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
label start() const noexcept
The patch start within the polyMesh face list.
Definition fvPatch.H:226
const fvPatch & patch() const noexcept
Return the patch.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
wordList sortedNames() const
The sorted names of all objects.
A class for managing temporary objects.
Definition tmp.H:75
@ BEGIN_BLOCK
Begin block [isseparator].
Definition token.H:178
@ END_BLOCK
End block [isseparator].
Definition token.H:179
@ NL
Newline [isspace].
Definition token.H:143
Mesh data needed to do the Finite Volume discretisation.
Definition volMesh.H:47
A class for handling words, derived from Foam::string.
Definition word.H:66
rDeltaTY field()
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for handling debugging switches.
Definition debug.C:45
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
List< word > wordList
List of word.
Definition fileName.H:60
List< label > labelList
A List of labels.
Definition List.H:62
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
errorManip< error > abort(error &err)
Definition errorManip.H:139
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Functor to negate primitives. Dummy for most other types.
Definition flipOp.H:67
Dispatch tag: Construct 'one-sided' from the non-master local sizes using gather but no broadcast.