All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
sampledSurfacesTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "sampledSurfaces.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "ListListOps.H"
30 #include "stringListOps.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 template<class Type>
35 void Foam::sampledSurfaces::writeSurface
36 (
37  const Field<Type>& values,
38  const label surfI,
39  const word& fieldName,
40  const fileName& outputDir
41 )
42 {
43  const sampledSurface& s = operator[](surfI);
44 
45  if (Pstream::parRun())
46  {
47  // Collect values from all processors
48  List<Field<Type>> gatheredValues(Pstream::nProcs());
49  gatheredValues[Pstream::myProcNo()] = values;
50  Pstream::gatherList(gatheredValues);
51 
52  if (Pstream::master())
53  {
54  // Combine values into single field
55  Field<Type> allValues
56  (
57  ListListOps::combine<Field<Type>>
58  (
59  gatheredValues,
60  accessOp<Field<Type>>()
61  )
62  );
63 
64  // Renumber (point data) to correspond to merged points
65  if (mergeList_[surfI].pointsMap.size() == allValues.size())
66  {
67  inplaceReorder(mergeList_[surfI].pointsMap, allValues);
68  allValues.setSize(mergeList_[surfI].points.size());
69  }
70 
71  // Write to time directory under outputPath_
72  // skip surface without faces (eg, a failed cut-plane)
73  if (mergeList_[surfI].faces.size())
74  {
75  formatter_->write
76  (
77  outputDir,
78  s.name(),
79  mergeList_[surfI].points,
80  mergeList_[surfI].faces,
81  fieldName,
82  allValues,
83  s.interpolate()
84  );
85  }
86  }
87  }
88  else
89  {
90  // Write to time directory under outputPath_
91  // skip surface without faces (eg, a failed cut-plane)
92  if (s.faces().size())
93  {
94  formatter_->write
95  (
96  outputDir,
97  s.name(),
98  s.points(),
99  s.faces(),
100  fieldName,
101  values,
102  s.interpolate()
103  );
104  }
105  }
106 }
107 
108 
109 template<class Type>
110 void Foam::sampledSurfaces::sampleAndWrite
111 (
112  const GeometricField<Type, fvPatchField, volMesh>& vField
113 )
114 {
115  // interpolator for this field
116  autoPtr<interpolation<Type>> interpolatorPtr;
117 
118  const word& fieldName = vField.name();
119  const fileName outputDir = outputPath_/vField.time().timeName();
120 
121  forAll(*this, surfI)
122  {
123  const sampledSurface& s = operator[](surfI);
124 
125  Field<Type> values;
126 
127  if (s.interpolate())
128  {
129  if (interpolatorPtr.empty())
130  {
131  interpolatorPtr = interpolation<Type>::New
132  (
133  interpolationScheme_,
134  vField
135  );
136  }
137 
138  values = s.interpolate(interpolatorPtr());
139  }
140  else
141  {
142  values = s.sample(vField);
143  }
144 
145  writeSurface<Type>(values, surfI, fieldName, outputDir);
146  }
147 }
148 
149 
150 template<class Type>
151 void Foam::sampledSurfaces::sampleAndWrite
152 (
153  const GeometricField<Type, fvsPatchField, surfaceMesh>& sField
154 )
155 {
156  const word& fieldName = sField.name();
157  const fileName outputDir = outputPath_/sField.time().timeName();
158 
159  forAll(*this, surfI)
160  {
161  const sampledSurface& s = operator[](surfI);
162  Field<Type> values(s.sample(sField));
163  writeSurface<Type>(values, surfI, fieldName, outputDir);
164  }
165 }
166 
167 
168 template<class GeoField>
169 void Foam::sampledSurfaces::sampleAndWrite(const IOobjectList& objects)
170 {
171  wordList names;
172  if (loadFromFiles_)
173  {
174  IOobjectList fieldObjects(objects.lookupClass(GeoField::typeName));
175  names = fieldObjects.names();
176  }
177  else
178  {
179  names = mesh_.thisDb().names<GeoField>();
180  }
181 
182  labelList nameIDs(findStrings(fieldSelection_, names));
183 
184  wordHashSet fieldNames(wordList(names, nameIDs));
185 
187  {
188  const word& fieldName = iter.key();
189 
190  if ((Pstream::master()) && verbose_)
191  {
192  Pout<< "sampleAndWrite: " << fieldName << endl;
193  }
194 
195  if (loadFromFiles_)
196  {
197  const GeoField fld
198  (
199  IOobject
200  (
201  fieldName,
202  mesh_.time().timeName(),
203  mesh_,
205  ),
206  mesh_
207  );
208 
209  sampleAndWrite(fld);
210  }
211  else
212  {
213  sampleAndWrite
214  (
215  mesh_.thisDb().lookupObject<GeoField>(fieldName)
216  );
217  }
218  }
219 }
220 
221 
222 // ************************************************************************* //
Foam::surfaceFields.
void inplaceReorder(const labelUList &oldToNew, ListType &)
Inplace reorder the elements of a list.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
wordList names() const
Return the list of names of the IOobjects.
const T & operator[](const label) const
Return element const reference.
Definition: UPtrListI.H:96
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:622
Operations on lists of strings.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
AccessType combine(const List< T > &, AccessOp aop=accessOp< T >())
Combines sublists into one big list.
Definition: ListListOps.C:34
static List< word > fieldNames
Definition: globalFoam.H:46
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:245
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const pointField & points
List< label > labelList
A List of labels.
Definition: labelList.H:56
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:208
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
List< word > wordList
A List of words.
Definition: fileName.H:54
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.