singleCellMesh.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  singleCellMesh
26 
27 Description
28  Reads all fields and maps them to a mesh with all internal faces removed
29  (singleCellFvMesh) which gets written to region "singleCell".
30 
31  Used to generate mesh and fields that can be used for boundary-only data.
32  Might easily result in illegal mesh though so only look at boundaries
33  in paraview.
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "argList.H"
38 #include "fvMesh.H"
39 #include "volFields.H"
40 #include "Time.H"
41 #include "ReadFields.H"
42 #include "singleCellFvMesh.H"
43 #include "timeSelector.H"
44 
45 using namespace Foam;
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 // Name of region to create
50 const string singleCellName = "singleCell";
51 
52 
53 template<class GeoField>
54 void interpolateFields
55 (
56  const singleCellFvMesh& scMesh,
57  const PtrList<GeoField>& flds
58 )
59 {
60  forAll(flds, i)
61  {
62  tmp<GeoField> scFld = scMesh.interpolate(flds[i]);
63  GeoField* scFldPtr = scFld.ptr();
64  scFldPtr->writeOpt() = IOobject::AUTO_WRITE;
65  scFldPtr->store();
66  }
67 }
68 
69 
70 
71 int main(int argc, char *argv[])
72 {
73  // constant, not false
74  timeSelector::addOptions(true, false);
76  (
77  "noFields",
78  "do not update fields"
79  );
80 
81  #include "setRootCase.H"
82  #include "createTime.H"
83 
84  const bool fields = !args.optionFound("noFields");
85 
87 
88  #include "createNamedMesh.H"
89 
90  if (regionName == singleCellName)
91  {
93  << "Cannot convert region " << singleCellName
94  << " since result would overwrite it. Please rename your region."
95  << exit(FatalError);
96  }
97 
98  // Create the mesh
99  Info<< "Creating singleCell mesh" << nl << endl;
101  (
102  new singleCellFvMesh
103  (
104  IOobject
105  (
106  singleCellName,
107  mesh.polyMesh::instance(),
108  runTime,
111  ),
112  mesh
113  )
114  );
115 
116  forAll(timeDirs, timeI)
117  {
118  runTime.setTime(timeDirs[timeI], timeI);
119 
120  Info<< nl << "Time = " << runTime.name() << endl;
121 
122 
123  // Check for new mesh
124  if (mesh.readUpdate() != polyMesh::UNCHANGED)
125  {
126  Info<< "Detected changed mesh. Recreating singleCell mesh." << endl;
127  scMesh.clear(); // remove any registered objects
128  scMesh.reset
129  (
130  new singleCellFvMesh
131  (
132  IOobject
133  (
134  singleCellName,
135  mesh.polyMesh::instance(),
136  runTime,
139  ),
140  mesh
141  )
142  );
143  }
144 
145 
146  // Read objects in time directory
147  IOobjectList objects(mesh, runTime.name());
148 
149  if (fields) Info<< "Reading geometric fields" << nl << endl;
150 
151  #include "readVolFields.H"
152 
153  // Map and store the fields on the scMesh.
154  if (fields) interpolateFields(scMesh(), vsFlds);
155  if (fields) interpolateFields(scMesh(), vvFlds);
156  if (fields) interpolateFields(scMesh(), vstFlds);
157  if (fields) interpolateFields(scMesh(), vsymtFlds);
158  if (fields) interpolateFields(scMesh(), vtFlds);
159 
160 
161  // Write
162  Info<< "Writing mesh to time " << runTime.name() << endl;
163  scMesh().write();
164  }
165 
166 
167  Info<< "End\n" << endl;
168 
169  return 0;
170 }
171 
172 
173 // ************************************************************************* //
Field reading functions for post-processing utilities.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
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
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1736
void clear()
Remove all regIOobject owned by the registry.
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
fvMesh as subset of other mesh. Consists of one cell and all original boundary faces....
tmp< VolField< Type > > interpolate(const VolField< Type > &) const
Map volField. Internal field set to average, patch fields straight.
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
A class for managing temporary objects.
Definition: tmp.H:55
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:205
Foam::word regionName
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
int main(int argc, char *argv[])
Definition: financialFoam.C:44
static instantList timeDirs
Definition: globalFoam.H:44
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
error FatalError
static const char nl
Definition: Ostream.H:260
objects
PtrList< volSphericalTensorField > vstFlds
Definition: readVolFields.H:9
PtrList< volScalarField > vsFlds
Definition: readVolFields.H:3
PtrList< volTensorField > vtFlds
Definition: readVolFields.H:15
PtrList< volSymmTensorField > vsymtFlds
Definition: readVolFields.H:12
PtrList< volVectorField > vvFlds
Definition: readVolFields.H:6
Foam::argList args(argc, argv)