sampledSets.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-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2020 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 \*---------------------------------------------------------------------------*/
28 
29 #include "sampledSets.H"
30 #include "dictionary.H"
31 #include "Time.H"
32 #include "volFields.H"
33 #include "volPointInterpolation.H"
34 #include "mapPolyMesh.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(sampledSets, 0);
42 
44  (
45  functionObject,
46  sampledSets,
47  dictionary
48  );
49 }
50 
51 bool Foam::sampledSets::verbose_ = false;
52 
53 
54 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55 
56 void Foam::sampledSets::combineSampledSets
57 (
58  PtrList<coordSet>& masterSampledSets,
59  labelListList& indexSets
60 )
61 {
62  // Combine sampleSets from processors. Sort by curveDist. Return
63  // ordering in indexSets.
64  // Note: only master results are valid
65 
66  masterSampledSets_.clear();
67  masterSampledSets_.setSize(size());
68  indexSets_.setSize(size());
69 
70  const PtrList<sampledSet>& sampledSets = *this;
71 
72  forAll(sampledSets, setI)
73  {
74  labelList segments;
75  masterSampledSets.set
76  (
77  setI,
78  sampledSets[setI].gather(indexSets[setI], segments)
79  );
80  }
81 }
82 
83 
84 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
85 
86 Foam::sampledSets::sampledSets
87 (
88  const word& name,
89  const Time& runTime,
90  const dictionary& dict
91 )
92 :
95  mesh_(refCast<const fvMesh>(obr_)),
96  loadFromFiles_(false),
97  outputPath_(fileName::null),
98  searchEngine_(mesh_),
99  interpolationScheme_(word::null),
100  writeFormat_(word::null),
101  writeFormatOptions_(dict.subOrEmptyDict("formatOptions"))
102 {
103  outputPath_ =
104  (
105  mesh_.time().globalPath()/functionObject::outputPrefix/name
106  );
107 
108  if (mesh_.name() != polyMesh::defaultRegion)
109  {
110  outputPath_ /= mesh_.name();
111  }
112 
113  outputPath_.clean(); // Remove unneeded ".."
114 
115  read(dict);
116 }
117 
118 
119 Foam::sampledSets::sampledSets
120 (
121  const word& name,
122  const objectRegistry& obr,
123  const dictionary& dict,
124  const bool loadFromFiles
125 )
126 :
129  mesh_(refCast<const fvMesh>(obr)),
130  loadFromFiles_(loadFromFiles),
131  outputPath_(fileName::null),
132  searchEngine_(mesh_),
133  interpolationScheme_(word::null),
134  writeFormat_(word::null),
135  writeFormatOptions_(dict.subOrEmptyDict("formatOptions"))
136 {
137  outputPath_ =
138  (
139  mesh_.time().globalPath()/functionObject::outputPrefix/name
140  );
141 
142  if (mesh_.name() != polyMesh::defaultRegion)
143  {
144  outputPath_ /= mesh_.name();
145  }
146 
147  outputPath_.clean(); // Remove unneeded ".."
148 
149  read(dict);
150 }
151 
152 
153 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
154 
155 bool Foam::sampledSets::verbose(const bool on)
156 {
157  bool old(verbose_);
158  verbose_ = on;
159  return old;
160 }
161 
162 
164 {
165  return true;
166 }
167 
168 
170 {
171  if (size())
172  {
173  const label nFields = classifyFields();
174 
175  if (Pstream::master())
176  {
177  if (debug)
178  {
179  Pout<< "timeName = " << mesh_.time().timeName() << nl
180  << "scalarFields " << scalarFields_ << nl
181  << "vectorFields " << vectorFields_ << nl
182  << "sphTensorFields " << sphericalTensorFields_ << nl
183  << "symTensorFields " << symmTensorFields_ <<nl
184  << "tensorFields " << tensorFields_ <<nl;
185  }
186 
187  if (nFields)
188  {
189  if (debug)
190  {
191  Pout<< "Creating directory "
192  << outputPath_/mesh_.time().timeName()
193  << nl << endl;
194  }
195 
196  mkDir(outputPath_/mesh_.time().timeName());
197  }
198  else
199  {
200  Info<< "No fields to sample" << endl;
201  }
202  }
203 
204  if (nFields)
205  {
206  sampleAndWrite(scalarFields_);
207  sampleAndWrite(vectorFields_);
208  sampleAndWrite(sphericalTensorFields_);
209  sampleAndWrite(symmTensorFields_);
210  sampleAndWrite(tensorFields_);
211  }
212  }
213 
214  return true;
215 }
216 
217 
219 {
220  dict_ = dict;
221 
222  if (dict_.found("sets"))
223  {
224  dict_.readEntry("fields", fieldSelection_);
225  clearFieldGroups();
226 
227  dict.readEntry("interpolationScheme", interpolationScheme_);
228  dict.readEntry("setFormat", writeFormat_);
229 
230  PtrList<sampledSet> newList
231  (
232  dict_.lookup("sets"),
233  sampledSet::iNew(mesh_, searchEngine_)
234  );
235  transfer(newList);
236  combineSampledSets(masterSampledSets_, indexSets_);
237 
238  if (this->size())
239  {
240  Info<< "Reading set description:" << nl;
241  forAll(*this, setI)
242  {
243  Info<< " " << operator[](setI).name() << nl;
244  }
245  Info<< endl;
246  }
247  }
248 
249  if (Pstream::master() && debug)
250  {
251  Pout<< "sample fields:" << fieldSelection_ << nl
252  << "sample sets:" << nl << "(" << nl;
253 
254  forAll(*this, setI)
255  {
256  Pout<< " " << operator[](setI) << endl;
257  }
258  Pout<< ")" << endl;
259  }
260 
261  return true;
262 }
263 
264 
266 {
267  if (dict_.found("sets"))
268  {
269  searchEngine_.correct();
270 
271  PtrList<sampledSet> newList
272  (
273  dict_.lookup("sets"),
274  sampledSet::iNew(mesh_, searchEngine_)
275  );
276  transfer(newList);
277  combineSampledSets(masterSampledSets_, indexSets_);
278  }
279 }
280 
281 
283 {
284  if (&mpm.mesh() == &mesh_)
285  {
286  correct();
287  }
288 }
289 
290 
292 {
293  if (&mesh == &mesh_)
294  {
295  correct();
296  }
297 }
298 
299 
301 {
302  if (state != polyMesh::UNCHANGED)
303  {
304  correct();
305  }
306 }
307 
308 
309 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
volFields.H
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::sampledSets::readUpdate
virtual void readUpdate(const polyMesh::readUpdateState state)
Update for changes of mesh due to readUpdate.
Definition: sampledSets.C:300
Foam::sampledSets::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: sampledSets.C:282
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::polyMesh::defaultRegion
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:318
mapPolyMesh.H
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::sampledSets::write
virtual bool write()
Sample and write.
Definition: sampledSets.C:169
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
sampledSets.H
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
correct
fvOptions correct(rho)
Foam::sampledSets::verbose
bool verbose(const bool on)
Enable/disable verbose output.
Definition: sampledSets.C:155
Foam::polyMesh::UNCHANGED
Definition: polyMesh.H:92
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:302
Foam::PtrList< sampledSet >
Foam::functionObject::outputPrefix
static word outputPrefix
Directory prefix.
Definition: functionObject.H:376
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::sampledSets::correct
void correct()
Correct for mesh changes.
Definition: sampledSets.C:265
Foam::dictionary::subOrEmptyDict
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:540
Foam::fileName::null
static const fileName null
An empty fileName.
Definition: fileName.H:101
Foam::polyMesh::readUpdateState
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:90
volPointInterpolation.H
Time.H
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::sampledSets::execute
virtual bool execute()
Execute, currently does nothing.
Definition: sampledSets.C:163
Foam::nl
constexpr char nl
Definition: Ostream.H:404
dictionary.H
Foam::sampledSets::movePoints
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: sampledSets.C:291
Foam::word::null
static const word null
An empty word.
Definition: word.H:80
Foam::List::set
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:341
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::functionObjects::regionFunctionObject
Specialization of Foam::functionObject for a region and providing a reference to the region Foam::obj...
Definition: regionFunctionObject.H:90
Foam::mapPolyMesh::mesh
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:363
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:507
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::sampledSet::iNew
Class used for the read-construction of.
Definition: sampledSet.H:204
Foam::sampledSets::read
virtual bool read(const dictionary &)
Read the sampledSets.
Definition: sampledSets.C:218