All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
foamToTetDualMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 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(),
89  Zero
90  )
91  )
92  );
93 
94  Field<Type>& fld = tetFields[i].primitiveFieldRef();
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 "createMeshNoChangers.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  false
175  );
176  // From tet vertices to poly cells/faces
177  const labelIOList pointDualAddressing
178  (
179  IOobject
180  (
181  "pointDualAddressing",
182  tetDualMesh.facesInstance(),
183  tetDualMesh.meshSubDir,
184  tetDualMesh,
186  )
187  );
188 
189  if (pointDualAddressing.size() != tetDualMesh.nPoints())
190  {
192  << "Size " << pointDualAddressing.size()
193  << " of addressing map "
194  << pointDualAddressing.relativeObjectPath()
195  << " differs from number of points in mesh "
196  << tetDualMesh.nPoints()
197  << exit(FatalError);
198  }
199 
200 
201  // Some stats on addressing
202  label nCells = 0;
203  label nPatchFaces = 0;
204  label nUnmapped = 0;
205  forAll(pointDualAddressing, pointi)
206  {
207  label index = pointDualAddressing[pointi];
208 
209  if (index > 0)
210  {
211  nCells++;
212  }
213  else if (index == 0)
214  {
215  nUnmapped++;
216  }
217  else
218  {
219  label facei = -index-1;
220  if (facei < mesh.nInternalFaces())
221  {
223  << "Face " << facei << " from index " << index
224  << " is not a boundary face."
225  << " nInternalFaces:" << mesh.nInternalFaces()
226  << exit(FatalError);
227  }
228  else
229  {
230  nPatchFaces++;
231  }
232  }
233  }
234 
235  reduce(nCells, sumOp<label>());
236  reduce(nPatchFaces, sumOp<label>());
237  reduce(nUnmapped, sumOp<label>());
238  Info<< "tetDualMesh points : " << tetDualMesh.nPoints()
239  << " of which mapped to" << nl
240  << " cells : " << nCells << nl
241  << " patch faces : " << nPatchFaces << nl
242  << " not mapped : " << nUnmapped << nl
243  << endl;
244 
245 
246  // Read objects in time directory
247  IOobjectList objects(mesh, runTime.timeName());
248 
249  // Read vol fields, interpolate onto tet points
251  ReadAndMapFields<volScalarField, pointScalarField>
252  (
253  mesh,
254  objects,
255  tetDualMesh,
256  pointDualAddressing,
257  Zero, // nullValue
258  psFlds
259  );
260 
262  ReadAndMapFields<volVectorField, pointVectorField>
263  (
264  mesh,
265  objects,
266  tetDualMesh,
267  pointDualAddressing,
268  Zero, // nullValue
269  pvFlds
270  );
271 
273  ReadAndMapFields<volSphericalTensorField, pointSphericalTensorField>
274  (
275  mesh,
276  objects,
277  tetDualMesh,
278  pointDualAddressing,
279  Zero, // nullValue
280  pstFlds
281  );
282 
284  ReadAndMapFields<volSymmTensorField, pointSymmTensorField>
285  (
286  mesh,
287  objects,
288  tetDualMesh,
289  pointDualAddressing,
290  Zero, // nullValue
291  psymmtFlds
292  );
293 
295  ReadAndMapFields<volTensorField, pointTensorField>
296  (
297  mesh,
298  objects,
299  tetDualMesh,
300  pointDualAddressing,
301  Zero, // nullValue
302  ptFlds
303  );
304 
305  tetDualMesh.objectRegistry::write();
306 
307  Info<< "End\n" << endl;
308 
309  return 0;
310 }
311 
312 
313 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:888
if(fields) ReadFields(pointMesh PtrList< pointTensorField > ptFlds
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
label nInternalFaces() const
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:328
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const labelList & patchID() const
Per boundary face label the patch index.
Generic dimensioned Type class.
fvMesh & mesh
Pre-declare SubField and related Field type.
Definition: Field.H:56
static const zero Zero
Definition: zero.H:97
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
static const char nl
Definition: Ostream.H:260
objects
if(fields) ReadFields(pointMesh PtrList< pointVectorField > pvFlds
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
IOobjectList lookupClass(const word &className) const
Return the list for all IOobjects of a given class.
Definition: IOobjectList.C:190
void setSize(const label)
Reset size of List.
Definition: List.C:281
label patchi
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
messageStream Info
PtrList< pointScalarField > psFlds
label nPoints() const
if(fields) ReadFields(pointMesh PtrList< pointSphericalTensorField > pstFlds
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
Namespace for OpenFOAM.