ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / sampling / sampledSurface / sampledSurfaces / sampledSurfacesTemplates.C
blob78d7808d7af2cd09c5d80ddae203420eaa4a55e6
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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 "sampledSurfaces.H"
27 #include "volFields.H"
28 #include "ListListOps.H"
30 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
32 template<class Type>
33 void Foam::sampledSurfaces::sampleAndWrite
35     const GeometricField<Type, fvPatchField, volMesh>& vField
38     // interpolator for this field
39     autoPtr< interpolation<Type> > interpolator;
41     const word& fieldName   = vField.name();
42     const fileName outputDir = outputPath_/vField.time().timeName();
44     forAll(*this, surfI)
45     {
46         const sampledSurface& s = operator[](surfI);
48         Field<Type> values;
50         if (s.interpolate())
51         {
52             if (interpolator.empty())
53             {
54                 interpolator = interpolation<Type>::New
55                 (
56                     interpolationScheme_,
57                     vField
58                 );
59             }
61             values = s.interpolate(interpolator());
62         }
63         else
64         {
65             values = s.sample(vField);
66         }
68         if (Pstream::parRun())
69         {
70             // Collect values from all processors
71             List<Field<Type> > gatheredValues(Pstream::nProcs());
72             gatheredValues[Pstream::myProcNo()] = values;
73             Pstream::gatherList(gatheredValues);
75             if (Pstream::master())
76             {
77                 // Combine values into single field
78                 Field<Type> allValues
79                 (
80                     ListListOps::combine<Field<Type> >
81                     (
82                         gatheredValues,
83                         accessOp<Field<Type> >()
84                     )
85                 );
87                 // Renumber (point data) to correspond to merged points
88                 if (mergeList_[surfI].pointsMap.size() == allValues.size())
89                 {
90                     inplaceReorder(mergeList_[surfI].pointsMap, allValues);
91                     allValues.setSize(mergeList_[surfI].points.size());
92                 }
94                 // Write to time directory under outputPath_
95                 // skip surface without faces (eg, a failed cut-plane)
96                 if (mergeList_[surfI].faces.size())
97                 {
98                     formatter_->write
99                     (
100                         outputDir,
101                         s.name(),
102                         mergeList_[surfI].points,
103                         mergeList_[surfI].faces,
104                         fieldName,
105                         allValues,
106                         s.interpolate()
107                     );
108                 }
109             }
110         }
111         else
112         {
113             // Write to time directory under outputPath_
114             // skip surface without faces (eg, a failed cut-plane)
115             if (s.faces().size())
116             {
117                 formatter_->write
118                 (
119                     outputDir,
120                     s.name(),
121                     s.points(),
122                     s.faces(),
123                     fieldName,
124                     values,
125                     s.interpolate()
126                 );
127             }
128         }
129     }
133 template<class Type>
134 void Foam::sampledSurfaces::sampleAndWrite
136     fieldGroup<Type>& fields
139     if (fields.size())
140     {
141         forAll(fields, fieldI)
142         {
143             if (Pstream::master() && verbose_)
144             {
145                 Pout<< "sampleAndWrite: " << fields[fieldI] << endl;
146             }
148             if (loadFromFiles_)
149             {
150                 sampleAndWrite
151                 (
152                     GeometricField<Type, fvPatchField, volMesh>
153                     (
154                         IOobject
155                         (
156                             fields[fieldI],
157                             mesh_.time().timeName(),
158                             mesh_,
159                             IOobject::MUST_READ,
160                             IOobject::NO_WRITE,
161                             false
162                         ),
163                         mesh_
164                     )
165                 );
166             }
167             else
168             {
169                 objectRegistry::const_iterator iter =
170                     mesh_.find(fields[fieldI]);
172                 if
173                 (
174                     iter != objectRegistry::end()
175                  && iter()->type()
176                  == GeometricField<Type, fvPatchField, volMesh>::typeName
177                 )
178                 {
179                    sampleAndWrite
180                    (
181                        mesh_.lookupObject
182                        <GeometricField<Type, fvPatchField, volMesh> >
183                        (
184                            fields[fieldI]
185                        )
186                    );
187                 }
188             }
189         }
190     }
194 // ************************************************************************* //