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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  foamToTetDualMesh
29 
30 Group
31  grpPostProcessingUtilities
32 
33 Description
34  Convert polyMesh results to tetDualMesh.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "argList.H"
39 #include "fvMesh.H"
40 #include "volFields.H"
41 #include "pointFields.H"
42 #include "Time.H"
43 #include "IOobjectList.H"
44 
45 using namespace Foam;
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 template<class ReadGeoField, class MappedGeoField>
50 void ReadAndMapFields
51 (
52  const fvMesh& mesh,
53  const IOobjectList& objects,
54  const fvMesh& tetDualMesh,
55  const labelList& map,
56  const typename MappedGeoField::value_type& nullValue,
57  PtrList<MappedGeoField>& tetFields
58 )
59 {
60  typedef typename MappedGeoField::value_type Type;
61 
62  // Objects of wanted type
63  const UPtrList<const IOobject> fieldObjects
64  (
65  objects.csorted<ReadGeoField>()
66  );
67 
68  tetFields.resize(fieldObjects.size());
69 
70  label i = 0;
71  for (const IOobject& io : fieldObjects)
72  {
73  Info<< "Converting "
74  << ReadGeoField::typeName << ' ' << io.name() << endl;
75 
76  ReadGeoField readField(io, mesh);
77 
78  tetFields.set
79  (
80  i,
81  new MappedGeoField
82  (
83  IOobject
84  (
85  readField.name(),
86  readField.instance(),
87  readField.local(),
88  tetDualMesh,
91  readField.registerObject()
92  ),
93  pointMesh::New(tetDualMesh),
94  dimensioned<Type>(readField.dimensions(), Zero)
95  )
96  );
97 
98  Field<Type>& fld = tetFields[i].primitiveFieldRef();
99 
100  // Map from read field. Set unmapped entries to nullValue.
101  fld.setSize(map.size(), nullValue);
102  forAll(map, pointi)
103  {
104  label index = map[pointi];
105 
106  if (index > 0)
107  {
108  label celli = index-1;
109  fld[pointi] = readField[celli];
110  }
111  else if (index < 0)
112  {
113  const label facei = -index-1;
114  const label patchi = mesh.boundaryMesh().patchID(facei);
115  if (patchi >= 0)
116  {
117  label localFacei = mesh.boundaryMesh()[patchi].whichFace
118  (
119  facei
120  );
121  fld[pointi] = readField.boundaryField()[patchi][localFacei];
122  }
123  //else
124  //{
125  // FatalErrorInFunction
126  // << "Face " << facei << " from index " << index
127  // << " is not a boundary face." << abort(FatalError);
128  //}
129 
130  }
131  //else
132  //{
133  // WarningInFunction
134  // << "Point " << pointi << " at "
135  // << tetDualMesh.points()[pointi]
136  // << " has no dual correspondence." << endl;
137  //}
138  }
139  tetFields[i].correctBoundaryConditions();
140 
141  i++;
142  }
143 }
144 
145 
146 int main(int argc, char *argv[])
147 {
149  (
150  "Convert polyMesh results to tetDualMesh"
151  );
152 
153  #include "addOverwriteOption.H"
154  #include "addTimeOptions.H"
155 
156  #include "setRootCase.H"
157  #include "createTime.H"
158  // Get times list
159  instantList Times = runTime.times();
160  #include "checkTimeOptions.H"
162 
163  // Read the mesh
164  #include "createNamedMesh.H"
165 
166  // Read the tetDualMesh
167  Info<< "Create tetDualMesh for time = "
168  << runTime.timeName() << nl << endl;
169 
170  fvMesh tetDualMesh
171  (
172  IOobject
173  (
174  "tetDualMesh",
175  runTime.timeName(),
176  runTime,
178  )
179  );
180  // From tet vertices to poly cells/faces
181  const labelIOList pointDualAddressing
182  (
183  IOobject
184  (
185  "pointDualAddressing",
186  tetDualMesh.facesInstance(),
188  tetDualMesh,
190  )
191  );
192 
193  if (pointDualAddressing.size() != tetDualMesh.nPoints())
194  {
196  << "Size " << pointDualAddressing.size()
197  << " of addressing map " << pointDualAddressing.objectPath()
198  << " differs from number of points in mesh "
199  << tetDualMesh.nPoints()
200  << exit(FatalError);
201  }
202 
203 
204  // Some stats on addressing
205  label nCells = 0;
206  label nPatchFaces = 0;
207  label nUnmapped = 0;
208  forAll(pointDualAddressing, pointi)
209  {
210  label index = pointDualAddressing[pointi];
211 
212  if (index > 0)
213  {
214  nCells++;
215  }
216  else if (index == 0)
217  {
218  nUnmapped++;
219  }
220  else
221  {
222  label facei = -index-1;
223  if (facei < mesh.nInternalFaces())
224  {
226  << "Face " << facei << " from index " << index
227  << " is not a boundary face."
228  << " nInternalFaces:" << mesh.nInternalFaces()
229  << exit(FatalError);
230  }
231  else
232  {
233  nPatchFaces++;
234  }
235  }
236  }
237 
238  reduce(nCells, sumOp<label>());
239  reduce(nPatchFaces, sumOp<label>());
240  reduce(nUnmapped, sumOp<label>());
241  Info<< "tetDualMesh points : " << tetDualMesh.nPoints()
242  << " of which mapped to" << nl
243  << " cells : " << nCells << nl
244  << " patch faces : " << nPatchFaces << nl
245  << " not mapped : " << nUnmapped << nl
246  << endl;
247 
248 
249  // Read objects in time directory
250  IOobjectList objects(mesh, runTime.timeName());
251 
252  // Read vol fields, interpolate onto tet points
254  ReadAndMapFields<volScalarField, pointScalarField>
255  (
256  mesh,
257  objects,
258  tetDualMesh,
259  pointDualAddressing,
260  Zero, // nullValue
261  psFlds
262  );
263 
265  ReadAndMapFields<volVectorField, pointVectorField>
266  (
267  mesh,
268  objects,
269  tetDualMesh,
270  pointDualAddressing,
271  Zero, // nullValue
272  pvFlds
273  );
274 
276  ReadAndMapFields<volSphericalTensorField, pointSphericalTensorField>
277  (
278  mesh,
279  objects,
280  tetDualMesh,
281  pointDualAddressing,
282  Zero, // nullValue
283  pstFlds
284  );
285 
287  ReadAndMapFields<volSymmTensorField, pointSymmTensorField>
288  (
289  mesh,
290  objects,
291  tetDualMesh,
292  pointDualAddressing,
293  Zero, // nullValue
294  psymmtFlds
295  );
296 
298  ReadAndMapFields<volTensorField, pointTensorField>
299  (
300  mesh,
301  objects,
302  tetDualMesh,
303  pointDualAddressing,
304  Zero, // nullValue
305  ptFlds
306  );
307 
308  tetDualMesh.objectRegistry::write();
309 
310  Info<< "End\n" << endl;
311 
312  return 0;
313 }
314 
315 
316 // ************************************************************************* //
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:859
label nPoints() const noexcept
Number of mesh points.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:410
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const labelList & patchID() const
Per boundary face label the patch index.
Required Classes.
Generic dimensioned Type class.
Required Classes.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
dynamicFvMesh & mesh
Generic templated field type.
Definition: Field.H:62
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:935
label nInternalFaces() const noexcept
Number of internal faces.
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: HashTable.H:106
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
void resize(const label newLen)
Adjust size of PtrList.
Definition: PtrList.C:95
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
instantList times() const
Search the case for valid time directories.
Definition: TimePaths.C:118
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
Automatically write from objectRegistry::writeObject()
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
messageStream Info
Information stream (stdout output on master, null elsewhere)
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Foam::label startTime
UPtrList< const IOobject > csorted() const
The sorted list of IOobjects with headerClassName == Type::typeName.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127