Loading...
Searching...
No Matches
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-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
36template<class Type>
37void 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
49 propsDict.add
50 (
51 "file",
52 time_.relativePath(outputName, true)
53 );
54 setProperty(fieldName, propsDict);
55}
56
57
58template<class Type, class GeoMeshType>
59void Foam::sampledSurfaces::storeRegistryField
60(
61 const sampledSurface& s,
62 const word& fieldName,
63 const dimensionSet& dims,
64 Field<Type>&& values
65)
66{
67 s.storeRegistryField<Type, GeoMeshType>
68 (
69 storedObjects(),
70 fieldName,
71 dims,
72 std::move(values),
73 IOobject::groupName(name(), s.name()) // surfaceName
74 );
75}
76
77
78template<class Type>
79void Foam::sampledSurfaces::performAction
80(
82 unsigned request
83)
84{
85 // The sampler for this field
87
88 // The interpolator for this field
90
91 const word& fieldName = fld.name();
92
93 const dimensionSet& dims = fld.dimensions();
94
95 forAll(*this, surfi)
96 {
97 // Skip empty surfaces (eg, failed cut-plane)
98 if (!hasFaces_[surfi]) continue;
99
100 const sampledSurface& s = operator[](surfi);
101
103
104 if (s.isPointData())
105 {
106 if (!interpPtr)
107 {
108 interpPtr = interpolation<Type>::New
109 (
110 sampleNodeScheme_,
111 fld
112 );
113 }
114
115 values = s.interpolate(*interpPtr);
116 }
117 else
118 {
119 if (!samplePtr)
120 {
121 samplePtr = interpolation<Type>::New
122 (
123 sampleFaceScheme_,
124 fld
125 );
126 }
127
128 values = s.sample(*samplePtr);
129 }
130
131 if ((request & actions_[surfi]) & ACTION_WRITE)
132 {
133 writeSurface<Type>(writers_[surfi], values, fieldName);
134 }
135
136 if ((request & actions_[surfi]) & ACTION_STORE)
137 {
138 if (s.isPointData())
139 {
140 storeRegistryField<Type, polySurfacePointGeoMesh>
141 (
142 s, fieldName, dims, std::move(values)
143 );
144 }
145 else
146 {
147 storeRegistryField<Type, polySurfaceGeoMesh>
148 (
149 s, fieldName, dims, std::move(values)
150 );
151 }
152 }
153 }
154}
155
156
157template<class Type>
158void Foam::sampledSurfaces::performAction
159(
160 const SurfaceField<Type>& fld,
161 unsigned request
162)
163{
164 const word& fieldName = fld.name();
165
166 const dimensionSet& dims = fld.dimensions();
167
168 forAll(*this, surfi)
169 {
170 // Skip empty surfaces (eg, failed cut-plane)
171 if (!hasFaces_[surfi]) continue;
172
173 const sampledSurface& s = (*this)[surfi];
174
175 Field<Type> values(s.sample(fld));
176
177 if ((request & actions_[surfi]) & ACTION_WRITE)
178 {
179 writeSurface<Type>(writers_[surfi], values, fieldName);
180 }
181
182 if ((request & actions_[surfi]) & ACTION_STORE)
183 {
184 storeRegistryField<Type, polySurfaceGeoMesh>
185 (
186 s, fieldName, dims, std::move(values)
187 );
188 }
189 }
190}
191
192
193template<class GeoField>
194void Foam::sampledSurfaces::performAction
195(
196 const IOobjectList& objects,
197 unsigned request
198)
199{
200 wordList fieldNames;
201 if (loadFromFiles_)
202 {
203 fieldNames = objects.sortedNames<GeoField>(fieldSelection_);
204 }
205 else
206 {
207 fieldNames = mesh_.thisDb().sortedNames<GeoField>(fieldSelection_);
208 }
209
210 for (const word& fieldName : fieldNames)
211 {
212 if (verbose_)
213 {
214 Info<< "sampleWrite: " << fieldName << endl;
215 }
216
217 refPtr<GeoField> tfield;
218
219 if (loadFromFiles_)
220 {
221 tfield.emplace
222 (
224 (
225 fieldName,
226 time_.timeName(),
227 mesh_.thisDb(),
231 ),
232 mesh_
233 );
234 }
235 else
236 {
237 tfield.cref(mesh_.thisDb().cfindObject<GeoField>(fieldName));
238 }
239
240 if (tfield)
241 {
242 performAction(tfield(), request);
243 }
244 }
245}
246
247
248template<class Container, class Predicate>
249bool Foam::sampledSurfaces::testAny
250(
251 const Container& items,
252 const Predicate& pred
253)
254{
255 for (const auto& item : items)
256 {
257 if (pred(item))
258 {
259 return true;
260 }
261 }
262
263 return false;
264}
265
266
267// ************************************************************************* //
IOdictionary propsDict(dictIO)
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edgesCentres")))
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
@ NO_REGISTER
Do not request registration (bool: false).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
void setProperty(const word &entryName, const Type &value)
Add generic property.
const Time & time_
Reference to the time database.
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
A class for managing references or pointers (no reference counting).
Definition refPtr.H:54
An abstract class for surfaces with sampling.
A class for handling words, derived from Foam::string.
Definition word.H:66
word outputName("finiteArea-edges.obj")
auto & name
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))
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition HashOps.H:164
List< word > wordList
List of word.
Definition fileName.H:60
GeometricField< Type, fvPatchField, volMesh > VolumeField
A volume field for a given type.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< Type, fvsPatchField, surfaceMesh > SurfaceField
A (volume) surface field for a given type.
Fields (face and point) for polySurface.
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.