cellSourceTemplates.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-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 "cellSource.H"
27 #include "volFields.H"
28 
29 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
30 
31 template<class Type>
32 bool Foam::fieldValues::cellSource::validField(const word& fieldName) const
33 {
35 
36  if (obr_.foundObject<vf>(fieldName))
37  {
38  return true;
39  }
40 
41  return false;
42 }
43 
44 
45 template<class Type>
47 (
48  const word& fieldName,
49  const bool mustGet
50 ) const
51 {
53 
54  if (obr_.foundObject<vf>(fieldName))
55  {
56  return filterField(obr_.lookupObject<vf>(fieldName));
57  }
58 
59  if (mustGet)
60  {
62  << "Field " << fieldName << " not found in database"
63  << abort(FatalError);
64  }
65 
66  return tmp<Field<Type> >(new Field<Type>(0.0));
67 }
68 
69 
70 template<class Type>
72 (
73  const Field<Type>& values,
74  const scalarField& V,
75  const scalarField& weightField
76 ) const
77 {
78  Type result = pTraits<Type>::zero;
79  switch (operation_)
80  {
81  case opSum:
82  {
83  result = gSum(values);
84  break;
85  }
86  case opSumMag:
87  {
88  result = gSum(cmptMag(values));
89  break;
90  }
91  case opAverage:
92  {
93  label n = returnReduce(values.size(), sumOp<label>());
94  result = gSum(values)/(scalar(n) + ROOTVSMALL);
95  break;
96  }
97  case opWeightedAverage:
98  {
99  label wSize = returnReduce(weightField.size(), sumOp<label>());
100 
101  if (wSize > 0)
102  {
103  result =
104  gSum(weightField*values)/(gSum(weightField) + ROOTVSMALL);
105  }
106  else
107  {
108  label n = returnReduce(values.size(), sumOp<label>());
109  result = gSum(values)/(scalar(n) + ROOTVSMALL);
110  }
111  break;
112  }
113  case opVolAverage:
114  {
115  result = gSum(values*V)/(gSum(V) + ROOTVSMALL);
116  break;
117  }
118  case opWeightedVolAverage:
119  {
120  result = gSum(weightField*V*values)/gSum(weightField*V);
121  break;
122  }
123  case opVolIntegrate:
124  {
125  result = gSum(V*values);
126  break;
127  }
128  case opMin:
129  {
130  result = gMin(values);
131  break;
132  }
133  case opMax:
134  {
135  result = gMax(values);
136  break;
137  }
138  case opCoV:
139  {
140  const scalar sumV = gSum(V);
141 
142  Type meanValue = gSum(V*values)/sumV;
143 
144  const label nComp = pTraits<Type>::nComponents;
145 
146  for (direction d=0; d<nComp; ++d)
147  {
148  scalarField vals(values.component(d));
149  scalar mean = component(meanValue, d);
150  scalar& res = setComponent(result, d);
151 
152  res = sqrt(gSum(V*sqr(vals - mean))/sumV)/(mean + ROOTVSMALL);
153  }
154 
155  break;
156  }
157  default:
158  {
159  // Do nothing
160  }
161  }
162 
163  return result;
164 }
165 
166 
167 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
168 
169 template<class Type>
171 (
172  const word& fieldName,
173  const scalarField& weightField
174 )
175 {
176  const bool ok = validField<Type>(fieldName);
177 
178  if (ok)
179  {
180  Field<Type> values(setFieldValues<Type>(fieldName));
181  scalarField V(filterField(mesh().V()));
182 
183  if (valueOutput_)
184  {
185  Field<Type> allValues(values);
186  combineFields(allValues);
187 
188  if (Pstream::master())
189  {
191  (
192  IOobject
193  (
194  fieldName + "_" + sourceTypeNames_[source_] + "-"
195  + sourceName_,
196  obr_.time().timeName(),
197  obr_,
200  ),
201  allValues
202  ).write();
203  }
204  }
205 
206  // Apply scale factor
207  values *= scaleFactor_;
208 
209  Type result = processValues(values, V, weightField);
210 
211  file()<< tab << result;
212 
213  if (log_) Info
214  << " " << operationTypeNames_[operation_]
215  << "(" << sourceName_ << ") of " << fieldName
216  << " = " << result << endl;
217 
218  // write state/results information
219  const word& opName = operationTypeNames_[operation_];
220  word resultName = opName + '(' + sourceName_ + ',' + fieldName + ')';
221  this->setResult(resultName, result);
222  }
223 
224  return ok;
225 }
226 
227 
228 template<class Type>
230 (
231  const Field<Type>& field
232 ) const
233 {
234  return tmp<Field<Type> >(new Field<Type>(field, cellId_));
235 }
236 
237 
238 // ************************************************************************* //
volFields.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::setComponent
label & setComponent(label &l, const direction)
Definition: label.H:79
Foam::component
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Definition: FieldFieldFunctions.C:41
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::fieldValues::cellSource::setFieldValues
tmp< Field< Type > > setFieldValues(const word &fieldName, const bool mustGet=false) const
Insert field values into values list.
Foam::regIOobject::write
virtual bool write() const
Write using setting from DB.
Definition: regIOobjectWrite.C:126
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::fieldValues::cellSource::processValues
Type processValues(const Field< Type > &values, const scalarField &V, const scalarField &weightField) const
Apply the 'operation' to the values.
Definition: cellSourceTemplates.C:72
Foam::cmptMag
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:359
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::label
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
Foam::Field< Type >
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::Info
messageStream Info
Foam::fieldValues::cellSource::writeValues
bool writeValues(const word &fieldName, const scalarField &weightField)
Templated helper function to output field values.
Definition: cellSourceTemplates.C:171
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:399
Foam::objectRegistry::foundObject
bool foundObject(const word &name) const
Is the named Type found?
Definition: objectRegistryTemplates.C:142
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::sumOp
Definition: ops.H:162
Foam::tab
static const char tab
Definition: Ostream.H:259
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::fieldValues::cellSource::filterField
tmp< Field< Type > > filterField(const Field< Type > &field) const
Filter a field according to cellIds.
Foam::fieldValues::cellSource::validField
bool validField(const word &fieldName) const
Return true if the field name is valid.
Definition: cellSourceTemplates.C:32
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
cellSource.H
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562
Foam::fieldValue::obr_
const objectRegistry & obr_
Database this class is registered to.
Definition: fieldValue.H:74