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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Application
27  foamToTetDualMesh
28 
29 Group
30  grpPostProcessingUtilities
31 
32 Description
33  Convert polyMesh results to tetDualMesh.
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "argList.H"
38 #include "fvMesh.H"
39 #include "volFields.H"
40 #include "pointFields.H"
41 #include "Time.H"
42 #include "IOobjectList.H"
43 
44 using namespace Foam;
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 template<class ReadGeoField, class MappedGeoField>
49 void ReadAndMapFields
50 (
51  const fvMesh& mesh,
52  const IOobjectList& objects,
53  const fvMesh& tetDualMesh,
54  const labelList& map,
55  const typename MappedGeoField::value_type& nullValue,
56  PtrList<MappedGeoField>& tetFields
57 )
58 {
59  typedef typename MappedGeoField::value_type Type;
60 
61  // Search list of objects for wanted type
62  IOobjectList fieldObjects(objects.lookupClass(ReadGeoField::typeName));
63 
64  tetFields.setSize(fieldObjects.size());
65 
66  label i = 0;
67  forAllConstIters(fieldObjects, iter)
68  {
69  Info<< "Converting " << ReadGeoField::typeName << ' ' << iter.key()
70  << endl;
71 
72  ReadGeoField readField(*iter(), mesh);
73 
74  tetFields.set
75  (
76  i,
77  new MappedGeoField
78  (
79  IOobject
80  (
81  readField.name(),
82  readField.instance(),
83  readField.local(),
84  tetDualMesh,
87  readField.registerObject()
88  ),
89  pointMesh::New(tetDualMesh),
90  dimensioned<Type>(readField.dimensions(), Zero)
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 int main(int argc, char *argv[])
144 {
146  (
147  "Convert polyMesh results to tetDualMesh"
148  );
149 
150  #include "addOverwriteOption.H"
151  #include "addTimeOptions.H"
152 
153  #include "setRootCase.H"
154  #include "createTime.H"
155  // Get times list
156  instantList Times = runTime.times();
157  #include "checkTimeOptions.H"
159 
160  // Read the mesh
161  #include "createNamedMesh.H"
162 
163  // Read the tetDualMesh
164  Info<< "Create tetDualMesh for time = "
165  << runTime.timeName() << nl << endl;
166 
167  fvMesh tetDualMesh
168  (
169  IOobject
170  (
171  "tetDualMesh",
172  runTime.timeName(),
173  runTime,
175  )
176  );
177  // From tet vertices to poly cells/faces
178  const labelIOList pointDualAddressing
179  (
180  IOobject
181  (
182  "pointDualAddressing",
183  tetDualMesh.facesInstance(),
184  tetDualMesh.meshSubDir,
185  tetDualMesh,
187  )
188  );
189 
190  if (pointDualAddressing.size() != tetDualMesh.nPoints())
191  {
193  << "Size " << pointDualAddressing.size()
194  << " of addressing map " << pointDualAddressing.objectPath()
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 // ************************************************************************* //
volFields.H
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:165
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:190
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
addOverwriteOption.H
Foam::polyMesh::meshSubDir
static word meshSubDir
Definition: polyMesh.H:317
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Definition: MeshObject.C:41
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:773
Foam::argList::addNote
static void addNote(const string &note)
Definition: argList.C:405
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Definition: polyMesh.C:845
IOobjectList.H
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Definition: polyMesh.H:440
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::PtrList::set
const T * set(const label i) const
Definition: PtrList.H:134
addTimeOptions.H
Foam::sumOp
Definition: ops.H:207
Foam::primitiveMesh::nPoints
label nPoints() const noexcept
Definition: primitiveMeshI.H:30
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::polyBoundaryMesh::patchID
const labelList & patchID() const
Definition: polyBoundaryMesh.C:449
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:45
Foam::Field< Type >
Foam::Info
messageStream Info
argList.H
Foam::TimePaths::times
instantList times() const
Definition: TimePaths.C:142
Foam::PtrList::setSize
void setSize(const label newLen)
Definition: PtrList.H:147
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:55
fld
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))
Definition: gmvOutputLagrangian.H:23
Foam::FatalError
error FatalError
createNamedMesh.H
Required Variables.
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned< Type >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:81
fvMesh.H
Foam
Definition: atmBoundaryLayer.C:26
checkTimeOptions.H
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Time.H
setRootCase.H
FatalErrorInFunction
#define FatalErrorInFunction
Definition: error.H:465
Foam::nl
constexpr char nl
Definition: Ostream.H:424
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const noexcept
Definition: primitiveMeshI.H:71
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::Time::setTime
virtual void setTime(const Time &t)
Definition: Time.C:996
Foam::IOList< label >
createTime.H
startTime
Foam::label startTime
Definition: checkTimeOptions.H:1
Foam::IOobjectList::lookupClass
IOobjectList lookupClass(const char *clsName) const
Definition: IOobjectList.C:283
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:184
pointFields.H
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:181