sampledSetsTemplates.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-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2021 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 "volFields.H"
31 #include "globalIndex.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 template<class Type>
36 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
37 (
38  const word& interpolationScheme,
39  const GeometricField<Type, fvPatchField, volMesh>& field,
40  const PtrList<sampledSet>& samplers
41 )
42 :
43  List<Field<Type>>(samplers.size()),
44  name_(field.name())
45 {
46  autoPtr<interpolation<Type>> interpolator
47  (
48  interpolation<Type>::New(interpolationScheme, field)
49  );
50 
51  forAll(samplers, setI)
52  {
53  Field<Type>& values = this->operator[](setI);
54  const sampledSet& samples = samplers[setI];
55 
56  values.setSize(samples.size());
57  forAll(samples, sampleI)
58  {
59  const point& samplePt = samples[sampleI];
60  label celli = samples.cells()[sampleI];
61  label facei = samples.faces()[sampleI];
62 
63  if (celli == -1 && facei == -1)
64  {
65  // Special condition for illegal sampling points
66  values[sampleI] = pTraits<Type>::max;
67  }
68  else
69  {
70  values[sampleI] = interpolator().interpolate
71  (
72  samplePt,
73  celli,
74  facei
75  );
76  }
77  }
78  }
79 }
80 
81 
82 template<class Type>
83 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
84 (
85  const GeometricField<Type, fvPatchField, volMesh>& field,
86  const PtrList<sampledSet>& samplers
87 )
88 :
89  List<Field<Type>>(samplers.size()),
90  name_(field.name())
91 {
92  forAll(samplers, setI)
93  {
94  Field<Type>& values = this->operator[](setI);
95  const sampledSet& samples = samplers[setI];
96 
97  values.setSize(samples.size());
98  forAll(samples, sampleI)
99  {
100  label celli = samples.cells()[sampleI];
101 
102  if (celli ==-1)
103  {
104  values[sampleI] = pTraits<Type>::max;
105  }
106  else
107  {
108  values[sampleI] = field[celli];
109  }
110  }
111  }
112 }
113 
114 
115 template<class Type>
116 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
117 (
118  const List<Field<Type>>& values,
119  const word& name
120 )
121 :
122  List<Field<Type>>(values),
123  name_(name)
124 {}
125 
126 
127 template<class Type>
128 Foam::fileName Foam::sampledSets::writeSampleFile
129 (
130  const coordSet& masterSampleSet,
131  const PtrList<volFieldSampler<Type>>& masterFields,
132  const label setI,
133  const fileName& timeDir,
134  const writer<Type>& formatter
135 )
136 {
137  wordList valueSetNames(masterFields.size());
138  List<const Field<Type>*> valueSets(masterFields.size());
139  forAll(masterFields, fieldi)
140  {
141  const word& fieldName = masterFields[fieldi].name();
142 
143  valueSetNames[fieldi] = fieldName;
144 
145  // Values only available on master
146  Type averageValue, minValue, maxValue;
147  label sizeValue;
148  if (Pstream::master())
149  {
150  valueSets[fieldi] = &masterFields[fieldi][setI];
151  averageValue = average(*valueSets[fieldi]);
152  minValue = min(*valueSets[fieldi]);
153  maxValue = max(*valueSets[fieldi]);
154  sizeValue = valueSets[fieldi]->size();
155  }
156  Pstream::scatter(averageValue);
159  Pstream::scatter(sizeValue);
160 
161  // Set results
162 
163  setResult("average(" + fieldName + ")", averageValue);
164  setResult("min(" + fieldName + ")", minValue);
165  setResult("max(" + fieldName + ")", maxValue);
166  setResult("size(" + fieldName + ")", sizeValue);
167  }
168 
169  fileName fName;
170  if (Pstream::master())
171  {
172  fName = timeDir/formatter.getFileName(masterSampleSet, valueSetNames);
173 
174  OFstream ofs(fName);
175  if (ofs.opened())
176  {
177  formatter.write
178  (
179  masterSampleSet,
180  valueSetNames,
181  valueSets,
182  ofs
183  );
184  }
185  else
186  {
188  << "File " << ofs.name() << " could not be opened. "
189  << "No data will be written" << endl;
190  }
191  }
192 
193  Pstream::scatter(fName);
194 
195  return fName;
196 }
197 
198 
199 template<class T>
200 void Foam::sampledSets::combineSampledValues
201 (
202  const PtrList<volFieldSampler<T>>& sampledFields,
203  const labelListList& indexSets,
204  PtrList<volFieldSampler<T>>& masterFields
205 )
206 {
207  forAll(sampledFields, fieldi)
208  {
209  List<Field<T>> masterValues(indexSets.size());
210 
211  forAll(indexSets, setI)
212  {
213  // Collect data from all processors
214 
215  Field<T> allData;
216  globalIndex::gatherOp(sampledFields[fieldi][setI], allData);
217 
218  if (Pstream::master())
219  {
220  masterValues[setI] = UIndirectList<T>
221  (
222  allData,
223  indexSets[setI]
224  )();
225  }
226  }
227 
228  masterFields.set
229  (
230  fieldi,
231  new volFieldSampler<T>
232  (
233  masterValues,
234  sampledFields[fieldi].name()
235  )
236  );
237  }
238 }
239 
240 
241 template<class Type>
242 void Foam::sampledSets::sampleAndWrite(fieldGroup<Type>& fields)
243 {
244  if (fields.size())
245  {
246  const bool interpolate = interpolationScheme_ != "cell";
247 
248  // Create or use existing writer
249  if (!fields.formatter)
250  {
251  fields.setFormatter(writeFormat_, writeFormatOptions_);
252  }
253 
254  // Storage for interpolated values
255  PtrList<volFieldSampler<Type>> sampledFields(fields.size());
256 
257  forAll(fields, fieldi)
258  {
259  if (Pstream::master() && verbose_)
260  {
261  Pout<< "sampledSets::sampleAndWrite: "
262  << fields[fieldi] << endl;
263  }
264 
265  if (loadFromFiles_)
266  {
267  GeometricField<Type, fvPatchField, volMesh> vf
268  (
269  IOobject
270  (
271  fields[fieldi],
272  mesh_.time().timeName(),
273  mesh_,
276  false
277  ),
278  mesh_
279  );
280 
281  if (interpolate)
282  {
283  sampledFields.set
284  (
285  fieldi,
286  new volFieldSampler<Type>
287  (
288  interpolationScheme_,
289  vf,
290  *this
291  )
292  );
293  }
294  else
295  {
296  sampledFields.set
297  (
298  fieldi,
299  new volFieldSampler<Type>(vf, *this)
300  );
301  }
302  }
303  else
304  {
305  if (interpolate)
306  {
307  sampledFields.set
308  (
309  fieldi,
310  new volFieldSampler<Type>
311  (
312  interpolationScheme_,
313  mesh_.lookupObject
314  <GeometricField<Type, fvPatchField, volMesh>>
315  (fields[fieldi]),
316  *this
317  )
318  );
319  }
320  else
321  {
322  sampledFields.set
323  (
324  fieldi,
325  new volFieldSampler<Type>
326  (
327  mesh_.lookupObject
328  <GeometricField<Type, fvPatchField, volMesh>>
329  (fields[fieldi]),
330  *this
331  )
332  );
333  }
334  }
335  }
336 
337  // Combine sampled fields from processors.
338  // Note: only master results are valid
339 
340  PtrList<volFieldSampler<Type>> masterFields(sampledFields.size());
341  combineSampledValues(sampledFields, indexSets_, masterFields);
342 
343  forAll(masterSampledSets_, setI)
344  {
345  fileName sampleFile = writeSampleFile
346  (
347  masterSampledSets_[setI],
348  masterFields,
349  setI,
350  outputPath_/mesh_.time().timeName(),
351  fields.formatter()
352  );
353 
354  if (sampleFile.size())
355  {
356  // Case-local file name with "<case>" to make relocatable
357 
358  forAll(masterFields, fieldi)
359  {
360  dictionary propsDict;
361  propsDict.add
362  (
363  "file",
364  time_.relativePath(sampleFile, true)
365  );
366 
367  const word& fieldName = masterFields[fieldi].name();
368  setProperty(fieldName, propsDict);
369  }
370  }
371  }
372  }
373 }
374 
375 
376 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
volFields.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::OFstream::name
virtual const fileName & name() const
Read/write access to the name of the stream.
Definition: OSstream.H:107
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
globalIndex.H
Foam::writer::write
virtual void write(const coordSet &, const wordList &, const List< const Field< Type > * > &, Ostream &) const =0
General entry point for writing.
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::Pstream::scatter
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
Definition: gatherScatter.C:150
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
sampledSets.H
Foam::writer::getFileName
virtual fileName getFileName(const coordSet &, const wordList &) const =0
Generate file name with correct extension.
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Foam::Field< T >
Foam::interpolation::New
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
Definition: interpolationNew.C:36
Foam::IOstream::opened
bool opened() const noexcept
True if stream has been opened.
Definition: IOstream.H:221
propsDict
IOdictionary propsDict(IOobject("particleTrackProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED))
samples
scalarField samples(nIntervals, Zero)
field
rDeltaTY field()
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::globalIndex::gatherOp
static void gatherOp(const UList< Type > &fld, List< Type > &allFld, const int tag=UPstream::msgType(), const Pstream::commsTypes=Pstream::commsTypes::nonBlocking)
Collect data in processor order on master.
Definition: globalIndexTemplates.C:472
Foam::writer
Base class for graphics format writing. Entry points are.
Definition: writer.H:81
Foam::coordSet
Holds list of sampling positions.
Definition: coordSet.H:53
Foam::entry::name
virtual const fileName & name() const =0
Return the entry name.
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::List< word >
minValue
scalar minValue
Definition: LISASMDCalcMethod2.H:12
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::dictionary::add
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:640
Foam::point
vector point
Point is a vector.
Definition: point.H:43
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
Foam::average
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:328
maxValue
scalar maxValue
Definition: LISASMDCalcMethod1.H:5
Foam::IOobject::MUST_READ
Definition: IOobject.H:185