foamToTetDualMesh.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 |
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 Application
25  foamToTetDualMesh
26 
27 Description
28  Converts polyMesh results to tetDualMesh.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "argList.H"
33 #include "fvMesh.H"
34 #include "volFields.H"
35 #include "pointFields.H"
36 #include "Time.H"
37 #include "IOobjectList.H"
38 
39 using namespace Foam;
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 template<class ReadGeoField, class MappedGeoField>
44 void ReadAndMapFields
45 (
46  const fvMesh& mesh,
47  const IOobjectList& objects,
48  const fvMesh& tetDualMesh,
49  const labelList& map,
50  const typename MappedGeoField::value_type& nullValue,
51  PtrList<MappedGeoField>& tetFields
52 )
53 {
54  typedef typename MappedGeoField::value_type Type;
55 
56  // Search list of objects for wanted type
57  IOobjectList fieldObjects(objects.lookupClass(ReadGeoField::typeName));
58 
59  tetFields.setSize(fieldObjects.size());
60 
61  label i = 0;
62  forAllConstIter(IOobjectList, fieldObjects, iter)
63  {
64  Info<< "Converting " << ReadGeoField::typeName << ' ' << iter.key()
65  << endl;
66 
67  ReadGeoField readField(*iter(), mesh);
68 
69  tetFields.set
70  (
71  i,
72  new MappedGeoField
73  (
74  IOobject
75  (
76  readField.name(),
77  readField.instance(),
78  readField.local(),
79  tetDualMesh,
82  readField.registerObject()
83  ),
84  pointMesh::New(tetDualMesh),
86  (
87  "zero",
88  readField.dimensions(),
90  )
91  )
92  );
93 
94  Field<Type>& fld = tetFields[i].internalField();
95 
96  // Map from read field. Set unmapped entries to nullValue.
97  fld.setSize(map.size(), nullValue);
98  forAll(map, pointI)
99  {
100  label index = map[pointI];
101 
102  if (index > 0)
103  {
104  label cellI = index-1;
105  fld[pointI] = readField[cellI];
106  }
107  else if (index < 0)
108  {
109  label faceI = -index-1;
110  label bFaceI = faceI - mesh.nInternalFaces();
111  if (bFaceI >= 0)
112  {
113  label patchI = mesh.boundaryMesh().patchID()[bFaceI];
114  label localFaceI = mesh.boundaryMesh()[patchI].whichFace
115  (
116  faceI
117  );
118  fld[pointI] = readField.boundaryField()[patchI][localFaceI];
119  }
120  //else
121  //{
122  // FatalErrorInFunction
123  // << "Face " << faceI << " from index " << index
124  // << " is not a boundary face." << abort(FatalError);
125  //}
126 
127  }
128  //else
129  //{
130  // WarningInFunction
131  // << "Point " << pointI << " at "
132  // << tetDualMesh.points()[pointI]
133  // << " has no dual correspondence." << endl;
134  //}
135  }
136  tetFields[i].correctBoundaryConditions();
137 
138  i++;
139  }
140 }
141 
142 
143 
144 
145 int main(int argc, char *argv[])
146 {
147  #include "addOverwriteOption.H"
148  #include "addTimeOptions.H"
149 
150  #include "setRootCase.H"
151  #include "createTime.H"
152  // Get times list
153  instantList Times = runTime.times();
154  #include "checkTimeOptions.H"
155  runTime.setTime(Times[startTime], startTime);
156 
157 
158  // Read the mesh
159  #include "createMesh.H"
160 
161  // Read the tetDualMesh
162  Info<< "Create tetDualMesh for time = "
163  << runTime.timeName() << nl << endl;
164 
165  fvMesh tetDualMesh
166  (
167  IOobject
168  (
169  "tetDualMesh",
170  runTime.timeName(),
171  runTime,
173  )
174  );
175  // From tet vertices to poly cells/faces
176  const labelIOList pointDualAddressing
177  (
178  IOobject
179  (
180  "pointDualAddressing",
181  tetDualMesh.facesInstance(),
182  tetDualMesh.meshSubDir,
183  tetDualMesh,
185  )
186  );
187 
188  if (pointDualAddressing.size() != tetDualMesh.nPoints())
189  {
191  << "Size " << pointDualAddressing.size()
192  << " of addressing map " << pointDualAddressing.objectPath()
193  << " differs from number of points in mesh "
194  << tetDualMesh.nPoints()
195  << exit(FatalError);
196  }
197 
198 
199  // Some stats on addressing
200  label nCells = 0;
201  label nPatchFaces = 0;
202  label nUnmapped = 0;
203  forAll(pointDualAddressing, pointI)
204  {
205  label index = pointDualAddressing[pointI];
206 
207  if (index > 0)
208  {
209  nCells++;
210  }
211  else if (index == 0)
212  {
213  nUnmapped++;
214  }
215  else
216  {
217  label faceI = -index-1;
218  if (faceI < mesh.nInternalFaces())
219  {
221  << "Face " << faceI << " from index " << index
222  << " is not a boundary face."
223  << " nInternalFaces:" << mesh.nInternalFaces()
224  << exit(FatalError);
225  }
226  else
227  {
228  nPatchFaces++;
229  }
230  }
231  }
232 
233  reduce(nCells, sumOp<label>());
234  reduce(nPatchFaces, sumOp<label>());
235  reduce(nUnmapped, sumOp<label>());
236  Info<< "tetDualMesh points : " << tetDualMesh.nPoints()
237  << " of which mapped to" << nl
238  << " cells : " << nCells << nl
239  << " patch faces : " << nPatchFaces << nl
240  << " not mapped : " << nUnmapped << nl
241  << endl;
242 
243 
244  // Read objects in time directory
245  IOobjectList objects(mesh, runTime.timeName());
246 
247  // Read vol fields, interpolate onto tet points
249  ReadAndMapFields<volScalarField, pointScalarField>
250  (
251  mesh,
252  objects,
253  tetDualMesh,
254  pointDualAddressing,
255  pTraits<scalar>::zero, // nullValue
256  psFlds
257  );
258 
260  ReadAndMapFields<volVectorField, pointVectorField>
261  (
262  mesh,
263  objects,
264  tetDualMesh,
265  pointDualAddressing,
266  pTraits<vector>::zero, // nullValue
267  pvFlds
268  );
269 
271  ReadAndMapFields<volSphericalTensorField, pointSphericalTensorField>
272  (
273  mesh,
274  objects,
275  tetDualMesh,
276  pointDualAddressing,
277  pTraits<sphericalTensor>::zero, // nullValue
278  pstFlds
279  );
280 
282  ReadAndMapFields<volSymmTensorField, pointSymmTensorField>
283  (
284  mesh,
285  objects,
286  tetDualMesh,
287  pointDualAddressing,
288  pTraits<symmTensor>::zero, // nullValue
289  psymmtFlds
290  );
291 
293  ReadAndMapFields<volTensorField, pointTensorField>
294  (
295  mesh,
296  objects,
297  tetDualMesh,
298  pointDualAddressing,
299  pTraits<tensor>::zero, // nullValue
300  ptFlds
301  );
302 
303  tetDualMesh.objectRegistry::write();
304 
305  Info<< "End\n" << endl;
306 
307  return 0;
308 }
309 
310 
311 // ************************************************************************* //
volFields.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
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
readField
void readField(const IOobject &io, const fvMesh &mesh, const label i, PtrList< GeoField > &fields)
Definition: redistributePar.C:557
addOverwriteOption.H
Foam::polyMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:309
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:778
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
IOobjectList.H
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
addTimeOptions.H
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::polyBoundaryMesh::patchID
const labelList & patchID() const
Per boundary face label the patch index.
Definition: polyBoundaryMesh.C:399
Foam::primitiveMesh::nPoints
label nPoints() const
Definition: primitiveMeshI.H:35
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::PtrList::set
bool set(const label) const
Is element set.
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< Type >
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
argList.H
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
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::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
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::dimensioned< Type >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::IOobjectList::lookupClass
IOobjectList lookupClass(const word &className) const
Return the list for all IOobjects of a given class.
Definition: IOobjectList.C:199
checkTimeOptions.H
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
setRootCase.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sumOp
Definition: ops.H:162
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::pTraits
Traits class for primitives.
Definition: pTraits.H:50
createMesh.H
Foam::IOList< label >
Foam::PtrList::setSize
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:142
createTime.H
startTime
Foam::label startTime
Definition: checkTimeOptions.H:5
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
pointFields.H