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-2014 OpenFOAM Foundation
6  \\/ M anipulation |
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>
60 void Foam::fvMeshDistribute::saveBoundaryFields
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>
89 void Foam::fvMeshDistribute::mapBoundaryFields
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  {
107  FatalErrorIn("fvMeshDistribute::mapBoundaryFields(..)") << "problem"
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 =
117  fld.boundaryField();
118 
119  const FieldField<fvsPatchField, T>& oldBfld = oldBflds[fieldI++];
120 
121  // Pull from old boundary field into bfld.
122 
123  forAll(bfld, patchI)
124  {
125  fvsPatchField<T>& patchFld = bfld[patchI];
126  label faceI = patchFld.patch().start();
127 
128  forAll(patchFld, i)
129  {
130  label oldFaceI = faceMap[faceI++];
131 
132  // Find patch and local patch face oldFaceI was in.
133  forAll(oldPatchStarts, oldPatchI)
134  {
135  label oldLocalI = oldFaceI - oldPatchStarts[oldPatchI];
136 
137  if (oldLocalI >= 0 && oldLocalI < oldBfld[oldPatchI].size())
138  {
139  patchFld[i] = oldBfld[oldPatchI][oldLocalI];
140  }
141  }
142  }
143  }
144  }
145 }
146 
147 
148 // Init patch fields of certain type
149 template<class GeoField, class PatchFieldType>
150 void Foam::fvMeshDistribute::initPatchFields
151 (
152  const typename GeoField::value_type& initVal
153 )
154 {
156  (
157  mesh_.objectRegistry::lookupClass<GeoField>()
158  );
159 
160  forAllIter(typename HashTable<GeoField*>, flds, iter)
161  {
162  GeoField& fld = *iter();
163 
164  typename GeoField::GeometricBoundaryField& bfld =
165  fld.boundaryField();
166 
167  forAll(bfld, patchI)
168  {
169  if (isA<PatchFieldType>(bfld[patchI]))
170  {
171  bfld[patchI] == initVal;
172  }
173  }
174  }
175 }
176 
177 
178 // correctBoundaryConditions patch fields of certain type
179 template<class GeoField>
180 void Foam::fvMeshDistribute::correctBoundaryConditions()
181 {
183  (
184  mesh_.objectRegistry::lookupClass<GeoField>()
185  );
186 
187  forAllIter(typename HashTable<GeoField*>, flds, iter)
188  {
189  const GeoField& fld = *iter();
190  fld.correctBoundaryConditions();
191  }
192 }
193 
194 
195 // Send fields. Note order supplied so we can receive in exactly the same order.
196 // Note that field gets written as entry in dictionary so we
197 // can construct from subdictionary.
198 // (since otherwise the reading as-a-dictionary mixes up entries from
199 // consecutive fields)
200 // The dictionary constructed is:
201 // volScalarField
202 // {
203 // p {internalField ..; boundaryField ..;}
204 // k {internalField ..; boundaryField ..;}
205 // }
206 // volVectorField
207 // {
208 // U {internalField ... }
209 // }
210 
211 // volVectorField {U {internalField ..; boundaryField ..;}}
212 //
213 template<class GeoField>
214 void Foam::fvMeshDistribute::sendFields
215 (
216  const label domain,
217  const wordList& fieldNames,
218  const fvMeshSubset& subsetter,
219  Ostream& toNbr
220 )
221 {
222  toNbr << GeoField::typeName << token::NL << token::BEGIN_BLOCK << token::NL;
223  forAll(fieldNames, i)
224  {
225  if (debug)
226  {
227  Pout<< "Subsetting field " << fieldNames[i]
228  << " for domain:" << domain << endl;
229  }
230 
231  // Send all fieldNames. This has to be exactly the same set as is
232  // being received!
233  const GeoField& fld =
234  subsetter.baseMesh().lookupObject<GeoField>(fieldNames[i]);
235 
236  tmp<GeoField> tsubfld = subsetter.interpolate(fld);
237 
238  toNbr
239  << fieldNames[i] << token::NL << token::BEGIN_BLOCK
240  << tsubfld
241  << token::NL << token::END_BLOCK << token::NL;
242  }
243  toNbr << token::END_BLOCK << token::NL;
244 }
245 
246 
247 // Opposite of sendFields
248 template<class GeoField>
249 void Foam::fvMeshDistribute::receiveFields
250 (
251  const label domain,
252  const wordList& fieldNames,
253  fvMesh& mesh,
254  PtrList<GeoField>& fields,
255  const dictionary& fieldDicts
256 )
257 {
258  if (debug)
259  {
260  Pout<< "Receiving fields " << fieldNames
261  << " from domain:" << domain << endl;
262  }
263 
264  fields.setSize(fieldNames.size());
265 
266  forAll(fieldNames, i)
267  {
268  if (debug)
269  {
270  Pout<< "Constructing field " << fieldNames[i]
271  << " from domain:" << domain << endl;
272  }
273 
274  fields.set
275  (
276  i,
277  new GeoField
278  (
279  IOobject
280  (
281  fieldNames[i],
282  mesh.time().timeName(),
283  mesh,
286  ),
287  mesh,
288  fieldDicts.subDict(fieldNames[i])
289  )
290  );
291  }
292 }
293 
294 
295 // ************************************************************************* //
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.
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:142
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Generic field type.
Definition: FieldField.H:51
bool set(const label) const
Is element set.
Definition: PtrListI.H:107
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
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
const fvPatch & patch() const
Return patch.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
An STL-conforming hash table.
Definition: HashTable.H:61
#define forAllIter(Container, container, iter)
Definition: UList.H:440
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
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
static void printFieldInfo(const fvMesh &)
Print some field info.
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:638
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
const fvMesh & baseMesh() const
Original mesh.
Definition: fvMeshSubset.H:179
Type type() const
Return the file type: FILE, DIRECTORY or UNDEFINED.
Definition: fileName.C:51
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
#define forAll(list, i)
Definition: UList.H:421
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
const labelList & faceMap() const
Old face map.
Definition: mapPolyMesh.H:406
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
errorManip< error > abort(error &err)
Definition: errorManip.H:131
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
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 ))
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
Generic GeometricField class.
error FatalError
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const labelList & oldPatchStarts() const
Return list of the old patch start labels.
Definition: mapPolyMesh.H:628
A class for managing temporary objects.
Definition: PtrList.H:118
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Post-processing mesh subset tool. Given the original mesh and the list of selected cells...
Definition: fvMeshSubset.H:71
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53