Loading...
Searching...
No Matches
scaledFixedValueFvPatchField.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-------------------------------------------------------------------------------
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
32// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33
34template<class Type>
36(
37 const fvPatch& p,
39)
40:
42 scalePtr_(),
43 refValuePtr_(fvPatchField<Type>::New("refValue", p, iF))
44{}
45
46
47template<class Type>
49(
50 const fvPatch& p,
52 const dictionary& dict
53)
54:
55 fixedValueFvPatchField<Type>(p, iF, dict, IOobjectOption::NO_READ),
56 scalePtr_(PatchFunction1<scalar>::New(p.patch(), "scale", dict)),
57 refValuePtr_(fvPatchField<Type>::New(p, iF, dict.subDict("refValue")))
58{
60 {
62 << typeName << " condition can only be applied to fixed value "
63 << "conditions"
65 }
67 const scalarField s(scalePtr_->value(this->db().time().timeOutputValue()));
68 this->operator==(s*refValuePtr_());
69}
70
71
72template<class Type>
74(
76 const fvPatch& p,
78 const fvPatchFieldMapper& mapper
79)
80:
81 fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
83 refValuePtr_(fvPatchField<Type>::New(ptf.refValue(), p, iF, mapper))
84{}
85
86
87template<class Type>
89(
91)
92:
104)
105:
106 fixedValueFvPatchField<Type>(spf, iF),
107 scalePtr_(spf.scalePtr_.clone(spf.patch().patch())),
109{}
110
111
112// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113
114template<class Type>
116(
117 const fvPatchFieldMapper& m
118)
119{
121 refValuePtr_->autoMap(m);
122
123 scalePtr_().autoMap(m);
124
125 if (scalePtr_().constant())
126 {
127 // If mapper is not dependent on time we're ok to evaluate
128 this->evaluate();
129 }
130}
131
132
133template<class Type>
135(
136 const fvPatchField<Type>& ptf,
137 const labelList& addr
138)
139{
141
142 const scaledFixedValueFvPatchField& sptf =
144
145 refValuePtr_->rmap(sptf.refValue(), addr);
146
147 scalePtr_().rmap(sptf.scalePtr_(), addr);
148}
149
150
151template<class Type>
153{
154 if (this->updated())
155 {
156 return;
157 }
158
159 refValuePtr_->evaluate();
160
161 const scalarField s(scalePtr_->value(this->db().time().timeOutputValue()));
162
163 // Note: setting this field value using = operator (not ==)
174
175 scalePtr_->writeData(os);
176
177 os.beginBlock("refValue");
178 refValuePtr_->write(os);
179 os.endBlock();
180}
181
182
183template<class Type>
185(
186 const fvPatchField<Type>& ptf
187)
188{
189 const scalarField s(scalePtr_->value(this->db().time().timeOutputValue()));
190
191 forAll(s, facei)
192 {
193 const scalar si = s[facei];
194 if (mag(si) > ROOTVSMALL)
195 {
196 refValuePtr_->operator[](facei) = ptf[facei]/si;
197 }
199
201}
202
203
204template<class Type>
206{
207 const scalarField s(scalePtr_->value(this->db().time().timeOutputValue()));
208
209 forAll(s, facei)
210 {
211 const scalar si = s[facei];
212 if (mag(si) > ROOTVSMALL)
213 {
214 refValuePtr_->operator[](facei) = tf[facei]/si;
215 }
217
219}
220
221
222template<class Type>
224{
225 const scalarField s(scalePtr_->value(this->db().time().timeOutputValue()));
226
227 forAll(s, facei)
228 {
229 const scalar si = s[facei];
230 if (mag(si) > ROOTVSMALL)
231 {
232 refValuePtr_->operator[](facei) = t/si;
233 }
234 }
235
237}
238
239
240// ************************************************************************* //
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...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
void operator=(const Field< Type > &)
Copy assignment.
Definition Field.C:781
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
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
fixedValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
static tmp< fvPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field, sets updated() to false.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
This condition applies a scalar multiplier to the value of another boundary condition.
tmp< fvPatchField< Type > > refValuePtr_
Condition to supply the reference value.
virtual void write(Ostream &) const
Write.
virtual void operator==(const fvPatchField< Type > &)
scaledFixedValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
autoPtr< PatchFunction1< scalar > > scalePtr_
Scalar scale factor.
const fvPatchField< Type > & refValue() const
Return the reference value condition.
virtual void rmap(const fvPatchField< Type > &ptf, const labelList &addr)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void autoMap(const fvPatchFieldMapper &m)
Map (and resize as needed) from self given a mapping object.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual tmp< fvPatchField< Type > > clone() const
Return a clone.
volScalarField & p
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
OBJstream os(runTime.globalPath()/outputName)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Different types of constants.
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299