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 | Copyright (C) 2011-2015 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);
75 
76  #include "setRootCase.H"
77  #include "createTime.H"
78 
79  instantList timeDirs = timeSelector::select0(runTime, args);
80 
81  #include "createNamedMesh.H"
82 
83  if (regionName == singleCellName)
84  {
86  << "Cannot convert region " << singleCellName
87  << " since result would overwrite it. Please rename your region."
88  << exit(FatalError);
89  }
90 
91  // Create the mesh
92  Info<< "Creating singleCell mesh" << nl << endl;
94  (
96  (
97  IOobject
98  (
99  singleCellName,
100  mesh.polyMesh::instance(),
101  runTime,
104  ),
105  mesh
106  )
107  );
108  // For convenience create any fv* files
109  if (!exists(scMesh().fvSolution::objectPath()))
110  {
111  mkDir(scMesh().fvSolution::path());
112  ln("../fvSolution", scMesh().fvSolution::objectPath());
113  }
114  if (!exists(scMesh().fvSchemes::objectPath()))
115  {
116  mkDir(scMesh().fvSolution::path());
117  ln("../fvSchemes", scMesh().fvSchemes::objectPath());
118  }
119 
120 
121  forAll(timeDirs, timeI)
122  {
123  runTime.setTime(timeDirs[timeI], timeI);
124 
125  Info<< nl << "Time = " << runTime.timeName() << endl;
126 
127 
128  // Check for new mesh
130  {
131  Info<< "Detected changed mesh. Recreating singleCell mesh." << endl;
132  scMesh.clear(); // remove any registered objects
133  scMesh.reset
134  (
135  new singleCellFvMesh
136  (
137  IOobject
138  (
139  singleCellName,
140  mesh.polyMesh::instance(),
141  runTime,
144  ),
145  mesh
146  )
147  );
148  }
149 
150 
151  // Read objects in time directory
152  IOobjectList objects(mesh, runTime.timeName());
153 
154  // Read vol fields.
156  ReadFields(mesh, objects, vsFlds);
157 
159  ReadFields(mesh, objects, vvFlds);
160 
162  ReadFields(mesh, objects, vstFlds);
163 
164  PtrList<volSymmTensorField> vsymtFlds;
165  ReadFields(mesh, objects, vsymtFlds);
166 
168  ReadFields(mesh, objects, vtFlds);
169 
170  // Map and store the fields on the scMesh.
171  interpolateFields(scMesh(), vsFlds);
172  interpolateFields(scMesh(), vvFlds);
173  interpolateFields(scMesh(), vstFlds);
174  interpolateFields(scMesh(), vsymtFlds);
175  interpolateFields(scMesh(), vtFlds);
176 
177 
178  // Write
179  Info<< "Writing mesh to time " << runTime.timeName() << endl;
180  scMesh().write();
181  }
182 
183 
184  Info<< "End\n" << endl;
185 
186  return 0;
187 }
188 
189 
190 // ************************************************************************* //
fileName objectPath() const
Return complete path + object name.
Definition: IOobject.H:363
wordList ReadFields(const Mesh &mesh, const IOobjectList &objects, PtrList< GeoField > &fields, const bool syncPar=true)
Read all fields of the specified type.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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.
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Foam::word regionName
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Field reading functions for post-processing utilities.
tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Map volField. Internal field set to average, patch fields straight.
fileName path() const
Return complete path.
Definition: IOobject.C:275
dynamicFvMesh & mesh
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:494
void clear()
Clear all entries from table.
Definition: HashTable.C:464
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:720
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:198
static const char nl
Definition: Ostream.H:262
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Definition: timeSelector.C:253
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:295
bool exists(const fileName &, const bool checkGzip=true)
Does the name exist (as DIRECTORY or FILE) in the file system?
Definition: POSIX.C:480
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
messageStream Info
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
A class for managing temporary objects.
Definition: PtrList.H:54
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:870
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
Namespace for OpenFOAM.