basicSymmetryFvPatchField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "symmTransformField.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const fvPatch& p,
36 )
37 :
38  transformFvPatchField<Type>(p, iF)
39 {}
40 
41 
42 template<class Type>
44 (
45  const fvPatch& p,
47  const dictionary& dict
48 )
49 :
50  transformFvPatchField<Type>(p, iF, dict)
51 {
52  this->evaluate();
53 }
54 
55 
56 template<class Type>
58 (
60  const fvPatch& p,
62  const fvPatchFieldMapper& mapper
63 )
64 :
65  transformFvPatchField<Type>(ptf, p, iF, mapper)
66 {}
67 
68 
69 template<class Type>
71 (
74 )
75 :
76  transformFvPatchField<Type>(ptf, iF)
77 {}
78 
79 
80 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
81 
82 template<class Type>
85 {
86  tmp<vectorField> nHat = this->patch().nf();
87 
88  const Field<Type> iF(this->patchInternalField());
89 
90  return
91  (transform(I - 2.0*sqr(nHat), iF) - iF)
92  *(this->patch().deltaCoeffs()/2.0);
93 }
94 
95 
96 template<class Type>
98 {
99  if (!this->updated())
100  {
101  this->updateCoeffs();
102  }
103 
104  tmp<vectorField> nHat = this->patch().nf();
105 
106  const Field<Type> iF(this->patchInternalField());
107 
109  (
110  (iF + transform(I - 2.0*sqr(nHat), iF))/2.0
111  );
112 
114 }
115 
116 
117 template<class Type>
120 {
121  const vectorField nHat(this->patch().nf());
122 
123  vectorField diag(nHat.size());
124 
125  diag.replace(vector::X, mag(nHat.component(vector::X)));
126  diag.replace(vector::Y, mag(nHat.component(vector::Y)));
127  diag.replace(vector::Z, mag(nHat.component(vector::Z)));
128 
129  return transformFieldMask<Type>(pow<vector, pTraits<Type>::rank>(diag));
130 }
131 
132 
133 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Pre-declare SubField and related Field type.
Definition: Field.H:82
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:455
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
commsTypes
Types of communications.
Definition: UPstream.H:65
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
basicSymmetryFvPatchField(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.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Foam::fvPatchFieldMapper.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:211
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Traits class for primitives.
Definition: pTraits.H:53
A class for managing temporary objects.
Definition: tmp.H:55
Foam::transformFvPatchField.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static const Identity< scalar > I
Definition: Identity.H:93
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:483
dimensioned< scalar > mag(const dimensioned< Type > &)
dictionary dict
volScalarField & p