ENH: RASModel.C: clipping input to log
[OpenFOAM-1.7.x.git] / src / sampling / sampledSet / sampledSets / sampledSetsTemplates.C
blobb3221d68fe3921de8d3b45b8a96b0dc88bfd80bc
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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  * * * * * * * * * * * //
32 template <class Type>
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()),
41     name_(field.name())
43     autoPtr<interpolation<Type> > interpolator
44     (
45         interpolation<Type>::New(interpolationScheme, field)
46     );
48     forAll(samplers, seti)
49     {
50         Field<Type>& values = this->operator[](seti);
51         const sampledSet& samples = samplers[seti];
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];
60             values[samplei] = interpolator().interpolate
61             (
62                 samplePt,
63                 celli,
64                 facei
65             );
66         }
67     }
71 template <class Type>
72 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
74     const GeometricField<Type, fvPatchField, volMesh>& field,
75     const PtrList<sampledSet>& samplers
78     List<Field<Type> >(samplers.size()),
79     name_(field.name())
81     forAll(samplers, seti)
82     {
83         Field<Type>& values = this->operator[](seti);
84         const sampledSet& samples = samplers[seti];
86         values.setSize(samples.size());
87         forAll(samples, samplei)
88         {
89             values[samplei] = field[samples.cells()[samplei]];
90         }
91     }
95 template <class Type>
96 Foam::sampledSets::volFieldSampler<Type>::volFieldSampler
98     const List<Field<Type> >& values,
99     const word& name
102     List<Field<Type> >(values),
103     name_(name)
107 template<class Type>
108 Foam::label Foam::sampledSets::grep
110     fieldGroup<Type>& fieldList,
111     const wordList& fieldTypes
112 ) const
114     fieldList.setSize(fieldNames_.size());
115     label nFields = 0;
117     forAll(fieldNames_, fieldi)
118     {
119         if
120         (
121             fieldTypes[fieldi]
122          == GeometricField<Type, fvPatchField, volMesh>::typeName
123         )
124         {
125             fieldList[nFields] = fieldNames_[fieldi];
126             nFields++;
127         }
128     }
130     fieldList.setSize(nFields);
132     return nFields;
136 template<class Type>
137 void Foam::sampledSets::writeSampleFile
139     const coordSet& masterSampleSet,
140     const PtrList<volFieldSampler<Type> >& masterFields,
141     const label seti,
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)
150     {
151         valueSetNames[fieldi] = masterFields[fieldi].name();
152         valueSets[fieldi] = &masterFields[fieldi][seti];
153     }
155     fileName fName
156     (
157         timeDir/formatter.getFileName(masterSampleSet, valueSetNames)
158     );
160     OFstream ofs(fName);
161     if (ofs.opened())
162     {
163         formatter.write
164         (
165             masterSampleSet,
166             valueSetNames,
167             valueSets,
168             ofs
169         );
170     }
171     else
172     {
173         WarningIn
174         (
175             "void Foam::sampledSets::writeSampleFile"
176             "("
177                 "const coordSet&, "
178                 "const PtrList<volFieldSampler<Type> >&, "
179                 "const label, "
180                 "const fileName&, "
181                 "const writer<Type>&"
182             ")"
183         )   << "File " << ofs.name() << " could not be opened. "
184             << "No data will be written" << endl;
185     }
189 template<class T>
190 void Foam::sampledSets::combineSampledValues
192     const PtrList<volFieldSampler<T> >& sampledFields,
193     const labelListList& indexSets,
194     PtrList<volFieldSampler<T> >& masterFields
197     forAll(sampledFields, fieldi)
198     {
199         List<Field<T> > masterValues(indexSets.size());
201         forAll(indexSets, seti)
202         {
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())
209             {
210                 Field<T> allData
211                 (
212                     ListListOps::combine<Field<T> >
213                     (
214                         gatheredData,
215                         Foam::accessOp<Field<T> >()
216                     )
217                 );
219                 masterValues[seti] = UIndirectList<T>
220                 (
221                     allData,
222                     indexSets[seti]
223                 )();
224             }
225         }
227         masterFields.set
228         (
229             fieldi,
230             new volFieldSampler<T>
231             (
232                 masterValues,
233                 sampledFields[fieldi].name()
234             )
235         );
236     }
240 template<class Type>
241 void Foam::sampledSets::sampleAndWrite
243     fieldGroup<Type>& fields
246     if (fields.size())
247     {
248         bool interpolate = interpolationScheme_ != "cell";
250         // Create or use existing writer
251         if (fields.formatter.empty())
252         {
253             fields.formatter = writer<Type>::New(writeFormat_);
254         }
256         // Storage for interpolated values
257         PtrList<volFieldSampler<Type> > sampledFields(fields.size());
259         forAll(fields, fieldi)
260         {
261             if (Pstream::master() && verbose_)
262             {
263                 Pout<< "sampledSets::sampleAndWrite: "
264                     << fields[fieldi] << endl;
265             }
267             if (loadFromFiles_)
268             {
269                 GeometricField<Type, fvPatchField, volMesh> vf
270                 (
271                     IOobject
272                     (
273                         fields[fieldi],
274                         mesh_.time().timeName(),
275                         mesh_,
276                         IOobject::MUST_READ,
277                         IOobject::NO_WRITE,
278                         false
279                     ),
280                     mesh_
281                 );
283                 if (interpolate)
284                 {
285                     sampledFields.set
286                     (
287                         fieldi,
288                         new volFieldSampler<Type>
289                         (
290                             interpolationScheme_,
291                             vf,
292                             *this
293                         )
294                     );
295                 }
296                 else
297                 {
298                     sampledFields.set
299                     (
300                         fieldi,
301                         new volFieldSampler<Type>(vf, *this)
302                     );
303                 }
304             }
305             else
306             {
307                 if (interpolate)
308                 {
309                     sampledFields.set
310                     (
311                         fieldi,
312                         new volFieldSampler<Type>
313                         (
314                             interpolationScheme_,
315                             mesh_.lookupObject
316                             <GeometricField<Type, fvPatchField, volMesh> >
317                             (fields[fieldi]),
318                             *this
319                         )
320                     );
321                 }
322                 else
323                 {
324                     sampledFields.set
325                     (
326                         fieldi,
327                         new volFieldSampler<Type>
328                         (
329                             mesh_.lookupObject
330                             <GeometricField<Type, fvPatchField, volMesh> >
331                             (fields[fieldi]),
332                             *this
333                         )
334                     );
335                 }
336             }
337         }
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())
346         {
347             forAll(masterSampledSets_, seti)
348             {
349                 writeSampleFile
350                 (
351                     masterSampledSets_[seti],
352                     masterFields,
353                     seti,
354                     outputPath_/mesh_.time().timeName(),
355                     fields.formatter()
356                 );
357             }
358         }
359     }
363 // ************************************************************************* //