1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
26 Simple central-difference snGrad scheme with non-orthogonal correction.
28 \*---------------------------------------------------------------------------*/
30 #include "patchCorrectedSnGrad.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
35 #include "gaussGrad.H"
36 #include "fixedValueCorrectedFvPatchFields.H"
37 #include "fixedGradientCorrectedFvPatchFields.H"
38 #include "zeroGradientCorrectedFvPatchFields.H"
39 #include "fixedGradientFvPatchFields.H"
40 #include "zeroGradientFvPatchFields.H"
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
55 patchCorrectedSnGrad<Type>::~patchCorrectedSnGrad()
59 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
63 patchCorrectedSnGrad<Type>::correction
65 const GeometricField<Type, fvPatchField, volMesh>& vf
68 const fvMesh& mesh = this->mesh();
70 // construct GeometricField<Type, fvsPatchField, surfaceMesh>
71 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tssf
73 new GeometricField<Type, fvsPatchField, surfaceMesh>
77 "snGradCorr("+vf.name()+')',
84 vf.dimensions()*mesh.deltaCoeffs().dimensions()
87 GeometricField<Type, fvsPatchField, surfaceMesh>& ssf = tssf();
90 for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; cmpt++)
95 mesh.correctionVectors()
99 outerProduct<vector, typename pTraits<Type>::cmptType>::type
102 gradScheme<typename pTraits<Type>::cmptType>::New
105 mesh.schemesDict().gradScheme(ssf.name())
107 //gaussGrad<typename pTraits<Type>::cmptType>(mesh)
108 .grad(vf.component(cmpt))
113 forAll(ssf.boundaryField(), patchI)
117 vf.boundaryField()[patchI].type()
118 == fixedValueCorrectedFvPatchField<Type>::typeName
121 const fixedValueCorrectedFvPatchField<Type>& pField =
122 refCast<const fixedValueCorrectedFvPatchField<Type> >
124 vf.boundaryField()[patchI]
127 ssf.boundaryField()[patchI] = pField.snGradCorrection();
135 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 } // End namespace fv
139 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141 } // End namespace Foam
143 // ************************************************************************* //