Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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)
42 fixedValueFvPatchVectorField(p, iF),
43 refValueFunc_(nullptr),
44 ramp_(nullptr)
45{}
46
47
50(
51 const fvPatch& p,
53 const dictionary& dict
54)
55:
56 fixedValueFvPatchVectorField(p, iF, dict, IOobjectOption::NO_READ),
57 refValueFunc_(PatchFunction1<scalar>::New(p.patch(), "uniformValue", dict)),
58 ramp_(Function1<scalar>::NewIfPresent("ramp", dict, word::null, &db()))
59{
60 if (!this->readValueEntry(dict))
61 {
62 this->evaluate();
63 }
64}
65
66
69(
70 const uniformNormalFixedValueFvPatchVectorField& ptf,
71 const fvPatch& p,
72 const DimensionedField<vector, volMesh>& iF,
73 const fvPatchFieldMapper& mapper
74)
75:
76 fixedValueFvPatchVectorField(p, iF), // Don't map
77 refValueFunc_(ptf.refValueFunc_.clone(p.patch())),
78 ramp_(ptf.ramp_.clone())
79{
80 if (mapper.direct() && !mapper.hasUnmapped())
81 {
82 // Use mapping instead of re-evaluation
83 this->map(ptf, mapper);
84 }
85 else
86 {
87 // Evaluate since value not mapped
88 this->evaluate();
89 }
90}
91
92
95(
97)
99 fixedValueFvPatchVectorField(ptf),
100 refValueFunc_(ptf.refValueFunc_.clone(this->patch().patch())),
101 ramp_(ptf.ramp_.clone())
102{}
103
104
107(
110)
111:
112 fixedValueFvPatchVectorField(ptf, iF),
113 refValueFunc_(ptf.refValueFunc_.clone(this->patch().patch())),
114 ramp_(ptf.ramp_.clone())
115{}
116
117
118// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119
121(
122 const fvPatchFieldMapper& mapper
123)
124{
125 fixedValueFvPatchVectorField::autoMap(mapper);
126 refValueFunc_().autoMap(mapper);
127
128 if (refValueFunc_().constant())
130 // If mapper is not dependent on time we're ok to evaluate
131 this->evaluate();
132 }
133}
134
135
137(
138 const fvPatchVectorField& ptf,
139 const labelList& addr
140)
141{
142 fixedValueFvPatchVectorField::rmap(ptf, addr);
143
146
147 refValueFunc_().rmap(tiptf.refValueFunc_(), addr);
148}
149
150
152{
153 if (updated())
154 {
155 return;
156 }
157
158 const scalar t = this->db().time().timeOutputValue();
159
160 tmp<vectorField> tvalues(refValueFunc_->value(t)*patch().nf());
161
162 if (ramp_)
163 {
164 tvalues.ref() *= ramp_->value(t);
175 refValueFunc_->writeData(os);
176 if (ramp_)
177 {
178 ramp_->writeData(os);
179 }
182
183
184// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185
186namespace Foam
187{
189 (
191 uniformNormalFixedValueFvPatchVectorField
192 );
193}
194
195// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual bool direct() const =0
Is it a direct (non-interpolating) mapper?
virtual bool hasUnmapped() const =0
Any unmapped values?
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition Function1.H:92
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual void operator=(const UList< vector > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
This boundary condition provides a uniform surface-normal vector boundary condition by its magnitude.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual tmp< fvPatchField< vector > > clone() const
Return a clone.
uniformNormalFixedValueFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Different types of constants.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
List< label > labelList
A List of labels.
Definition List.H:62
fvPatchField< vector > fvPatchVectorField
dictionary dict