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-2023 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 "timeSelector.H"
34 #include "fvMesh.H"
35 #include "volFields.H"
36 #include "pointFields.H"
37 #include "Time.H"
38 #include "IOobjectList.H"
39 
40 using namespace Foam;
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 template<class ReadGeoField, class MappedGeoField>
45 void ReadAndMapFields
46 (
47  const fvMesh& mesh,
48  const IOobjectList& objects,
49  const fvMesh& tetDualMesh,
50  const labelList& map,
51  const typename MappedGeoField::value_type& nullValue,
52  PtrList<MappedGeoField>& tetFields
53 )
54 {
55  typedef typename MappedGeoField::value_type Type;
56 
57  // Search list of objects for wanted type
58  IOobjectList fieldObjects(objects.lookupClass(ReadGeoField::typeName));
59 
60  tetFields.setSize(fieldObjects.size());
61 
62  label i = 0;
63  forAllConstIter(IOobjectList, fieldObjects, iter)
64  {
65  Info<< "Converting " << ReadGeoField::typeName << ' ' << iter.key()
66  << endl;
67 
68  ReadGeoField readField(*iter(), mesh);
69 
70  tetFields.set
71  (
72  i,
73  new MappedGeoField
74  (
75  IOobject
76  (
77  readField.name(),
78  readField.instance(),
79  readField.local(),
80  tetDualMesh,
83  readField.registerObject()
84  ),
85  pointMesh::New(tetDualMesh),
87  (
88  "zero",
89  readField.dimensions(),
90  Zero
91  )
92  )
93  );
94 
95  Field<Type>& fld = tetFields[i].primitiveFieldRef();
96 
97  // Map from read field. Set unmapped entries to nullValue.
98  fld.setSize(map.size(), nullValue);
99  forAll(map, pointi)
100  {
101  label index = map[pointi];
102 
103  if (index > 0)
104  {
105  label celli = index-1;
106  fld[pointi] = readField[celli];
107  }
108  else if (index < 0)
109  {
110  label facei = -index-1;
111  label bFacei = facei - mesh.nInternalFaces();
112  if (bFacei >= 0)
113  {
114  label patchi = mesh.boundaryMesh().patchID()[bFacei];
115  label localFacei = mesh.boundaryMesh()[patchi].whichFace
116  (
117  facei
118  );
119  fld[pointi] = readField.boundaryField()[patchi][localFacei];
120  }
121  // else
122  //{
123  // FatalErrorInFunction
124  // << "Face " << facei << " from index " << index
125  // << " is not a boundary face." << abort(FatalError);
126  //}
127 
128  }
129  // else
130  //{
131  // WarningInFunction
132  // << "Point " << pointi << " at "
133  // << tetDualMesh.points()[pointi]
134  // << " has no dual correspondence." << endl;
135  //}
136  }
137  tetFields[i].correctBoundaryConditions();
138 
139  i++;
140  }
141 }
142 
143 
144 
145 
146 int main(int argc, char *argv[])
147 {
149  #include "addOverwriteOption.H"
150 
151  #include "setRootCase.H"
152  #include "createTime.H"
153  timeSelector::select0(runTime, args);
154  #include "createMeshNoChangers.H"
155 
156  // Read the tetDualMesh
157  Info<< "Create tetDualMesh for time = "
158  << runTime.name() << nl << endl;
159 
160  fvMesh tetDualMesh
161  (
162  IOobject
163  (
164  "tetDualMesh",
165  runTime.name(),
166  runTime,
168  ),
169  false
170  );
171  // From tet vertices to poly cells/faces
172  const labelIOList pointDualAddressing
173  (
174  IOobject
175  (
176  "pointDualAddressing",
177  tetDualMesh.facesInstance(),
178  tetDualMesh.meshSubDir,
179  tetDualMesh,
181  )
182  );
183 
184  if (pointDualAddressing.size() != tetDualMesh.nPoints())
185  {
187  << "Size " << pointDualAddressing.size()
188  << " of addressing map "
189  << pointDualAddressing.relativeObjectPath()
190  << " differs from number of points in mesh "
191  << tetDualMesh.nPoints()
192  << exit(FatalError);
193  }
194 
195 
196  // Some stats on addressing
197  label nCells = 0;
198  label nPatchFaces = 0;
199  label nUnmapped = 0;
200  forAll(pointDualAddressing, pointi)
201  {
202  label index = pointDualAddressing[pointi];
203 
204  if (index > 0)
205  {
206  nCells++;
207  }
208  else if (index == 0)
209  {
210  nUnmapped++;
211  }
212  else
213  {
214  label facei = -index-1;
215  if (facei < mesh.nInternalFaces())
216  {
218  << "Face " << facei << " from index " << index
219  << " is not a boundary face."
220  << " nInternalFaces:" << mesh.nInternalFaces()
221  << exit(FatalError);
222  }
223  else
224  {
225  nPatchFaces++;
226  }
227  }
228  }
229 
230  reduce(nCells, sumOp<label>());
231  reduce(nPatchFaces, sumOp<label>());
232  reduce(nUnmapped, sumOp<label>());
233  Info<< "tetDualMesh points : " << tetDualMesh.nPoints()
234  << " of which mapped to" << nl
235  << " cells : " << nCells << nl
236  << " patch faces : " << nPatchFaces << nl
237  << " not mapped : " << nUnmapped << nl
238  << endl;
239 
240 
241  // Read objects in time directory
242  IOobjectList objects(mesh, runTime.name());
243 
244  // Read vol fields, interpolate onto tet points
246  ReadAndMapFields<volScalarField, pointScalarField>
247  (
248  mesh,
249  objects,
250  tetDualMesh,
251  pointDualAddressing,
252  Zero, // nullValue
253  psFlds
254  );
255 
257  ReadAndMapFields<volVectorField, pointVectorField>
258  (
259  mesh,
260  objects,
261  tetDualMesh,
262  pointDualAddressing,
263  Zero, // nullValue
264  pvFlds
265  );
266 
268  ReadAndMapFields<volSphericalTensorField, pointSphericalTensorField>
269  (
270  mesh,
271  objects,
272  tetDualMesh,
273  pointDualAddressing,
274  Zero, // nullValue
275  pstFlds
276  );
277 
279  ReadAndMapFields<volSymmTensorField, pointSymmTensorField>
280  (
281  mesh,
282  objects,
283  tetDualMesh,
284  pointDualAddressing,
285  Zero, // nullValue
286  psymmtFlds
287  );
288 
290  ReadAndMapFields<volTensorField, pointTensorField>
291  (
292  mesh,
293  objects,
294  tetDualMesh,
295  pointDualAddressing,
296  Zero, // nullValue
297  ptFlds
298  );
299 
300  tetDualMesh.objectRegistry::write();
301 
302  Info<< "End\n" << endl;
303 
304  return 0;
305 }
306 
307 
308 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Pre-declare SubField and related Field type.
Definition: Field.H:82
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const labelList & patchID() const
Per boundary face label the patch index.
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:1025
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:405
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:271
label nPoints() const
label nInternalFaces() const
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Definition: timeSelector.C:252
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
error FatalError
static const char nl
Definition: Ostream.H:260
if(fields) ReadFields(pointMesh PtrList< pointVectorField > pvFlds
PtrList< pointScalarField > psFlds
if(fields) ReadFields(pointMesh PtrList< pointSphericalTensorField > pstFlds
if(fields) ReadFields(pointMesh PtrList< pointTensorField > ptFlds
objects
Foam::argList args(argc, argv)