symmetryPlaneFvPatchField.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) 2013-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
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 \*---------------------------------------------------------------------------*/
27 
29 
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  const fvPatch& p,
38 )
39 :
41  symmetryPlanePatch_(refCast<const symmetryPlaneFvPatch>(p))
42 {}
43 
44 
45 template<class Type>
47 (
49  const fvPatch& p,
51  const fvPatchFieldMapper& mapper
52 )
53 :
54  basicSymmetryFvPatchField<Type>(ptf, p, iF, mapper),
55  symmetryPlanePatch_(refCast<const symmetryPlaneFvPatch>(p))
56 {
57  if (!isType<symmetryPlaneFvPatch>(this->patch()))
58  {
60  << "\n patch type '" << p.type()
61  << "' not constraint type '" << typeName << "'"
62  << "\n for patch " << p.name()
63  << " of field " << this->internalField().name()
64  << " in file " << this->internalField().objectPath()
66  }
67 }
68 
69 
70 template<class Type>
72 (
73  const fvPatch& p,
74  const DimensionedField<Type, volMesh>& iF,
75  const dictionary& dict
76 )
77 :
78  basicSymmetryFvPatchField<Type>(p, iF, dict),
79  symmetryPlanePatch_(refCast<const symmetryPlaneFvPatch>(p, dict))
80 {
81  if (!isType<symmetryPlaneFvPatch>(p))
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 
93 
94 template<class Type>
96 (
97  const symmetryPlaneFvPatchField<Type>& ptf
98 )
99 :
101  symmetryPlanePatch_(ptf.symmetryPlanePatch_)
102 {}
103 
104 
105 template<class Type>
107 (
110 )
111 :
112  basicSymmetryFvPatchField<Type>(ptf, iF),
113  symmetryPlanePatch_(ptf.symmetryPlanePatch_)
114 {}
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
119 template<class Type>
122 {
123  vector nHat(symmetryPlanePatch_.n());
124 
125  const Field<Type> iF(this->patchInternalField());
126 
127  return
128  (transform(I - 2.0*sqr(nHat), iF) - iF)
129  *(this->patch().deltaCoeffs()/2.0);
130 }
131 
132 
133 template<class Type>
135 {
136  if (!this->updated())
137  {
138  this->updateCoeffs();
139  }
140 
141  vector nHat(symmetryPlanePatch_.n());
142 
143  const Field<Type> iF(this->patchInternalField());
144 
146  (
147  (iF + transform(I - 2.0*sqr(nHat), iF))/2.0
148  );
151 }
152 
153 
154 template<class Type>
157 {
158  vector nHat(symmetryPlanePatch_.n());
159 
160  const vector diag
161  (
162  mag(nHat.component(vector::X)),
163  mag(nHat.component(vector::Y)),
164  mag(nHat.component(vector::Z))
165  );
166 
167  return tmp<Field<Type>>
168  (
169  new Field<Type>
170  (
171  this->size(),
172  transformMask<Type>
173  (
174  //pow<vector, pTraits<Type>::rank>(diag)
175  pow
176  (
177  diag,
179  ::type>::zero
180  )
181  )
182  )
183  );
184 }
185 
186 
187 // ************************************************************************* //
dictionary dict
Symmetry-plane patch.
virtual tmp< Field< Type > > snGradTransformDiag() const
Return face-gradient transform diagonal.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
commsTypes
Communications types.
Definition: UPstream.H:72
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:269
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Type & refCast(U &obj)
A dynamic_cast (for references). Generates a FatalError on failed casts and uses the virtual type() m...
Definition: typeInfo.H:159
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:75
symmetryPlaneFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field.
static const char *const typeName
Typename for Field.
Definition: Field.H:86
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
friend Ostream & operator(Ostream &, const Field< Type > &)
static const Identity< scalar > I
Definition: Identity.H:100
Generic templated field type.
Definition: Field.H:62
A FieldMapper for finite-volume patch fields.
Vector< scalar > vector
Definition: vector.H:57
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
virtual tmp< Field< Type > > snGrad() const
Return gradient at boundary.
This boundary condition enforces a symmetryPlane constraint.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets updated() to false.
Definition: fvPatchField.C:329
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
PtrList< volScalarField > & Y
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const DimensionedField< Type, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
Definition: fvPatchField.H:662
const std::string patch
OpenFOAM patch number as a std::string.
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...