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-2024 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 "timeSelector.H"
40 #include "fvMesh.H"
41 #include "volFields.H"
42 #include "pointFields.H"
43 #include "Time.H"
44 #include "IOobjectList.H"
45 
46 using namespace Foam;
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 template<class ReadGeoField, class MappedGeoField>
51 void ReadAndMapFields
52 (
53  const fvMesh& mesh,
54  const IOobjectList& objects,
55  const fvMesh& tetDualMesh,
56  const labelList& map,
57  const typename MappedGeoField::value_type& nullValue,
58  PtrList<MappedGeoField>& tetFields
59 )
60 {
61  typedef typename MappedGeoField::value_type Type;
62 
63  // Objects of wanted type
64  const UPtrList<const IOobject> fieldObjects
65  (
66  objects.csorted<ReadGeoField>()
67  );
68 
69  tetFields.resize(fieldObjects.size());
70 
71  label i = 0;
72  for (const IOobject& io : fieldObjects)
73  {
74  Info<< "Converting "
75  << ReadGeoField::typeName << ' ' << io.name() << endl;
76 
77  ReadGeoField readField(io, mesh);
78 
79  tetFields.set
80  (
81  i,
82  new MappedGeoField
83  (
84  IOobject
85  (
86  readField.name(),
87  readField.instance(),
88  readField.local(),
89  tetDualMesh,
92  readField.registerObject()
93  ),
94  pointMesh::New(tetDualMesh),
95  dimensioned<Type>(readField.dimensions(), Zero)
96  )
97  );
98 
99  Field<Type>& fld = tetFields[i].primitiveFieldRef();
100 
101  // Map from read field. Set unmapped entries to nullValue.
102  fld.setSize(map.size(), nullValue);
103  forAll(map, pointi)
104  {
105  label index = map[pointi];
106 
107  if (index > 0)
108  {
109  label celli = index-1;
110  fld[pointi] = readField[celli];
111  }
112  else if (index < 0)
113  {
114  const label facei = -index-1;
115  const label patchi = mesh.boundaryMesh().patchID(facei);
116  if (patchi >= 0)
117  {
118  label localFacei = mesh.boundaryMesh()[patchi].whichFace
119  (
120  facei
121  );
122  fld[pointi] = readField.boundaryField()[patchi][localFacei];
123  }
124  //else
125  //{
126  // FatalErrorInFunction
127  // << "Face " << facei << " from index " << index
128  // << " is not a boundary face." << abort(FatalError);
129  //}
130 
131  }
132  //else
133  //{
134  // WarningInFunction
135  // << "Point " << pointi << " at "
136  // << tetDualMesh.points()[pointi]
137  // << " has no dual correspondence." << endl;
138  //}
139  }
140  tetFields[i].correctBoundaryConditions();
141 
142  i++;
143  }
144 }
145 
146 
147 int main(int argc, char *argv[])
148 {
149  argList::noFunctionObjects(); // Never use function objects
150 
151  timeSelector::addOptions_singleTime(); // Single-time options
152 
154  (
155  "Convert polyMesh results to tetDualMesh"
156  );
157 
158  #include "addOverwriteOption.H"
159 
160  #include "setRootCase.H"
161  #include "createTime.H"
162 
163  // Set time from specified time options, or force start from Time=0
165 
166  // Read the mesh
167  #include "createNamedMesh.H"
168 
169  // Read the tetDualMesh
170  Info<< "Create tetDualMesh for time = "
171  << runTime.timeName() << nl << endl;
172 
173  fvMesh tetDualMesh
174  (
175  IOobject
176  (
177  "tetDualMesh",
178  runTime.timeName(),
179  runTime,
181  )
182  );
183  // From tet vertices to poly cells/faces
184  const labelIOList pointDualAddressing
185  (
186  IOobject
187  (
188  "pointDualAddressing",
189  tetDualMesh.facesInstance(),
191  tetDualMesh,
193  )
194  );
195 
196  if (pointDualAddressing.size() != tetDualMesh.nPoints())
197  {
199  << "Size " << pointDualAddressing.size()
200  << " of addressing map " << pointDualAddressing.objectPath()
201  << " differs from number of points in mesh "
202  << tetDualMesh.nPoints()
203  << exit(FatalError);
204  }
205 
206 
207  // Some stats on addressing
208  label nCells = 0;
209  label nPatchFaces = 0;
210  label nUnmapped = 0;
211  forAll(pointDualAddressing, pointi)
212  {
213  label index = pointDualAddressing[pointi];
214 
215  if (index > 0)
216  {
217  nCells++;
218  }
219  else if (index == 0)
220  {
221  nUnmapped++;
222  }
223  else
224  {
225  label facei = -index-1;
226  if (facei < mesh.nInternalFaces())
227  {
229  << "Face " << facei << " from index " << index
230  << " is not a boundary face."
231  << " nInternalFaces:" << mesh.nInternalFaces()
232  << exit(FatalError);
233  }
234  else
235  {
236  nPatchFaces++;
237  }
238  }
239  }
240 
241  reduce(nCells, sumOp<label>());
242  reduce(nPatchFaces, sumOp<label>());
243  reduce(nUnmapped, sumOp<label>());
244  Info<< "tetDualMesh points : " << tetDualMesh.nPoints()
245  << " of which mapped to" << nl
246  << " cells : " << nCells << nl
247  << " patch faces : " << nPatchFaces << nl
248  << " not mapped : " << nUnmapped << nl
249  << endl;
250 
251 
252  // Read objects in time directory
253  IOobjectList objects(mesh, runTime.timeName());
254 
255  // Read vol fields, interpolate onto tet points
257  ReadAndMapFields<volScalarField, pointScalarField>
258  (
259  mesh,
260  objects,
261  tetDualMesh,
262  pointDualAddressing,
263  Zero, // nullValue
264  psFlds
265  );
266 
268  ReadAndMapFields<volVectorField, pointVectorField>
269  (
270  mesh,
271  objects,
272  tetDualMesh,
273  pointDualAddressing,
274  Zero, // nullValue
275  pvFlds
276  );
277 
279  ReadAndMapFields<volSphericalTensorField, pointSphericalTensorField>
280  (
281  mesh,
282  objects,
283  tetDualMesh,
284  pointDualAddressing,
285  Zero, // nullValue
286  pstFlds
287  );
288 
290  ReadAndMapFields<volSymmTensorField, pointSymmTensorField>
291  (
292  mesh,
293  objects,
294  tetDualMesh,
295  pointDualAddressing,
296  Zero, // nullValue
297  psymmtFlds
298  );
299 
301  ReadAndMapFields<volTensorField, pointTensorField>
302  (
303  mesh,
304  objects,
305  tetDualMesh,
306  pointDualAddressing,
307  Zero, // nullValue
308  ptFlds
309  );
310 
311  tetDualMesh.objectRegistry::write();
312 
313  Info<< "End\n" << endl;
314 
315  return 0;
316 }
317 
318 
319 // ************************************************************************* //
static void noFunctionObjects(bool addWithOption=false)
Remove &#39;-noFunctionObjects&#39; option and ignore any occurrences.
Definition: argList.C:547
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:608
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 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:411
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.
static bool setTimeIfPresent(Time &runTime, const argList &args, const bool forceInitial=false)
Set the runTime based on -constant (if present), -time (value), or -latestTime.
Definition: timeSelector.C:314
Generic dimensioned Type class.
#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:609
static void addOptions_singleTime()
Add single-time timeSelector options to argList::validOptions()
Definition: timeSelector.C:146
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))
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()
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)
Foam::argList args(argc, argv)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
UPtrList< const IOobject > csorted() const
The sorted list of IOobjects with headerClassName == Type::typeName.
void reduce(T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) using linear/tree communication schedule.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127