uniformFixedGradientFvPatchField.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) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  const fvPatch& p,
38 )
39 :
41  refGradFunc_(nullptr)
42 {}
43 
44 
45 template<class Type>
47 (
48  const fvPatch& p,
50  const Field<Type>& fld
51 )
52 :
54  refGradFunc_(nullptr)
55 {}
56 
57 
58 template<class Type>
60 (
61  const fvPatch& p,
63  const dictionary& dict
64 )
65 :
66  fixedGradientFvPatchField<Type>(p, iF), // Bypass dictionary constructor
67  refGradFunc_
68  (
69  PatchFunction1<Type>::New(p.patch(), "uniformGradient", dict)
70  )
71 {
73 
74  if (!this->readValueEntry(dict))
75  {
76  // Ensure field has initialised values
77  this->extrapolateInternal();
78 
79  // Evaluate to assign a value
80  this->evaluate();
81  }
82 }
83 
84 
85 template<class Type>
87 (
88  const uniformFixedGradientFvPatchField<Type>& ptf,
89  const fvPatch& p,
90  const DimensionedField<Type, volMesh>& iF,
91  const fvPatchFieldMapper& mapper
92 )
93 :
94  fixedGradientFvPatchField<Type>(ptf, p, iF, mapper),
95  refGradFunc_(ptf.refGradFunc_.clone())
96 {}
97 
98 
99 template<class Type>
101 (
103 )
104 :
106  refGradFunc_(ptf.refGradFunc_.clone())
107 {}
108 
109 
110 template<class Type>
112 (
115 )
116 :
117  fixedGradientFvPatchField<Type>(ptf, iF),
118  refGradFunc_(ptf.refGradFunc_.clone())
119 {
120  // Evaluate the profile if defined
121  if (refGradFunc_)
122  {
123  this->evaluate();
124  }
125 }
126 
127 
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129 
130 template<class Type>
132 {
133  if (this->updated())
134  {
135  return;
136  }
137 
138  const scalar t = this->db().time().timeOutputValue();
139  // Extra safety on the evaluation
140  if (refGradFunc_)
141  {
142  this->gradient() = refGradFunc_->value(t);
143  }
144  else
145  {
146  this->gradient() = Zero;
147  }
148 
150 }
151 
152 
153 template<class Type>
155 {
157  if (refGradFunc_)
158  {
159  refGradFunc_->writeData(os);
160  }
162 }
163 
164 
165 // ************************************************************************* //
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
Definition: fvPatchField.C:32
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
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.
void extrapolateInternal()
Assign the patch field from the internal field.
Definition: fvPatchField.C:62
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void readDict(const dictionary &dict)
Read dictionary entries.
virtual void write(Ostream &os) const
Write.
uniformFixedGradientFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual void write(Ostream &) const
Write.
Generic templated field type.
Definition: Field.H:62
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
This boundary condition supplies a fixed gradient condition, such that the patch values are calculate...
OBJstream os(runTime.globalPath()/outputName)
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
This boundary condition provides a uniform fixed gradient condition.
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;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
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.
volScalarField & p
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127