Loading...
Searching...
No Matches
fixedGradientFvPatchField.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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\*---------------------------------------------------------------------------*/
30#include "dictionary.H"
31
32// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
33
34template<class Type>
36(
37 const dictionary& dict,
39)
40{
41 if (!IOobjectOption::isAnyRead(readOpt)) return false;
42 const auto& p = fvPatchFieldBase::patch();
43
44
45 const auto* eptr = dict.findEntry("gradient", keyType::LITERAL);
46
47 if (eptr)
48 {
49 gradient_.assign(*eptr, p.size());
50 return true;
51 }
52
54 {
56 << "Required entry 'gradient' : missing for patch " << p.name()
57 << " in dictionary " << dict.relativeName() << nl
59 }
60
61 return false;
62}
63
64
65// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66
67template<class Type>
69(
70 const fvPatch& p,
72)
74 fvPatchField<Type>(p, iF),
75 gradient_(p.size(), Zero)
76{}
77
78
79template<class Type>
81(
82 const fvPatch& p,
84 const dictionary& dict,
86)
87:
88 fvPatchField<Type>(p, iF, dict, IOobjectOption::NO_READ),
89 gradient_(p.size())
90{
91 if (readGradientEntry(dict, requireGrad))
92 {
93 evaluate();
94 }
95 else
96 {
97 // Not read (eg, optional and missing):
98 // - treat as zero-gradient, do not evaluate
100 gradient_ = Zero;
101 }
102}
103
104
105template<class Type>
107(
108 const fixedGradientFvPatchField<Type>& ptf,
109 const fvPatch& p,
110 const DimensionedField<Type, volMesh>& iF,
111 const fvPatchFieldMapper& mapper
112)
113:
114 fvPatchField<Type>(ptf, p, iF, mapper),
115 gradient_(ptf.gradient_, mapper)
116{
117 if (notNull(iF) && mapper.hasUnmapped())
118 {
119 WarningInFunction
120 << "On field " << iF.name() << " patch " << p.name()
121 << " patchField " << this->type()
122 << " : mapper does not map all values." << nl
123 << " To avoid this warning fully specify the mapping in derived"
124 << " patch fields." << endl;
125 }
126}
127
128
129template<class Type>
131(
133)
135 fvPatchField<Type>(ptf),
136 gradient_(ptf.gradient_)
137{}
138
139
140template<class Type>
142(
145)
146:
147 fvPatchField<Type>(ptf, iF),
148 gradient_(ptf.gradient_)
149{}
150
151
152// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
153
154template<class Type>
156(
157 const fvPatchFieldMapper& m
158)
161 gradient_.autoMap(m);
162}
163
164
165template<class Type>
167(
168 const fvPatchField<Type>& ptf,
169 const labelList& addr
170)
171{
172 fvPatchField<Type>::rmap(ptf, addr);
173
176
177 gradient_.rmap(fgptf.gradient_, addr);
178}
179
180
181template<class Type>
183{
184 if (!this->updated())
185 {
186 this->updateCoeffs();
187 }
188
189 Field<Type>::operator=
190 (
191 this->patchInternalField() + gradient_/this->patch().deltaCoeffs()
192 );
195}
196
197
198template<class Type>
201(
202 const tmp<scalarField>&
203) const
206}
207
208
209template<class Type>
212(
213 const tmp<scalarField>&
214) const
216 return gradient()/this->patch().deltaCoeffs();
217}
218
219
220template<class Type>
224 return tmp<Field<Type>>::New(this->size(), Zero);
225}
226
227
228template<class Type>
235
236template<class Type>
238{
240 gradient_.writeEntry("gradient", os);
241}
242
243
244// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual bool hasUnmapped() const =0
Any unmapped values?
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
static bool isReadRequired(readOption opt) noexcept
True if (MUST_READ | READ_MODIFIED) bits are set.
readOption
Enumeration defining read preferences.
static bool isAnyRead(readOption opt) noexcept
True if any reading may be required (ie, != NO_READ).
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
void size(const label n)
Definition UList.H:118
commsTypes
Communications types.
Definition UPstream.H:81
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 gradient condition, such that the patch values are calculate...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual Field< Type > & gradient()
Return gradient at boundary.
virtual void write(Ostream &) const
Write.
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the.
bool readGradientEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "gradient" entry into corresponding member.
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
fixedGradientFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal 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.
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.
virtual void write(Ostream &) const
Write.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
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.
void extrapolateInternal()
Assign the patch field from the internal field.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field, sets updated() to false.
fvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
@ LITERAL
String literal.
Definition keyType.H:82
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
A class for managing temporary objects.
Definition tmp.H:75
volScalarField & p
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
OBJstream os(runTime.globalPath()/outputName)
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
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
bool notNull(const T *ptr) noexcept
True if ptr is not a pointer (of type T) to the nullObject.
Definition nullObject.H:267
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict