fieldExtentsTemplates.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2018 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "fieldExtents.H"
29 #include "volFields.H"
30 #include "boundBox.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 template<class Type>
36 (
38 ) const
39 {
40  return
41  pos
42  (
43  mag(field)
44  - dimensionedScalar("t", field.dimensions(), threshold_)
45  );
46 }
47 
48 
49 template<class Type>
51 (
52  const word& fieldName,
53  const bool calcMag
54 )
55 {
57 
58  const VolFieldType* fieldPtr =
59  obr_.findObject<VolFieldType>(fieldName);
60 
61  if (!fieldPtr)
62  {
63  return;
64  }
65 
66  auto extents = [this](const scalarField& mask, const vectorField& C)
67  {
68  boundBox extents(boundBox::invertedBox);
69  forAll(mask, i)
70  {
71  if (mask[i] > 0.5)
72  {
73  extents.add(C[i] - C0_);
74  }
75  };
76 
77  extents.reduce();
78 
79  if (extents.empty())
80  {
81  extents.add(point::zero);
82  }
83 
84  return extents;
85  };
86 
87  Log << "field: " << fieldName << nl;
88 
89  writeCurrentTime(file());
90 
91  tmp<volScalarField> tmask = calcMask<Type>(*fieldPtr);
92  const volScalarField& mask = tmask();
93 
94  // Internal field
95  if (internalField_)
96  {
97  boundBox bb(extents(mask, mesh_.C()));
98  Log << " internal field: " << bb << nl;
99  file() << bb;
100 
101  this->setResult(fieldName + "_internal_min" , bb.min());
102  this->setResult(fieldName + "_internal_max", bb.max());
103  }
104 
105  // Patches
106  for (const label patchi : patchIDs_)
107  {
108  const fvPatchScalarField& maskp = mask.boundaryField()[patchi];
109  boundBox bb(extents(maskp, maskp.patch().Cf()));
110  const word& patchName = maskp.patch().name();
111  Log << " patch " << patchName << ": " << bb << nl;
112  file() << bb;
113  this->setResult(fieldName + "_" + patchName + "_min", bb.min());
114  this->setResult(fieldName + "_" + patchName + "_max", bb.max());
115  }
116 
117  Log << endl;
118  file() << endl;
119 }
120 
121 
122 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
Foam::boundBox::reduce
void reduce()
Parallel reduction of min/max values.
Definition: boundBox.C:184
Log
#define Log
Definition: PDRblock.C:35
Foam::functionObjects::fieldExtents::calcMask
tmp< volScalarField > calcMask(const GeometricField< Type, fvPatchField, volMesh > &field) const
Return the field mask.
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
fieldExtents.H
Foam::functionObjects::fieldExtents::calcFieldExtents
void calcFieldExtents(const word &fieldName, const bool calcMag=false)
Main calculation.
Definition: fieldExtentsTemplates.C:51
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
Foam::fvPatch::name
virtual const word & name() const
Return name.
Definition: fvPatch.H:167
C
volScalarField & C
Definition: readThermalProperties.H:102
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::Field< scalar >
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
field
rDeltaTY field()
Foam::fvPatch::Cf
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:113
Foam::boundBox::empty
bool empty() const
Bounding box is inverted, contains no points.
Definition: boundBoxI.H:62
boundBox.H
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:343
Foam::C
Graphite solid properties.
Definition: C.H:50
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::boundBox::add
void add(const boundBox &bb)
Extend to include the second box.
Definition: boundBoxI.H:191
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177