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-2020 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 GeoField>
35 {
37  (
38  mesh.objectRegistry::lookupClass<GeoField>()
39  );
40 
41  forAllConstIters(flds, iter)
42  {
43  const GeoField& fld = *iter();
44 
45  Pout<< "Field:" << iter.key() << " internalsize:" << fld.size()
46  //<< " value:" << fld
47  << endl;
48 
49  for (const auto& patchFld : fld.boundaryField())
50  {
51  Pout<< " " << patchFld.patch().index()
52  << ' ' << patchFld.patch().name()
53  << ' ' << patchFld.type()
54  << ' ' << patchFld.size()
55  << nl;
56  }
57  }
58 }
59 
60 
61 template<class T, class Mesh>
62 void Foam::fvMeshDistribute::saveBoundaryFields
63 (
65 ) const
66 {
67  // Save whole boundary field
68 
70 
72  (
73  mesh_.objectRegistry::lookupClass<const fldType>()
74  );
75 
76  bflds.setSize(flds.size());
77 
78  label i = 0;
79  forAllConstIters(flds, iter)
80  {
81  const fldType& fld = *iter();
82 
83  bflds.set(i, fld.boundaryField().clone().ptr());
84 
85  ++i;
86  }
87 }
88 
89 
90 template<class T, class Mesh>
91 void Foam::fvMeshDistribute::mapBoundaryFields
92 (
93  const mapPolyMesh& map,
94  const PtrList<FieldField<fvsPatchField, T>>& oldBflds
95 )
96 {
97  // Map boundary field
98 
99  const labelList& oldPatchStarts = map.oldPatchStarts();
100  const labelList& faceMap = map.faceMap();
101 
102  typedef GeometricField<T, fvsPatchField, Mesh> fldType;
103 
104  HashTable<fldType*> flds
105  (
106  mesh_.objectRegistry::lookupClass<fldType>()
107  );
108 
109  if (flds.size() != oldBflds.size())
110  {
112  << abort(FatalError);
113  }
114 
115  label fieldi = 0;
116 
117  forAllIters(flds, iter)
118  {
119  fldType& fld = *iter();
120  auto& bfld = fld.boundaryFieldRef();
121 
122  const FieldField<fvsPatchField, T>& oldBfld = oldBflds[fieldi++];
123 
124  // Pull from old boundary field into bfld.
125 
126  forAll(bfld, patchi)
127  {
128  fvsPatchField<T>& patchFld = bfld[patchi];
129  label facei = patchFld.patch().start();
130 
131  forAll(patchFld, i)
132  {
133  label oldFacei = faceMap[facei++];
134 
135  // Find patch and local patch face oldFacei was in.
136  forAll(oldPatchStarts, oldPatchi)
137  {
138  label oldLocalI = oldFacei - oldPatchStarts[oldPatchi];
139 
140  if (oldLocalI >= 0 && oldLocalI < oldBfld[oldPatchi].size())
141  {
142  patchFld[i] = oldBfld[oldPatchi][oldLocalI];
143  }
144  }
145  }
146  }
147  }
148 }
149 
150 
151 template<class T>
152 void Foam::fvMeshDistribute::saveInternalFields
153 (
154  PtrList<Field<T>>& iflds
155 ) const
156 {
157  typedef GeometricField<T, fvsPatchField, surfaceMesh> fldType;
158 
159  HashTable<const fldType*> flds
160  (
161  mesh_.objectRegistry::lookupClass<const fldType>()
162  );
163 
164  iflds.setSize(flds.size());
165 
166  label i = 0;
167 
168  forAllConstIters(flds, iter)
169  {
170  const fldType& fld = *iter();
171 
172  iflds.set(i, fld.primitiveField().clone());
173 
174  ++i;
175  }
176 }
177 
178 
179 template<class T>
180 void Foam::fvMeshDistribute::mapExposedFaces
181 (
182  const mapPolyMesh& map,
183  const PtrList<Field<T>>& oldFlds
184 )
185 {
186  // Set boundary values of exposed internal faces
187 
188  const labelList& faceMap = map.faceMap();
189 
190  typedef GeometricField<T, fvsPatchField, surfaceMesh> fldType;
191 
192  HashTable<fldType*> flds
193  (
194  mesh_.objectRegistry::lookupClass<fldType>()
195  );
196 
197  if (flds.size() != oldFlds.size())
198  {
200  << "problem"
201  << abort(FatalError);
202  }
203 
204 
205  label fieldI = 0;
206 
207  forAllIters(flds, iter)
208  {
209  fldType& fld = *iter();
210  const bool oriented = fld.oriented()();
211 
212  typename fldType::Boundary& bfld = fld.boundaryFieldRef();
213 
214  const Field<T>& oldInternal = oldFlds[fieldI++];
215 
216  // Pull from old internal field into bfld.
217 
218  forAll(bfld, patchi)
219  {
220  fvsPatchField<T>& patchFld = bfld[patchi];
221 
222  forAll(patchFld, i)
223  {
224  const label faceI = patchFld.patch().start()+i;
225 
226  label oldFaceI = faceMap[faceI];
227 
228  if (oldFaceI < oldInternal.size())
229  {
230  patchFld[i] = oldInternal[oldFaceI];
231 
232  if (oriented && map.flipFaceFlux().found(faceI))
233  {
234  patchFld[i] = flipOp()(patchFld[i]);
235  }
236  }
237  }
238  }
239  }
240 }
241 
242 
243 template<class GeoField, class PatchFieldType>
244 void Foam::fvMeshDistribute::initPatchFields
245 (
246  const typename GeoField::value_type& initVal
247 )
248 {
249  // Init patch fields of certain type
250 
251  HashTable<GeoField*> flds
252  (
253  mesh_.objectRegistry::lookupClass<GeoField>()
254  );
255 
256  forAllIters(flds, iter)
257  {
258  GeoField& fld = *iter();
259 
260  auto& bfld = fld.boundaryFieldRef();
261 
262  forAll(bfld, patchi)
263  {
264  if (isA<PatchFieldType>(bfld[patchi]))
265  {
266  bfld[patchi] == initVal;
267  }
268  }
269  }
270 }
271 
272 
273 template<class GeoField>
274 void Foam::fvMeshDistribute::correctBoundaryConditions()
275 {
276  // CorrectBoundaryConditions patch fields of certain type
277 
278  HashTable<GeoField*> flds
279  (
280  mesh_.objectRegistry::lookupClass<GeoField>()
281  );
282 
283  forAllIters(flds, iter)
284  {
285  GeoField& fld = *iter();
286  fld.correctBoundaryConditions();
287  }
288 }
289 
290 
291 template<class GeoField>
292 void Foam::fvMeshDistribute::getFieldNames
293 (
294  const fvMesh& mesh,
295  HashTable<wordList>& allFieldNames,
296  const bool syncPar
297 )
298 {
299  wordList& list = allFieldNames(GeoField::typeName);
300  list = mesh.sortedNames<GeoField>();
301 
302  // Check all procs have same names
303  if (syncPar)
304  {
305  List<wordList> allNames(Pstream::nProcs());
306  allNames[Pstream::myProcNo()] = list;
307  Pstream::gatherList(allNames);
308  Pstream::scatterList(allNames);
309 
310  for (const int proci : Pstream::subProcs())
311  {
312  if (allNames[proci] != allNames[0])
313  {
315  << "When checking for equal "
316  << GeoField::typeName
317  << " :" << nl
318  << "processor0 has:" << allNames[0] << endl
319  << "processor" << proci << " has:" << allNames[proci] << nl
320  << GeoField::typeName
321  << " need to be synchronised on all processors."
322  << exit(FatalError);
323  }
324  }
325  }
326 }
327 
328 
329 template<class GeoField>
330 void Foam::fvMeshDistribute::sendFields
331 (
332  const label domain,
333  const HashTable<wordList>& allFieldNames,
334  const fvMeshSubset& subsetter,
335  Ostream& toNbr
336 )
337 {
338  // Send fields. Note order supplied so we can receive in exactly the same
339  // order.
340  // Note that field gets written as entry in dictionary so we
341  // can construct from subdictionary.
342  // (since otherwise the reading as-a-dictionary mixes up entries from
343  // consecutive fields)
344  // The dictionary constructed is:
345  // volScalarField
346  // {
347  // p {internalField ..; boundaryField ..;}
348  // k {internalField ..; boundaryField ..;}
349  // }
350  // volVectorField
351  // {
352  // U {internalField ... }
353  // }
354 
355  // volVectorField {U {internalField ..; boundaryField ..;}}
356 
357  const wordList& fieldNames =
358  allFieldNames.lookup(GeoField::typeName, wordList::null());
359 
360  toNbr << GeoField::typeName << token::NL << token::BEGIN_BLOCK << token::NL;
361 
362  for (const word& fieldName : fieldNames)
363  {
364  if (debug)
365  {
366  Pout<< "Subsetting field " << fieldName
367  << " for domain:" << domain << endl;
368  }
369 
370  // Send all fieldNames. This has to be exactly the same set as is
371  // being received!
372  const GeoField& fld =
373  subsetter.baseMesh().lookupObject<GeoField>(fieldName);
374 
375  // Note: use subsetter to get sub field. Override default behaviour
376  // to warn for unset fields since they will be reset later on
377  tmp<GeoField> tsubfld = subsetter.interpolate(fld, true);
378 
379  toNbr
380  << fieldName << token::NL << token::BEGIN_BLOCK
381  << tsubfld
383  }
384  toNbr << token::END_BLOCK << token::NL;
385 }
386 
387 
388 template<class GeoField>
389 void Foam::fvMeshDistribute::receiveFields
390 (
391  const label domain,
392  const HashTable<wordList>& allFieldNames,
393  fvMesh& mesh,
394  PtrList<GeoField>& fields,
395  const dictionary& allFieldsDict
396 )
397 {
398  // Opposite of sendFields
399 
400  const wordList& fieldNames =
401  allFieldNames.lookup(GeoField::typeName, wordList::null());
402 
403  const dictionary& fieldDicts =
404  allFieldsDict.subDict(GeoField::typeName);
405 
406 
407  if (debug)
408  {
409  Pout<< "Receiving fields " << fieldNames
410  << " from domain:" << domain << endl;
411  }
412 
413  fields.resize(fieldNames.size());
414 
415  label fieldi = 0;
416  for (const word& fieldName : fieldNames)
417  {
418  if (debug)
419  {
420  Pout<< "Constructing field " << fieldName
421  << " from domain:" << domain << endl;
422  }
423 
424  fields.set
425  (
426  fieldi++,
427  new GeoField
428  (
429  IOobject
430  (
431  fieldName,
432  mesh.time().timeName(),
433  mesh,
436  ),
437  mesh,
438  fieldDicts.subDict(fieldName)
439  )
440  );
441  }
442 }
443 
444 
445 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
Foam::objectRegistry::sortedNames
wordList sortedNames() const
The sorted names of all objects.
Definition: objectRegistry.C:170
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::List::null
static const List< T > & null()
Return a null List.
Definition: ListI.H:108
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeTopological.C:94
Foam::FieldField
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:53
Foam::Pstream::scatterList
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
Definition: gatherScatterList.C:215
Foam::fvMeshDistribute::printFieldInfo
static void printFieldInfo(const fvMesh &)
Print some field info.
Definition: fvMeshDistributeTemplates.C:34
mapPolyMesh.H
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:59
Foam::UPstream::subProcs
static rangeType subProcs(const label communicator=worldComm)
Range of process indices for sub-processes.
Definition: UPstream.H:516
Foam::OSstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:107
Foam::fileName::type
Type type(bool followLink=true, bool checkGzip=false) const
Definition: fileName.C:268
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:528
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:62
fld
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))
Definition: gmvOutputLagrangian.H:23
Foam::token::END_BLOCK
End block [isseparator].
Definition: token.H:127
forAllIters
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:223
Foam::token::NL
Newline [isspace].
Definition: token.H:119
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::token::BEGIN_BLOCK
Begin block [isseparator].
Definition: token.H:126
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::HashTable
A HashTable similar to std::unordered_map.
Definition: HashTable.H:105
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:52
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:464
Foam::nl
constexpr char nl
Definition: Ostream.H:385
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:275
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53
Foam::IOobject::NO_READ
Definition: IOobject.H:123
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:446