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-2023 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 "limitFields.H"
29 #include "volFields.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 template<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 (withBounds_ & limitType::CLAMP_MIN)
51  {
52  auto tmField = volScalarField::New
53  (
54  IOobject::scopedName(typeName, "mag" + field.name()),
56  mag(field)
57  );
58  auto& mField = tmField.ref();
59 
60  Log << " min(|" << gMin(mField) << "|)";
61  //field.normalise();
62  field /= mag(field) + eps;
63  mField.clamp_min(min_);
64  field *= tmField;
65  }
66 
67  if (withBounds_ & limitType::CLAMP_MAX)
68  {
69  auto tmField = volScalarField::New
70  (
71  IOobject::scopedName(typeName, "mag" + field.name()),
73  mag(field)
74  );
75  auto& mField = tmField.ref();
76 
77  Log << " max(|" << gMax(mField) << "|)";
78  //field.normalise();
79  field /= mag(field) + eps;
80  mField.clamp_max(max_);
81  field *= tmField;
82  }
83 
84  return true;
85 }
86 
87 
88 // ************************************************************************* //
scalar min_
Minimum limit.
Definition: limitFields.H:212
rDeltaTY field()
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Type gMin(const FieldField< Field, Type > &f)
Generic GeometricField class.
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
bool limitField(const word &fieldName)
Limit a field.
A class for handling words, derived from Foam::string.
Definition: word.H:63
scalar max_
Maximum limit.
Definition: limitFields.H:217
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
Type gMax(const FieldField< Field, Type > &f)
limitType withBounds_
Limiting type.
Definition: limitFields.H:207
#define Log
Definition: PDRblock.C:28
Do not request registration (bool: false)