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 | Copyright (C) 2013-2016 OpenFOAM Foundation
6  \\/ M anipulation |
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 "AveragingMethod.H"
27 #include "runTimeSelectionTables.H"
28 #include "pointMesh.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class Type>
34 (
35  const IOobject& io,
36  const dictionary& dict,
37  const fvMesh& mesh,
38  const labelList& size
39 )
40 :
41  regIOobject(io),
43  dict_(dict),
44  mesh_(mesh)
45 {
46  forAll(size, i)
47  {
49  (
50  new Field<Type>(size[i], Zero)
51  );
52  }
53 }
54 
55 
56 template<class Type>
58 (
59  const AveragingMethod<Type>& am
60 )
61 :
62  regIOobject(am),
64  dict_(am.dict_),
65  mesh_(am.mesh_)
66 {}
67 
68 
69 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
70 
71 template<class Type>
74 (
75  const IOobject& io,
76  const dictionary& dict,
77  const fvMesh& mesh
78 )
79 {
80  word averageType(dict.lookup(typeName));
81 
82  //Info<< "Selecting averaging method "
83  // << averageType << endl;
84 
85  typename dictionaryConstructorTable::iterator cstrIter =
86  dictionaryConstructorTablePtr_->find(averageType);
87 
88  if (cstrIter == dictionaryConstructorTablePtr_->end())
89  {
91  << "Unknown averaging method " << averageType
92  << ", constructor not in hash table" << nl << nl
93  << " Valid averaging methods are:" << nl
94  << dictionaryConstructorTablePtr_->sortedToc()
95  << abort(FatalError);
96  }
97 
98  return autoPtr<AveragingMethod<Type>>(cstrIter()(io, dict, mesh));
99 }
100 
101 
102 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
103 
104 template<class Type>
106 {}
107 
108 
109 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
110 
111 template<class Type>
113 {
114  // do nothing
115 }
116 
117 
118 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119 
120 template<class Type>
122 {
123  updateGrad();
124 }
125 
126 
127 template<class Type>
129 (
130  const AveragingMethod<scalar>& weight
131 )
132 {
133  updateGrad();
134 
135  *this /= max(weight, SMALL);
136 }
137 
138 
139 template<class Type>
141 {
142  return os.good();
143 }
144 
145 
146 template<class Type>
148 {
149  const pointMesh pointMesh_(mesh_);
150 
151  // point volumes
152  Field<scalar> pointVolume(mesh_.nPoints(), 0);
153 
154  // output fields
156  (
157  IOobject
158  (
159  this->name() + ":cellValue",
160  this->time().timeName(),
161  mesh_
162  ),
163  mesh_,
164  dimensioned<Type>("zero", dimless, Zero)
165  );
167  (
168  IOobject
169  (
170  this->name() + ":cellGrad",
171  this->time().timeName(),
172  mesh_
173  ),
174  mesh_,
176  );
178  (
179  IOobject
180  (
181  this->name() + ":pointValue",
182  this->time().timeName(),
183  mesh_
184  ),
185  pointMesh_,
186  dimensioned<Type>("zero", dimless, Zero)
187  );
189  (
190  IOobject
191  (
192  this->name() + ":pointGrad",
193  this->time().timeName(),
194  mesh_
195  ),
196  pointMesh_,
198  );
199 
200  // tet-volume weighted sums
201  forAll(mesh_.C(), celli)
202  {
203  const List<tetIndices> cellTets =
204  polyMeshTetDecomposition::cellTetIndices(mesh_, celli);
205 
206  forAll(cellTets, tetI)
207  {
208  const tetIndices& tetIs = cellTets[tetI];
209  const scalar v = tetIs.tet(mesh_).mag();
210 
211  cellValue[celli] += v*interpolate(mesh_.C()[celli], tetIs);
212  cellGrad[celli] += v*interpolateGrad(mesh_.C()[celli], tetIs);
213 
214  const face& f = mesh_.faces()[tetIs.face()];
215  labelList vertices(3);
216  vertices[0] = f[tetIs.faceBasePt()];
217  vertices[1] = f[tetIs.facePtA()];
218  vertices[2] = f[tetIs.facePtB()];
219 
220  forAll(vertices, vertexI)
221  {
222  const label pointi = vertices[vertexI];
223 
224  pointVolume[pointi] += v;
225  pointValue[pointi] +=
226  v*interpolate(mesh_.points()[pointi], tetIs);
227  pointGrad[pointi] +=
228  v*interpolateGrad(mesh_.points()[pointi], tetIs);
229  }
230  }
231  }
232 
233  // average
234  cellValue.primitiveFieldRef() /= mesh_.V();
235  cellGrad.primitiveFieldRef() /= mesh_.V();
236  pointValue.primitiveFieldRef() /= pointVolume;
237  pointGrad.primitiveFieldRef() /= pointVolume;
238 
239  // write
240  if (!cellValue.write()) return false;
241  if (!cellGrad.write()) return false;
242  if (!pointValue.write()) return false;
243  if (!pointGrad.write()) return false;
244 
245  return true;
246 }
247 
248 
249 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
label facePtA() const
Return face point A.
Definition: tetIndicesI.H:48
static autoPtr< AveragingMethod< Type > > New(const IOobject &io, const dictionary &dict, const fvMesh &mesh)
Selector.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
Generic GeometricField class.
Generic dimensioned Type class.
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Generic field type.
Definition: FieldField.H:51
label faceBasePt() const
Return the face base point.
Definition: tetIndicesI.H:42
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
virtual bool write() const
Write using setting from DB.
dynamicFvMesh & mesh
label facePtB() const
Return face point B.
Definition: tetIndicesI.H:54
Pre-declare SubField and related Field type.
Definition: Field.H:57
A class for handling words, derived from string.
Definition: word.H:59
Base class for lagrangian averaging methods.
Definition: MPPICParcel.H:55
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:91
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
label face() const
Return the face.
Definition: tetIndicesI.H:36
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const fvMesh & mesh_
The mesh on which the averaging is to be done.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
static const char nl
Definition: Ostream.H:262
const dictionary & dict_
Protected data.
scalar mag() const
Return volume.
Definition: tetrahedronI.H:171
labelList f(nPoints)
virtual void updateGrad()
Protected member functions.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
virtual void average()
Calculate the average.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual bool write() const
Write using setting from DB.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:60
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
Macros to ease declaration of run-time selection tables.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
AveragingMethod(const IOobject &io, const dictionary &dict, const fvMesh &mesh, const labelList &size)
Constructors.
tetPointRef tet(const polyMesh &mesh) const
Return the geometry corresponding to this tet from the.
Definition: tetIndicesI.H:66
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451
virtual bool writeData(Ostream &) const
Dummy write.
virtual ~AveragingMethod()
Destructor.