limitFieldsTemplates.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) 2019 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 "limitFields.H"
29#include "volFields.H"
30
31// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32
33template<class Type>
35{
37
38 auto* fieldPtr = getObjectPtr<VolFieldType>(fieldName);
39 if (!fieldPtr)
40 {
41 return false;
42 }
43
44 auto& field = *fieldPtr;
45
46 Log << " Limiting field " << fieldName << ":";
47
48 const dimensionedScalar eps("eps", field.dimensions(), ROOTVSMALL);
49
50 if (limit_ & MIN)
51 {
52 volScalarField mField(typeName + ":mag" + field.name(), mag(field));
53 Log << " min(|" << gMin(mField) << "|)";
54 //field.normalise();
55 field /= mag(field) + eps;
56 mField.max(dimensionedScalar("min", field.dimensions(), min_));
57 field *= mField;
58 }
59
60 if (limit_ & MAX)
61 {
62 volScalarField mField(typeName + ":mag" + field.name(), mag(field));
63 Log << " max(|" << gMax(mField) << "|)";
64 //field.normalise();
65 field /= mag(field) + eps;
66 mField.min(dimensionedScalar("max", field.dimensions(), max_));
67 field *= mField;
68 }
69
70 return true;
71}
72
73
74// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Generic GeometricField class.
void max(const dimensioned< Type > &dt)
Use the maximum of the field and specified value.
void min(const dimensioned< Type > &dt)
Use the minimum of the field and specified value.
@ MIN
limit by minimum value
Definition: limitFields.H:187
@ MAX
limit by maximum value
Definition: limitFields.H:188
limitType limit_
Limiting type.
Definition: limitFields.H:201
bool limitField(const word &fieldName)
Limit a field.
Computes the magnitude of an input field.
Definition: mag.H:153
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.
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)