uniformNormalFixedValueFvPatchVectorField.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-2021 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 
30 #include "volFields.H"
31 #include "fvPatchFieldMapper.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const fvPatch& p,
40 )
41 :
42  fixedValueFvPatchVectorField(p, iF),
43  uniformValue_(nullptr),
44  ramp_(nullptr)
45 {}
46 
47 
50 (
51  const fvPatch& p,
53  const dictionary& dict
54 )
55 :
56  fixedValueFvPatchVectorField(p, iF, dict, false),
57  uniformValue_(PatchFunction1<scalar>::New(p.patch(), "uniformValue", dict)),
58  ramp_(Function1<scalar>::NewIfPresent("ramp", dict, word::null, &db()))
59 {
60  if (dict.found("value"))
61  {
63  (
64  vectorField("value", dict, p.size())
65  );
66  }
67  else
68  {
69  this->evaluate();
70  }
71 }
72 
73 
76 (
77  const uniformNormalFixedValueFvPatchVectorField& ptf,
78  const fvPatch& p,
79  const DimensionedField<vector, volMesh>& iF,
80  const fvPatchFieldMapper& mapper
81 )
82 :
83  fixedValueFvPatchVectorField(p, iF), // Don't map
84  uniformValue_(ptf.uniformValue_.clone(p.patch())),
85  ramp_(ptf.ramp_.clone())
86 {
87  if (mapper.direct() && !mapper.hasUnmapped())
88  {
89  // Use mapping instead of re-evaluation
90  this->map(ptf, mapper);
91  }
92  else
93  {
94  // Evaluate since value not mapped
95  this->evaluate();
96  }
97 }
98 
99 
102 (
103  const uniformNormalFixedValueFvPatchVectorField& ptf
104 )
105 :
106  fixedValueFvPatchVectorField(ptf),
107  uniformValue_(ptf.uniformValue_.clone(this->patch().patch())),
108  ramp_(ptf.ramp_.clone())
109 {}
110 
111 
114 (
117 )
118 :
119  fixedValueFvPatchVectorField(ptf, iF),
120  uniformValue_(ptf.uniformValue_.clone(this->patch().patch())),
121  ramp_(ptf.ramp_.clone())
122 {}
123 
124 
125 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
126 
128 (
129  const fvPatchFieldMapper& mapper
130 )
131 {
132  fixedValueFvPatchVectorField::autoMap(mapper);
133  uniformValue_().autoMap(mapper);
134 
135  if (uniformValue_().constant())
136  {
137  // If mapper is not dependent on time we're ok to evaluate
138  this->evaluate();
139  }
140 }
141 
142 
144 (
145  const fvPatchVectorField& ptf,
146  const labelList& addr
147 )
148 {
149  fixedValueFvPatchVectorField::rmap(ptf, addr);
150 
152  refCast<const uniformNormalFixedValueFvPatchVectorField>(ptf);
153 
154  uniformValue_().rmap(tiptf.uniformValue_(), addr);
155 }
156 
157 
159 {
160  if (updated())
161  {
162  return;
163  }
164 
165  const scalar t = this->db().time().timeOutputValue();
166 
167  tmp<vectorField> tvalues(uniformValue_->value(t)*patch().nf());
168 
169  if (ramp_)
170  {
171  tvalues.ref() *= ramp_->value(t);
172  }
173 
176 }
177 
178 
180 {
182  uniformValue_->writeData(os);
183  if (ramp_)
184  {
185  ramp_->writeData(os);
186  }
187  this->writeEntry("value", os);
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 namespace Foam
194 {
196  (
198  uniformNormalFixedValueFvPatchVectorField
199  );
200 }
201 
202 // ************************************************************************* //
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
dictionary dict
fvPatchField< vector > fvPatchVectorField
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:68
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:333
This boundary condition provides a uniform surface-normal vector boundary condition by its magnitude...
Macros for easy insertion into run-time selection tables.
A class for handling words, derived from Foam::string.
Definition: word.H:63
A FieldMapper for finite-volume patch fields.
string evaluate(label fieldWidth, const std::string &s, size_t pos=0, size_t len=std::string::npos)
String evaluation with specified (positive, non-zero) field width.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
OBJstream os(runTime.globalPath()/outputName)
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:270
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:352
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
Field< vector > vectorField
Specialisation of Field<T> for vector.
uniformNormalFixedValueFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.