correctedSnGrad.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 -------------------------------------------------------------------------------
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 
28 #include "correctedSnGrad.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 #include "linear.H"
32 #include "fvcGrad.H"
33 #include "gaussGrad.H"
34 
35 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
36 
37 template<class Type>
39 {}
40 
41 
42 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
43 
44 template<class Type>
47 (
49 ) const
50 {
51  const fvMesh& mesh = this->mesh();
52 
53  // construct GeometricField<Type, fvsPatchField, surfaceMesh>
56  (
57  mesh.nonOrthCorrectionVectors(),
59  (
60  mesh,
61  mesh.gradScheme("grad(" + vf.name() + ')')
62  )().grad(vf, "grad(" + vf.name() + ')')
63  );
64  tssf.ref().rename("snGradCorr(" + vf.name() + ')');
65 
66  return tssf;
67 }
68 
69 
70 template<class Type>
73 (
75 ) const
76 {
77  const fvMesh& mesh = this->mesh();
78 
79  // construct GeometricField<Type, fvsPatchField, surfaceMesh>
81  (
83  (
84  IOobject
85  (
86  "snGradCorr("+vf.name()+')',
87  vf.instance(),
88  mesh,
89  IOobject::NO_READ,
90  IOobject::NO_WRITE
91  ),
92  mesh,
93  vf.dimensions()*mesh.nonOrthDeltaCoeffs().dimensions()
94  )
95  );
97  ssf.setOriented();
98 
99  for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; cmpt++)
100  {
101  ssf.replace
102  (
103  cmpt,
105  .fullGradCorrection(vf.component(cmpt))
106  );
107  }
108 
109  return tssf;
110 }
111 
112 
113 // ************************************************************************* //
volFields.H
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::GeometricField::component
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
Foam::GeometricField::replace
void replace(const direction d, const GeometricField< cmptType, PatchField, GeoMesh > &gcf)
Replace specified field component with content from another field.
Foam::fvc::dotInterpolate
static tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &tvf)
Interpolate field onto faces.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::fac::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:56
Foam::fv::correctedSnGrad::~correctedSnGrad
virtual ~correctedSnGrad()
Destructor.
Definition: correctedSnGrad.C:38
surfaceFields.H
Foam::surfaceFields.
gaussGrad.H
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
correctedSnGrad.H
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:54
fvcGrad.H
Calculate the gradient of the given field.
Foam::fv::correctedSnGrad
Simple central-difference snGrad scheme with non-orthogonal correction.
Definition: correctedSnGrad.H:59
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::linear
Central-differencing interpolation scheme class.
Definition: linear.H:55
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::fv::correctedSnGrad::fullGradCorrection
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > fullGradCorrection(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the correctedSnGrad.
Definition: correctedSnGrad.C:47
Foam::fv::correctedSnGrad::correction
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the correctedSnGrad.
Definition: correctedSnGrad.C:73