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  Copyright (C) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 \*---------------------------------------------------------------------------*/
28 
29 #include "correctedSnGrad.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "linear.H"
33 #include "fvcGrad.H"
34 #include "gaussGrad.H"
35 
36 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
37 
38 template<class Type>
41 (
43 ) const
44 {
45  const fvMesh& mesh = this->mesh();
46 
47  // construct GeometricField<Type, fvsPatchField, surfaceMesh>
50  (
51  mesh.nonOrthCorrectionVectors(),
53  (
54  mesh,
55  mesh.gradScheme("grad(" + vf.name() + ')')
56  )().grad(vf, "grad(" + vf.name() + ')')
57  );
58  tssf.ref().rename("snGradCorr(" + vf.name() + ')');
59 
60  return tssf;
61 }
62 
63 
64 template<class Type>
67 (
69 ) const
70 {
71  const fvMesh& mesh = this->mesh();
72 
73  // construct GeometricField<Type, fvsPatchField, surfaceMesh>
75  (
77  (
78  IOobject
79  (
80  "snGradCorr("+vf.name()+')',
81  vf.instance(),
82  mesh,
83  IOobject::NO_READ,
84  IOobject::NO_WRITE
85  ),
86  mesh,
87  vf.dimensions()*mesh.nonOrthDeltaCoeffs().dimensions()
88  )
89  );
91  ssf.setOriented();
92 
93  for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; ++cmpt)
94  {
95  ssf.replace
96  (
97  cmpt,
99  .fullGradCorrection(vf.component(cmpt))
100  );
101  }
102 
103  return tssf;
104 }
105 
106 
107 // ************************************************************************* //
volFields.H
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
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
surfaceFields.H
Foam::surfaceFields.
gaussGrad.H
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
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
A traits class, which is primarily used for primitives.
Definition: pTraits.H:56
fvcGrad.H
Calculate the gradient of the given field.
Foam::fv::correctedSnGrad
Surface gradient scheme with full explicit non-orthogonal correction.
Definition: correctedSnGrad.H:69
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
Definition: correctedSnGrad.C:41
Foam::fv::correctedSnGrad::correction
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Definition: correctedSnGrad.C:67