surfaceRegionTemplates.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 "surfaceRegion.H"
27 #include "surfaceFields.H"
28 #include "volFields.H"
29 #include "sampledSurface.H"
30 #include "surfaceWriter.H"
31 #include "interpolationCellPoint.H"
32 
33 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  const word& fieldName
39 ) const
40 {
43 
44  if (regionType_ != stSampledSurface && obr_.foundObject<sf>(fieldName))
45  {
46  return true;
47  }
48  else if (obr_.foundObject<vf>(fieldName))
49  {
50  return true;
51  }
52 
53  return false;
54 }
55 
56 
57 template<class Type>
60 (
61  const word& fieldName,
62  const bool mustGet,
63  const bool applyOrientation
64 ) const
65 {
68 
69  if (regionType_ != stSampledSurface && obr_.foundObject<sf>(fieldName))
70  {
71  return filterField(obr_.lookupObject<sf>(fieldName), applyOrientation);
72  }
73  else if (obr_.foundObject<vf>(fieldName))
74  {
75  const vf& fld = obr_.lookupObject<vf>(fieldName);
76 
77  if (surfacePtr_.valid())
78  {
79  if (surfacePtr_().interpolate())
80  {
81  const interpolationCellPoint<Type> interp(fld);
82  tmp<Field<Type>> tintFld(surfacePtr_().interpolate(interp));
83  const Field<Type>& intFld = tintFld();
84 
85  // Average
86  const faceList& faces = surfacePtr_().faces();
87  tmp<Field<Type>> tavg
88  (
89  new Field<Type>(faces.size(), Zero)
90  );
91  Field<Type>& avg = tavg.ref();
92 
93  forAll(faces, facei)
94  {
95  const face& f = faces[facei];
96  forAll(f, fp)
97  {
98  avg[facei] += intFld[f[fp]];
99  }
100  avg[facei] /= f.size();
101  }
102 
103  return tavg;
104  }
105  else
106  {
107  return surfacePtr_().sample(fld);
108  }
109  }
110  else
111  {
112  return filterField(fld, applyOrientation);
113  }
114  }
115 
116  if (mustGet)
117  {
119  << "Field " << fieldName << " not found in database"
120  << abort(FatalError);
121  }
122 
123  return tmp<Field<Type>>(new Field<Type>(0));
124 }
125 
126 
127 template<class Type>
129 (
130  const Field<Type>& values,
131  const vectorField& Sf,
132  const scalarField& weightField
133 ) const
134 {
135  Type result = Zero;
136  switch (operation_)
137  {
138  case opSum:
139  {
140  result = sum(values);
141  break;
142  }
143  case opSumMag:
144  {
145  result = sum(cmptMag(values));
146  break;
147  }
148  case opSumDirection:
149  {
151  << "Operation " << operationTypeNames_[operation_]
152  << " not available for values of type "
154  << exit(FatalError);
155 
156  result = Zero;
157  break;
158  }
159  case opSumDirectionBalance:
160  {
162  << "Operation " << operationTypeNames_[operation_]
163  << " not available for values of type "
165  << exit(FatalError);
166 
167  result = Zero;
168  break;
169  }
170  case opAverage:
171  {
172  result = sum(values)/values.size();
173  break;
174  }
175  case opWeightedAverage:
176  {
177  if (weightField.size())
178  {
179  result = sum(weightField*values)/sum(weightField);
180  }
181  else
182  {
183  result = sum(values)/values.size();
184  }
185  break;
186  }
187  case opAreaAverage:
188  {
189  const scalarField magSf(mag(Sf));
190 
191  result = sum(magSf*values)/sum(magSf);
192  break;
193  }
194  case opWeightedAreaAverage:
195  {
196  const scalarField magSf(mag(Sf));
197 
198  if (weightField.size())
199  {
200  result = sum(weightField*magSf*values)/sum(magSf*weightField);
201  }
202  else
203  {
204  result = sum(magSf*values)/sum(magSf);
205  }
206  break;
207  }
208  case opAreaIntegrate:
209  {
210  const scalarField magSf(mag(Sf));
211 
212  result = sum(magSf*values);
213  break;
214  }
215  case opMin:
216  {
217  result = min(values);
218  break;
219  }
220  case opMax:
221  {
222  result = max(values);
223  break;
224  }
225  case opCoV:
226  {
227  const scalarField magSf(mag(Sf));
228 
229  Type meanValue = sum(values*magSf)/sum(magSf);
230 
231  const label nComp = pTraits<Type>::nComponents;
232 
233  for (direction d=0; d<nComp; ++d)
234  {
235  scalarField vals(values.component(d));
236  scalar mean = component(meanValue, d);
237  scalar& res = setComponent(result, d);
238 
239  res = sqrt(sum(magSf*sqr(vals - mean))/sum(magSf))/mean;
240  }
241 
242  break;
243  }
244  case opAreaNormalAverage:
245  {}
246  case opAreaNormalIntegrate:
247  {}
248  case opNone:
249  {}
250  }
251 
252  return result;
253 }
254 
255 
256 template<class Type>
258 (
259  const Field<Type>& values,
260  const vectorField& Sf,
261  const scalarField& weightField
262 ) const
263 {
264  return processSameTypeValues(values, Sf, weightField);
265 }
266 
267 
268 
269 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
270 
271 template<class Type>
273 (
274  const word& fieldName,
275  const scalarField& weightField,
276  const bool orient
277 )
278 {
279  const bool ok = validField<Type>(fieldName);
280 
281  if (ok)
282  {
283  Field<Type> values(getFieldValues<Type>(fieldName, true, orient));
284 
285  vectorField Sf;
286  if (surfacePtr_.valid())
287  {
288  // Get oriented Sf
289  Sf = surfacePtr_().Sf();
290  }
291  else
292  {
293  // Get oriented Sf
294  Sf = filterField(mesh().Sf(), true);
295  }
296 
297  // Combine onto master
298  combineFields(values);
299  combineFields(Sf);
300 
301  // Write raw values on surface if specified
302  if (surfaceWriterPtr_.valid())
303  {
304  faceList faces;
306 
307  if (surfacePtr_.valid())
308  {
309  combineSurfaceGeometry(faces, points);
310  }
311  else
312  {
313  combineMeshGeometry(faces, points);
314  }
315 
316  if (Pstream::master())
317  {
318  fileName outputDir =
319  baseFileDir()/name()/"surface"/obr_.time().timeName();
320 
321  surfaceWriterPtr_->write
322  (
323  outputDir,
324  word(regionTypeNames_[regionType_]) + "_" + regionName_,
325  points,
326  faces,
327  fieldName,
328  values,
329  false
330  );
331  }
332  }
333 
334 
335  // Apply scale factor
336  values *= scaleFactor_;
337 
338  if (Pstream::master())
339  {
340  Type result = processValues(values, Sf, weightField);
341 
342  // Add to result dictionary, over-writing any previous entry
343  resultDict_.add(fieldName, result, true);
344 
345  file()<< tab << result;
346 
347  Log << " " << operationTypeNames_[operation_]
348  << "(" << regionName_ << ") of " << fieldName
349  << " = " << result << endl;
350  }
351  }
352 
353  return ok;
354 }
355 
356 
357 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.ref();
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  << regionTypeNames_[regionType_] << "(" << regionName_ << "):"
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>
402 (
404  const bool applyOrientation
405 ) const
406 {
407  tmp<Field<Type>> tvalues(new Field<Type>(faceId_.size()));
408  Field<Type>& values = tvalues.ref();
409 
410  forAll(values, i)
411  {
412  label facei = faceId_[i];
413  label patchi = facePatchId_[i];
414  if (patchi >= 0)
415  {
416  values[i] = field.boundaryField()[patchi][facei];
417  }
418  else
419  {
420  values[i] = field[facei];
421  }
422  }
423 
424  if (applyOrientation)
425  {
426  forAll(values, i)
427  {
428  values[i] *= faceSign_[i];
429  }
430  }
431 
432  return tvalues;
433 }
434 
435 
436 // ************************************************************************* //
Foam::surfaceFields.
bool validField(const word &fieldName) const
Return true if the field name is valid.
#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
uint8_t direction
Definition: direction.H:46
A class for handling file names.
Definition: fileName.H:69
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
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 > > getFieldValues(const word &fieldName, const bool mustGet=false, const bool applyOrientation=false) const
Return field values by looking up field name.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
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.
tmp< surfaceScalarField > interpolate(const RhoType &rho)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
dynamicFvMesh & mesh
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))
const pointField & points
const Boundary & boundaryField() const
Return const-reference to the boundary field.
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
Type processValues(const Field< Type > &values, const vectorField &Sf, const scalarField &weightField) const
Apply the &#39;operation&#39; to the values. Wrapper around.
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
static const char nl
Definition: Ostream.H:262
tmp< Field< Type > > filterField(const GeometricField< Type, fvsPatchField, surfaceMesh > &field, const bool applyOrientation) const
Filter a surface field according to faceIds.
bool writeValues(const word &fieldName, const scalarField &weightField, const bool orient)
Templated helper function to output field values.
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:650
labelList f(nPoints)
Type processSameTypeValues(const Field< Type > &values, const vectorField &Sf, const scalarField &weightField) const
Apply the &#39;operation&#39; to the values. Operation has to.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
volScalarField sf(fieldObject, mesh)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label patchi
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
#define Log
Report write to Foam::Info if the local log switch is true.
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:54
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
label & setComponent(label &l, const direction)
Definition: label.H:79
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
const word & name() const
Return name.
Definition: IOobject.H:260