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-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 "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 
140  forAll(masterFields, fieldi)
141  {
142  valueSetNames[fieldi] = masterFields[fieldi].name();
143  valueSets[fieldi] = &masterFields[fieldi][setI];
144  }
145 
146  fileName fName
147  (
148  timeDir/formatter.getFileName(masterSampleSet, valueSetNames)
149  );
150 
151  OFstream ofs(fName);
152  if (ofs.opened())
153  {
154  formatter.write
155  (
156  masterSampleSet,
157  valueSetNames,
158  valueSets,
159  ofs
160  );
161  return fName;
162  }
163  else
164  {
166  << "File " << ofs.name() << " could not be opened. "
167  << "No data will be written" << endl;
168  return fileName::null;
169  }
170 }
171 
172 
173 template<class T>
174 void Foam::sampledSets::combineSampledValues
175 (
176  const PtrList<volFieldSampler<T>>& sampledFields,
177  const labelListList& indexSets,
178  PtrList<volFieldSampler<T>>& masterFields
179 )
180 {
181  forAll(sampledFields, fieldi)
182  {
183  List<Field<T>> masterValues(indexSets.size());
184 
185  forAll(indexSets, setI)
186  {
187  // Collect data from all processors
188 
189  Field<T> allData;
190  globalIndex::gatherOp(sampledFields[fieldi][setI], allData);
191 
192  if (Pstream::master())
193  {
194  masterValues[setI] = UIndirectList<T>
195  (
196  allData,
197  indexSets[setI]
198  )();
199  }
200  }
201 
202  masterFields.set
203  (
204  fieldi,
205  new volFieldSampler<T>
206  (
207  masterValues,
208  sampledFields[fieldi].name()
209  )
210  );
211  }
212 }
213 
214 
215 template<class Type>
216 void Foam::sampledSets::sampleAndWrite(fieldGroup<Type>& fields)
217 {
218  if (fields.size())
219  {
220  const bool interpolate = interpolationScheme_ != "cell";
221 
222  // Create or use existing writer
223  if (!fields.formatter)
224  {
225  fields = writeFormat_;
226  }
227 
228  // Storage for interpolated values
229  PtrList<volFieldSampler<Type>> sampledFields(fields.size());
230 
231  forAll(fields, fieldi)
232  {
233  if (Pstream::master() && verbose_)
234  {
235  Pout<< "sampledSets::sampleAndWrite: "
236  << fields[fieldi] << endl;
237  }
238 
239  if (loadFromFiles_)
240  {
241  GeometricField<Type, fvPatchField, volMesh> vf
242  (
243  IOobject
244  (
245  fields[fieldi],
246  mesh_.time().timeName(),
247  mesh_,
250  false
251  ),
252  mesh_
253  );
254 
255  if (interpolate)
256  {
257  sampledFields.set
258  (
259  fieldi,
260  new volFieldSampler<Type>
261  (
262  interpolationScheme_,
263  vf,
264  *this
265  )
266  );
267  }
268  else
269  {
270  sampledFields.set
271  (
272  fieldi,
273  new volFieldSampler<Type>(vf, *this)
274  );
275  }
276  }
277  else
278  {
279  if (interpolate)
280  {
281  sampledFields.set
282  (
283  fieldi,
284  new volFieldSampler<Type>
285  (
286  interpolationScheme_,
287  mesh_.lookupObject
288  <GeometricField<Type, fvPatchField, volMesh>>
289  (fields[fieldi]),
290  *this
291  )
292  );
293  }
294  else
295  {
296  sampledFields.set
297  (
298  fieldi,
299  new volFieldSampler<Type>
300  (
301  mesh_.lookupObject
302  <GeometricField<Type, fvPatchField, volMesh>>
303  (fields[fieldi]),
304  *this
305  )
306  );
307  }
308  }
309  }
310 
311  // Combine sampled fields from processors.
312  // Note: only master results are valid
313 
314  PtrList<volFieldSampler<Type>> masterFields(sampledFields.size());
315  combineSampledValues(sampledFields, indexSets_, masterFields);
316 
317  forAll(masterSampledSets_, setI)
318  {
319  fileName sampleFile;
320  if (Pstream::master())
321  {
322  sampleFile = writeSampleFile
323  (
324  masterSampledSets_[setI],
325  masterFields,
326  setI,
327  outputPath_/mesh_.time().timeName(),
328  fields.formatter()
329  );
330  }
331 
332  Pstream::scatter(sampleFile);
333  if (sampleFile.size())
334  {
335  // Case-local file name with "<case>" to make relocatable
336 
337  forAll(masterFields, fieldi)
338  {
339  dictionary propsDict;
340  propsDict.add
341  (
342  "file",
343  time_.relativePath(sampleFile, true)
344  );
345 
346  const word& fieldName = masterFields[fieldi].name();
347  setProperty(fieldName, propsDict);
348  }
349  }
350  }
351  }
352 }
353 
354 
355 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
volFields.H
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
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:458
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:350
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
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::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
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:62
Foam::globalIndex::gatherOp
static void gatherOp(const UList< Type > &fld, List< Type > &allFld, const int tag=UPstream::msgType(), const Pstream::commsTypes commsType=Pstream::commsTypes::nonBlocking)
Collect data in processor order on master.
Definition: globalIndexTemplates.C:208
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::writer
Base class for graphics format writing. Entry points are.
Definition: writer.H:80
Foam::coordSet
Holds list of sampling positions.
Definition: coordSet.H:53
Foam::entry::name
virtual const fileName & name() const =0
Return the dictionary name.
Foam::fileName::null
static const fileName null
An empty fileName.
Definition: fileName.H:97
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::List< word >
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:708
Foam::point
vector point
Point is a vector.
Definition: point.H:43
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
Foam::IOobject::MUST_READ
Definition: IOobject.H:120