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-2021 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\*---------------------------------------------------------------------------*/
28
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(), 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(), Zero),
74 valueFraction_("valueFraction", dict, p.size()),
75 writeValue_(dict.getOrDefault("writeValue", false))
76{
77 this->patchType() = dict.getOrDefault<word>("patchType", word::null);
78
79 // Backwards compatibility - leave refValue as zero unless specified
80 if (dict.found("refValue"))
81 {
82 refValue_ = Field<Type>("refValue", dict, p.size());
83 }
84
85 evaluate();
86}
87
88
89template<class Type>
91(
93)
94:
95 parent_bctype(ptf),
96 refValue_(ptf.refValue_),
97 valueFraction_(ptf.valueFraction_),
98 writeValue_(ptf.writeValue_)
99{}
100
101
102template<class Type>
104(
107)
108:
109 parent_bctype(ptf, iF),
110 refValue_(ptf.refValue_),
111 valueFraction_(ptf.valueFraction_),
112 writeValue_(ptf.writeValue_)
113{}
114
115
116// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117
118template<class Type>
120(
121 const fvPatchFieldMapper& m
122)
123{
124 parent_bctype::autoMap(m);
125 refValue_.autoMap(m);
126 valueFraction_.autoMap(m);
127}
128
129
130template<class Type>
132(
133 const fvPatchField<Type>& ptf,
134 const labelList& addr
135)
136{
137 parent_bctype::rmap(ptf, addr);
138
139 const auto& dmptf =
140 refCast<const partialSlipFvPatchField<Type>>(ptf);
141
142 refValue_.rmap(dmptf.refValue_, addr);
143 valueFraction_.rmap(dmptf.valueFraction_, addr);
144}
145
146
147template<class Type>
150{
151 tmp<vectorField> nHat = this->patch().nf();
152 const Field<Type> pif(this->patchInternalField());
153
154 return
155 (
156 valueFraction_*refValue_
157 + (1.0 - valueFraction_)*transform(I - sqr(nHat), pif) - pif
158 )*this->patch().deltaCoeffs();
159}
160
161
162template<class Type>
164(
166)
167{
168 if (!this->updated())
169 {
170 this->updateCoeffs();
171 }
172
173 tmp<vectorField> nHat = this->patch().nf();
174
176 (
177 valueFraction_*refValue_
178 +
179 (1.0 - valueFraction_)
180 *transform(I - sqr(nHat), this->patchInternalField())
181 );
182
183 parent_bctype::evaluate();
184}
185
186
187template<class Type>
190{
191 const vectorField nHat(this->patch().nf());
192 vectorField diag(nHat.size());
193
194 diag.replace(vector::X, mag(nHat.component(vector::X)));
195 diag.replace(vector::Y, mag(nHat.component(vector::Y)));
196 diag.replace(vector::Z, mag(nHat.component(vector::Z)));
197
198 return
199 valueFraction_*pTraits<Type>::one
200 + (1.0 - valueFraction_)
201 *transformFieldMask<Type>(pow<vector, pTraits<Type>::rank>(diag));
202}
203
204
205template<class Type>
207{
208 this->parent_bctype::write(os);
209 refValue_.writeEntry("refValue", os);
210 valueFraction_.writeEntry("valueFraction", os);
211
212 if (writeValue_)
213 {
214 os.writeEntry("writeValue", "true");
215 this->writeEntry("value", os);
216 }
217}
218
219
220// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type.
Definition: Field.H:82
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:545
void evaluate()
Evaluate boundary conditions.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
commsTypes
Types of communications.
Definition: UPstream.H:67
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
virtual bool write()
Write the output fields.
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...
Definition: fvPatchField.H:82
const word & patchType() const
Optional patch type.
Definition: fvPatchField.H:380
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:62
A traits class, which is primarily used for primitives.
Definition: pTraits.H:59
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 evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch 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.
A class for managing temporary objects.
Definition: tmp.H:65
Foam::transformFvPatchField.
A class for handling words, derived from Foam::string.
Definition: word.H:68
static const word null
An empty word.
Definition: word.H:80
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static const Identity< scalar > I
Definition: Identity.H:94
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
dictionary dict