sampledSurfacesTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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>
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  (
58  (
59  gatheredValues,
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  fileName fName = 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 
87  propsDict.add("file", fName);
88  setProperty(fieldName, propsDict);
89  }
90  }
91  }
92  else
93  {
94  // Write to time directory under outputPath_
95  // skip surface without faces (eg, a failed cut-plane)
96  if (s.faces().size())
97  {
98  fileName fName = formatter_->write
99  (
100  outputDir,
101  s.name(),
102  s.points(),
103  s.faces(),
104  fieldName,
105  values,
106  s.interpolate()
107  );
108 
110  propsDict.add("file", fName);
111  setProperty(fieldName, propsDict);
112  }
113  }
114 }
115 
116 
117 template<class Type>
119 (
121 )
122 {
123  // interpolator for this field
124  autoPtr<interpolation<Type> > interpolatorPtr;
125 
126  const word& fieldName = vField.name();
127  const fileName outputDir = outputPath_/vField.time().timeName();
128 
129  forAll(*this, surfI)
130  {
131  const sampledSurface& s = operator[](surfI);
132 
133  Field<Type> values;
134 
135  if (s.interpolate())
136  {
137  if (interpolatorPtr.empty())
138  {
139  interpolatorPtr = interpolation<Type>::New
140  (
141  interpolationScheme_,
142  vField
143  );
144  }
145 
146  values = s.interpolate(interpolatorPtr());
147  }
148  else
149  {
150  values = s.sample(vField);
151  }
152 
153  writeSurface<Type>(values, surfI, fieldName, outputDir);
154  }
155 }
156 
157 
158 template<class Type>
160 (
162 )
163 {
164  const word& fieldName = sField.name();
165  const fileName outputDir = outputPath_/sField.time().timeName();
166 
167  forAll(*this, surfI)
168  {
169  const sampledSurface& s = operator[](surfI);
170  Field<Type> values(s.sample(sField));
171  writeSurface<Type>(values, surfI, fieldName, outputDir);
172  }
173 }
174 
175 
176 template<class GeoField>
178 {
179  wordList names;
180  const fvMesh& mesh = refCast<const fvMesh>(obr_);
181 
182  if (loadFromFiles_)
183  {
184  IOobjectList fieldObjects(objects.lookupClass(GeoField::typeName));
185  names = fieldObjects.names();
186  }
187  else
188  {
189  names = mesh.thisDb().names<GeoField>();
190  }
191 
192  labelList nameIDs(findStrings(fieldSelection_, names));
193 
194  wordHashSet fieldNames(wordList(names, nameIDs));
195 
197  {
198  const word& fieldName = iter.key();
199 
200  if ((Pstream::master()) && verbose_)
201  {
202  Pout<< "sampleAndWrite: " << fieldName << endl;
203  }
204 
205  if (loadFromFiles_)
206  {
207  const GeoField fld
208  (
209  IOobject
210  (
211  fieldName,
212  mesh.time().timeName(),
213  mesh,
215  ),
216  mesh
217  );
218 
220  }
221  else
222  {
224  (
225  mesh.thisDb().lookupObject<GeoField>(fieldName)
226  );
227  }
228  }
229 }
230 
231 
232 // ************************************************************************* //
volFields.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::sampledSurfaces::obr_
const objectRegistry & obr_
Const reference to database.
Definition: sampledSurfaces.H:97
ListListOps.H
Foam::autoPtr::empty
bool empty() const
Return true if the autoPtr is empty (ie, no pointer set).
Definition: autoPtrI.H:76
Foam::ListListOps::combine
AccessType combine(const List< T > &, AccessOp aop=accessOp< T >())
Combines sublists into one big list.
Definition: ListListOps.C:34
Foam::accessOp
Definition: ListListOps.H:98
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::IOobjectList::names
wordList names() const
Return the list of names of the IOobjects.
Definition: IOobjectList.C:221
Foam::HashSet
A HashTable with keys but without contents.
Definition: HashSet.H:59
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
sampledSurfaces.H
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:54
Foam::inplaceReorder
void inplaceReorder(const labelUList &oldToNew, ListType &)
Inplace reorder the elements of a list.
Definition: ListOpsTemplates.C:102
Foam::label
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
Foam::Field< Type >
propsDict
IOdictionary propsDict(IOobject("particleTrackProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED))
Foam::sampledSurface
An abstract class for surfaces with sampling.
Definition: sampledSurface.H:77
Foam::sampledSurfaces::writeSurface
void writeSurface(const Field< Type > &values, const label surfI, const word &fieldName, const fileName &outputDir)
Write sampled fieldName on surface and on outputDir path.
Definition: sampledSurfacesTemplates.C:36
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fld
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))
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::IOobjectList::lookupClass
IOobjectList lookupClass(const word &className) const
Return the list for all IOobjects of a given class.
Definition: IOobjectList.C:199
s
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){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:399
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::findStrings
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
Foam::sampledSurfaces::loadFromFiles_
const bool loadFromFiles_
Load fields from files (not from objectRegistry)
Definition: sampledSurfaces.H:100
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::fieldNames
wordList fieldNames(const IOobjectList &objects, const bool syncPar)
Get sorted names of fields of type. If syncPar and running in parallel.
Definition: ReadFields.C:36
Foam::sampledSurfaces::verbose_
static bool verbose_
Output verbosity.
Definition: sampledSurfaces.H:88
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
points
const pointField & points
Definition: gmvOutputHeader.H:1
stringListOps.H
Operations on lists of strings.
Foam::sampledSurfaces::sampleAndWrite
void sampleAndWrite(const GeometricField< Type, fvPatchField, volMesh > &)
Sample and write a particular volume field.
Definition: sampledSurfacesTemplates.C:119
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::sampledSurfaces::fieldSelection_
wordReList fieldSelection_
Names of fields to sample.
Definition: sampledSurfaces.H:109