limitMag.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2016-2022 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 "limitMag.H"
27 #include "volFields.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace fv
35 {
36  defineTypeNameAndDebug(limitMag, 0);
38  (
39  fvConstraint,
40  limitMag,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::fv::limitMag::readCoeffs()
50 {
51  fieldName_ = coeffs().lookup<word>("field");
52  max_ = coeffs().lookup<scalar>("max");
53 }
54 
55 
56 template<class Type>
57 inline bool Foam::fv::limitMag::constrainType
58 (
59  GeometricField<Type, fvPatchField, volMesh>& psi
60 ) const
61 {
62  const scalar maxSqrPsi = sqr(max_);
63 
64  Field<Type>& psiif = psi.primitiveFieldRef();
65 
66  const labelList& cells = set_.cells();
67 
68  forAll(cells, i)
69  {
70  const label celli = cells[i];
71 
72  const scalar magSqrPsii = magSqr(psiif[celli]);
73 
74  if (magSqrPsii > maxSqrPsi)
75  {
76  psiif[celli] *= sqrt(maxSqrPsi/magSqrPsii);
77  }
78  }
79 
80  // handle boundaries in the case of 'all'
81  if (set_.selectionMode() == fvCellSet::selectionModeType::all)
82  {
84  psi.boundaryFieldRef();
85 
86  forAll(psibf, patchi)
87  {
88  fvPatchField<Type>& psip = psibf[patchi];
89 
90  if (!psip.fixesValue())
91  {
92  forAll(psip, facei)
93  {
94  const scalar magSqrPsii = magSqr(psip[facei]);
95 
96  if (magSqrPsii > maxSqrPsi)
97  {
98  psip[facei] *= sqrt(maxSqrPsi/magSqrPsii);
99  }
100  }
101  }
102  }
103  }
104 
105  return cells.size();
106 }
107 
108 
109 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
110 
112 (
113  const word& name,
114  const word& modelType,
115  const dictionary& dict,
116  const fvMesh& mesh
117 )
118 :
119  fvConstraint(name, modelType, dict, mesh),
120  set_(coeffs(), mesh),
121  fieldName_(word::null),
122  max_(vGreat)
123 {
124  readCoeffs();
125 }
126 
127 
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129 
131 {
132  return wordList(1, fieldName_);
133 }
134 
135 
137 (
140 );
141 
142 
144 {
145  set_.movePoints();
146  return true;
147 }
148 
149 
151 {
152  set_.topoChange(map);
153 }
154 
155 
157 {
158  set_.mapMesh(map);
159 }
160 
161 
163 {
164  set_.distribute(map);
165 }
166 
167 
169 {
170  if (fvConstraint::read(dict))
171  {
172  set_.read(coeffs());
173  readCoeffs();
174  return true;
175  }
176  else
177  {
178  return false;
179  }
180 }
181 
182 
183 // ************************************************************************* //
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: limitMag.C:162
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: limitMag.C:168
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual wordList constrainedFields() const
Return the list of fields constrained by the fvConstraint.
Definition: limitMag.C:130
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
virtual bool movePoints()
Update for mesh motion.
Definition: limitMag.C:143
Macros for easy insertion into run-time selection tables.
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
#define IMPLEMENT_FV_CONSTRAINT_CONSTRAIN_FIELD(Type, constraintType)
Definition: fvConstraintM.H:49
static const word null
An empty word.
Definition: word.H:77
Limits the magnitude of the specified field to the specified max value.
Definition: limitMag.H:66
List< label > labelList
A List of labels.
Definition: labelList.H:56
dimensioned< scalar > magSqr(const dimensioned< Type > &)
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: limitMag.C:150
List< word > wordList
A List of words.
Definition: fileName.H:54
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: limitMag.C:156
label patchi
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
FOR_ALL_FIELD_TYPES(DefineFvWallInfoType)
limitMag(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition: limitMag.C:112
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvConstraint.C:143
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
Finite volume options abstract base class.
Definition: fvConstraint.H:56
Namespace for OpenFOAM.