1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "sampledSets.H"
28 #include "volFields.H"
29 #include "ListListOps.H"
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
36 const word& interpolationScheme,
37 const GeometricField<Type, fvPatchField, volMesh>& field,
38 const PtrList<sampledSet>& samplers
41 List<Field<Type> >(samplers.size()),
44 autoPtr<interpolation<Type> > interpolator
46 interpolation<Type>::New(interpolationScheme, field)
49 forAll(samplers, seti)
51 Field<Type>& values = this->operator[](seti);
52 const sampledSet& samples = samplers[seti];
54 values.setSize(samples.size());
55 forAll(samples, samplei)
57 const point& samplePt = samples[samplei];
58 label celli = samples.cells()[samplei];
59 label facei = samples.faces()[samplei];
61 values[samplei] = interpolator().interpolate
73 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
75 const GeometricField<Type, fvPatchField, volMesh>& field,
76 const PtrList<sampledSet>& samplers
79 List<Field<Type> >(samplers.size()),
82 forAll(samplers, seti)
84 Field<Type>& values = this->operator[](seti);
85 const sampledSet& samples = samplers[seti];
87 values.setSize(samples.size());
88 forAll(samples, samplei)
90 values[samplei] = field[samples.cells()[samplei]];
97 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
99 const List<Field<Type> >& values,
103 List<Field<Type> >(values),
109 Foam::label Foam::sampledSets::grep
111 fieldGroup<Type>& fieldList,
112 const wordList& fieldTypes
115 fieldList.setSize(fieldNames_.size());
118 forAll(fieldNames_, fieldi)
123 == GeometricField<Type, fvPatchField, volMesh>::typeName
126 fieldList[nFields] = fieldNames_[fieldi];
131 fieldList.setSize(nFields);
138 void Foam::sampledSets::writeSampleFile
140 const coordSet& masterSampleSet,
141 const PtrList<volFieldSampler<Type> >& masterFields,
143 const fileName& timeDir,
144 const writer<Type>& formatter
147 wordList valueSetNames(masterFields.size());
148 List<const Field<Type>*> valueSets(masterFields.size());
150 forAll(masterFields, fieldi)
152 valueSetNames[fieldi] = masterFields[fieldi].name();
153 valueSets[fieldi] = &masterFields[fieldi][seti];
158 timeDir/formatter.getFileName(masterSampleSet, valueSetNames)
172 void Foam::sampledSets::combineSampledValues
174 const PtrList<volFieldSampler<T> >& sampledFields,
175 const labelListList& indexSets,
176 PtrList<volFieldSampler<T> >& masterFields
179 forAll(sampledFields, fieldi)
181 List<Field<T> > masterValues(indexSets.size());
183 forAll(indexSets, seti)
185 // Collect data from all processors
186 List<Field<T> > gatheredData(Pstream::nProcs());
187 gatheredData[Pstream::myProcNo()] = sampledFields[fieldi][seti];
188 Pstream::gatherList(gatheredData);
190 if (Pstream::master())
194 ListListOps::combine<Field<T> >
197 Foam::accessOp<Field<T> >()
201 masterValues[seti] = UIndirectList<T>
212 new volFieldSampler<T>
215 sampledFields[fieldi].name()
223 void Foam::sampledSets::sampleAndWrite
225 fieldGroup<Type>& fields
230 bool interpolate = interpolationScheme_ != "cell";
232 // Create or use existing writer
233 if (fields.formatter.empty())
235 fields.formatter = writer<Type>::New(writeFormat_);
238 // Storage for interpolated values
239 PtrList<volFieldSampler<Type> > sampledFields(fields.size());
241 forAll(fields, fieldi)
243 if (Pstream::master() && verbose_)
245 Pout<< "sampledSets::sampleAndWrite: "
246 << fields[fieldi] << endl;
251 GeometricField<Type, fvPatchField, volMesh> vf
256 mesh_.time().timeName(),
270 new volFieldSampler<Type>
272 interpolationScheme_,
283 new volFieldSampler<Type>(vf, *this)
294 new volFieldSampler<Type>
296 interpolationScheme_,
298 <GeometricField<Type, fvPatchField, volMesh> >
309 new volFieldSampler<Type>
312 <GeometricField<Type, fvPatchField, volMesh> >
321 // Combine sampled fields from processors.
322 // Note: only master results are valid
324 PtrList<volFieldSampler<Type> > masterFields(sampledFields.size());
325 combineSampledValues(sampledFields, indexSets_, masterFields);
327 if (Pstream::master())
329 forAll(masterSampledSets_, seti)
333 masterSampledSets_[seti],
336 outputPath_/mesh_.time().timeName(),
345 // ************************************************************************* //