Loading...
Searching...
No Matches
directionMixedFvPatchField.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-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
29#include "symmTransformField.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class Type>
35(
36 const fvPatch& p,
38)
39:
40 transformFvPatchField<Type>(p, iF),
41 refValue_(p.size()),
42 refGrad_(p.size()),
43 valueFraction_(p.size())
44{}
45
46
47template<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 refGrad_(ptf.refGrad_, mapper),
59 valueFraction_(ptf.valueFraction_, mapper)
60{
61 if (notNull(iF) && mapper.hasUnmapped())
62 {
63 WarningInFunction
64 << "On field " << iF.name() << " patch " << p.name()
65 << " patchField " << this->type()
66 << " : mapper does not map all values." << nl
67 << " To avoid this warning fully specify the mapping in derived"
68 << " patch fields." << endl;
69 }
70}
71
72
73template<class Type>
75(
76 const fvPatch& p,
78 const dictionary& dict
79)
80:
81 transformFvPatchField<Type>(p, iF, dict),
82 refValue_("refValue", dict, p.size()),
83 refGrad_("refGradient", dict, p.size()),
84 valueFraction_("valueFraction", dict, p.size())
85{
86 evaluate();
87}
88
89
90template<class Type>
92(
95)
96:
97 transformFvPatchField<Type>(ptf, iF),
98 refValue_(ptf.refValue_),
99 refGrad_(ptf.refGrad_),
100 valueFraction_(ptf.valueFraction_)
101{}
102
103
104// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105
106template<class Type>
108(
109 const fvPatchFieldMapper& m
110)
111{
113 refValue_.autoMap(m);
114 refGrad_.autoMap(m);
115 valueFraction_.autoMap(m);
116}
117
118
119template<class Type>
121(
122 const fvPatchField<Type>& ptf,
123 const labelList& addr
124)
125{
127
130
131 refValue_.rmap(dmptf.refValue_, addr);
132 refGrad_.rmap(dmptf.refGrad_, addr);
133 valueFraction_.rmap(dmptf.valueFraction_, addr);
134}
135
136
137template<class Type>
140{
141 const Field<Type> pif(this->patchInternalField());
142
143 tmp<Field<Type>> normalValue = transform(valueFraction_, refValue_);
144
145 tmp<Field<Type>> gradValue = pif + refGrad_/this->patch().deltaCoeffs();
146
147 tmp<Field<Type>> transformGradValue =
148 transform(I - valueFraction_, gradValue);
149
150 return
151 (normalValue + transformGradValue - pif)*
152 this->patch().deltaCoeffs();
153}
154
155
156template<class Type>
158{
159 if (!this->updated())
160 {
161 this->updateCoeffs();
162 }
163
164 tmp<Field<Type>> normalValue = transform(valueFraction_, refValue_);
165
166 tmp<Field<Type>> gradValue =
167 this->patchInternalField() + refGrad_/this->patch().deltaCoeffs();
168
169 tmp<Field<Type>> transformGradValue =
170 transform(I - valueFraction_, gradValue);
171
172 Field<Type>::operator=(normalValue + transformGradValue);
175}
176
177
178template<class Type>
181{
182 // static_assert(!std::is_arithmetic_v<T>, "direction-mixed with scalar??");
183 // if constexpr (!is_rotational_vectorspace_v<Type>)
184 // {
185 // // Rotational-invariant type
186 // return tmp<Field<Type>>::New(this->size(), Foam::zero{});
187 // }
188
189 {
190 vectorField diag(valueFraction_.size());
191
192 std::transform
193 (
194 valueFraction_.cbegin(),
195 valueFraction_.cend(),
196 diag.begin(),
197
198 [](const symmTensor& t)
199 {
200 return vector
201 (
202 sqrt(mag(t.xx())), sqrt(mag(t.yy())), sqrt(mag(t.zz()))
203 );
204 }
205 );
208 }
209}
210
211
212template<class Type>
214{
216 refValue_.writeEntry("refValue", os);
217 refGrad_.writeEntry("refGradient", os);
218 valueFraction_.writeEntry("valueFraction", os);
220}
221
222
223// ************************************************************************* //
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
void operator=(const Field< Type > &)
Copy assignment.
Definition Field.C:781
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
Base class for direction-mixed boundary conditions.
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 > > 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.
directionMixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal 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.
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
Intermediate layer (not used directly as a user boundary condition). The "value" entry is NO_READ,...
transformFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
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
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.
bool notNull(const T *ptr) noexcept
True if ptr is not a pointer (of type T) to the nullObject.
Definition nullObject.H:267
Vector< scalar > vector
Definition vector.H:57
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition symmTensor.H:55
dictionary dict