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-2018 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"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(sampledSets, 0);
41 
43  (
44  functionObject,
45  sampledSets,
46  dictionary
47  );
48 }
49 
50 bool Foam::sampledSets::verbose_ = false;
51 
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
55 void Foam::sampledSets::combineSampledSets
56 (
57  PtrList<coordSet>& masterSampledSets,
58  labelListList& indexSets
59 )
60 {
61  // Combine sampleSets from processors. Sort by curveDist. Return
62  // ordering in indexSets.
63  // Note: only master results are valid
64 
65  masterSampledSets_.clear();
66  masterSampledSets_.setSize(size());
67  indexSets_.setSize(size());
68 
69  const PtrList<sampledSet>& sampledSets = *this;
70 
71  forAll(sampledSets, setI)
72  {
73  const sampledSet& samplePts = sampledSets[setI];
74 
75  // Collect data from all processors
76  List<List<point>> gatheredPts(Pstream::nProcs());
77  gatheredPts[Pstream::myProcNo()] = samplePts;
78  Pstream::gatherList(gatheredPts);
79 
80  List<labelList> gatheredSegments(Pstream::nProcs());
81  gatheredSegments[Pstream::myProcNo()] = samplePts.segments();
82  Pstream::gatherList(gatheredSegments);
83 
84  List<scalarList> gatheredDist(Pstream::nProcs());
85  gatheredDist[Pstream::myProcNo()] = samplePts.curveDist();
86  Pstream::gatherList(gatheredDist);
87 
88 
89  // Combine processor lists into one big list.
90  List<point> allPts
91  (
92  ListListOps::combine<List<point>>
93  (
94  gatheredPts, accessOp<List<point>>()
95  )
96  );
97  labelList allSegments
98  (
99  ListListOps::combine<labelList>
100  (
101  gatheredSegments, accessOp<labelList>()
102  )
103  );
104  scalarList allCurveDist
105  (
106  ListListOps::combine<scalarList>
107  (
108  gatheredDist, accessOp<scalarList>()
109  )
110  );
111 
112 
113  if (Pstream::master() && allCurveDist.size() == 0)
114  {
116  << "Sample set " << samplePts.name()
117  << " has zero points." << endl;
118  }
119 
120  // Sort curveDist and use to fill masterSamplePts
121  SortableList<scalar> sortedDist(allCurveDist);
122  indexSets[setI] = sortedDist.indices();
123 
124  masterSampledSets.set
125  (
126  setI,
127  new coordSet
128  (
129  samplePts.name(),
130  samplePts.axis(),
131  List<point>(UIndirectList<point>(allPts, indexSets[setI])),
132  scalarList(UIndirectList<scalar>(allCurveDist, indexSets[setI]))
133  )
134  );
135  }
136 }
137 
138 
139 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
140 
142 (
143  const word& name,
144  const Time& t,
145  const dictionary& dict
146 )
147 :
148  functionObject(name),
150  mesh_
151  (
152  refCast<const fvMesh>
153  (
155  (
157  )
158  )
159  ),
160  loadFromFiles_(false),
161  outputPath_(fileName::null),
162  searchEngine_(mesh_),
163  interpolationScheme_(word::null),
164  writeFormat_(word::null)
165 {
166  if (Pstream::parRun())
167  {
168  outputPath_ = mesh_.time().path()/".."/"postProcessing"/name;
169  }
170  else
171  {
172  outputPath_ = mesh_.time().path()/"postProcessing"/name;
173  }
174  if (mesh_.name() != fvMesh::defaultRegion)
175  {
176  outputPath_ = outputPath_/mesh_.name();
177  }
178  // Remove ".."
179  outputPath_.clean();
180 
181  read(dict);
182 }
183 
184 
186 (
187  const word& name,
188  const objectRegistry& obr,
189  const dictionary& dict,
190  const bool loadFromFiles
191 )
192 :
193  functionObject(name),
195  mesh_(refCast<const fvMesh>(obr)),
196  loadFromFiles_(loadFromFiles),
197  outputPath_(fileName::null),
198  searchEngine_(mesh_),
199  interpolationScheme_(word::null),
200  writeFormat_(word::null)
201 {
202  if (Pstream::parRun())
203  {
204  outputPath_ = mesh_.time().path()/".."/"postProcessing"/name;
205  }
206  else
207  {
208  outputPath_ = mesh_.time().path()/"postProcessing"/name;
209  }
210  if (mesh_.name() != fvMesh::defaultRegion)
211  {
212  outputPath_ = outputPath_/mesh_.name();
213  }
214  // Remove ".."
215  outputPath_.clean();
216 
217  read(dict);
218 }
219 
220 
221 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
222 
224 {}
225 
226 
227 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
228 
229 void Foam::sampledSets::verbose(const bool verbosity)
230 {
231  verbose_ = verbosity;
232 }
233 
234 
236 {
237  return true;
238 }
239 
240 
242 {
243  if (size())
244  {
245  const label nFields = classifyFields();
246 
247  if (Pstream::master())
248  {
249  if (debug)
250  {
251  Pout<< "timeName = " << mesh_.time().timeName() << nl
252  << "scalarFields " << scalarFields_ << nl
253  << "vectorFields " << vectorFields_ << nl
254  << "sphTensorFields " << sphericalTensorFields_ << nl
255  << "symTensorFields " << symmTensorFields_ <<nl
256  << "tensorFields " << tensorFields_ <<nl;
257  }
258 
259  if (nFields)
260  {
261  if (debug)
262  {
263  Pout<< "Creating directory "
264  << outputPath_/mesh_.time().timeName()
265  << nl << endl;
266  }
267 
268  mkDir(outputPath_/mesh_.time().timeName());
269  }
270  }
271 
272  if (nFields)
273  {
274  sampleAndWrite(scalarFields_);
275  sampleAndWrite(vectorFields_);
276  sampleAndWrite(sphericalTensorFields_);
277  sampleAndWrite(symmTensorFields_);
278  sampleAndWrite(tensorFields_);
279  }
280  }
281 
282  return true;
283 }
284 
285 
287 {
288  dict_ = dict;
289 
290  bool setsFound = dict_.found("sets");
291  if (setsFound)
292  {
293  dict_.lookup("fields") >> fieldSelection_;
294  clearFieldGroups();
295 
296  dict.lookup("interpolationScheme") >> interpolationScheme_;
297  dict.lookup("setFormat") >> writeFormat_;
298 
299  PtrList<sampledSet> newList
300  (
301  dict_.lookup("sets"),
302  sampledSet::iNew(mesh_, searchEngine_)
303  );
304  transfer(newList);
305  combineSampledSets(masterSampledSets_, indexSets_);
306 
307  if (this->size())
308  {
309  Info<< "Reading set description:" << nl;
310  forAll(*this, setI)
311  {
312  Info<< " " << operator[](setI).name() << nl;
313  }
314  Info<< endl;
315  }
316  }
317 
318  if (Pstream::master() && debug)
319  {
320  Pout<< "sample fields:" << fieldSelection_ << nl
321  << "sample sets:" << nl << "(" << nl;
322 
323  forAll(*this, setI)
324  {
325  Pout<< " " << operator[](setI) << endl;
326  }
327  Pout<< ")" << endl;
328  }
329 
330  return true;
331 }
332 
333 
335 {
336  bool setsFound = dict_.found("sets");
337  if (setsFound)
338  {
339  searchEngine_.correct();
340 
341  PtrList<sampledSet> newList
342  (
343  dict_.lookup("sets"),
344  sampledSet::iNew(mesh_, searchEngine_)
345  );
346  transfer(newList);
347  combineSampledSets(masterSampledSets_, indexSets_);
348  }
349 }
350 
351 
353 {
354  if (&mpm.mesh() == &mesh_)
355  {
356  correct();
357  }
358 }
359 
360 
362 {
363  if (&mesh == &mesh_)
364  {
365  correct();
366  }
367 }
368 
369 
371 {
372  if (state != polyMesh::UNCHANGED)
373  {
374  correct();
375  }
376 }
377 
378 
379 // ************************************************************************* //
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:438
#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:361
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
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:256
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:286
Abstract base-class for Time/database function objects.
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
A class for handling words, derived from string.
Definition: word.H:59
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:334
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: sampledSets.C:352
List< label > labelList
A List of labels.
Definition: labelList.H:56
virtual bool execute()
Execute, currently does nothing.
Definition: sampledSets.C:235
Class used for the read-construction of.
Definition: sampledSet.H:126
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:265
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-7/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
sampledSets(const word &name, const Time &time, const dictionary &dict)
Construct from Time and dictionary.
Definition: sampledSets.C:142
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:229
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:370
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
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:241
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:223
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583