Loading...
Searching...
No Matches
wedgeFvPatchField.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) 2024-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\*---------------------------------------------------------------------------*/
28
29#include "wedgeFvPatch.H"
31#include "transformField.H"
32#include "symmTransform.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36template<class Type>
38(
39 const fvPatch& p,
42:
43 parent_bctype(p, iF)
44{}
45
46
47template<class Type>
49(
50 const wedgeFvPatchField<Type>& ptf,
51 const fvPatch& p,
53 const fvPatchFieldMapper& mapper
54)
55:
56 parent_bctype(ptf, p, iF, mapper)
57{
58 if (!isType<wedgeFvPatch>(this->patch()))
59 {
61 << "\n patch type '" << p.type()
62 << "' not constraint type '" << typeName << "'"
63 << "\n for patch " << p.name()
64 << " of field " << this->internalField().name()
65 << " in file " << this->internalField().objectPath()
66 << exit(FatalError);
67 }
68}
69
70
71template<class Type>
73(
74 const fvPatch& p,
75 const DimensionedField<Type, volMesh>& iF,
76 const dictionary& dict
77)
78:
79 parent_bctype(p, iF, dict) // "value" is NO_READ
80{
82 {
84 << "\n patch type '" << p.type()
85 << "' not constraint type '" << typeName << "'"
86 << "\n for patch " << p.name()
87 << " of field " << this->internalField().name()
88 << " in file " << this->internalField().objectPath()
90 }
91
92 evaluate();
93}
94
95
96template<class Type>
98(
99 const wedgeFvPatchField<Type>& ptf,
100 const DimensionedField<Type, volMesh>& iF
102:
103 parent_bctype(ptf, iF)
104{}
105
106
107template<class Type>
109(
110 const wedgeFvPatchField<Type>& ptf
111)
112:
114{}
115
116
117// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118
119template<class Type>
121{
123 {
124 // Rotational-invariant type : treat like zero-gradient
125 return tmp<Field<Type>>::New(this->size(), Foam::zero{});
126 }
127 else
128 {
129 const auto& rot = refCast<const wedgeFvPatch>(this->patch()).cellT();
130
131 const Field<Type> pif(this->patchInternalField());
132
133 const auto& dc = this->patch().deltaCoeffs();
134
135 return
136 (
137 (0.5*dc)
138 * (transform(rot, pif) - pif)
139 );
140 }
141}
142
143
144template<class Type>
146{
147 if (!this->updated())
148 {
149 this->updateCoeffs();
150 }
151
153 {
154 // Rotational-invariant type : treat like zero-gradient
155 this->extrapolateInternal();
156 }
157 else
158 {
159 const auto& rot = refCast<const wedgeFvPatch>(this->patch()).faceT();
160
162 (
163 transform(rot, this->patchInternalField())
164 );
165 }
166}
167
168
169template<class Type>
172{
174 {
175 // Rotational-invariant type
176 // FatalErrorInFunction
177 // << "Should not be called for this type"
178 // << ::Foam::abort(FatalError);
179 return tmp<Field<Type>>::New(this->size(), Foam::zero{});
180 }
181 else
182 {
183 const auto& rot = refCast<const wedgeFvPatch>(this->patch()).cellT();
184
185 const vector diag = 0.5*(I - rot).diag();
186
188 (
189 this->size(),
191 (
192 pow
193 (
194 diag,
195 pTraits
196 <
198 >::zero
199 )
200 )
201 );
202 }
203}
204
205
206// ************************************************************************* //
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
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
fileName objectPath() const
The complete path + object name.
Definition IOobjectI.H:313
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
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 tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
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.
void extrapolateInternal()
Assign the patch field from the internal field.
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
This boundary condition is similar to the cyclic condition, except that it is applied to 2-D geometri...
wedgeFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
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 representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
volScalarField & p
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const std::string patch
OpenFOAM patch number as a std::string.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
static const Identity< scalar > I
Definition Identity.H:100
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Type1 transformMask(const Type2 &t)
Definition transform.H:435
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
bool isType(const U &obj)
Check if typeid of the object and Type are identical.
Definition typeInfo.H:112
constexpr bool is_rotational_vectorspace_v
The is_rotational_vectorspace value of Type.
Definition pTraits.H:251
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dictionary dict
A non-counting (dummy) refCount.
Definition refCount.H:55
3D symmetric tensor transformation operations.
Spatial transformation functions for primitive fields.