volumeExprDriverFields.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) 2019-2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "volumeExprDriver.H"
29#include "fvPatch.H"
30#include "error.H"
31
32// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
33
36(
37 const word& name,
38 enum topoSetSource::sourceType setType
39) const
40{
41 auto tresult = volScalarField::New
42 (
43 "selected",
44 mesh(),
46 );
47
48 refPtr<labelList> tselected;
49 switch (setType)
50 {
53 {
54 tselected = getTopoSetLabels(name, setType);
55 break;
56 }
57
58 default:
59 {
61 << "Unexpected sourceType: " << int(setType) << nl
62 << exit(FatalError);
63 break;
64 }
65 }
66 const auto& selected = tselected();
67
68 auto& fld = tresult.ref().primitiveFieldRef();
69 UIndirectList<scalar>(fld, selected) = scalar(1);
70
71 return tresult;
72}
73
74
77(
78 const word& name,
79 enum topoSetSource::sourceType setType
80) const
81{
82 auto tresult = surfaceScalarField::New
83 (
84 "selected",
85 mesh(),
87 );
88
89 refPtr<labelList> tselected;
90 switch (setType)
91 {
94 {
95 tselected = getTopoSetLabels(name, setType);
96 break;
97 }
98
99 default:
100 {
102 << "Unexpected sourceType: " << int(setType) << nl
103 << exit(FatalError);
104 break;
105 }
106 }
107 const auto& selected = tselected();
108
109 const auto& bmesh = mesh().boundaryMesh();
110
111 auto& result = tresult.ref();
112 auto& fld = result.primitiveFieldRef();
113 auto& bfld = result.boundaryFieldRef();
114
115 label nErrors = 0;
116
117 for (const label facei : selected)
118 {
119 if (facei < mesh().nInternalFaces())
120 {
121 fld[facei] = scalar(1);
122 }
123 else
124 {
125 const label patchi = bmesh.whichPatch(facei);
126
127 if (patchi < 0)
128 {
129 ++nErrors;
130 }
131 else
132 {
133 bfld[patchi][facei-bmesh[patchi].start()] = scalar(1);
134 }
135 }
136 }
137
138 if (nErrors)
139 {
141 << "The faceSet/faceZone " << name << " contained "
142 << nErrors << " faces outside of the addressing range" << nl
143 << nl;
144 }
145
146
147 return tresult;
148}
149
150
153(
154 const word& name,
155 enum topoSetSource::sourceType setType
156) const
157{
158 auto tresult = pointScalarField::New
159 (
160 "selected",
163 );
164
165 refPtr<labelList> tselected;
166 switch (setType)
167 {
170 {
171 tselected = getTopoSetLabels(name, setType);
172 break;
173 }
174
175 default:
176 {
178 << "Unexpected sourceType: " << int(setType) << nl
179 << exit(FatalError);
180 break;
181 }
182 }
183 const auto& selected = tselected();
184
185 auto& fld = tresult.ref().primitiveFieldRef();
186 UIndirectList<scalar>(fld, selected) = scalar(1);
187
188 return tresult;
189}
190
191
192
195{
197 (
198 "vol",
199 mesh(),
200 dimVol,
201 mesh().V()
202 );
203}
204
205
208{
209 return tmp<volVectorField>::New(mesh().C());
210}
211
212
215{
217 (
218 "face",
219 mesh(),
220 dimless,
221 mesh().magSf()
222 );
223}
224
225
228{
230 (
231 "fpos",
232 mesh(),
233 dimless,
234 mesh().Cf()
235 );
236}
237
238
241{
243 (
244 "face",
245 mesh(),
246 dimless,
247 mesh().Sf()
248 );
249}
250
251
254{
256 (
257 "pts",
259 dimless,
260 mesh().points()
261 );
262}
263
264
267(
268 label seed,
269 bool gaussian
270) const
271{
272 auto tfld = volScalarField::New("rand", mesh(), dimless);
273
274 exprDriver::fill_random(tfld.ref().primitiveFieldRef(), seed, gaussian);
275
276 return tfld;
277}
278
279
280// ************************************************************************* //
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Graphite solid properties.
Definition: C.H:53
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A List with indirect addressing. Like IndirectList but does not store addressing.
Definition: IndirectList.H:79
void fill_random(scalarField &field, label seed=0, const bool gaussian=false) const
Fill a random field.
refPtr< labelList > getTopoSetLabels(const word &name, enum topoSetSource::sourceType setType) const
Get the labels associated with the topo set.
tmp< volScalarField > field_cellVolume() const
The cell volumes - (swak = vol)
tmp< surfaceScalarField > field_faceArea() const
The face area magnitudes [magSf] - (swak = area)
tmp< pointVectorField > field_pointField() const
The mesh point locations - (swak = pts)
tmp< surfaceScalarField > field_faceSelection(const word &name, enum topoSetSource::sourceType setType) const
Face selections (as logical)
tmp< surfaceVectorField > field_areaNormal() const
The face areas with their vector direction [Sf] - (swak = face)
tmp< surfaceVectorField > field_faceCentre() const
The face centres - (swak = fpos)
tmp< volScalarField > field_cellSelection(const word &name, enum topoSetSource::sourceType setType) const
Cell selections (as logical)
tmp< volScalarField > field_rand(label seed=0, bool gaussian=false) const
A uniform random field.
tmp< volVectorField > field_cellCentre() const
The cell centres - (swak = pos)
virtual const fvMesh & mesh() const
The mesh we are attached to.
tmp< pointScalarField > field_pointSelection(const word &name, enum topoSetSource::sourceType setType) const
Point selections (as logical)
A class for managing references or pointers (no reference counting)
Definition: refPtr.H:58
A class for managing temporary objects.
Definition: tmp.H:65
sourceType
Enumeration defining the types of sources.
Definition: topoSetSource.H:75
@ POINTSET_SOURCE
Points as set.
Definition: topoSetSource.H:84
@ FACESET_SOURCE
Faces as set.
Definition: topoSetSource.H:83
@ FACEZONE_SOURCE
Faces as zone.
Definition: topoSetSource.H:88
@ POINTZONE_SOURCE
Points as zone.
Definition: topoSetSource.H:89
@ CELLSET_SOURCE
Cells as set.
Definition: topoSetSource.H:82
@ CELLZONE_SOURCE
Cells as zone.
Definition: topoSetSource.H:87
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
Definition: dimensionSets.H:64
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53