BUGFIX: Uninitialised member variables
[foam-extend-3.2.git] / applications / solvers / surfaceTracking / freeSurface / patchCorrectedSnGrad / patchCorrectedSnGrad.C
blob0e85c1a98cb55759001c5b7ac3d39066ec3a1a82
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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
25 Description
26     Simple central-difference snGrad scheme with non-orthogonal correction.
28 \*---------------------------------------------------------------------------*/
30 #include "patchCorrectedSnGrad.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 #include "linear.H"
34 #include "fvcGrad.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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 namespace Foam
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 namespace fv
52 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
54 template<class Type>
55 patchCorrectedSnGrad<Type>::~patchCorrectedSnGrad()
59 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
61 template<class Type>
62 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
63 patchCorrectedSnGrad<Type>::correction
65     const GeometricField<Type, fvPatchField, volMesh>& vf
66 ) const
68     const fvMesh& mesh = this->mesh();
70     // construct GeometricField<Type, fvsPatchField, surfaceMesh>
71     tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tssf
72     (
73         new GeometricField<Type, fvsPatchField, surfaceMesh>
74         (
75             IOobject
76             (
77                 "snGradCorr("+vf.name()+')',
78                 vf.instance(),
79                 mesh,
80                 IOobject::NO_READ,
81                 IOobject::NO_WRITE
82             ),
83             mesh,
84             vf.dimensions()*mesh.deltaCoeffs().dimensions()
85         )
86     );
87     GeometricField<Type, fvsPatchField, surfaceMesh>& ssf = tssf();
90     for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; cmpt++)
91     {
92         ssf.replace
93         (
94             cmpt,
95             mesh.correctionVectors()
96           & linear
97             <
98                 typename
99                 outerProduct<vector, typename pTraits<Type>::cmptType>::type
100             >(mesh).interpolate
101             (
102                 gradScheme<typename pTraits<Type>::cmptType>::New
103                 (
104                     mesh,
105                     mesh.schemesDict().gradScheme(ssf.name())
106                 )()
107                 //gaussGrad<typename pTraits<Type>::cmptType>(mesh)
108                .grad(vf.component(cmpt))
109             )
110         );
111     }
113     forAll(ssf.boundaryField(), patchI)
114     {
115         if
116         (
117             vf.boundaryField()[patchI].type()
118          == fixedValueCorrectedFvPatchField<Type>::typeName
119         )
120         {
121             const fixedValueCorrectedFvPatchField<Type>& pField =
122                 refCast<const fixedValueCorrectedFvPatchField<Type> >
123                 (
124                     vf.boundaryField()[patchI]
125                 );
127             ssf.boundaryField()[patchI] = pField.snGradCorrection();
128         }
129     }
131     return tssf;
135 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 } // End namespace fv
139 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141 } // End namespace Foam
143 // ************************************************************************* //