Loading...
Searching...
No Matches
AveragingMethod.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) 2013-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2023 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
31#include "pointMesh.H"
32
33// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
34
35template<class Type>
37{}
38
39
40// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41
42template<class Type>
44(
45 const IOobject& io,
46 const dictionary& dict,
47 const fvMesh& mesh,
48 const labelList& size
49)
50:
52 FieldField<Field, Type>(),
53 dict_(dict),
54 mesh_(mesh)
55{
56 forAll(size, i)
57 {
59 (
61 );
62 }
63}
64
65
66template<class Type>
68(
69 const AveragingMethod<Type>& am
70)
71:
72 regIOobject(am),
73 FieldField<Field, Type>(am),
74 dict_(am.dict_),
75 mesh_(am.mesh_)
76{}
77
78
79// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
80
81template<class Type>
84(
85 const IOobject& io,
86 const dictionary& dict,
87 const fvMesh& mesh
88)
89{
90 const word modelType
91 (
92 dict.template getOrDefault<word>(typeName, "basic")
93 );
94
95 auto* ctorPtr = dictionaryConstructorTable(modelType);
96
97 if (!ctorPtr)
98 {
100 (
101 dict,
102 "averaging limiter",
103 modelType,
104 *dictionaryConstructorTablePtr_
105 ) << abort(FatalIOError);
106 }
107
109}
110
111
112// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113
114template<class Type>
116{
117 updateGrad();
118}
119
120
121template<class Type>
123(
124 const AveragingMethod<scalar>& weight
125)
126{
128
129 *this /= max(weight, SMALL);
130}
131
132
133template<class Type>
135{
136 return os.good();
137}
138
139
140template<class Type>
141bool Foam::AveragingMethod<Type>::write(const bool writeOnProc) const
142{
143 const pointMesh pointMesh_(mesh_);
144
145 // point volumes
146 Field<scalar> pointVolume(mesh_.nPoints(), Zero);
147
148 // output fields
150 (
151 IOobject::scopedName(this->name(), "cellValue"),
153 mesh_,
155 );
156 auto& cellValue = tcellValue.ref();
157
159 (
160 IOobject::scopedName(this->name(), "cellGrad"),
162 mesh_,
164 );
165 auto& cellGrad = tcellGrad.ref();
166
168 (
169 IOobject::scopedName(this->name(), "pointValue"),
171 pointMesh_,
173 );
174 auto& pointValue = tpointValue.ref();
175
177 (
178 IOobject::scopedName(this->name(), "pointGrad"),
180 pointMesh_,
182 );
183 auto& pointGrad = tpointGrad.ref();
184
185 // Barycentric coordinates of the tet vertices
187 tetCrds
188 ({
189 barycentric(1, 0, 0, 0),
190 barycentric(0, 1, 0, 0),
191 barycentric(0, 0, 1, 0),
192 barycentric(0, 0, 0, 1)
193 });
194
195 // tet-volume weighted sums
196 forAll(mesh_.C(), celli)
197 {
198 const List<tetIndices> cellTets =
200
201 forAll(cellTets, tetI)
202 {
203 const tetIndices& tetIs = cellTets[tetI];
204 const triFace triIs = tetIs.faceTriIs(mesh_);
205 const scalar v = tetIs.tet(mesh_).mag();
206
207 cellValue[celli] += v*interpolate(tetCrds[0], tetIs);
208 cellGrad[celli] += v*interpolateGrad(tetCrds[0], tetIs);
210 forAll(triIs, vertexI)
211 {
212 const label pointi = triIs[vertexI];
213
214 pointVolume[pointi] += v;
215 pointValue[pointi] += v*interpolate(tetCrds[vertexI], tetIs);
216 pointGrad[pointi] += v*interpolateGrad(tetCrds[vertexI], tetIs);
217 }
218 }
219 }
220
221 // average
222 cellValue.primitiveFieldRef() /= mesh_.V();
223 cellGrad.primitiveFieldRef() /= mesh_.V();
224 pointValue.primitiveFieldRef() /= pointVolume;
225 pointGrad.primitiveFieldRef() /= pointVolume;
226
227 // write
228 if (!cellValue.write(writeOnProc)) return false;
229 if (!cellGrad.write(writeOnProc)) return false;
230 if (!pointValue.write(writeOnProc)) return false;
231 if (!pointGrad.write(writeOnProc)) return false;
232
233 return true;
234}
235
236
237// ************************************************************************* //
Base class for lagrangian averaging methods.
const fvMesh & mesh_
The mesh on which the averaging is to be done.
virtual bool writeData(Ostream &) const
Dummy write.
virtual TypeGrad interpolateGrad(const barycentric &coordinates, const tetIndices &tetIs) const=0
AveragingMethod(const IOobject &io, const dictionary &dict, const fvMesh &mesh, const labelList &size)
Constructors.
static autoPtr< AveragingMethod< Type > > New(const IOobject &io, const dictionary &dict, const fvMesh &mesh)
Selector.
const dictionary & dict_
Protected data.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
virtual void updateGrad()
Protected member functions.
virtual void average()
Calculate the average.
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
constexpr FieldField() noexcept
Construct null.
Definition FieldField.C:101
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=PatchField< Type >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
void append(Field< Type > *ptr)
Definition UPtrList.H:877
label size() const noexcept
Definition UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
static List< tetIndices > cellTetIndices(const polyMesh &mesh, label cI)
Return the tet decomposition of the given cell, see.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition regIOobject.H:71
regIOobject(const IOobject &io, const bool isTimeObject=false)
Construct from IOobject. The optional flag adds special handling if the object is the top-level regIO...
Definition regIOobject.C:43
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition tetIndices.H:79
triFace faceTriIs(const polyMesh &mesh, const bool warn=true) const
The indices corresponding to the tri on the face for this tet. The normal of the tri points out of th...
Definition tetIndicesI.H:48
tetPointRef tet(const polyMesh &mesh) const
The tet geometry for this tet, where point0 is the cell centre.
scalar mag() const
Return volume.
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition triFace.H:68
A class for handling words, derived from Foam::string.
Definition word.H:66
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
OBJstream os(runTime.globalPath()/outputName)
const auto & io
auto & name
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
Definition List.H:62
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
errorManip< error > abort(error &err)
Definition errorManip.H:139
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition curveTools.C:75
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition barycentric.H:45
Macros to ease declaration of run-time selection tables.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299