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 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_(nullptr)
59 {
60  if (dict.found("ramp"))
61  {
62  ramp_ = Function1<scalar>::New("ramp", dict);
63  }
64 
65  if (dict.found("value"))
66  {
67  fvPatchVectorField::operator=
68  (
69  vectorField("value", dict, p.size())
70  );
71  }
72  else
73  {
74  this->evaluate();
75  }
76 }
77 
78 
81 (
83  const fvPatch& p,
85  const fvPatchFieldMapper& mapper
86 )
87 :
88  fixedValueFvPatchVectorField(p, iF), // Don't map
89  uniformValue_(ptf.uniformValue_.clone(p.patch())),
90  ramp_(ptf.ramp_.clone())
91 {
92  if (mapper.direct() && !mapper.hasUnmapped())
93  {
94  // Use mapping instead of re-evaluation
95  this->map(ptf, mapper);
96  }
97  else
98  {
99  // Evaluate since value not mapped
100  this->evaluate();
101  }
102 }
103 
104 
107 (
109 )
110 :
111  fixedValueFvPatchVectorField(ptf),
112  uniformValue_(ptf.uniformValue_.clone(this->patch().patch())),
113  ramp_(ptf.ramp_.clone())
114 {}
115 
116 
119 (
122 )
123 :
124  fixedValueFvPatchVectorField(ptf, iF),
125  uniformValue_(ptf.uniformValue_.clone(this->patch().patch())),
126  ramp_(ptf.ramp_.clone())
127 {}
128 
129 
130 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
131 
133 (
134  const fvPatchFieldMapper& mapper
135 )
136 {
137  fixedValueFvPatchVectorField::autoMap(mapper);
138  uniformValue_().autoMap(mapper);
139 
140  if (uniformValue_().constant())
141  {
142  // If mapper is not dependent on time we're ok to evaluate
143  this->evaluate();
144  }
145 }
146 
147 
149 (
150  const fvPatchVectorField& ptf,
151  const labelList& addr
152 )
153 {
154  fixedValueFvPatchVectorField::rmap(ptf, addr);
155 
157  refCast<const uniformNormalFixedValueFvPatchVectorField>(ptf);
158 
159  uniformValue_().rmap(tiptf.uniformValue_(), addr);
160 }
161 
162 
164 {
165  if (updated())
166  {
167  return;
168  }
169 
170  const scalar t = this->db().time().timeOutputValue();
171 
172  tmp<vectorField> tvalues(uniformValue_->value(t)*patch().nf());
173 
174  if (ramp_)
175  {
176  tvalues.ref() *= ramp_->value(t);
177  }
178 
181 }
182 
183 
185 {
187  uniformValue_->writeData(os);
188  if (ramp_)
189  {
190  ramp_->writeData(os);
191  }
192  this->writeEntry("value", os);
193 }
194 
195 
196 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
197 
198 namespace Foam
199 {
201  (
204  );
205 }
206 
207 // ************************************************************************* //
Foam::fvPatchField< vector >
volFields.H
Foam::fvPatchField< vector >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::uniformNormalFixedValueFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: uniformNormalFixedValueFvPatchVectorField.C:163
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::uniformNormalFixedValueFvPatchVectorField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: uniformNormalFixedValueFvPatchVectorField.C:133
fvPatchFieldMapper.H
Foam::FieldMapper::direct
virtual bool direct() const =0
Foam::uniformNormalFixedValueFvPatchVectorField
This boundary condition provides a uniform surface-normal vector boundary condition by its magnitude.
Definition: uniformNormalFixedValueFvPatchVectorField.H:97
Foam::FieldMapper::hasUnmapped
virtual bool hasUnmapped() const =0
Are there unmapped values? I.e. do all size() elements get.
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::uniformNormalFixedValueFvPatchVectorField::uniformNormalFixedValueFvPatchVectorField
uniformNormalFixedValueFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: uniformNormalFixedValueFvPatchVectorField.C:37
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::PatchFunction1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: PatchFunction1.H:63
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fvPatchField< vector >::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:313
Foam::uniformNormalFixedValueFvPatchVectorField::write
virtual void write(Ostream &os) const
Write.
Definition: uniformNormalFixedValueFvPatchVectorField.C:184
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
uniformNormalFixedValueFvPatchVectorField.H
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
constant
constant condensation/saturation model.
Foam::fvPatchField< vector >::operator=
virtual void operator=(const UList< vector > &)
Definition: fvPatchField.C:379
Foam::stringOps::evaluate
string evaluate(const std::string &s, size_t pos=0, size_t len=std::string::npos)
Definition: stringOpsEvaluate.C:37
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::uniformNormalFixedValueFvPatchVectorField::rmap
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: uniformNormalFixedValueFvPatchVectorField.C:149