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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "sampledSurfaces.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "polySurfaceFields.H"
33 
34 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35 
36 template<class Type>
37 void Foam::sampledSurfaces::writeSurface
38 (
39  surfaceWriter& writer,
40  const Field<Type>& values,
41  const word& fieldName
42 )
43 {
44  fileName outputName = writer.write(fieldName, values);
45 
46  // Case-local file name with "<case>" to make relocatable
47 
48  dictionary propsDict;
49  propsDict.add
50  (
51  "file",
53  );
54  setProperty(fieldName, propsDict);
55 }
56 
57 
58 template<class Type, class GeoMeshType>
59 bool Foam::sampledSurfaces::storeRegistryField
60 (
61  const sampledSurface& s,
62  const word& fieldName,
63  const dimensionSet& dims,
64  Field<Type>&& values
65 )
66 {
67  return s.storeRegistryField<Type, GeoMeshType>
68  (
69  storedObjects(),
70  fieldName,
71  dims,
72  std::move(values),
73  IOobject::groupName(name(), s.name())
74  );
75 }
76 
77 
78 template<class Type>
79 void Foam::sampledSurfaces::performAction
80 (
81  const VolumeField<Type>& fld,
82  unsigned request
83 )
84 {
85  // The sampler for this field
86  autoPtr<interpolation<Type>> samplePtr;
87 
88  // The interpolator for this field
89  autoPtr<interpolation<Type>> interpPtr;
90 
91  const word& fieldName = fld.name();
92 
93  const dimensionSet& dims = fld.dimensions();
94 
95  forAll(*this, surfi)
96  {
97  const sampledSurface& s = operator[](surfi);
98 
99  // Skip surface without faces (eg, failed cut-plane)
100  if (!nFaces_[surfi])
101  {
102  continue;
103  }
104 
105  Field<Type> values;
106 
107  if (s.isPointData())
108  {
109  if (!interpPtr)
110  {
111  interpPtr = interpolation<Type>::New
112  (
113  sampleNodeScheme_,
114  fld
115  );
116  }
117 
118  values = s.interpolate(*interpPtr);
119  }
120  else
121  {
122  if (!samplePtr)
123  {
124  samplePtr = interpolation<Type>::New
125  (
126  sampleFaceScheme_,
127  fld
128  );
129  }
130 
131  values = s.sample(*samplePtr);
132  }
133 
134  if ((request & actions_[surfi]) & ACTION_WRITE)
135  {
136  writeSurface<Type>(writers_[surfi], values, fieldName);
137  }
138 
139  if ((request & actions_[surfi]) & ACTION_STORE)
140  {
141  if (s.isPointData())
142  {
143  storeRegistryField<Type, polySurfacePointGeoMesh>
144  (
145  s, fieldName, dims, std::move(values)
146  );
147  }
148  else
149  {
150  storeRegistryField<Type, polySurfaceGeoMesh>
151  (
152  s, fieldName, dims, std::move(values)
153  );
154  }
155  }
156  }
157 }
158 
159 
160 template<class Type>
161 void Foam::sampledSurfaces::performAction
162 (
163  const SurfaceField<Type>& fld,
164  unsigned request
165 )
166 {
167  const word& fieldName = fld.name();
168 
169  const dimensionSet& dims = fld.dimensions();
170 
171  forAll(*this, surfi)
172  {
173  const sampledSurface& s = (*this)[surfi];
174 
175  // Skip surface without faces (eg, failed cut-plane)
176  if (!nFaces_[surfi])
177  {
178  continue;
179  }
180 
181  Field<Type> values(s.sample(fld));
182 
183  if ((request & actions_[surfi]) & ACTION_WRITE)
184  {
185  writeSurface<Type>(writers_[surfi], values, fieldName);
186  }
187 
188  if ((request & actions_[surfi]) & ACTION_STORE)
189  {
190  storeRegistryField<Type, polySurfaceGeoMesh>
191  (
192  s, fieldName, dims, std::move(values)
193  );
194  }
195  }
196 }
197 
198 
199 template<class GeoField>
200 void Foam::sampledSurfaces::performAction
201 (
202  const IOobjectList& objects,
203  unsigned request
204 )
205 {
206  wordList fieldNames;
207  if (loadFromFiles_)
208  {
209  fieldNames = objects.sortedNames<GeoField>(fieldSelection_);
210  }
211  else
212  {
213  fieldNames = mesh_.thisDb().sortedNames<GeoField>(fieldSelection_);
214  }
215 
216  for (const word& fieldName : fieldNames)
217  {
218  if (verbose_)
219  {
220  Info<< "sampleWrite: " << fieldName << endl;
221  }
222 
223  if (loadFromFiles_)
224  {
225  const GeoField fld
226  (
227  IOobject
228  (
229  fieldName,
230  time_.timeName(),
231  mesh_,
233  ),
234  mesh_
235  );
236 
237  performAction(fld, request);
238  }
239  else
240  {
241  performAction
242  (
243  mesh_.thisDb().lookupObject<GeoField>(fieldName),
244  request
245  );
246  }
247  }
248 }
249 
250 
251 template<class Container, class Predicate>
252 bool Foam::sampledSurfaces::testAny
253 (
254  const Container& items,
255  const Predicate& pred
256 )
257 {
258  for (const auto& item : items)
259  {
260  if (pred(item))
261  {
262  return true;
263  }
264  }
265 
266  return false;
267 }
268 
269 
270 // ************************************************************************* //
Foam::surfaceFields.
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
Fields (face and point) for polySurface.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
fileName relativePath(const fileName &input, const bool caseTag=false) const
Return the input relative to the globalPath by stripping off a leading value of the globalPath...
Definition: TimePathsI.H:80
IOdictionary propsDict(dictIO)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
word outputName("finiteArea-edges.obj")
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
vtk::internalMeshWriter writer(topoMesh, topoCells, vtk::formatType::INLINE_ASCII, runTime.path()/"blockTopology")
void setProperty(const word &entryName, const Type &value)
Add generic property.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
List< word > wordList
A List of words.
Definition: fileName.H:58
messageStream Info
Information stream (stdout output on master, null elsewhere)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
const Time & time_
Reference to the time database.