boundConstraint.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) 2023 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 "boundConstraint.H"
27 #include "volFields.H"
28 #include "bound.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
39  (
41  bound,
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 void Foam::fv::bound::readCoeffs()
51 {
52  fieldName_ = coeffs().lookup<word>("field");
53  min_ = coeffs().lookup<scalar>("min");
54 }
55 
56 
57 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 
60 (
61  const word& name,
62  const word& modelType,
63  const fvMesh& mesh,
64  const dictionary& dict
65 )
66 :
67  fvConstraint(name, modelType, mesh, dict),
68  fieldName_(word::null),
69  min_(0)
70 {
71  readCoeffs();
72 }
73 
74 
75 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76 
78 {
79  return wordList(1, fieldName_);
80 }
81 
82 
84 {
85  return Foam::bound(f, dimensionedScalar(f.dimensions(), min_));
86 }
87 
88 
90 {
91  return true;
92 }
93 
94 
96 {}
97 
98 
100 {}
101 
102 
104 {}
105 
106 
108 {
110  {
111  readCoeffs();
112  return true;
113  }
114  else
115  {
116  return false;
117  }
118 }
119 
120 
121 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Bound the given scalar field where it is below the specified minimum.
Generic GeometricField class.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
Finite volume options abstract base class.
Definition: fvConstraint.H:57
const dictionary & coeffs() const
Return dictionary.
Definition: fvConstraintI.H:40
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvConstraint.C:166
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Bound the specified scalar field where it is below the specified minimum.
virtual bool movePoints()
Update for mesh motion.
virtual bool constrain(volScalarField &p) const
Constrain the pressure field.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read dictionary.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
bound(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
virtual wordList constrainedFields() const
Return the list of fields constrained by the fvConstraint.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList f(nPoints)
labelList fv(nPoints)
dictionary dict