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-2021 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 "mapPolyMesh.H"
34 #include "writeFile.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  const sampledSet& samplePts = sampledSets[setI];
75 
76  // Collect data from all processors
77  List<List<point>> gatheredPts(Pstream::nProcs());
78  gatheredPts[Pstream::myProcNo()] = samplePts;
79  Pstream::gatherList(gatheredPts);
80 
81  List<labelList> gatheredSegments(Pstream::nProcs());
82  gatheredSegments[Pstream::myProcNo()] = samplePts.segments();
83  Pstream::gatherList(gatheredSegments);
84 
85  List<scalarList> gatheredDist(Pstream::nProcs());
86  gatheredDist[Pstream::myProcNo()] = samplePts.curveDist();
87  Pstream::gatherList(gatheredDist);
88 
89 
90  // Combine processor lists into one big list.
91  List<point> allPts
92  (
93  ListListOps::combine<List<point>>
94  (
95  gatheredPts, accessOp<List<point>>()
96  )
97  );
98  labelList allSegments
99  (
100  ListListOps::combine<labelList>
101  (
102  gatheredSegments, accessOp<labelList>()
103  )
104  );
105  scalarList allCurveDist
106  (
107  ListListOps::combine<scalarList>
108  (
109  gatheredDist, accessOp<scalarList>()
110  )
111  );
112 
113 
114  if (Pstream::master() && allCurveDist.size() == 0)
115  {
117  << "Sample set " << samplePts.name()
118  << " has zero points." << endl;
119  }
120 
121  // Sort curveDist and use to fill masterSamplePts
122  SortableList<scalar> sortedDist(allCurveDist);
123  indexSets[setI] = sortedDist.indices();
124 
125  masterSampledSets.set
126  (
127  setI,
128  new coordSet
129  (
130  samplePts.name(),
131  samplePts.axis(),
132  List<point>(UIndirectList<point>(allPts, indexSets[setI])),
133  scalarList(UIndirectList<scalar>(allCurveDist, indexSets[setI]))
134  )
135  );
136  }
137 }
138 
139 
140 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
141 
143 (
144  const word& name,
145  const Time& t,
146  const dictionary& dict
147 )
148 :
149  functionObject(name),
151  mesh_
152  (
153  refCast<const fvMesh>
154  (
156  (
158  )
159  )
160  ),
161  loadFromFiles_(false),
162  outputPath_(fileName::null),
163  searchEngine_(mesh_),
164  interpolationScheme_(word::null),
165  writeFormat_(word::null)
166 {
167  outputPath_ =
168  mesh_.time().globalPath()/functionObjects::writeFile::outputPrefix/name;
169 
170  if (mesh_.name() != fvMesh::defaultRegion)
171  {
172  outputPath_ = outputPath_/mesh_.name();
173  }
174 
175  read(dict);
176 }
177 
178 
180 (
181  const word& name,
182  const objectRegistry& obr,
183  const dictionary& dict,
184  const bool loadFromFiles
185 )
186 :
187  functionObject(name),
189  mesh_(refCast<const fvMesh>(obr)),
190  loadFromFiles_(loadFromFiles),
191  outputPath_(fileName::null),
192  searchEngine_(mesh_),
193  interpolationScheme_(word::null),
194  writeFormat_(word::null)
195 {
196  outputPath_ =
197  mesh_.time().globalPath()/functionObjects::writeFile::outputPrefix/name;
198 
199  if (mesh_.name() != fvMesh::defaultRegion)
200  {
201  outputPath_ = outputPath_/mesh_.name();
202  }
203 
204  read(dict);
205 }
206 
207 
208 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
209 
211 {}
212 
213 
214 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
215 
216 void Foam::sampledSets::verbose(const bool verbosity)
217 {
218  verbose_ = verbosity;
219 }
220 
221 
223 {
224  return true;
225 }
226 
227 
229 {
230  if (size())
231  {
232  const label nFields = classifyFields();
233 
234  if (Pstream::master())
235  {
236  if (debug)
237  {
238  Pout<< "timeName = " << mesh_.time().timeName() << nl
239  << "scalarFields " << scalarFields_ << nl
240  << "vectorFields " << vectorFields_ << nl
241  << "sphTensorFields " << sphericalTensorFields_ << nl
242  << "symTensorFields " << symmTensorFields_ <<nl
243  << "tensorFields " << tensorFields_ <<nl;
244  }
245 
246  if (nFields)
247  {
248  if (debug)
249  {
250  Pout<< "Creating directory "
251  << outputPath_/mesh_.time().timeName()
252  << nl << endl;
253  }
254 
255  mkDir(outputPath_/mesh_.time().timeName());
256  }
257  }
258 
259  if (nFields)
260  {
261  sampleAndWrite(scalarFields_);
262  sampleAndWrite(vectorFields_);
263  sampleAndWrite(sphericalTensorFields_);
264  sampleAndWrite(symmTensorFields_);
265  sampleAndWrite(tensorFields_);
266  }
267  }
268 
269  return true;
270 }
271 
272 
274 {
275  dict_ = dict;
276 
277  bool setsFound = dict_.found("sets");
278  if (setsFound)
279  {
280  dict_.lookup("fields") >> fieldSelection_;
281  clearFieldGroups();
282 
283  dict.lookup("interpolationScheme") >> interpolationScheme_;
284  dict.lookup("setFormat") >> writeFormat_;
285 
286  PtrList<sampledSet> newList
287  (
288  dict_.lookup("sets"),
289  sampledSet::iNew(mesh_, searchEngine_)
290  );
291  transfer(newList);
292  combineSampledSets(masterSampledSets_, indexSets_);
293 
294  if (this->size())
295  {
296  Info<< "Reading set description:" << nl;
297  forAll(*this, setI)
298  {
299  Info<< " " << operator[](setI).name() << nl;
300  }
301  Info<< endl;
302  }
303  }
304 
305  if (Pstream::master() && debug)
306  {
307  Pout<< "sample fields:" << fieldSelection_ << nl
308  << "sample sets:" << nl << "(" << nl;
309 
310  forAll(*this, setI)
311  {
312  Pout<< " " << operator[](setI) << endl;
313  }
314  Pout<< ")" << endl;
315  }
316 
317  return true;
318 }
319 
320 
322 {
323  bool setsFound = dict_.found("sets");
324  if (setsFound)
325  {
326  searchEngine_.correct();
327 
328  PtrList<sampledSet> newList
329  (
330  dict_.lookup("sets"),
331  sampledSet::iNew(mesh_, searchEngine_)
332  );
333  transfer(newList);
334  combineSampledSets(masterSampledSets_, indexSets_);
335  }
336 }
337 
338 
340 {
341  if (&mpm.mesh() == &mesh_)
342  {
343  correct();
344  }
345 }
346 
347 
349 {
350  if (&mesh == &mesh_)
351  {
352  correct();
353  }
354 }
355 
356 
358 {
359  if (state != polyMesh::UNCHANGED)
360  {
361  correct();
362  }
363 }
364 
365 
366 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: sampledSets.C:348
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
static const fileName null
An empty fileName.
Definition: fileName.H:97
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
virtual bool read(const dictionary &)
Read the sampledSets.
Definition: sampledSets.C:273
Abstract base-class for Time/database functionObjects.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
AccessType combine(const List< T > &, AccessOp aop=accessOp< T >())
Combines sublists into one big list.
Definition: ListListOps.C:34
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Macros for easy insertion into run-time selection tables.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dynamicFvMesh & mesh
static const word outputPrefix
Directory prefix.
Definition: writeFile.H:72
A class for handling words, derived from string.
Definition: word.H:59
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-9/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
static const word null
An empty word.
Definition: word.H:77
void correct()
Correct for mesh changes.
Definition: sampledSets.C:321
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: sampledSets.C:339
List< label > labelList
A List of labels.
Definition: labelList.H:56
virtual bool execute()
Execute, currently does nothing.
Definition: sampledSets.C:222
Class used for the read-construction of.
Definition: sampledSet.H:126
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:260
sampledSets(const word &name, const Time &time, const dictionary &dict)
Construct from Time and dictionary.
Definition: sampledSets.C:143
defineTypeNameAndDebug(combustionModel, 0)
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
void verbose(const bool verbosity=true)
Set verbosity level.
Definition: sampledSets.C:216
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void readUpdate(const polyMesh::readUpdateState state)
Update for changes of mesh due to readUpdate.
Definition: sampledSets.C:357
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
virtual bool write()
Sample and write.
Definition: sampledSets.C:228
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:352
messageStream Info
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:88
Registry of regIOobjects.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
virtual ~sampledSets()
Destructor.
Definition: sampledSets.C:210
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844