Loading...
Searching...
No Matches
fixedNormalSlipFvPatchField.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\*---------------------------------------------------------------------------*/
30#include "symmTransformField.H"
31
32// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33
34template<class Type>
36(
37 const fvPatch& p,
39)
40:
41 parent_bctype(p, iF),
42 fixedValue_(p.size(), Foam::zero{}),
43 writeValue_(false)
44{}
45
46
47template<class Type>
49(
51 const fvPatch& p,
53 const fvPatchFieldMapper& mapper
54)
55:
56 parent_bctype(ptf, p, iF, mapper),
57 fixedValue_(ptf.fixedValue_, mapper),
58 writeValue_(ptf.writeValue_)
59{}
60
61
62template<class Type>
64(
65 const fvPatch& p,
67 const dictionary& dict
68)
69:
70 parent_bctype(p, iF),
71 fixedValue_("fixedValue", dict, p.size()),
72 writeValue_(dict.getOrDefault("writeValue", false))
73{
84)
85:
86 parent_bctype(ptf),
87 fixedValue_(ptf.fixedValue_),
88 writeValue_(ptf.writeValue_)
89{}
90
91
92template<class Type>
94(
97)
98:
99 parent_bctype(ptf, iF),
100 fixedValue_(ptf.fixedValue_),
101 writeValue_(ptf.writeValue_)
102{}
103
104
105// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106
107template<class Type>
109(
110 const fvPatchFieldMapper& m
111)
114 fixedValue_.autoMap(m);
115}
116
117
118template<class Type>
120(
121 const fvPatchField<Type>& ptf,
122 const labelList& addr
123)
124{
125 parent_bctype::rmap(ptf, addr);
126
127 const auto& dmptf =
130 fixedValue_.rmap(dmptf.fixedValue_, addr);
131}
132
133
134template<class Type>
137{
138 static_assert
139 (
141 "normal-slip with vector, tensor only!"
142 );
143
144 const vectorField nHat(this->patch().nf());
145 const Field<Type> pif(this->patchInternalField());
146
147 return
149 (nHat*(nHat & fixedValue_) + transform(I - sqr(nHat), pif)) - pif
150 )*this->patch().deltaCoeffs();
151}
152
153
154template<class Type>
156(
158)
159{
160 if (!this->updated())
161 {
162 this->updateCoeffs();
163 }
164
165 const vectorField nHat(this->patch().nf());
166
167 Field<Type>::operator=
168 (
169 nHat*(nHat & fixedValue_)
170 + transform(I - sqr(nHat), this->patchInternalField())
171 );
174}
175
176
177template<class Type>
180{
181 static_assert
182 (
184 "normal-slip with vector, tensor only!"
185 );
186 // if constexpr (!is_rotational_vectorspace_v<Type>)
187 // {
188 // // Rotational-invariant type
189 // return tmp<Field<Type>>::New(this->size(), Foam::zero{});
190 // }
191
192 {
193 tmp<vectorField> diag(cmptMag(this->patch().nf()));
196 }
197}
198
199
200template<class Type>
202{
203 this->parent_bctype::write(os);
204 if (writeValue_)
205 {
206 os.writeEntry("writeValue", "true");
207 }
208 fixedValue_.writeEntry("fixedValue", os);
209 if (writeValue_)
210 {
212 }
213}
214
215
216// ************************************************************************* //
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
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 sets the patch-normal component to the field (vector or tensor) to the patch-...
virtual void write(Ostream &) const
Write.
fixedNormalSlipFvPatchField(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 > &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 tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field.
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.
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
A class for managing temporary objects.
Definition tmp.H:75
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
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)
Field< vector > vectorField
Specialisation of Field<T> for vector.
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 > &)
dictionary dict
A non-counting (dummy) refCount.
Definition refCount.H:55