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-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  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 word singleCellName = "singleCell";
57 
58 
59 template<class GeoField>
60 wordList interpolateFields
61 (
62  const singleCellFvMesh& subsetter,
63  const PtrList<GeoField>& flds
64 )
65 {
66  wordList names(flds.size());
67 
68  forAll(flds, i)
69  {
70  GeoField* subFld = subsetter.interpolate(flds[i]).ptr();
71 
72  subFld->writeOpt(IOobject::AUTO_WRITE);
73  subFld->store();
74 
75  names[i] = subFld->name();
76  }
77 
78  return names;
79 }
80 
81 
82 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
83 
84 int main(int argc, char *argv[])
85 {
87  (
88  "Map fields to a mesh with all internal faces removed"
89  " (singleCellFvMesh) which gets written to region 'singleCell'"
90  );
91 
92  timeSelector::addOptions(true, false); // constant(true), zero(false)
93 
94  #include "setRootCase.H"
95  #include "createTime.H"
96 
98 
99  #include "createNamedMesh.H"
100 
101  if (regionName == singleCellName)
102  {
104  << "Cannot convert region " << regionName
105  << " since result would overwrite it. Please rename your region."
106  << exit(FatalError);
107  }
108 
109  // Create the mesh
110  Info<< "Creating singleCell mesh" << nl << endl;
112  (
113  new singleCellFvMesh
114  (
115  IOobject
116  (
117  singleCellName,
119  runTime,
123  ),
124  mesh
125  )
126  );
127  scMesh().setInstance(mesh.pointsInstance());
128 
129  // For convenience create any fv* files
130  if (!exists(scMesh().fvSolution::objectPath()))
131  {
132  mkDir(scMesh().fvSolution::path());
133  ln("../fvSolution", scMesh().fvSolution::objectPath());
134  }
135  if (!exists(scMesh().fvSchemes::objectPath()))
136  {
137  mkDir(scMesh().fvSolution::path());
138  ln("../fvSchemes", scMesh().fvSchemes::objectPath());
139  }
140 
141 
142  // List of stored objects to clear prior to reading
143  DynamicList<word> storedObjects;
144 
145  forAll(timeDirs, timei)
146  {
147  runTime.setTime(timeDirs[timei], timei);
148 
149  Info<< nl << "Time = " << runTime.timeName() << endl;
150 
151  // Purge any previously interpolated fields
152  if (!storedObjects.empty())
153  {
154  static_cast<objectRegistry&>(scMesh()).erase(storedObjects);
155  storedObjects.clear();
156  }
157 
158  // Check for new mesh
160  {
161  Info<< "Detected changed mesh. Recreating singleCell mesh." << endl;
162  scMesh.reset(nullptr); // first free any stored objects
163  scMesh.reset
164  (
165  new singleCellFvMesh
166  (
167  IOobject
168  (
169  singleCellName,
171  runTime,
175  ),
176  mesh
177  )
178  );
179  }
180 
181  // Read objects in time directory
182  IOobjectList objects(mesh, runTime.timeName());
183  storedObjects.reserve(objects.size());
184 
185  // Read fvMesh fields, map and store the interpolated fields
186  #define doLocalCode(GeoField) \
187  { \
188  PtrList<GeoField> flds; \
189  ReadFields(mesh, objects, flds); \
190  storedObjects.push_back(interpolateFields(scMesh(), flds)); \
191  }
192 
198 
199  #undef doLocalCode
200 
201  // Write
202  Info<< "Writing " << singleCellName
203  << " mesh/fields to time " << runTime.timeName() << endl;
204  scMesh().write();
205  }
206 
207 
208  Info<< "End\n" << endl;
209 
210  return 0;
211 }
212 
213 
214 // ************************************************************************* //
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 boundary faces. Useful when manipulating boundary data. Single internal cell only needed to be able to manipulate in a standard way.
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
List of names generated by calling name() for each list item and filtered for matches.
tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Map volField. Internal field set to average, patch fields straight copies.
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
srcOptions erase("case")
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
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:675
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Required Classes.
Field reading functions for post-processing utilities.
#define doLocalCode(FieldType, Variable)
fileName objectPath() const
The complete path + object name.
Definition: IOobjectI.H:284
fileName path() const
The complete path for the object (with instance, local,...).
Definition: IOobject.C:480
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:853
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:715
A class for handling words, derived from Foam::string.
Definition: word.H:63
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
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 bool write(const token &tok)=0
Write token to stream or otherwise handle it.
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:902
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:1237
void reserve(const label len)
Reserve allocation space for at least this size, allocating new space if required and retaining old c...
Definition: DynamicListI.H:333
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 clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:405
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:260
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
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()
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
Registry of regIOobjects.
Foam::argList args(argc, argv)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Do not request registration (bool: false)
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
Definition: timeSelector.C:101
Namespace for OpenFOAM.