sampledSets.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-2025 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 \*---------------------------------------------------------------------------*/
25 
26 #include "sampledSets.H"
27 #include "dictionary.H"
28 #include "Time.H"
29 #include "volFields.H"
30 #include "ListListOps.H"
31 #include "SortableList.H"
32 #include "volPointInterpolation.H"
33 #include "polyTopoChangeMap.H"
34 #include "polyMeshMap.H"
35 #include "polyDistributionMap.H"
36 #include "writeFile.H"
37 #include "OFstream.H"
38 #include "OSspecific.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace functionObjects
46 {
48 
50  (
54  );
55 }
56 }
57 
58 
59 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
60 
61 void Foam::functionObjects::sampledSets::combineSampledSets()
62 {
63  masterSets_.setSize(size());
64  masterSetOrders_.setSize(size());
65 
66  forAll(*this, seti)
67  {
68  const sampledSet& s = operator[](seti);
69 
70  Tuple2<coordSet, labelList> g = s.gather();
71 
72  masterSets_.set(seti, new coordSet(g.first()));
73  masterSetOrders_[seti] = g.second();
74 
75  if (Pstream::master() && masterSets_[seti].size() == 0)
76  {
78  << "Sample set " << s.name()
79  << " has zero points." << endl;
80  }
81  }
82 }
83 
84 
85 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
86 
88 (
89  const word& name,
90  const Time& t,
91  const dictionary& dict
92 )
93 :
96  outputPath_
97  (
98  mesh_.time().globalPath()
99  /writeFile::outputPrefix
100  /(mesh_.name() != polyMesh::defaultRegion ? mesh_.name() : word())
101  /name
102  ),
103  searchEngine_(mesh_),
104  interpolationScheme_(word::null),
105  formatter_(nullptr)
106 {
107  read(dict);
108 }
109 
110 
111 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
112 
114 {}
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
120 {
121  bool setsFound = dict.found("sets");
122  if (setsFound)
123  {
124  dict.lookup("fields") >> fields_;
125 
126  dict.lookup("interpolationScheme") >> interpolationScheme_;
127 
128  const word writeType(dict.lookup("setFormat"));
129 
130  // Define the set formatter
131  formatter_ = setWriter::New(writeType, dict);
132 
133  if (dict.isDict("sets"))
134  {
135  const dictionary& setsDict = dict.subDict("sets");
136 
137  setSize(setsDict.size());
138 
139  label i = 0;
140 
141  forAllConstIter(dictionary, setsDict, iter)
142  {
143  set
144  (
145  i++,
147  (
148  iter().keyword(),
149  mesh_,
150  searchEngine_,
151  iter().dict()
152  )
153  );
154  }
155  }
156  else
157  {
158  PtrList<sampledSet> newList
159  (
160  dict.lookup("sets"),
161  sampledSet::iNew(mesh_, searchEngine_)
162  );
163  transfer(newList);
164  }
165 
166  combineSampledSets();
167 
168  if (this->size())
169  {
170  Info<< "Reading set description:" << nl;
171  forAll(*this, seti)
172  {
173  Info<< " " << operator[](seti).name() << nl;
174  }
175  Info<< endl;
176  }
177  }
178 
179  return true;
180 }
181 
182 
184 {
185  return fields_;
186 }
187 
188 
190 {
191  return true;
192 }
193 
194 
196 {
197  if (size())
198  {
199  if (Pstream::master())
200  {
201  if (debug)
202  {
203  Pout<< "Creating directory "
204  << outputPath_/mesh_.time().name() << nl << endl;
205  }
206 
207  mkDir(outputPath_/mesh_.time().name());
208  }
209 
210  // Create a list of names of fields that are actually available
212  forAll(fields_, fieldi)
213  {
214  #define FoundFieldType(Type, nullArg) \
215  || foundObject<VolField<Type>>(fields_[fieldi])
217  {
218  fieldNames.append(fields_[fieldi]);
219  }
220  else
221  {
222  cannotFindObject(fields_[fieldi]);
223  }
224  #undef FoundFieldType
225  }
226 
227  // Create table of cached interpolations, to prevent unnecessary work
228  // when interpolating fields over multiple surfaces
229  #define DeclareInterpolations(Type, nullArg) \
230  HashPtrTable<interpolation<Type>> interpolation##Type##s;
232  #undef DeclareInterpolations
233 
234  // Sample and write the sets
235  forAll(*this, seti)
236  {
237  #define GenerateFieldTypeValues(Type, nullArg) \
238  PtrList<Field<Type>> field##Type##Values = \
239  sampleType<Type>(seti, fieldNames, interpolation##Type##s);
241  #undef GenerateFieldTypeValues
242 
243  if (Pstream::parRun())
244  {
245  if (Pstream::master() && masterSets_[seti].size())
246  {
247  formatter_->write
248  (
249  outputPath_/mesh_.time().name(),
250  operator[](seti).name(),
251  masterSets_[seti],
252  fieldNames
253  #define FieldTypeValuesParameter(Type, nullArg) \
254  , field##Type##Values
257  );
258 
259  }
260  }
261  else
262  {
263  if (operator[](seti).size())
264  {
265  formatter_->write
266  (
267  outputPath_/mesh_.time().name(),
268  operator[](seti).name(),
269  operator[](seti),
270  fieldNames
271  #define FieldTypeValuesParameter(Type, nullArg) \
272  , field##Type##Values
275  );
276  }
277  }
278  }
279  }
280 
281  return true;
282 }
283 
284 
286 {
287  if (&mesh == &mesh_)
288  {
289  if (this->size())
290  {
291  searchEngine_.correct();
292  }
293 
294  forAll(*this, seti)
295  {
296  operator[](seti).movePoints();
297  }
298 
299  if (this->size())
300  {
301  combineSampledSets();
302  }
303  }
304 }
305 
306 
308 (
309  const polyTopoChangeMap& map
310 )
311 {
312  if (&map.mesh() == &mesh_)
313  {
314  if (this->size())
315  {
316  searchEngine_.correct();
317  }
318 
319  forAll(*this, seti)
320  {
321  operator[](seti).topoChange(map);
322  }
323 
324  if (this->size())
325  {
326  combineSampledSets();
327  }
328  }
329 }
330 
331 
333 {
334  if (&map.mesh() == &mesh_)
335  {
336  if (this->size())
337  {
338  searchEngine_.correct();
339  }
340 
341  forAll(*this, seti)
342  {
343  operator[](seti).mapMesh(map);
344  }
345 
346  if (this->size())
347  {
348  combineSampledSets();
349  }
350  }
351 }
352 
353 
355 (
356  const polyDistributionMap& map
357 )
358 {
359  if (&map.mesh() == &mesh_)
360  {
361  if (this->size())
362  {
363  searchEngine_.correct();
364  }
365 
366  forAll(*this, seti)
367  {
368  operator[](seti).distribute(map);
369  }
370 
371  if (this->size())
372  {
373  combineSampledSets();
374  }
375  }
376 }
377 
378 
379 // ************************************************************************* //
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Macros for easy insertion into run-time selection tables.
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void setSize(const label)
Reset size of List.
Definition: List.C:281
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:85
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
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
const T & operator[](const label) const
Return element const reference.
Definition: UPtrListI.H:105
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Abstract base-class for Time/database functionObjects.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
sampledSets(const word &name, const Time &time, const dictionary &dict)
Construct from Time and dictionary.
Definition: sampledSets.C:88
virtual wordList fields() const
Return the list of fields required.
Definition: sampledSets.C:183
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: sampledSets.C:308
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: sampledSets.C:355
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: sampledSets.C:332
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: sampledSets.C:285
virtual ~sampledSets()
Destructor.
Definition: sampledSets.C:113
virtual bool execute()
Execute, currently does nothing.
Definition: sampledSets.C:189
virtual bool write()
Sample and write.
Definition: sampledSets.C:195
virtual bool read(const dictionary &)
Read the sampledSets.
Definition: sampledSets.C:119
Writes run time, CPU time and clock time and optionally the CPU and clock times per time step.
functionObject base class for writing single files
Definition: writeFile.H:56
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
const polyMesh & mesh() const
Return polyMesh.
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
const polyMesh & mesh() const
Return polyMesh.
Definition: polyMeshMap.H:75
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const polyMesh & mesh() const
Return polyMesh.
Class used for the read-construction of.
Definition: sampledSet.H:142
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:67
static autoPtr< sampledSet > New(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Return a reference to the selected sampledSet.
Definition: sampledSet.C:190
static autoPtr< setWriter > New(const word &writeType, const IOstream::streamFormat writeFormat=IOstream::ASCII, const IOstream::compressionType writeCompression=IOstream::UNCOMPRESSED)
Select given write options.
Definition: setWriter.C:286
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
static List< word > fieldNames
Definition: globalFoam.H:46
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define WarningInFunction
Report a warning using Foam::Warning.
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
static const char nl
Definition: Ostream.H:267
points setSize(newPointi)
#define FoundFieldType(Type, nullArg)
#define FieldTypeValuesParameter(Type, nullArg)
#define DeclareInterpolations(Type, nullArg)
#define GenerateFieldTypeValues(Type, nullArg)
dictionary dict