Loading...
Searching...
No Matches
partialSlipFvPatchField.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) 2017-2025 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 "symmTransformField.H"
31
32// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33
34template<class Type>
36(
37 const fvPatch& p,
39)
40:
41 parent_bctype(p, iF),
42 refValue_(p.size(), Foam::zero{}),
43 valueFraction_(p.size(), 1.0),
44 writeValue_(false)
45{}
46
47
48template<class Type>
50(
52 const fvPatch& p,
54 const fvPatchFieldMapper& mapper
55)
56:
57 parent_bctype(ptf, p, iF, mapper),
58 refValue_(ptf.refValue_, mapper),
59 valueFraction_(ptf.valueFraction_, mapper),
60 writeValue_(ptf.writeValue_)
61{}
62
63
64template<class Type>
66(
67 const fvPatch& p,
69 const dictionary& dict
70)
71:
72 parent_bctype(p, iF),
73 refValue_(p.size(), Foam::zero{}),
74 valueFraction_("valueFraction", dict, p.size()),
75 writeValue_(dict.getOrDefault("writeValue", false))
76{
78
79 // Backwards compatibility - leave refValue as zero unless specified
80 refValue_.assign("refValue", dict, p.size(), IOobjectOption::LAZY_READ);
81
82 evaluate();
83}
84
85
86template<class Type>
88(
91)
92:
93 parent_bctype(ptf, iF),
94 refValue_(ptf.refValue_),
95 valueFraction_(ptf.valueFraction_),
96 writeValue_(ptf.writeValue_)
97{}
98
99
100template<class Type>
102(
104)
105:
107{}
108
109
110// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111
112template<class Type>
114(
115 const fvPatchFieldMapper& m
116)
117{
119 refValue_.autoMap(m);
120 valueFraction_.autoMap(m);
121}
122
123
124template<class Type>
126(
127 const fvPatchField<Type>& ptf,
128 const labelList& addr
129)
130{
131 parent_bctype::rmap(ptf, addr);
132
133 const auto& dmptf =
135
136 refValue_.rmap(dmptf.refValue_, addr);
137 valueFraction_.rmap(dmptf.valueFraction_, addr);
138}
139
140
141template<class Type>
144{
145 const Field<Type> pif(this->patchInternalField());
146
147 tmp<Field<Type>> rotated;
148
150 {
151 // Rotational-invariant type
152 rotated.cref(pif);
153 }
154 else
155 {
156 tmp<vectorField> nHat(this->patch().nf());
157 rotated = transform(I - sqr(nHat), pif);
158 }
159
160 return
162 lerp(rotated, refValue_, valueFraction_) - pif
163 )*this->patch().deltaCoeffs();
164}
165
166
167template<class Type>
169(
171)
172{
173 if (!this->updated())
174 {
175 this->updateCoeffs();
176 }
177
178 tmp<Field<Type>> rotated;
179
181 {
182 // Rotational-invariant type
183 rotated = this->patchInternalField();
184 }
185 else
186 {
187 tmp<vectorField> nHat(this->patch().nf());
188 rotated = transform(I - sqr(nHat), this->patchInternalField());
189 }
190
191
193 (
194 lerp(rotated, refValue_, valueFraction_)
195 );
198}
199
200
201template<class Type>
204{
205 // static_assert(!std::is_arithmetic_v<T>, "partial-slip with scalar??");
206 // if constexpr (!is_rotational_vectorspace_v<Type>)
207 // {
208 // // Rotational-invariant type
209 // return tmp<Field<Type>>::New(this->size(), Foam::zero{});
210 // }
211
212 {
213 tmp<vectorField> diag(cmptMag(this->patch().nf()));
214
215 return
216 (
217 valueFraction_*pTraits<Type>::one
218 + (1.0 - valueFraction_)
220 );
221 }
222}
223
224
225template<class Type>
227{
228 this->parent_bctype::write(os);
229 if (writeValue_)
230 {
231 os.writeEntry("writeValue", "true");
232 }
233 refValue_.writeEntry("refValue", os);
234 valueFraction_.writeEntry("valueFraction", os);
235 if (writeValue_)
236 {
238 }
239}
240
241
242// ************************************************************************* //
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
@ LAZY_READ
Reading is optional [identical to READ_IF_PRESENT].
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
virtual void readDict(const dictionary &dict)
Read dictionary entries.
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.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
const DimensionedField< Type, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
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
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
This boundary condition provides a partial slip condition. The amount of slip is controlled by a user...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void write(Ostream &) const
Write.
partialSlipFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field.
A class for managing temporary objects.
Definition tmp.H:75
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition tmpI.H:221
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< Field< Type1 > > transformFieldMask(const Field< Type2 > &fld)
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
static const Identity< scalar > I
Definition Identity.H:100
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
constexpr bool is_rotational_vectorspace_v
The is_rotational_vectorspace value of Type.
Definition pTraits.H:251
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Vector< scalar > vector
Definition vector.H:57
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
dictionary dict
A non-counting (dummy) refCount.
Definition refCount.H:55