Loading...
Searching...
No Matches
transformFaPatchField.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) 2016-2017 Wikki Ltd
9 Copyright (C) 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\*---------------------------------------------------------------------------*/
30#include "transformField.H"
31
32// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
33
34template<class Type>
36(
37 const faPatch& p,
40:
41 faPatchField<Type>(p, iF)
42{}
43
44
45template<class Type>
47(
49 const faPatch& p,
51 const faPatchFieldMapper& mapper
53:
54 faPatchField<Type>(ptf, p, iF, mapper)
55{}
56
57
58template<class Type>
60(
61 const faPatch& p,
63 const dictionary& dict
65:
66 faPatchField<Type>(p, iF, dict, IOobjectOption::NO_READ)
67{}
68
69
70template<class Type>
72(
75)
76:
77 faPatchField<Type>(ptf, iF)
78{}
79
80
81// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82
83template<class Type>
86(
87 const tmp<scalarField>&
88) const
89{
91 {
92 // Rotational-invariant type
93 return tmp<Field<Type>>::New(this->size(), pTraits<Type>::one);
94 }
95 else
96 {
98 }
99}
100
101
102template<class Type>
105(
106 const tmp<scalarField>&
107) const
108{
109 return
110 *this
112 (
113 valueInternalCoeffs(this->patch().weights()),
115 );
116}
117
118
119template<class Type>
122{
124 {
125 // Rotational-invariant type
126 return tmp<Field<Type>>::New(this->size(), Foam::zero{});
127 }
128 else
129 {
130 return -this->patch().deltaCoeffs()*snGradTransformDiag();
131 }
132}
133
134
135template<class Type>
138{
139 return
140 snGrad()
142}
143
144
145// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
146
147template<class Type>
149(
150 const faPatchField<Type>& ptf
151)
152{
153 this->evaluate();
154}
155
156
157// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
void size(const label n)
Definition UList.H:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const faPatch & patch() const noexcept
Return the patch.
A FieldMapper for finite-area patch fields.
faPatchField<Type> abstract base class. This class gives a fat-interface to all derived classes cover...
faPatchField(const faPatch &, const DimensionedField< Type, areaMesh > &)
Construct from patch and internal field.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field, sets updated() to false.
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition faPatch.H:76
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
Author Zeljko Tukovic, FMENA Hrvoje Jasak, Wikki Ltd.
virtual tmp< Field< Type > > snGradTransformDiag() const =0
Return face-gradient transform diagonal.
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the.
virtual tmp< Field< Type > > snGrad() const =0
Return gradient at boundary.
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
transformFaPatchField(const faPatch &, const DimensionedField< Type, areaMesh > &)
Construct from patch and internal 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
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.
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
constexpr bool is_rotational_vectorspace_v
The is_rotational_vectorspace value of Type.
Definition pTraits.H:251
dictionary dict
Spatial transformation functions for primitive fields.