fieldsDistributorTemplates.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) 2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29
30template<class GeoField>
32(
33 const IOobject& io,
34 const typename GeoField::Mesh& mesh,
35 const label i,
37)
38{
39 fields.set(i, new GeoField(io, mesh));
40}
41
42template<class Type, template<class> class PatchField, class GeoMesh>
44(
45 const IOobject& io,
46 const typename GeoMesh::Mesh& mesh,
47 const label i,
49)
50{
51 fields.set
52 (
53 i,
55 );
56}
57
58
59template<class Type, template<class> class PatchField, class GeoMesh>
61(
62 const typename GeoMesh::Mesh& mesh,
63 const IOobjectList& objects,
65 const bool readOldTime
66)
67{
69
70 // GeoField fields - sorted for consistent order on all processors
71 UPtrList<const IOobject> fieldObjects(objects.sorted<GeoField>());
72
73 // Construct the fields
74 fields.resize(fieldObjects.size());
75
76 forAll(fieldObjects, i)
77 {
78 fields.set(i, new GeoField(fieldObjects[i], mesh, readOldTime));
79 }
80}
81
82
83template<class Mesh, class GeoField>
85(
86 const Mesh& mesh,
87 const IOobjectList& objects,
89)
90{
91 // GeoField fields - sorted for consistent order on all processors
92 UPtrList<const IOobject> fieldObjects(objects.sorted<GeoField>());
93
94 // Construct the fields
95 fields.resize(fieldObjects.size());
96
97 forAll(fieldObjects, i)
98 {
99 fields.set(i, new GeoField(fieldObjects[i], mesh));
100 }
101}
102
103
104template<class GeoField, class MeshSubsetter>
106(
107 const boolList& haveMeshOnProc,
108 const typename GeoField::Mesh& mesh,
109 const autoPtr<MeshSubsetter>& subsetterPtr,
110 IOobjectList& allObjects,
112 const bool deregister
113)
114{
115 // Get my objects of type
116 IOobjectList objects(allObjects.lookupClass<GeoField>());
117
118 // Check that we all have all objects
119 wordList objectNames = objects.sortedNames();
120
121 // Get master names
122 wordList masterNames(objectNames);
123 Pstream::broadcast(masterNames);
124
125 if (haveMeshOnProc[Pstream::myProcNo()] && objectNames != masterNames)
126 {
128 << "Objects not synchronised across processors." << nl
129 << "Master has " << flatOutput(masterNames) << nl
130 << "Processor " << Pstream::myProcNo()
131 << " has " << flatOutput(objectNames)
132 << exit(FatalError);
133 }
134
135 fields.clear();
136 fields.resize(masterNames.size());
137
138 if (fields.empty())
139 {
140 if (deregister)
141 {
142 // Extra safety - remove all such types
144 (
145 mesh.thisDb().objectRegistry::template lookupClass<GeoField>()
146 );
147
148 forAllConstIters(other, iter)
149 {
150 GeoField& fld = const_cast<GeoField&>(*iter.val());
151
152 if (!fld.ownedByRegistry())
153 {
154 fld.checkOut();
155 }
156 }
157 }
158
159 // Early exit
160 return;
161 }
162
163
164 // Have master send all fields to processors that don't have a mesh. The
165 // issue is if a patchField does any parallel operations inside its
166 // construct-from-dictionary. This will not work when going to more
167 // processors (e.g. decompose = 1 -> many) ! We could make a special
168 // exception for decomposePar but nicer would be to have read-communicator
169 // ... For now detect if decomposing & disable parRun
170 if (Pstream::master())
171 {
172 // Work out if we're decomposing - none of the subprocs has a mesh
173 bool decompose = true;
174 for (const int proci : Pstream::subProcs())
175 {
176 if (haveMeshOnProc[proci])
177 {
178 decompose = false;
179 }
180 }
181
182 const bool oldParRun = Pstream::parRun();
183 if (decompose)
184 {
185 Pstream::parRun(false);
186 }
187
188 forAll(masterNames, i)
189 {
190 const word& name = masterNames[i];
191 IOobject& io = *objects[name];
193
194 // Load field (but not oldTime)
195 readField(io, mesh, i, fields);
196 }
197
198 Pstream::parRun(oldParRun); // Restore any changes
199 }
200 else if (haveMeshOnProc[Pstream::myProcNo()])
201 {
202 // Have mesh so just try to load
203 forAll(masterNames, i)
204 {
205 const word& name = masterNames[i];
206 IOobject& io = *objects[name];
208
210
211 // Load field (but not oldTime)
212 readField(io, mesh, i, fields);
213 }
214 }
215
216
217 // Missing fields on any processors?
218 // - construct from dictionary
219
220 PtrList<dictionary> fieldDicts;
221
222 if (Pstream::master())
223 {
224 // Broadcast zero sized fields everywhere (if needed)
225 // Send like a list of dictionaries
226
227 OPBstream toProcs(UPstream::masterNo()); // worldComm
228
229 const label nDicts = (subsetterPtr ? fields.size() : label(0));
230
231 toProcs << nDicts << token::BEGIN_LIST; // Begin list
232
233 if (nDicts)
234 {
235 // Disable communication for interpolate() method
236 const bool oldParRun = Pstream::parRun(false);
237
238 const auto& subsetter = subsetterPtr();
239
240 forAll(fields, i)
241 {
242 tmp<GeoField> tsubfld = subsetter.interpolate(fields[i]);
243
244 // Surround each with {} as dictionary entry
245 toProcs.beginBlock();
246 toProcs << tsubfld();
247 toProcs.endBlock();
248 }
249
250 Pstream::parRun(oldParRun); // Restore state
251 }
252
253 toProcs << token::END_LIST << token::NL; // End list
254 }
255 else
256 {
257 // Receive the broadcast...
258 IPBstream fromMaster(UPstream::masterNo()); // worldComm
259
260 // But only consume where needed...
261 if (!haveMeshOnProc[Pstream::myProcNo()])
262 {
263 fromMaster >> fieldDicts;
264 }
265 }
266
267
268 // Use the received dictionaries to create fields
269 // (will be empty if we didn't require them)
270
271 // Disable communication when constructing from dictionary
272 const bool oldParRun = Pstream::parRun(false);
273
274 forAll(fieldDicts, i)
275 {
276 fields.set
277 (
278 i,
279 new GeoField
280 (
282 (
283 masterNames[i],
284 mesh.time().timeName(),
285 mesh.thisDb(),
288 ),
289 mesh,
290 fieldDicts[i]
291 )
292 );
293 }
294
295 Pstream::parRun(oldParRun); // Restore any changes
296
297
298 // Finally. Can checkOut of registry as required
299 if (deregister)
300 {
302 for (auto& fld : fields)
303 {
305
306 // Ensure it is not destroyed by polyMesh deletion
307 fld.checkOut();
308 }
310
311 // Extra safety - remove all such types
313 (
314 mesh.thisDb().objectRegistry::template lookupClass<GeoField>()
315 );
316
317 forAllConstIters(other, iter)
318 {
319 GeoField& fld = const_cast<GeoField&>(*iter.val());
320
321 if (!fld.ownedByRegistry())
322 {
323 fld.checkOut();
324 }
325 }
326 }
327}
328
329
330// ************************************************************************* //
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))
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:49
MESH Mesh
Definition: GeoMesh.H:62
Generic GeometricField class.
A HashTable similar to std::unordered_map.
Definition: HashTable.H:123
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:59
wordList sortedNames() const
The sorted names of the IOobjects.
Definition: IOobjectList.C:383
IOobjectList lookupClass(const char *clsName) const
The list of IOobjects with the given headerClassName.
Definition: IOobjectList.C:313
UPtrList< const IOobject > sorted() const
The sorted list of IOobjects.
Definition: IOobjectList.C:336
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
writeOption writeOpt() const noexcept
The write option.
Definition: IOobjectI.H:179
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:105
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:87
UPstream::rangeType subProcs() const noexcept
Range of sub-processes indices associated with PstreamBuffers.
static void broadcast(Type &value, const label comm=UPstream::worldComm)
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
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static constexpr int masterNo() noexcept
Process index of the master (always 0)
Definition: UPstream.H:451
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:71
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
faMesh Mesh
The mesh type.
Definition: faMesh.H:491
static void readField(const IOobject &io, const typename GeoField::Mesh &mesh, const label i, PtrList< GeoField > &fields)
Generic mesh-based field reading.
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:158
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:302
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
int myProcNo() const noexcept
Return processor number.
splitCell * master() const
Definition: splitCell.H:113
A class for managing temporary objects.
Definition: tmp.H:65
@ BEGIN_LIST
Begin list [isseparator].
Definition: token.H:155
@ END_LIST
End list [isseparator].
Definition: token.H:156
@ NL
Newline [isspace].
Definition: token.H:124
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:215
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
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 forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278