initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / sampling / sampledSurface / sampledSurfaces / sampledSurfacesTemplates.C
blob9e536ef53b31636e289a6192abdbadeda832385e
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 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 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
19     for more details.
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 "sampledSurfaces.H"
28 #include "volFields.H"
29 #include "ListListOps.H"
31 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
33 template<class Type>
34 Foam::label Foam::sampledSurfaces::grep
36     fieldGroup<Type>& fieldList,
37     const wordList& fieldTypes
38 ) const
40     fieldList.setSize(fieldNames_.size());
41     label nFields = 0;
43     forAll(fieldNames_, fieldI)
44     {
45         if
46         (
47             fieldTypes[fieldI]
48          == GeometricField<Type, fvPatchField, volMesh>::typeName
49         )
50         {
51             fieldList[nFields] = fieldNames_[fieldI];
52             nFields++;
53         }
54     }
56     fieldList.setSize(nFields);
58     return nFields;
62 template<class Type>
63 Foam::autoPtr<Foam::interpolation<Type> >
64 Foam::sampledSurfaces::setInterpolator
66     const GeometricField<Type, fvPatchField, volMesh>& vField
69     if (!pMeshPtr_.valid() || !pInterpPtr_.valid())
70     {
71         // set up interpolation
72         pMeshPtr_.reset(new pointMesh(mesh_));
73         pInterpPtr_.reset(new volPointInterpolation(mesh_, pMeshPtr_()));
74     }
76     // interpolator for this field
77     return interpolation<Type>::New
78     (
79         interpolationScheme_,
80         pInterpPtr_(),
81         vField
82     );
86 template<class Type>
87 void Foam::sampledSurfaces::sampleAndWrite
89     const GeometricField<Type, fvPatchField, volMesh>& vField,
90     const surfaceWriter<Type>& formatter
93     // interpolator for this field
94     autoPtr<interpolation<Type> > interpolator;
96     const word& fieldName   = vField.name();
97     const fileName& timeDir = vField.time().timeName();
99     forAll(*this, surfI)
100     {
101         const sampledSurface& s = operator[](surfI);
103         Field<Type> values;
105         if (s.interpolate())
106         {
107             if (!interpolator.valid())
108             {
109                 interpolator = setInterpolator(vField);
110             }
112             values = s.interpolate(interpolator());
113         }
114         else
115         {
116             values = s.sample(vField);
117         }
119         if (Pstream::parRun())
120         {
121             // Collect values from all processors
122             List<Field<Type> > gatheredValues(Pstream::nProcs());
123             gatheredValues[Pstream::myProcNo()] = values;
124             Pstream::gatherList(gatheredValues);
126             if (Pstream::master())
127             {
128                 // Combine values into single field
129                 Field<Type> allValues
130                 (
131                     ListListOps::combine<Field<Type> >
132                     (
133                         gatheredValues,
134                         accessOp<Field<Type> >()
135                     )
136                 );
138                 // Renumber (point data) to correspond to merged points
139                 if (mergeList_[surfI].pointsMap.size() == allValues.size())
140                 {
141                     inplaceReorder(mergeList_[surfI].pointsMap, allValues);
142                     allValues.setSize(mergeList_[surfI].points.size());
143                 }
145                 // Write to time directory under outputPath_
146                 // skip surface without faces (eg, a failed cut-plane)
147                 if (mergeList_[surfI].faces.size())
148                 {
149                     formatter.write
150                     (
151                         outputPath_,
152                         timeDir,
153                         s.name(),
154                         mergeList_[surfI].points,
155                         mergeList_[surfI].faces,
156                         fieldName,
157                         allValues
158                     );
159                 }
160             }
161         }
162         else
163         {
164             // Write to time directory under outputPath_
165             // skip surface without faces (eg, a failed cut-plane)
166             if (s.faces().size())
167             {
168                 formatter.write
169                 (
170                     outputPath_,
171                     timeDir,
172                     s.name(),
173                     s.points(),
174                     s.faces(),
175                     fieldName,
176                     values
177                 );
178             }
179         }
180     }
184 template <class Type>
185 void Foam::sampledSurfaces::sampleAndWrite
187     fieldGroup<Type>& fields
190     if (fields.size())
191     {
192         // create or use existing surfaceWriter
193         if (!fields.formatter.valid())
194         {
195             fields.formatter = surfaceWriter<Type>::New(writeFormat_);
196         }
198         forAll(fields, fieldI)
199         {
200             if (Pstream::master() && verbose_)
201             {
202                 Pout<< "sampleAndWrite: " << fields[fieldI] << endl;
203             }
205             if (loadFromFiles_)
206             {
207                 sampleAndWrite
208                 (
209                     GeometricField<Type, fvPatchField, volMesh>
210                     (
211                         IOobject
212                         (
213                             fields[fieldI],
214                             mesh_.time().timeName(),
215                             mesh_,
216                             IOobject::MUST_READ,
217                             IOobject::NO_WRITE,
218                             false
219                         ),
220                         mesh_
221                     ),
222                     fields.formatter()
223                 );
224             }
225             else
226             {
227                 objectRegistry::const_iterator iter =
228                     mesh_.find(fields[fieldI]);
230                 if
231                 (
232                     iter != mesh_.objectRegistry::end()
233                  && iter()->type()
234                  == GeometricField<Type, fvPatchField, volMesh>::typeName
235                 )
236                 {
237                    sampleAndWrite
238                    (
239                        mesh_.lookupObject
240                            <GeometricField<Type, fvPatchField, volMesh> >
241                        (
242                            fields[fieldI]
243                        ),
244                        fields.formatter()
245                    );
246                 }
247             }
248         }
249     }
253 // ************************************************************************* //