All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
sampledSetsTemplates.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 "volFields.H"
28 #include "ListListOps.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
32 template<class Type>
33 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
34 (
35  const word& interpolationScheme,
36  const GeometricField<Type, fvPatchField, volMesh>& field,
37  const PtrList<sampledSet>& samplers
38 )
39 :
40  List<Field<Type>>(samplers.size()),
41  name_(field.name())
42 {
43  autoPtr<interpolation<Type>> interpolator
44  (
45  interpolation<Type>::New(interpolationScheme, field)
46  );
47 
48  forAll(samplers, setI)
49  {
50  Field<Type>& values = this->operator[](setI);
51  const sampledSet& samples = samplers[setI];
52 
53  values.setSize(samples.size());
54  forAll(samples, sampleI)
55  {
56  const point& samplePt = samples[sampleI];
57  label celli = samples.cells()[sampleI];
58  label facei = samples.faces()[sampleI];
59 
60  if (celli == -1 && facei == -1)
61  {
62  // Special condition for illegal sampling points
63  values[sampleI] = pTraits<Type>::max;
64  }
65  else
66  {
67  values[sampleI] = interpolator().interpolate
68  (
69  samplePt,
70  celli,
71  facei
72  );
73  }
74  }
75  }
76 }
77 
78 
79 template<class Type>
80 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
81 (
82  const GeometricField<Type, fvPatchField, volMesh>& field,
83  const PtrList<sampledSet>& samplers
84 )
85 :
86  List<Field<Type>>(samplers.size()),
87  name_(field.name())
88 {
89  forAll(samplers, setI)
90  {
91  Field<Type>& values = this->operator[](setI);
92  const sampledSet& samples = samplers[setI];
93 
94  values.setSize(samples.size());
95  forAll(samples, sampleI)
96  {
97  label celli = samples.cells()[sampleI];
98 
99  if (celli ==-1)
100  {
101  values[sampleI] = pTraits<Type>::max;
102  }
103  else
104  {
105  values[sampleI] = field[celli];
106  }
107  }
108  }
109 }
110 
111 
112 template<class Type>
113 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
114 (
115  const List<Field<Type>>& values,
116  const word& name
117 )
118 :
119  List<Field<Type>>(values),
120  name_(name)
121 {}
122 
123 
124 template<class Type>
125 void Foam::sampledSets::writeSampleFile
126 (
127  const coordSet& masterSampleSet,
128  const PtrList<volFieldSampler<Type>>& masterFields,
129  const label setI,
130  const fileName& timeDir,
131  const setWriter<Type>& formatter
132 )
133 {
134  wordList valueSetNames(masterFields.size());
135  List<const Field<Type>*> valueSets(masterFields.size());
136 
137  forAll(masterFields, fieldi)
138  {
139  valueSetNames[fieldi] = masterFields[fieldi].name();
140  valueSets[fieldi] = &masterFields[fieldi][setI];
141  }
142 
143  fileName fName
144  (
145  timeDir/formatter.getFileName(masterSampleSet, valueSetNames)
146  );
147 
148  OFstream ofs(fName);
149  if (ofs.opened())
150  {
151  formatter.write
152  (
153  masterSampleSet,
154  valueSetNames,
155  valueSets,
156  ofs
157  );
158  }
159  else
160  {
162  << "File " << ofs.name() << " could not be opened. "
163  << "No data will be written" << endl;
164  }
165 }
166 
167 
168 template<class T>
169 void Foam::sampledSets::combineSampledValues
170 (
171  const PtrList<volFieldSampler<T>>& sampledFields,
172  const labelListList& indexSets,
173  PtrList<volFieldSampler<T>>& masterFields
174 )
175 {
176  forAll(sampledFields, fieldi)
177  {
178  List<Field<T>> masterValues(indexSets.size());
179 
180  forAll(indexSets, setI)
181  {
182  // Collect data from all processors
183  List<Field<T>> gatheredData(Pstream::nProcs());
184  gatheredData[Pstream::myProcNo()] = sampledFields[fieldi][setI];
185  Pstream::gatherList(gatheredData);
186 
187  if (Pstream::master())
188  {
189  Field<T> allData
190  (
191  ListListOps::combine<Field<T>>
192  (
193  gatheredData,
194  Foam::accessOp<Field<T>>()
195  )
196  );
197 
198  masterValues[setI] = UIndirectList<T>
199  (
200  allData,
201  indexSets[setI]
202  )();
203  }
204  }
205 
206  masterFields.set
207  (
208  fieldi,
209  new volFieldSampler<T>
210  (
211  masterValues,
212  sampledFields[fieldi].name()
213  )
214  );
215  }
216 }
217 
218 
219 template<class Type>
220 void Foam::sampledSets::sampleAndWrite
221 (
222  fieldGroup<Type>& fields
223 )
224 {
225  if (fields.size())
226  {
227  bool interpolate = interpolationScheme_ != "cell";
228 
229  // Create or use existing writer
230  if (fields.formatter.empty())
231  {
232  fields = writeFormat_;
233  }
234 
235  // Storage for interpolated values
236  PtrList<volFieldSampler<Type>> sampledFields(fields.size());
237 
238  forAll(fields, fieldi)
239  {
240  if (Pstream::master() && verbose_)
241  {
242  Pout<< "sampledSets::sampleAndWrite: "
243  << fields[fieldi] << endl;
244  }
245 
246  if (loadFromFiles_)
247  {
248  GeometricField<Type, fvPatchField, volMesh> vf
249  (
250  IOobject
251  (
252  fields[fieldi],
253  mesh_.time().timeName(),
254  mesh_,
257  false
258  ),
259  mesh_
260  );
261 
262  if (interpolate)
263  {
264  sampledFields.set
265  (
266  fieldi,
267  new volFieldSampler<Type>
268  (
269  interpolationScheme_,
270  vf,
271  *this
272  )
273  );
274  }
275  else
276  {
277  sampledFields.set
278  (
279  fieldi,
280  new volFieldSampler<Type>(vf, *this)
281  );
282  }
283  }
284  else
285  {
286  if (interpolate)
287  {
288  sampledFields.set
289  (
290  fieldi,
291  new volFieldSampler<Type>
292  (
293  interpolationScheme_,
294  mesh_.lookupObject
295  <GeometricField<Type, fvPatchField, volMesh>>
296  (fields[fieldi]),
297  *this
298  )
299  );
300  }
301  else
302  {
303  sampledFields.set
304  (
305  fieldi,
306  new volFieldSampler<Type>
307  (
308  mesh_.lookupObject
309  <GeometricField<Type, fvPatchField, volMesh>>
310  (fields[fieldi]),
311  *this
312  )
313  );
314  }
315  }
316  }
317 
318  // Combine sampled fields from processors.
319  // Note: only master results are valid
320 
321  PtrList<volFieldSampler<Type>> masterFields(sampledFields.size());
322  combineSampledValues(sampledFields, indexSets_, masterFields);
323 
324  if (Pstream::master())
325  {
326  forAll(masterSampledSets_, setI)
327  {
328  writeSampleFile
329  (
330  masterSampledSets_[setI],
331  masterFields,
332  setI,
333  outputPath_/mesh_.time().timeName(),
334  fields.formatter()
335  );
336  }
337  }
338  }
339 }
340 
341 
342 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#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
const word & name() const
Return the name of this functionObject.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
Definition: interpolation.C:48
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
const T & operator[](const label) const
Return element const reference.
Definition: UPtrListI.H:96
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
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
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
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
List< word > wordList
A List of words.
Definition: fileName.H:54
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
vector point
Point is a vector.
Definition: point.H:41
#define WarningInFunction
Report a warning using Foam::Warning.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
PtrList()
Null Constructor.
Definition: PtrList.C:32
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.