BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / dynamicMesh / fvMeshDistribute / fvMeshDistributeTemplates.C
blobfa88841eb06f56a02ba300f8f34a60ebfa4db423
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "mapPolyMesh.H"
28 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
30 template<class GeoField>
31 void Foam::fvMeshDistribute::printFieldInfo(const fvMesh& mesh)
33     HashTable<const GeoField*> flds
34     (
35         mesh.objectRegistry::lookupClass<GeoField>()
36     );
38     forAllConstIter(typename HashTable<const GeoField*>, flds, iter)
39     {
40         const GeoField& fld = *iter();
42         Pout<< "Field:" << iter.key() << " internalsize:" << fld.size()
43             //<< " value:" << fld
44             << endl;
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     }
58 template<class GeoField>
59 void Foam::fvMeshDistribute::addPatchFields(const word& patchFieldType)
61     HashTable<const GeoField*> flds
62     (
63         mesh_.objectRegistry::lookupClass<GeoField>()
64     );
66     forAllConstIter(typename HashTable<const GeoField*>, flds, iter)
67     {
68         const GeoField& fld = *iter();
70         typename GeoField::GeometricBoundaryField& bfld =
71             const_cast<typename GeoField::GeometricBoundaryField&>
72             (
73                 fld.boundaryField()
74             );
76         label sz = bfld.size();
77         bfld.setSize(sz + 1);
78         bfld.set
79         (
80             sz,
81             GeoField::PatchFieldType::New
82             (
83                 patchFieldType,
84                 mesh_.boundary()[sz],
85                 fld.dimensionedInternalField()
86             )
87         );
88     }
92 // Delete trailing patch fields
93 template<class GeoField>
94 void Foam::fvMeshDistribute::deleteTrailingPatchFields()
96     HashTable<const GeoField*> flds
97     (
98         mesh_.objectRegistry::lookupClass<GeoField>()
99     );
101     forAllConstIter(typename HashTable<const GeoField*>, flds, iter)
102     {
103         const GeoField& fld = *iter();
105         typename GeoField::GeometricBoundaryField& bfld =
106             const_cast<typename GeoField::GeometricBoundaryField&>
107             (
108                 fld.boundaryField()
109             );
111         // Shrink patchFields
112         bfld.setSize(bfld.size() - 1);
113     }
117 // Save whole boundary field
118 template <class T, class Mesh>
119 void Foam::fvMeshDistribute::saveBoundaryFields
121     PtrList<FieldField<fvsPatchField, T> >& bflds
122 ) const
124     typedef GeometricField<T, fvsPatchField, Mesh> fldType;
126     HashTable<const fldType*> flds
127     (
128         mesh_.objectRegistry::lookupClass<fldType>()
129     );
131     bflds.setSize(flds.size());
133     label i = 0;
135     forAllConstIter(typename HashTable<const fldType*>, flds, iter)
136     {
137         const fldType& fld = *iter();
139         bflds.set(i, fld.boundaryField().clone().ptr());
141         i++;
142     }
146 // Map boundary field
147 template <class T, class Mesh>
148 void Foam::fvMeshDistribute::mapBoundaryFields
150     const mapPolyMesh& map,
151     const PtrList<FieldField<fvsPatchField, T> >& oldBflds
154     const labelList& oldPatchStarts = map.oldPatchStarts();
155     const labelList& faceMap = map.faceMap();
157     typedef GeometricField<T, fvsPatchField, Mesh> fldType;
159     HashTable<const fldType*> flds
160     (
161         mesh_.objectRegistry::lookupClass<fldType>()
162     );
164     if (flds.size() != oldBflds.size())
165     {
166         FatalErrorIn("fvMeshDistribute::mapBoundaryFields(..)") << "problem"
167             << abort(FatalError);
168     }
170     label fieldI = 0;
172     forAllConstIter(typename HashTable<const fldType*>, flds, iter)
173     {
174         const fldType& fld = *iter();
175         typename fldType::GeometricBoundaryField& bfld =
176             const_cast<typename fldType::GeometricBoundaryField&>
177             (
178                 fld.boundaryField()
179             );
182         const FieldField<fvsPatchField, T>& oldBfld = oldBflds[fieldI++];
184         // Pull from old boundary field into bfld.
186         forAll(bfld, patchI)
187         {
188             fvsPatchField<T>& patchFld = bfld[patchI];
189             label faceI = patchFld.patch().start();
191             forAll(patchFld, i)
192             {
193                 label oldFaceI = faceMap[faceI++];
195                 // Find patch and local patch face oldFaceI was in.
196                 forAll(oldPatchStarts, oldPatchI)
197                 {
198                     label oldLocalI = oldFaceI - oldPatchStarts[oldPatchI];
200                     if (oldLocalI >= 0 && oldLocalI < oldBfld[oldPatchI].size())
201                     {
202                         patchFld[i] = oldBfld[oldPatchI][oldLocalI];
203                     }
204                 }
205             }
206         }
207     }
211 // Init patch fields of certain type
212 template<class GeoField, class PatchFieldType>
213 void Foam::fvMeshDistribute::initPatchFields
215     const typename GeoField::value_type& initVal
218     HashTable<const GeoField*> flds
219     (
220         mesh_.objectRegistry::lookupClass<GeoField>()
221     );
223     forAllConstIter(typename HashTable<const GeoField*>, flds, iter)
224     {
225         const GeoField& fld = *iter();
227         typename GeoField::GeometricBoundaryField& bfld =
228             const_cast<typename GeoField::GeometricBoundaryField&>
229             (
230                 fld.boundaryField()
231             );
233         forAll(bfld, patchI)
234         {
235             if (isA<PatchFieldType>(bfld[patchI]))
236             {
237                 bfld[patchI] == initVal;
238             }
239         }
240     }
244 // Send fields. Note order supplied so we can receive in exactly the same order.
245 // Note that field gets written as entry in dictionary so we
246 // can construct from subdictionary.
247 // (since otherwise the reading as-a-dictionary mixes up entries from
248 // consecutive fields)
249 // The dictionary constructed is:
250 //  volScalarField
251 //  {
252 //      p {internalField ..; boundaryField ..;}
253 //      k {internalField ..; boundaryField ..;}
254 //  }
255 //  volVectorField
256 //  {
257 //      U {internalField ...  }
258 //  }
260 // volVectorField {U {internalField ..; boundaryField ..;}}
262 template<class GeoField>
263 void Foam::fvMeshDistribute::sendFields
265     const label domain,
266     const wordList& fieldNames,
267     const fvMeshSubset& subsetter,
268     UOPstream& toNbr
271     toNbr << GeoField::typeName << token::NL << token::BEGIN_BLOCK << token::NL;
272     forAll(fieldNames, i)
273     {
274         if (debug)
275         {
276             Pout<< "Subsetting field " << fieldNames[i]
277                 << " for domain:" << domain << endl;
278         }
280         // Send all fieldNames. This has to be exactly the same set as is
281         // being received!
282         const GeoField& fld =
283             subsetter.baseMesh().lookupObject<GeoField>(fieldNames[i]);
285         tmp<GeoField> tsubfld = subsetter.interpolate(fld);
287         toNbr
288             << fieldNames[i] << token::NL << token::BEGIN_BLOCK
289             << tsubfld
290             << token::NL << token::END_BLOCK << token::NL;
291     }
292     toNbr << token::END_BLOCK << token::NL;
296 // Opposite of sendFields
297 template<class GeoField>
298 void Foam::fvMeshDistribute::receiveFields
300     const label domain,
301     const wordList& fieldNames,
302     fvMesh& mesh,
303     PtrList<GeoField>& fields,
304     const dictionary& fieldDicts
307     if (debug)
308     {
309         Pout<< "Receiving fields " << fieldNames
310             << " from domain:" << domain << endl;
311     }
313     fields.setSize(fieldNames.size());
315     forAll(fieldNames, i)
316     {
317         if (debug)
318         {
319             Pout<< "Constructing field " << fieldNames[i]
320                 << " from domain:" << domain << endl;
321         }
323         fields.set
324         (
325             i,
326             new GeoField
327             (
328                 IOobject
329                 (
330                     fieldNames[i],
331                     mesh.time().timeName(),
332                     mesh,
333                     IOobject::NO_READ,
334                     IOobject::AUTO_WRITE
335                 ),
336                 mesh,
337                 fieldDicts.subDict(fieldNames[i])
338             )
339         );
340     }
344 // ************************************************************************* //