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-2018 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 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class Type>
34 (
35  const fvPatch& p,
37 )
38 :
40 {}
41 
42 
43 template<class Type>
45 (
47  const fvPatch& p,
49  const fvPatchFieldMapper& mapper
50 )
51 :
52  transformFvPatchField<Type>(ptf, p, iF, mapper)
53 {}
54 
55 
56 template<class Type>
58 (
59  const fvPatch& p,
61  const dictionary& dict
62 )
63 :
65 {
66  this->evaluate();
67 }
68 
69 
70 template<class Type>
72 (
74 )
75 :
77 {}
78 
79 
80 template<class Type>
82 (
85 )
86 :
88 {}
89 
90 
91 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92 
93 template<class Type>
96 {
97  tmp<vectorField> nHat = this->patch().nf();
98 
99  const Field<Type> iF(this->patchInternalField());
100 
101  return
102  (transform(I - 2.0*sqr(nHat), iF) - iF)
103  *(this->patch().deltaCoeffs()/2.0);
104 }
105 
106 
107 template<class Type>
109 {
110  if (!this->updated())
111  {
112  this->updateCoeffs();
113  }
114 
115  tmp<vectorField> nHat = this->patch().nf();
116 
117  const Field<Type> iF(this->patchInternalField());
118 
120  (
121  (iF + transform(I - 2.0*sqr(nHat), iF))/2.0
122  );
123 
125 }
126 
127 
128 template<class Type>
131 {
132  const vectorField nHat(this->patch().nf());
133 
134  vectorField diag(nHat.size());
135 
136  diag.replace(vector::X, mag(nHat.component(vector::X)));
137  diag.replace(vector::Y, mag(nHat.component(vector::Y)));
138  diag.replace(vector::Z, mag(nHat.component(vector::Z)));
139 
140  return transformFieldMask<Type>(pow<vector, pTraits<Type>::rank>(diag));
141 }
142 
143 
144 // ************************************************************************* //
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
commsTypes
Types of communications.
Definition: UPstream.H:64
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
static const Identity< scalar > I
Definition: Identity.H:93
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::fvPatchFieldMapper.
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:657
Foam::transformFvPatchField.
PtrList< volScalarField > & Y
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477
basicSymmetryFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.