1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd.
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
13 the Free Software Foundation, either version 3 of the License, or
14 (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "sampledSets.H"
27 #include "volFields.H"
28 #include "ListListOps.H"
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
35 const word& interpolationScheme,
36 const GeometricField<Type, fvPatchField, volMesh>& field,
37 const PtrList<sampledSet>& samplers
40 List<Field<Type> >(samplers.size()),
43 autoPtr<interpolation<Type> > interpolator
45 interpolation<Type>::New(interpolationScheme, field)
48 forAll(samplers, seti)
50 Field<Type>& values = this->operator[](seti);
51 const sampledSet& samples = samplers[seti];
53 values.setSize(samples.size());
54 forAll(samples, samplei)
56 const point& samplePt = samples[samplei];
57 label celli = samples.cells()[samplei];
58 label facei = samples.faces()[samplei];
60 values[samplei] = interpolator().interpolate
72 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
74 const GeometricField<Type, fvPatchField, volMesh>& field,
75 const PtrList<sampledSet>& samplers
78 List<Field<Type> >(samplers.size()),
81 forAll(samplers, seti)
83 Field<Type>& values = this->operator[](seti);
84 const sampledSet& samples = samplers[seti];
86 values.setSize(samples.size());
87 forAll(samples, samplei)
89 values[samplei] = field[samples.cells()[samplei]];
96 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
98 const List<Field<Type> >& values,
102 List<Field<Type> >(values),
108 Foam::label Foam::sampledSets::grep
110 fieldGroup<Type>& fieldList,
111 const wordList& fieldTypes
114 fieldList.setSize(fieldNames_.size());
117 forAll(fieldNames_, fieldi)
122 == GeometricField<Type, fvPatchField, volMesh>::typeName
125 fieldList[nFields] = fieldNames_[fieldi];
130 fieldList.setSize(nFields);
137 void Foam::sampledSets::writeSampleFile
139 const coordSet& masterSampleSet,
140 const PtrList<volFieldSampler<Type> >& masterFields,
142 const fileName& timeDir,
143 const writer<Type>& formatter
146 wordList valueSetNames(masterFields.size());
147 List<const Field<Type>*> valueSets(masterFields.size());
149 forAll(masterFields, fieldi)
151 valueSetNames[fieldi] = masterFields[fieldi].name();
152 valueSets[fieldi] = &masterFields[fieldi][seti];
157 timeDir/formatter.getFileName(masterSampleSet, valueSetNames)
175 "void Foam::sampledSets::writeSampleFile"
178 "const PtrList<volFieldSampler<Type> >&, "
181 "const writer<Type>&"
183 ) << "File " << ofs.name() << " could not be opened. "
184 << "No data will be written" << endl;
190 void Foam::sampledSets::combineSampledValues
192 const PtrList<volFieldSampler<T> >& sampledFields,
193 const labelListList& indexSets,
194 PtrList<volFieldSampler<T> >& masterFields
197 forAll(sampledFields, fieldi)
199 List<Field<T> > masterValues(indexSets.size());
201 forAll(indexSets, seti)
203 // Collect data from all processors
204 List<Field<T> > gatheredData(Pstream::nProcs());
205 gatheredData[Pstream::myProcNo()] = sampledFields[fieldi][seti];
206 Pstream::gatherList(gatheredData);
208 if (Pstream::master())
212 ListListOps::combine<Field<T> >
215 Foam::accessOp<Field<T> >()
219 masterValues[seti] = UIndirectList<T>
230 new volFieldSampler<T>
233 sampledFields[fieldi].name()
241 void Foam::sampledSets::sampleAndWrite
243 fieldGroup<Type>& fields
248 bool interpolate = interpolationScheme_ != "cell";
250 // Create or use existing writer
251 if (fields.formatter.empty())
253 fields.formatter = writer<Type>::New(writeFormat_);
256 // Storage for interpolated values
257 PtrList<volFieldSampler<Type> > sampledFields(fields.size());
259 forAll(fields, fieldi)
261 if (Pstream::master() && verbose_)
263 Pout<< "sampledSets::sampleAndWrite: "
264 << fields[fieldi] << endl;
269 GeometricField<Type, fvPatchField, volMesh> vf
274 mesh_.time().timeName(),
288 new volFieldSampler<Type>
290 interpolationScheme_,
301 new volFieldSampler<Type>(vf, *this)
312 new volFieldSampler<Type>
314 interpolationScheme_,
316 <GeometricField<Type, fvPatchField, volMesh> >
327 new volFieldSampler<Type>
330 <GeometricField<Type, fvPatchField, volMesh> >
339 // Combine sampled fields from processors.
340 // Note: only master results are valid
342 PtrList<volFieldSampler<Type> > masterFields(sampledFields.size());
343 combineSampledValues(sampledFields, indexSets_, masterFields);
345 if (Pstream::master())
347 forAll(masterSampledSets_, seti)
351 masterSampledSets_[seti],
354 outputPath_/mesh_.time().timeName(),
363 // ************************************************************************* //