volRegionTemplates.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) 2011-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 "volRegion.H"
27 #include "volFields.H"
28 
29 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const word& fieldName
35 ) const
36 {
38 
39  if (obr_.foundObject<vf>(fieldName))
40  {
41  return true;
42  }
43 
44  return false;
45 }
46 
47 
48 template<class Type>
51 (
52  const word& fieldName,
53  const bool mustGet
54 ) const
55 {
57 
58  if (obr_.foundObject<vf>(fieldName))
59  {
60  return filterField(obr_.lookupObject<vf>(fieldName));
61  }
62 
63  if (mustGet)
64  {
66  << "Field " << fieldName << " not found in database"
67  << abort(FatalError);
68  }
69 
70  return tmp<Field<Type>>(new Field<Type>(0.0));
71 }
72 
73 
74 template<class Type>
76 (
77  const Field<Type>& values,
78  const scalarField& V,
79  const scalarField& weightField
80 ) const
81 {
82  Type result = Zero;
83  switch (operation_)
84  {
85  case opSum:
86  {
87  result = sum(values);
88  break;
89  }
90  case opSumMag:
91  {
92  result = sum(cmptMag(values));
93  break;
94  }
95  case opAverage:
96  {
97  result = sum(values)/values.size();
98  break;
99  }
100  case opWeightedAverage:
101  {
102  result = sum(weightField*values)/sum(weightField);
103  break;
104  }
105  case opVolAverage:
106  {
107  result = sum(V*values)/sum(V);
108  break;
109  }
110  case opWeightedVolAverage:
111  {
112  result = sum(weightField*V*values)/sum(weightField*V);
113  break;
114  }
115  case opVolIntegrate:
116  {
117  result = sum(V*values);
118  break;
119  }
120  case opMin:
121  {
122  result = min(values);
123  break;
124  }
125  case opMax:
126  {
127  result = max(values);
128  break;
129  }
130  case opCoV:
131  {
132  Type meanValue = sum(values*V)/sum(V);
133 
134  const label nComp = pTraits<Type>::nComponents;
135 
136  for (direction d=0; d<nComp; ++d)
137  {
138  scalarField vals(values.component(d));
139  scalar mean = component(meanValue, d);
140  scalar& res = setComponent(result, d);
141 
142  res = sqrt(sum(V*sqr(vals - mean))/sum(V))/mean;
143  }
144 
145  break;
146  }
147  case opNone:
148  {}
149  }
150 
151  return result;
152 }
153 
154 
155 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
156 
157 template<class Type>
159 (
160  const word& fieldName
161 )
162 {
163  const bool ok = validField<Type>(fieldName);
164 
165  if (ok)
166  {
167  Field<Type> values(setFieldValues<Type>(fieldName));
168  scalarField V(filterField(mesh().V()));
169  scalarField weightField(values.size(), 1.0);
170 
171  if (weightFieldName_ != "none")
172  {
173  weightField = setFieldValues<scalar>(weightFieldName_, true);
174  }
175 
176  // Combine onto master
177  combineFields(values);
178  combineFields(V);
179  combineFields(weightField);
180 
181  if (Pstream::master())
182  {
183  Type result = processValues(values, V, weightField);
184 
185  // Add to result dictionary, over-writing any previous entry
186  resultDict_.add(fieldName, result, true);
187 
188  if (writeFields_)
189  {
191  (
192  IOobject
193  (
194  fieldName + "_" + regionTypeNames_[regionType_] + "-"
195  + regionName_,
196  obr_.time().timeName(),
197  obr_,
198  IOobject::NO_READ,
199  IOobject::NO_WRITE
200  ),
201  weightField*values
202  ).write();
203  }
204 
205 
206  file()<< tab << result;
207 
208  Log << " " << operationTypeNames_[operation_]
209  << "(" << regionName_ << ") of " << fieldName
210  << " = " << result << endl;
211  }
212  }
213 
214  return ok;
215 }
216 
217 
218 template<class Type>
221 (
222  const Field<Type>& field
223 ) const
224 {
225  return tmp<Field<Type>>(new Field<Type>(field, cellId_));
226 }
227 
228 
229 // ************************************************************************* //
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
uint8_t direction
Definition: direction.H:46
static const char tab
Definition: Ostream.H:261
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< Field< Type > > filterField(const Field< Type > &field) const
Filter a field according to cellIds.
void size(const label)
Override size to be inconsistent with allocated storage.
bool writeValues(const word &fieldName)
Templated helper function to output field values.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Traits class for primitives.
Definition: pTraits.H:50
Generic GeometricField class.
bool validField(const word &fieldName) const
Return true if the field name is valid.
tmp< Field< Type > > setFieldValues(const word &fieldName, const bool mustGet=false) const
Insert field values into values list.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
dynamicFvMesh & mesh
Pre-declare SubField and related Field type.
Definition: Field.H:57
A class for handling words, derived from string.
Definition: word.H:59
static const zero Zero
Definition: zero.H:91
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:650
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual bool write() const
Write using setting from DB.
#define Log
Report write to Foam::Info if the local log switch is true.
A class for managing temporary objects.
Definition: PtrList.H:54
label & setComponent(label &l, const direction)
Definition: label.H:79
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
Type processValues(const Field< Type > &values, const scalarField &V, const scalarField &weightField) const
Apply the &#39;operation&#39; to the values.