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-------------------------------------------------------------------------------
10License
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
34template<class Type>
36(
38) const
39{
40 return
41 pos
42 (
43 mag(field)
44 - dimensionedScalar("t", field.dimensions(), threshold_)
45 );
46}
47
48
49template<class Type>
51(
52 const word& fieldName,
53 const bool calcMag
54)
55{
57
58 const VolFieldType* fieldPtr = obr_.findObject<VolFieldType>(fieldName);
59
60 if (!fieldPtr)
61 {
62 return;
63 }
64
65 auto extents = [this](const scalarField& mask, const vectorField& C)
66 {
68 forAll(mask, i)
69 {
70 if (mask[i] > 0.5)
71 {
72 extents.add(C[i] - C0_);
73 }
74 };
75
76 extents.reduce();
77
78 if (extents.empty())
79 {
80 extents.add(point::zero);
81 }
82
83 return extents;
84 };
85
86 Log << "field: " << fieldName << nl;
87
88 writeCurrentTime(file());
89
90 tmp<volScalarField> tmask = calcMask<Type>(*fieldPtr);
91 const volScalarField& mask = tmask();
92
93 // Internal field
94 if (internalField_)
95 {
96 boundBox bb(extents(mask, mesh_.C()));
97 Log << " internal field: " << bb << nl;
98 file() << bb;
99
100 this->setResult(fieldName + "_internal_min" , bb.min());
101 this->setResult(fieldName + "_internal_max", bb.max());
102 }
103
104 // Patches
105 for (const label patchi : patchIDs_)
106 {
107 const fvPatchScalarField& maskp = mask.boundaryField()[patchi];
108 boundBox bb(extents(maskp, maskp.patch().Cf()));
109 const word& patchName = maskp.patch().name();
110 Log << " patch " << patchName << ": " << bb << nl;
111 file() << bb;
112 this->setResult(fieldName + "_" + patchName + "_min", bb.min());
113 this->setResult(fieldName + "_" + patchName + "_max", bb.max());
114 }
115
116 Log << endl;
117 file() << endl;
118}
119
120
121// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Graphite solid properties.
Definition: C.H:53
Generic GeometricField class.
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
void reduce()
Parallel reduction of min/max values.
Definition: boundBox.C:184
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:86
bool empty() const
Bounding box is inverted, contains no points.
Definition: boundBoxI.H:62
void add(const boundBox &bb)
Extend to include the second box.
Definition: boundBoxI.H:191
scalar threshold_
Threshold value.
Definition: fieldExtents.H:197
void calcFieldExtents(const word &fieldName, const bool calcMag=false)
Main calculation.
tmp< volScalarField > calcMask(const GeometricField< Type, fvPatchField, volMesh > &field) const
Return the field mask.
Computes the magnitude of an input field.
Definition: mag.H:153
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:362
virtual const word & name() const
Return name.
Definition: fvPatch.H:173
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:119
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
rDeltaTY field()
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pos(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
A non-counting (dummy) refCount.
Definition: refCount.H:59