faceSourceTemplates.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 "faceSource.H"
27 #include "surfaceFields.H"
28 #include "volFields.H"
29 #include "sampledSurface.H"
30 #include "interpolationCellPoint.H"
31 
32 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
33 
34 template<class Type>
35 bool Foam::fieldValues::faceSource::validField(const word& fieldName) const
36 {
39 
40  if (source_ != stSampledSurface && obr_.foundObject<sf>(fieldName))
41  {
42  return true;
43  }
44  else if (obr_.foundObject<vf>(fieldName))
45  {
46  return true;
47  }
48 
49  return false;
50 }
51 
52 
53 template<class Type>
55 (
56  const word& fieldName,
57  const bool mustGet,
58  const bool applyOrientation
59 ) const
60 {
63 
64  if (source_ != stSampledSurface && obr_.foundObject<sf>(fieldName))
65  {
66  return filterField(obr_.lookupObject<sf>(fieldName), applyOrientation);
67  }
68  else if (obr_.foundObject<vf>(fieldName))
69  {
70  const vf& fld = obr_.lookupObject<vf>(fieldName);
71 
72  if (surfacePtr_.valid())
73  {
74  if (surfacePtr_().interpolate())
75  {
76  const interpolationCellPoint<Type> interp(fld);
77  tmp<Field<Type> > tintFld(surfacePtr_().interpolate(interp));
78  const Field<Type>& intFld = tintFld();
79 
80  // Average
81  const faceList& faces = surfacePtr_().faces();
82  tmp<Field<Type> > tavg
83  (
85  );
86  Field<Type>& avg = tavg();
87 
88  forAll(faces, faceI)
89  {
90  const face& f = faces[faceI];
91  forAll(f, fp)
92  {
93  avg[faceI] += intFld[f[fp]];
94  }
95  avg[faceI] /= f.size();
96  }
97 
98  return tavg;
99  }
100  else
101  {
102  return surfacePtr_().sample(fld);
103  }
104  }
105  else
106  {
107  return filterField(fld, applyOrientation);
108  }
109  }
110 
111  if (mustGet)
112  {
114  << "Field " << fieldName << " not found in database"
115  << abort(FatalError);
116  }
117 
118  return tmp<Field<Type> >(new Field<Type>(0));
119 }
120 
121 
122 template<class Type>
124 (
125  const Field<Type>& values,
126  const vectorField& Sf,
127  const scalarField& weightField
128 ) const
129 {
130  Type result = pTraits<Type>::zero;
131  switch (operation_)
132  {
133  case opSum:
134  {
135  result = gSum(values);
136  break;
137  }
138  case opSumMag:
139  {
140  result = gSum(cmptMag(values));
141  break;
142  }
143  case opSumDirection:
144  {
146  << "Operation " << operationTypeNames_[operation_]
147  << " not available for values of type "
149  << exit(FatalError);
150 
151  result = pTraits<Type>::zero;
152  break;
153  }
154  case opSumDirectionBalance:
155  {
157  << "Operation " << operationTypeNames_[operation_]
158  << " not available for values of type "
160  << exit(FatalError);
161 
162  result = pTraits<Type>::zero;
163  break;
164  }
165  case opAverage:
166  {
167  label n = returnReduce(values.size(), sumOp<label>());
168  result = gSum(values)/(scalar(n) + ROOTVSMALL);
169  break;
170  }
171  case opWeightedAverage:
172  {
173  label wSize = returnReduce(weightField.size(), sumOp<label>());
174 
175  if (wSize > 0)
176  {
177  result =
178  gSum(weightField*values)/(gSum(weightField) + ROOTVSMALL);
179  }
180  else
181  {
182  label n = returnReduce(values.size(), sumOp<label>());
183  result = gSum(values)/(scalar(n) + ROOTVSMALL);
184  }
185  break;
186  }
187  case opAreaAverage:
188  {
189  const scalarField magSf(mag(Sf));
190 
191  result = gSum(magSf*values)/gSum(magSf);
192  break;
193  }
194  case opWeightedAreaAverage:
195  {
196  const scalarField magSf(mag(Sf));
197  label wSize = returnReduce(weightField.size(), sumOp<label>());
198 
199  if (wSize > 0)
200  {
201  result = gSum(weightField*magSf*values)/gSum(magSf*weightField);
202  }
203  else
204  {
205  result = gSum(magSf*values)/gSum(magSf);
206  }
207  break;
208  }
209  case opAreaIntegrate:
210  {
211  const scalarField magSf(mag(Sf));
212 
213  result = gSum(magSf*values);
214  break;
215  }
216  case opMin:
217  {
218  result = gMin(values);
219  break;
220  }
221  case opMax:
222  {
223  result = gMax(values);
224  break;
225  }
226  case opCoV:
227  {
228  const scalarField magSf(mag(Sf));
229 
230  const scalar gSumMagSf = gSum(magSf);
231 
232  Type meanValue = gSum(values*magSf)/gSumMagSf;
233 
234  const label nComp = pTraits<Type>::nComponents;
235 
236  for (direction d=0; d<nComp; ++d)
237  {
238  scalarField vals(values.component(d));
239  scalar mean = component(meanValue, d);
240  scalar& res = setComponent(result, d);
241 
242  res =
243  sqrt(gSum(magSf*sqr(vals - mean))/gSumMagSf)
244  /(mean + ROOTVSMALL);
245  }
246 
247  break;
248  }
249  default:
250  {
251  // Do nothing
252  }
253  }
254 
255  return result;
256 }
257 
258 
259 template<class Type>
261 (
262  const Field<Type>& values,
263  const vectorField& Sf,
264  const scalarField& weightField
265 ) const
266 {
267  return processSameTypeValues(values, Sf, weightField);
268 }
269 
270 
271 
272 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
273 
274 template<class Type>
276 (
277  const word& fieldName,
278  const scalarField& weightField,
279  const bool orient
280 )
281 {
282  const bool ok = validField<Type>(fieldName);
283 
284  if (ok)
285  {
286  Field<Type> values(setFieldValues<Type>(fieldName, true, orient));
287 
288  vectorField Sf;
289  if (surfacePtr_.valid())
290  {
291  // Get oriented Sf
292  Sf = surfacePtr_().Sf();
293  }
294  else
295  {
296  // Get oriented Sf
297  Sf = filterField(mesh().Sf(), true);
298  }
299 
300  // Write raw values on surface if specified
301  if (surfaceWriterPtr_.valid())
302  {
303  Field<Type> allValues(values);
304  combineFields(allValues);
305 
306  faceList faces;
308 
309  if (surfacePtr_.valid())
310  {
311  combineSurfaceGeometry(faces, points);
312  }
313  else
314  {
315  combineMeshGeometry(faces, points);
316  }
317 
318  if (Pstream::master())
319  {
320  fileName outputDir =
321  baseFileDir()/name_/"surface"/obr_.time().timeName();
322 
323  surfaceWriterPtr_->write
324  (
325  outputDir,
326  word(sourceTypeNames_[source_]) + "_" + sourceName_,
327  points,
328  faces,
329  fieldName,
330  allValues,
331  false
332  );
333  }
334  }
335 
336  // Apply scale factor
337  values *= scaleFactor_;
338 
339  Type result = processValues(values, Sf, weightField);
340 
341  file()<< tab << result;
342 
343  if (log_) Info
344  << " " << operationTypeNames_[operation_]
345  << "(" << sourceName_ << ") for " << fieldName
346  << " = " << result << endl;
347 
348  // Write state/results information
349  const word& opName = operationTypeNames_[operation_];
350  word resultName = opName + '(' + sourceName_ + ',' + fieldName + ')';
351  this->setResult(resultName, result);
352  }
353 
354  return ok;
355 }
356 
357 
358 template<class Type>
360 (
362  const bool applyOrientation
363 ) const
364 {
365  tmp<Field<Type> > tvalues(new Field<Type>(faceId_.size()));
366  Field<Type>& values = tvalues();
367 
368  forAll(values, i)
369  {
370  label faceI = faceId_[i];
371  label patchI = facePatchId_[i];
372  if (patchI >= 0)
373  {
374  values[i] = field.boundaryField()[patchI][faceI];
375  }
376  else
377  {
379  << type() << " " << name_ << ": "
380  << sourceTypeNames_[source_] << "(" << sourceName_ << "):"
381  << nl
382  << " Unable to process internal faces for volume field "
383  << field.name() << nl << abort(FatalError);
384  }
385  }
386 
387  if (applyOrientation)
388  {
389  forAll(values, i)
390  {
391  values[i] *= faceSign_[i];
392  }
393  }
394 
395  return tvalues;
396 }
397 
398 
399 template<class Type>
401 (
403  const bool applyOrientation
404 ) const
405 {
406  tmp<Field<Type> > tvalues(new Field<Type>(faceId_.size()));
407  Field<Type>& values = tvalues();
408 
409  forAll(values, i)
410  {
411  label faceI = faceId_[i];
412  label patchI = facePatchId_[i];
413  if (patchI >= 0)
414  {
415  values[i] = field.boundaryField()[patchI][faceI];
416  }
417  else
418  {
419  values[i] = field[faceI];
420  }
421  }
422 
423  if (applyOrientation)
424  {
425  forAll(values, i)
426  {
427  values[i] *= faceSign_[i];
428  }
429  }
430 
431  return tvalues;
432 }
433 
434 
435 // ************************************************************************* //
volFields.H
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::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::fieldValues::faceSource::validField
bool validField(const word &fieldName) const
Return true if the field name is valid.
Definition: faceSourceTemplates.C:35
interpolationCellPoint.H
Foam::fieldValues::faceSource::stSampledSurface
@ stSampledSurface
Definition: faceSource.H:317
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::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
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::fieldValues::faceSource::processValues
Type processValues(const Field< Type > &values, const vectorField &Sf, const scalarField &weightField) const
Apply the 'operation' to the values. Wrapper around.
Definition: faceSourceTemplates.C:261
Foam::interpolationCellPoint
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
Definition: interpolationCellPoint.H:48
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
sampledSurface.H
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::fieldValues::faceSource::setFieldValues
tmp< Field< Type > > setFieldValues(const word &fieldName, const bool mustGet=false, const bool applyOrientation=false) const
Return field values by looking up field name.
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
sf
volScalarField sf(fieldObject, mesh)
Foam::Field::component
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:644
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:399
Foam::fieldValues::faceSource::filterField
tmp< Field< Type > > filterField(const GeometricField< Type, fvsPatchField, surfaceMesh > &field, const bool applyOrientation) const
Filter a surface field according to faceIds.
Foam::fieldValues::faceSource::source_
sourceType source_
Source type.
Definition: faceSource.H:387
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::fieldValues::faceSource::processSameTypeValues
Type processSameTypeValues(const Field< Type > &values, const vectorField &Sf, const scalarField &weightField) const
Apply the 'operation' to the values. Operation has to.
Definition: faceSourceTemplates.C:124
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
f
labelList f(nPoints)
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fieldValues::faceSource::writeValues
bool writeValues(const word &fieldName, const scalarField &weightField, const bool orient)
Templated helper function to output field values.
Definition: faceSourceTemplates.C:276
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562
faceSource.H
Foam::fieldValue::obr_
const objectRegistry & obr_
Database this class is registered to.
Definition: fieldValue.H:74