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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "mapPolyMesh.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class GeoField>
32 {
34  (
35  mesh.objectRegistry::lookupClass<GeoField>()
36  );
37 
38  forAllConstIter(typename HashTable<const GeoField*>, flds, iter)
39  {
40  const GeoField& fld = *iter();
41 
42  Pout<< "Field:" << iter.key() << " internalsize:" << fld.size()
43  //<< " value:" << fld
44  << endl;
45 
46  forAll(fld.boundaryField(), patchI)
47  {
48  Pout<< " " << patchI
49  << ' ' << fld.boundaryField()[patchI].patch().name()
50  << ' ' << fld.boundaryField()[patchI].type()
51  << ' ' << fld.boundaryField()[patchI].size()
52  << endl;
53  }
54  }
55 }
56 
57 
58 // Save whole boundary field
59 template<class T, class Mesh>
61 (
63 ) const
64 {
66 
68  (
69  static_cast<const fvMesh&>(mesh_).objectRegistry::lookupClass<fldType>()
70  );
71 
72  bflds.setSize(flds.size());
73 
74  label i = 0;
75 
76  forAllConstIter(typename HashTable<const fldType*>, flds, iter)
77  {
78  const fldType& fld = *iter();
79 
80  bflds.set(i, fld.boundaryField().clone().ptr());
81 
82  i++;
83  }
84 }
85 
86 
87 // Map boundary field
88 template<class T, class Mesh>
90 (
91  const mapPolyMesh& map,
92  const PtrList<FieldField<fvsPatchField, T> >& oldBflds
93 )
94 {
95  const labelList& oldPatchStarts = map.oldPatchStarts();
96  const labelList& faceMap = map.faceMap();
97 
99 
101  (
102  mesh_.objectRegistry::lookupClass<fldType>()
103  );
104 
105  if (flds.size() != oldBflds.size())
106  {
108  << abort(FatalError);
109  }
110 
111  label fieldI = 0;
112 
113  forAllIter(typename HashTable<fldType*>, flds, iter)
114  {
115  fldType& fld = *iter();
116  typename fldType::GeometricBoundaryField& bfld = fld.boundaryField();
117 
118  const FieldField<fvsPatchField, T>& oldBfld = oldBflds[fieldI++];
119 
120  // Pull from old boundary field into bfld.
121 
122  forAll(bfld, patchI)
123  {
124  fvsPatchField<T>& patchFld = bfld[patchI];
125  label faceI = patchFld.patch().start();
126 
127  forAll(patchFld, i)
128  {
129  label oldFaceI = faceMap[faceI++];
130 
131  // Find patch and local patch face oldFaceI was in.
132  forAll(oldPatchStarts, oldPatchI)
133  {
134  label oldLocalI = oldFaceI - oldPatchStarts[oldPatchI];
135 
136  if (oldLocalI >= 0 && oldLocalI < oldBfld[oldPatchI].size())
137  {
138  patchFld[i] = oldBfld[oldPatchI][oldLocalI];
139  }
140  }
141  }
142  }
143  }
144 }
145 
146 
147 template<class T>
149 (
150  PtrList<Field<T> >& iflds
151 ) const
152 {
154 
156  (
157  static_cast<const fvMesh&>(mesh_).objectRegistry::lookupClass<fldType>()
158  );
159 
160  iflds.setSize(flds.size());
161 
162  label i = 0;
163 
164  forAllConstIter(typename HashTable<const fldType*>, flds, iter)
165  {
166  const fldType& fld = *iter();
167 
168  iflds.set(i, fld.internalField().clone());
169 
170  i++;
171  }
172 }
173 
174 
175 // Set boundary values of exposed internal faces
176 template<class T>
178 (
179  const mapPolyMesh& map,
180  const PtrList<Field<T> >& oldFlds
181 )
182 {
183  const labelList& faceMap = map.faceMap();
184 
186 
188  (
189  mesh_.objectRegistry::lookupClass<fldType>()
190  );
191 
192  if (flds.size() != oldFlds.size())
193  {
195  << "problem"
196  << abort(FatalError);
197  }
198 
199 
200  label fieldI = 0;
201 
202  forAllIter(typename HashTable<fldType*>, flds, iter)
203  {
204  fldType& fld = *iter();
205  typename fldType::GeometricBoundaryField& bfld = fld.boundaryField();
206 
207  const Field<T>& oldInternal = oldFlds[fieldI++];
208 
209  // Pull from old internal field into bfld.
210 
211  forAll(bfld, patchI)
212  {
213  fvsPatchField<T>& patchFld = bfld[patchI];
214 
215  forAll(patchFld, i)
216  {
217  const label faceI = patchFld.patch().start()+i;
218 
219  label oldFaceI = faceMap[faceI];
220 
221  if (oldFaceI < oldInternal.size())
222  {
223  patchFld[i] = oldInternal[oldFaceI];
224 
225  if (map.flipFaceFlux().found(faceI))
226  {
227  patchFld[i] = flipOp()(patchFld[i]);
228  }
229  }
230  }
231  }
232  }
233 }
234 
235 
236 // Init patch fields of certain type
237 template<class GeoField, class PatchFieldType>
239 (
240  const typename GeoField::value_type& initVal
241 )
242 {
244  (
245  mesh_.objectRegistry::lookupClass<GeoField>()
246  );
247 
248  forAllIter(typename HashTable<GeoField*>, flds, iter)
249  {
250  GeoField& fld = *iter();
251 
252  typename GeoField::GeometricBoundaryField& bfld =
253  fld.boundaryField();
254 
255  forAll(bfld, patchI)
256  {
257  if (isA<PatchFieldType>(bfld[patchI]))
258  {
259  bfld[patchI] == initVal;
260  }
261  }
262  }
263 }
264 
265 
266 // correctBoundaryConditions patch fields of certain type
267 template<class GeoField>
269 {
271  (
272  mesh_.objectRegistry::lookupClass<GeoField>()
273  );
274 
275  forAllIter(typename HashTable<GeoField*>, flds, iter)
276  {
277  const GeoField& fld = *iter();
278  fld.correctBoundaryConditions();
279  }
280 }
281 
282 
283 // Send fields. Note order supplied so we can receive in exactly the same order.
284 // Note that field gets written as entry in dictionary so we
285 // can construct from subdictionary.
286 // (since otherwise the reading as-a-dictionary mixes up entries from
287 // consecutive fields)
288 // The dictionary constructed is:
289 // volScalarField
290 // {
291 // p {internalField ..; boundaryField ..;}
292 // k {internalField ..; boundaryField ..;}
293 // }
294 // volVectorField
295 // {
296 // U {internalField ... }
297 // }
298 
299 // volVectorField {U {internalField ..; boundaryField ..;}}
300 //
301 template<class GeoField>
303 (
304  const label domain,
305  const wordList& fieldNames,
306  const fvMeshSubset& subsetter,
307  Ostream& toNbr
308 )
309 {
310  toNbr << GeoField::typeName << token::NL << token::BEGIN_BLOCK << token::NL;
311  forAll(fieldNames, i)
312  {
313  if (debug)
314  {
315  Pout<< "Subsetting field " << fieldNames[i]
316  << " for domain:" << domain << endl;
317  }
318 
319  // Send all fieldNames. This has to be exactly the same set as is
320  // being received!
321  const GeoField& fld =
322  subsetter.baseMesh().lookupObject<GeoField>(fieldNames[i]);
323 
324  tmp<GeoField> tsubfld = subsetter.interpolate(fld);
325 
326  toNbr
328  << tsubfld
330  }
331  toNbr << token::END_BLOCK << token::NL;
332 }
333 
334 
335 // Opposite of sendFields
336 template<class GeoField>
338 (
339  const label domain,
340  const wordList& fieldNames,
341  fvMesh& mesh,
343  const dictionary& fieldDicts
344 )
345 {
346  if (debug)
347  {
348  Pout<< "Receiving fields " << fieldNames
349  << " from domain:" << domain << endl;
350  }
351 
352  fields.setSize(fieldNames.size());
353 
354  forAll(fieldNames, i)
355  {
356  if (debug)
357  {
358  Pout<< "Constructing field " << fieldNames[i]
359  << " from domain:" << domain << endl;
360  }
361 
362  fields.set
363  (
364  i,
365  new GeoField
366  (
367  IOobject
368  (
369  fieldNames[i],
370  mesh.time().timeName(),
371  mesh,
374  ),
375  mesh,
376  fieldDicts.subDict(fieldNames[i])
377  )
378  );
379  }
380 }
381 
382 
383 // ************************************************************************* //
Foam::fvMeshSubset::interpolate
static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &cellMap, const labelList &faceMap)
Map volume field.
Definition: fvMeshSubsetInterpolate.C:43
Foam::mapPolyMesh::flipFaceFlux
const labelHashSet & flipFaceFlux() const
Map of flipped face flux faces.
Definition: mapPolyMesh.H:558
Foam::fvMeshDistribute::receiveFields
static void receiveFields(const label domain, const wordList &fieldNames, fvMesh &, PtrList< GeoField > &, const dictionary &fieldDicts)
Receive fields. Opposite of sendFields.
Definition: fvMeshDistributeTemplates.C:338
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeFast.C:90
Foam::fvMeshDistribute::saveBoundaryFields
void saveBoundaryFields(PtrList< FieldField< fvsPatchField, T > > &bflds) const
Save boundary fields.
Definition: fvMeshDistributeTemplates.C:61
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
Foam::FieldField
Generic field type.
Definition: FieldField.H:51
Foam::fvMeshDistribute::mapBoundaryFields
void mapBoundaryFields(const mapPolyMesh &map, const PtrList< FieldField< fvsPatchField, T > > &oldBflds)
Map boundary fields.
Definition: fvMeshDistributeTemplates.C:90
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::fvMeshSubset
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:73
fields
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar("dpdt", p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);p_rgh=p - rho *gh;mesh.setFluxRequired(p_rgh.name());multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:127
Foam::fvMeshDistribute::printFieldInfo
static void printFieldInfo(const fvMesh &)
Print some field info.
Definition: fvMeshDistributeTemplates.C:31
mapPolyMesh.H
Foam::flipOp
Class containing functor to negate primitives. Dummy for all other types.
Definition: flipOp.H:50
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::token::NL
@ NL
Definition: token.H:97
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::fvMeshDistribute::sendFields
static void sendFields(const label domain, const wordList &fieldNames, const fvMeshSubset &, Ostream &toNbr)
Send subset of fields.
Definition: fvMeshDistributeTemplates.C:303
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field< T >
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::HashTable::size
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::HashTable::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::fvMeshDistribute::saveInternalFields
void saveInternalFields(PtrList< Field< T > > &iflds) const
Save internal fields of surfaceFields.
Definition: fvMeshDistributeTemplates.C:149
Foam::HashTable
An STL-conforming hash table.
Definition: HashTable.H:61
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::fvsPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvsPatchField.H:278
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::fvMeshDistribute::correctBoundaryConditions
void correctBoundaryConditions()
Call correctBoundaryConditions on fields.
Definition: fvMeshDistributeTemplates.C:268
Foam::token::END_BLOCK
@ END_BLOCK
Definition: token.H:105
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::fieldNames
wordList fieldNames(const IOobjectList &objects, const bool syncPar)
Get sorted names of fields of type. If syncPar and running in parallel.
Definition: ReadFields.C:36
Foam::mapPolyMesh::oldPatchStarts
const labelList & oldPatchStarts() const
Return list of the old patch start labels.
Definition: mapPolyMesh.H:628
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::fvMeshDistribute::initPatchFields
void initPatchFields(const typename GeoField::value_type &initVal)
Init patch fields of certain type.
Definition: fvMeshDistributeTemplates.C:239
Foam::fvPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:155
Foam::mapPolyMesh::faceMap
const labelList & faceMap() const
Old face map.
Definition: mapPolyMesh.H:406
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::token::BEGIN_BLOCK
@ BEGIN_BLOCK
Definition: token.H:104
Foam::fvMeshDistribute::mapExposedFaces
void mapExposedFaces(const mapPolyMesh &map, const PtrList< Field< T > > &oldFlds)
Set value of patch faces resulting from internal faces.
Definition: fvMeshDistributeTemplates.C:178
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::fvMeshSubset::baseMesh
const fvMesh & baseMesh() const
Original mesh.
Definition: fvMeshSubset.H:216
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52