parFvFieldDistributorTemplates.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) 2015 OpenFOAM Foundation
9 Copyright (C) 2016-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
30#include "Time.H"
31#include "PtrList.H"
32#include "fvPatchFields.H"
33#include "emptyFvPatch.H"
34#include "emptyFvPatchField.H"
35#include "emptyFvsPatchField.H"
36#include "IOobjectList.H"
38#include "processorFvPatch.H"
39
42
43// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
44
45template<class Type>
48(
49 const DimensionedField<Type, volMesh>& fld
50) const
51{
52 // Create internalField by remote mapping
53
55 (
57 distMap_.cellMap()
58 );
59
60 auto tfield = tmp<DimensionedField<Type, volMesh>>::New
61 (
62 IOobject
63 (
64 fld.name(),
65 tgtMesh_.time().timeName(),
66 fld.local(),
67 tgtMesh_,
70 ),
71 tgtMesh_,
72 fld.dimensions(),
73 Field<Type>(fld, mapper)
74 );
75
76 tfield.ref().oriented() = fld.oriented();
77
78 return tfield;
79}
80
81
82template<class Type>
85(
86 const GeometricField<Type, fvPatchField, volMesh>& fld
87) const
88{
89 // Create internalField by remote mapping
91 (
93 distMap_.cellMap()
94 );
95
96 DimensionedField<Type, volMesh> internalField
97 (
98 IOobject
99 (
100 fld.name(),
101 tgtMesh_.time().timeName(),
102 fld.local(),
103 tgtMesh_,
106 ),
107 tgtMesh_,
108 fld.dimensions(),
109 Field<Type>(fld.internalField(), mapper)
110 );
111
112 internalField.oriented() = fld.oriented();
113
114
115 // Create patchFields by remote mapping
116 // Note: patchFields still on source mesh, not target mesh
117
118 PtrList<fvPatchField<Type>> oldPatchFields(fld.mesh().boundary().size());
119
120 const auto& bfld = fld.boundaryField();
121
122 forAll(bfld, patchi)
123 {
124 if (patchFaceMaps_.set(patchi))
125 {
126 // Clone local patch field
127 oldPatchFields.set(patchi, bfld[patchi].clone());
128
130 (
132 patchFaceMaps_[patchi]
133 );
134
135 // Map into local copy
136 oldPatchFields[patchi].autoMap(mapper);
137 }
138 }
139
140
141 // Clone the oldPatchFields onto the target patches. This is just to reset
142 // the reference to the patch, size and content stay the same.
143
144 PtrList<fvPatchField<Type>> newPatchFields(tgtMesh_.boundary().size());
145
146 forAll(oldPatchFields, patchi)
147 {
148 if (oldPatchFields.set(patchi))
149 {
150 const auto& pfld = oldPatchFields[patchi];
151
152 labelList dummyMap(identity(pfld.size()));
153 directFvPatchFieldMapper dummyMapper(dummyMap);
154
155 newPatchFields.set
156 (
157 patchi,
159 (
160 pfld,
161 tgtMesh_.boundary()[patchi],
163 dummyMapper
164 )
165 );
166 }
167 }
168
169 // Add some empty patches on remaining patches
170 // (... probably processor patches)
171
172 forAll(newPatchFields, patchi)
173 {
174 if (!newPatchFields.set(patchi))
175 {
176 newPatchFields.set
177 (
178 patchi,
180 (
182 tgtMesh_.boundary()[patchi],
184 )
185 );
186 }
187 }
188
189 // Return geometric field
190
191 return tmp<GeometricField<Type, fvPatchField, volMesh>>::New
192 (
193 std::move(internalField),
194 newPatchFields
195 );
196}
197
198
199template<class Type>
202(
203 const GeometricField<Type, fvsPatchField, surfaceMesh>& fld
204) const
205{
206 // Create internalField by remote mapping
208 (
210 distMap_.faceMap()
211 );
212
213
214 Field<Type> primitiveField;
215 {
216 // Create flat field of internalField + all patch fields
217 Field<Type> flatFld(fld.mesh().nFaces(), Type(Zero));
218 SubList<Type>(flatFld, fld.internalField().size())
219 = fld.internalField();
220
221 for (const fvsPatchField<Type>& fvp : fld.boundaryField())
222 {
223 SubList<Type>(flatFld, fvp.size(), fvp.patch().start()) = fvp;
224 }
225
226 // Map all faces
227 primitiveField = Field<Type>(flatFld, mapper, fld.oriented()());
228
229 // Trim to internal faces (note: could also have special mapper)
230 primitiveField.resize
231 (
232 min
233 (
234 primitiveField.size(),
235 tgtMesh_.nInternalFaces()
236 )
237 );
238 }
239
240
241 DimensionedField<Type, surfaceMesh> internalField
242 (
243 IOobject
244 (
245 fld.name(),
246 tgtMesh_.time().timeName(),
247 fld.local(),
248 tgtMesh_,
251 ),
252 tgtMesh_,
253 fld.dimensions(),
254 std::move(primitiveField)
255 );
256
257 internalField.oriented() = fld.oriented();
258
259
260 // Create patchFields by remote mapping
261 // Note: patchFields still on source mesh, not target mesh
262
263 PtrList<fvsPatchField<Type>> oldPatchFields(fld.mesh().boundary().size());
264
265 const auto& bfld = fld.boundaryField();
266
267 forAll(bfld, patchi)
268 {
269 if (patchFaceMaps_.set(patchi))
270 {
271 // Clone local patch field
272 oldPatchFields.set(patchi, bfld[patchi].clone());
273
275 (
277 patchFaceMaps_[patchi]
278 );
279
280 // Map into local copy
281 oldPatchFields[patchi].autoMap(mapper);
282 }
283 }
284
285
286 PtrList<fvsPatchField<Type>> newPatchFields(tgtMesh_.boundary().size());
287
288 // Clone the patchFields onto the base patches. This is just to reset
289 // the reference to the patch, size and content stay the same.
290 forAll(oldPatchFields, patchi)
291 {
292 if (oldPatchFields.set(patchi))
293 {
294 const fvsPatchField<Type>& pfld = oldPatchFields[patchi];
295
296 labelList dummyMap(identity(pfld.size()));
297 directFvPatchFieldMapper dummyMapper(dummyMap);
298
299 newPatchFields.set
300 (
301 patchi,
303 (
304 pfld,
305 tgtMesh_.boundary()[patchi],
307 dummyMapper
308 )
309 );
310 }
311 }
312
313 // Add some empty patches on remaining patches
314 // (... probably processor patches)
315 forAll(newPatchFields, patchi)
316 {
317 if (!newPatchFields.set(patchi))
318 {
319 newPatchFields.set
320 (
321 patchi,
323 (
325 tgtMesh_.boundary()[patchi],
327 )
328 );
329 }
330 }
331
332
333 // Return geometric field
334 return tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>::New
335 (
336 std::move(internalField),
337 newPatchFields
338 );
339}
340
341
342template<class Type>
345(
346 const IOobject& fieldObject
347) const
348{
349 // Read field
350 DimensionedField<Type, volMesh> fld
351 (
352 fieldObject,
353 srcMesh_
354 );
355
356 // Distribute
357 return distributeField(fld);
358}
359
360
361template<class Type>
364(
365 const IOobject& fieldObject
366) const
367{
368 // Read field
369 GeometricField<Type, fvPatchField, volMesh> fld
370 (
371 fieldObject,
372 srcMesh_
373 );
374
375 // Distribute
376 return distributeField(fld);
377}
378
379
380template<class Type>
383(
384 const IOobject& fieldObject
385) const
386{
387 // Read field
388 GeometricField<Type, fvsPatchField, surfaceMesh> fld
389 (
390 fieldObject,
391 srcMesh_
392 );
393
394 // Distribute
395 return distributeField(fld);
396}
397
398
399template<class Type>
401(
402 const IOobjectList& objects,
403 const wordRes& selectedFields
404) const
405{
406 typedef DimensionedField<Type, volMesh> fieldType;
407
408 // Available fields, sorted order
409 const wordList fieldNames =
410 (
411 selectedFields.empty()
412 ? objects.sortedNames<fieldType>()
413 : objects.sortedNames<fieldType>(selectedFields)
414 );
415
416 label nFields = 0;
417 for (const word& fieldName : fieldNames)
418 {
419 if ("cellDist" == fieldName)
420 {
421 // There is an odd chance this is an internal field
422 continue;
423 }
424 if (verbose_)
425 {
426 if (!nFields)
427 {
428 Info<< " Reconstructing "
429 << fieldType::typeName << "s\n" << nl;
430 }
431 Info<< " " << fieldName << nl;
432 }
433 ++nFields;
434
435 tmp<fieldType> tfld
436 (
437 distributeInternalField<Type>(*(objects[fieldName]))
438 );
439 if (isWriteProc_)
440 {
441 tfld().write();
442 }
443 }
444
445 if (nFields && verbose_) Info<< endl;
446 return nFields;
447}
448
449
450template<class Type>
452(
453 const IOobjectList& objects,
454 const wordRes& selectedFields
455) const
456{
457 typedef GeometricField<Type, fvPatchField, volMesh> fieldType;
458
459 // Available fields, sorted order
460 const wordList fieldNames =
461 (
462 selectedFields.empty()
463 ? objects.sortedNames<fieldType>()
464 : objects.sortedNames<fieldType>(selectedFields)
465 );
466
467 label nFields = 0;
468 for (const word& fieldName : fieldNames)
469 {
470 if ("cellDist" == fieldName)
471 {
472 continue;
473 }
474 if (verbose_)
475 {
476 if (!nFields)
477 {
478 Info<< " Reconstructing "
479 << fieldType::typeName << "s\n" << nl;
480 }
481 Info<< " " << fieldName << nl;
482 }
483 ++nFields;
484
485 tmp<fieldType> tfld
486 (
487 distributeVolumeField<Type>(*(objects[fieldName]))
488 );
489 if (isWriteProc_)
490 {
491 tfld().write();
492 }
493 }
494
495 if (nFields && verbose_) Info<< endl;
496 return nFields;
497}
498
499
500template<class Type>
502(
503 const IOobjectList& objects,
504 const wordRes& selectedFields
505) const
506{
507 typedef GeometricField<Type, fvsPatchField, surfaceMesh> fieldType;
508
509 // Available fields, sorted order
510 const wordList fieldNames =
511 (
512 selectedFields.empty()
513 ? objects.sortedNames<fieldType>()
514 : objects.sortedNames<fieldType>(selectedFields)
515 );
516
517 label nFields = 0;
518 for (const word& fieldName : fieldNames)
519 {
520 if (verbose_)
521 {
522 if (!nFields)
523 {
524 Info<< " Reconstructing "
525 << fieldType::typeName << "s\n" << nl;
526 }
527 Info<< " " << fieldName << nl;
528 }
529 ++nFields;
530
531 tmp<fieldType> tfld
532 (
533 distributeSurfaceField<Type>(*(objects[fieldName]))
534 );
535 if (isWriteProc_)
536 {
537 tfld().write();
538 }
539 }
540
541 if (nFields && verbose_) Info<< endl;
542 return nFields;
543}
544
545
546// ************************************************************************* //
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))
static const DimensionedField< Type, GeoMesh > & null()
Return a NullObjectRef DimensionedField.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
static const UList< label > & null()
Return a UList reference to a nullObject.
Definition: UListI.H:53
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
const mapDistribute & cellMap() const noexcept
Cell distribute map.
label distributeInternalFields(const IOobjectList &objects, const wordRes &selectedFields=wordRes()) const
Read, redistribute and write all/selected volume internal fields.
tmp< DimensionedField< Type, volMesh > > distributeField(const DimensionedField< Type, volMesh > &) const
Redistribute volume internal field.
label distributeSurfaceFields(const IOobjectList &objects, const wordRes &selectedFields=wordRes()) const
Read, reconstruct and write all/selected surface fields.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > distributeSurfaceField(const IOobject &fieldObject) const
Read and distribute surface field.
label distributeVolumeFields(const IOobjectList &objects, const wordRes &selectedFields=wordRes()) const
Read, redistribute and write all/selected volume fields.
tmp< DimensionedField< Type, volMesh > > distributeInternalField(const IOobject &fieldObject) const
Read and distribute volume internal field.
tmp< GeometricField< Type, fvPatchField, volMesh > > distributeVolumeField(const IOobject &fieldObject) const
Read and distribute volume field.
A class for managing temporary objects.
Definition: tmp.H:65
DirectFieldMapper< fvPatchFieldMapper > directFvPatchFieldMapper
A fvPatchFieldMapper with direct mapping.
DistributedFieldMapper< directFieldMapper > distributedFieldMapper
A directFieldMapper with distributed (with local or remote) quantities.
List< word > wordList
A List of words.
Definition: fileName.H:63
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
List< label > labelList
A List of labels.
Definition: List.H:66
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
DistributedFieldMapper< directFvPatchFieldMapper > distributedFvPatchFieldMapper
A directFvPatchFieldMapper with direct mapping from local or remote quantities.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
static const char *const typeName
The type name used in ensight case files.