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-2020 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 #include "symmTransformField.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class Type>
36 (
37  const fvPatch& p,
39 )
40 :
42  refValue_(p.size()),
43  valueFraction_(p.size(), 1.0)
44 {}
45 
46 
47 template<class Type>
49 (
51  const fvPatch& p,
53  const fvPatchFieldMapper& mapper
54 )
55 :
56  transformFvPatchField<Type>(ptf, p, iF, mapper),
57  refValue_(ptf.refValue_, mapper),
58  valueFraction_(ptf.valueFraction_, mapper)
59 {}
60 
61 
62 template<class Type>
64 (
65  const fvPatch& p,
67  const dictionary& dict
68 )
69 :
71  refValue_(p.size(), Zero),
72  valueFraction_("valueFraction", dict, p.size())
73 {
74  this->patchType() = dict.getOrDefault<word>("patchType", word::null);
75 
76  // Backwards compatibility - leave refValue as zero unless specified
77  if (dict.found("refValue"))
78  {
79  refValue_ = Field<Type>("refValue", dict, p.size());
80  }
81 
82  evaluate();
83 }
84 
85 
86 template<class Type>
88 (
90 )
91 :
93  refValue_(ptf.refValue_),
94  valueFraction_(ptf.valueFraction_)
95 {}
96 
97 
98 template<class Type>
100 (
103 )
104 :
106  refValue_(ptf.refValue_),
107  valueFraction_(ptf.valueFraction_)
108 {}
109 
110 
111 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112 
113 template<class Type>
115 (
116  const fvPatchFieldMapper& m
117 )
118 {
120  refValue_.autoMap(m);
121  valueFraction_.autoMap(m);
122 }
123 
124 
125 template<class Type>
127 (
128  const fvPatchField<Type>& ptf,
129  const labelList& addr
130 )
131 {
133 
134  const partialSlipFvPatchField<Type>& dmptf =
135  refCast<const partialSlipFvPatchField<Type>>(ptf);
136 
137  refValue_.rmap(dmptf.refValue_, addr);
138  valueFraction_.rmap(dmptf.valueFraction_, addr);
139 }
140 
141 
142 template<class Type>
145 {
146  tmp<vectorField> nHat = this->patch().nf();
147  const Field<Type> pif(this->patchInternalField());
148 
149  return
150  (
151  valueFraction_*refValue_
152  + (1.0 - valueFraction_)*transform(I - sqr(nHat), pif) - pif
153  )*this->patch().deltaCoeffs();
154 }
155 
156 
157 template<class Type>
159 (
160  const Pstream::commsTypes
161 )
162 {
163  if (!this->updated())
164  {
165  this->updateCoeffs();
166  }
167 
168  tmp<vectorField> nHat = this->patch().nf();
169 
171  (
172  valueFraction_*refValue_
173  +
174  (1.0 - valueFraction_)
175  *transform(I - sqr(nHat), this->patchInternalField())
176  );
177 
179 }
180 
181 
182 template<class Type>
185 {
186  const vectorField nHat(this->patch().nf());
187  vectorField diag(nHat.size());
188 
189  diag.replace(vector::X, mag(nHat.component(vector::X)));
190  diag.replace(vector::Y, mag(nHat.component(vector::Y)));
191  diag.replace(vector::Z, mag(nHat.component(vector::Z)));
192 
193  return
194  valueFraction_*pTraits<Type>::one
195  + (1.0 - valueFraction_)
196  *transformFieldMask<Type>(pow<vector, pTraits<Type>::rank>(diag));
197 }
198 
199 
200 template<class Type>
202 {
204  refValue_.writeEntry("refValue", os);
205  valueFraction_.writeEntry("valueFraction", os);
206 }
207 
208 
209 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:50
p
volScalarField & p
Definition: createFieldRefs.H:8
symmTransformField.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::diag
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
Definition: pointPatchFieldFunctions.H:287
Foam::partialSlipFvPatchField::partialSlipFvPatchField
partialSlipFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: partialSlipFvPatchField.C:36
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:519
partialSlipFvPatchField.H
Foam::partialSlipFvPatchField::snGradTransformDiag
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
Definition: partialSlipFvPatchField.C:184
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::partialSlipFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: partialSlipFvPatchField.C:201
Foam::partialSlipFvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
Definition: partialSlipFvPatchField.C:144
Foam::partialSlipFvPatchField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: partialSlipFvPatchField.C:115
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::Field::component
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:551
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::partialSlipFvPatchField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
Definition: partialSlipFvPatchField.C:159
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:69
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::List< label >
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:54
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::transformFvPatchField
Foam::transformFvPatchField.
Definition: transformFvPatchField.H:54
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:35
Foam::partialSlipFvPatchField::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: partialSlipFvPatchField.C:127
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::partialSlipFvPatchField
This boundary condition provides a partial slip condition. The amount of slip is controlled by a user...
Definition: partialSlipFvPatchField.H:91
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::stringOps::evaluate
string evaluate(const std::string &s, size_t pos=0, size_t len=std::string::npos)
Definition: stringOpsEvaluate.C:37
Foam::I
static const Identity< scalar > I
Definition: Identity.H:95
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54