singleCellMesh.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-2015 OpenFOAM Foundation
9  Copyright (C) 2016-2021 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  singleCellMesh
29 
30 Group
31  grpMeshManipulationUtilities
32 
33 Description
34  Reads all fields and maps them to a mesh with all internal faces removed
35  (singleCellFvMesh) which gets written to region "singleCell".
36 
37  Used to generate mesh and fields that can be used for boundary-only data.
38  Might easily result in illegal mesh though so only look at boundaries
39  in paraview.
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #include "argList.H"
44 #include "fvMesh.H"
45 #include "volFields.H"
46 #include "Time.H"
47 #include "ReadFields.H"
48 #include "singleCellFvMesh.H"
49 #include "timeSelector.H"
50 
51 using namespace Foam;
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 // Name of region to create
56 const string singleCellName = "singleCell";
57 
58 
59 template<class GeoField>
60 void interpolateFields
61 (
62  const singleCellFvMesh& scMesh,
63  const PtrList<GeoField>& flds
64 )
65 {
66  forAll(flds, i)
67  {
68  tmp<GeoField> scFld = scMesh.interpolate(flds[i]);
69  GeoField* scFldPtr = scFld.ptr();
70  scFldPtr->writeOpt(IOobject::AUTO_WRITE);
71  scFldPtr->store();
72  }
73 }
74 
75 
76 
77 int main(int argc, char *argv[])
78 {
80  (
81  "Map fields to a mesh with all internal faces removed"
82  " (singleCellFvMesh) which gets written to region 'singleCell'"
83  );
84 
85  timeSelector::addOptions(true, false); // constant(true), zero(false)
86 
87  #include "setRootCase.H"
88  #include "createTime.H"
89 
91 
92  #include "createNamedMesh.H"
93 
94  if (regionName == singleCellName)
95  {
97  << "Cannot convert region " << singleCellName
98  << " since result would overwrite it. Please rename your region."
99  << exit(FatalError);
100  }
101 
102  // Create the mesh
103  Info<< "Creating singleCell mesh" << nl << endl;
105  (
106  new singleCellFvMesh
107  (
108  IOobject
109  (
110  singleCellName,
112  runTime,
115  ),
116  mesh
117  )
118  );
119  scMesh().setInstance(mesh.pointsInstance());
120 
121  // For convenience create any fv* files
122  if (!exists(scMesh().fvSolution::objectPath()))
123  {
124  mkDir(scMesh().fvSolution::path());
125  ln("../fvSolution", scMesh().fvSolution::objectPath());
126  }
127  if (!exists(scMesh().fvSchemes::objectPath()))
128  {
129  mkDir(scMesh().fvSolution::path());
130  ln("../fvSchemes", scMesh().fvSchemes::objectPath());
131  }
132 
133 
134  forAll(timeDirs, timeI)
135  {
136  runTime.setTime(timeDirs[timeI], timeI);
137 
138  Info<< nl << "Time = " << runTime.timeName() << endl;
139 
140 
141  // Check for new mesh
143  {
144  Info<< "Detected changed mesh. Recreating singleCell mesh." << endl;
145  scMesh.clear(); // remove any registered objects
146  scMesh.reset
147  (
148  new singleCellFvMesh
149  (
150  IOobject
151  (
152  singleCellName,
154  runTime,
157  ),
158  mesh
159  )
160  );
161  }
162 
163 
164  // Read objects in time directory
165  IOobjectList objects(mesh, runTime.timeName());
166 
167  // Read vol fields.
169  ReadFields(mesh, objects, vsFlds);
170 
172  ReadFields(mesh, objects, vvFlds);
173 
175  ReadFields(mesh, objects, vstFlds);
176 
177  PtrList<volSymmTensorField> vsymtFlds;
178  ReadFields(mesh, objects, vsymtFlds);
179 
181  ReadFields(mesh, objects, vtFlds);
182 
183  // Map and store the fields on the scMesh.
184  interpolateFields(scMesh(), vsFlds);
185  interpolateFields(scMesh(), vvFlds);
186  interpolateFields(scMesh(), vstFlds);
187  interpolateFields(scMesh(), vsymtFlds);
188  interpolateFields(scMesh(), vtFlds);
189 
190 
191  // Write
192  Info<< "Writing mesh to time " << runTime.timeName() << endl;
193  scMesh().write();
194  }
195 
196 
197  Info<< "End\n" << endl;
198 
199  return 0;
200 }
201 
202 
203 // ************************************************************************* //
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
fvMesh as subset of other mesh. Consists of one cell and all original bounday faces. Useful when manipulating boundary data. Single internal cell only needed to be able to manipulate in a standard way.
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh >> &fields, const bool syncPar=true, const bool readOldTime=false)
Read Geometric fields of templated type.
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
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:578
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
Required Variables.
tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Map volField. Internal field set to average, patch fields straight copies.
Field reading functions for post-processing utilities.
fileName objectPath() const
The complete path + object name.
Definition: IOobjectI.H:251
fileName path() const
The complete path.
Definition: IOobject.C:449
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
Foam::word regionName(Foam::polyMesh::defaultRegion)
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:848
dynamicFvMesh & mesh
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:614
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:671
bool exists(const fileName &name, const bool checkGzip=true, const bool followLink=true)
Does the name exist (as DIRECTORY or FILE) in the file system?
Definition: POSIX.C:835
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:997
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1069
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:1237
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:770
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options and also set the runTime to the first i...
Definition: timeSelector.C:234
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition: polyMeshIO.C:29
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Nothing to be read.
Automatically write from objectRegistry::writeObject()
void clear()
Clear all entries from the registry.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
T * ptr() const
Return managed pointer for reuse, or clone() the object reference.
Definition: tmpI.H:265
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Foam::argList args(argc, argv)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:171
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
Definition: timeSelector.C:101
Namespace for OpenFOAM.