FreeFOAM The Cross-Platform CFD Toolkit
gaussLaplacianScheme.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 | Copyright (C) 1991-2010 OpenCFD Ltd.
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "gaussLaplacianScheme.H"
28 #include <finiteVolume/fvcDiv.H>
29 #include <finiteVolume/fvcGrad.H>
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace fv
40 {
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 template<class Type, class GType>
45 tmp<fvMatrix<Type> >
46 gaussLaplacianScheme<Type, GType>::fvmLaplacianUncorrected
47 (
48  const surfaceScalarField& gammaMagSf,
49  GeometricField<Type, fvPatchField, volMesh>& vf
50 )
51 {
52  tmp<surfaceScalarField> tdeltaCoeffs =
53  this->tsnGradScheme_().deltaCoeffs(vf);
54  const surfaceScalarField& deltaCoeffs = tdeltaCoeffs();
55 
56  tmp<fvMatrix<Type> > tfvm
57  (
58  new fvMatrix<Type>
59  (
60  vf,
61  deltaCoeffs.dimensions()*gammaMagSf.dimensions()*vf.dimensions()
62  )
63  );
64  fvMatrix<Type>& fvm = tfvm();
65 
66  fvm.upper() = deltaCoeffs.internalField()*gammaMagSf.internalField();
67  fvm.negSumDiag();
68 
69  forAll(fvm.psi().boundaryField(), patchI)
70  {
71  const fvPatchField<Type>& psf = fvm.psi().boundaryField()[patchI];
72  const fvsPatchScalarField& patchGamma =
73  gammaMagSf.boundaryField()[patchI];
74 
75  fvm.internalCoeffs()[patchI] = patchGamma*psf.gradientInternalCoeffs();
76  fvm.boundaryCoeffs()[patchI] = -patchGamma*psf.gradientBoundaryCoeffs();
77  }
78 
79  return tfvm;
80 }
81 
82 
83 template<class Type, class GType>
84 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
85 gaussLaplacianScheme<Type, GType>::gammaSnGradCorr
86 (
87  const surfaceVectorField& SfGammaCorr,
88  const GeometricField<Type, fvPatchField, volMesh>& vf
89 )
90 {
91  const fvMesh& mesh = this->mesh();
92 
93  tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tgammaSnGradCorr
94  (
95  new GeometricField<Type, fvsPatchField, surfaceMesh>
96  (
97  IOobject
98  (
99  "gammaSnGradCorr("+vf.name()+')',
100  vf.instance(),
101  mesh,
104  ),
105  mesh,
106  SfGammaCorr.dimensions()
107  *vf.dimensions()*mesh.deltaCoeffs().dimensions()
108  )
109  );
110 
111  for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; cmpt++)
112  {
113  tgammaSnGradCorr().replace
114  (
115  cmpt,
116  SfGammaCorr & fvc::interpolate(fvc::grad(vf.component(cmpt)))
117  );
118  }
119 
120  return tgammaSnGradCorr;
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126 template<class Type, class GType>
127 tmp<GeometricField<Type, fvPatchField, volMesh> >
129 (
131 )
132 {
133  const fvMesh& mesh = this->mesh();
134 
136  (
137  fvc::div(this->tsnGradScheme_().snGrad(vf)*mesh.magSf())
138  );
139 
140  tLaplacian().rename("laplacian(" + vf.name() + ')');
141 
142  return tLaplacian;
143 }
144 
145 
146 template<class Type, class GType>
149 (
152 )
153 {
154  const fvMesh& mesh = this->mesh();
155 
156  surfaceVectorField Sn = mesh.Sf()/mesh.magSf();
157 
158  surfaceVectorField SfGamma = mesh.Sf() & gamma;
160  surfaceVectorField SfGammaCorr = SfGamma - SfGammaSn*Sn;
161 
162  tmp<fvMatrix<Type> > tfvm = fvmLaplacianUncorrected(SfGammaSn, vf);
163  fvMatrix<Type>& fvm = tfvm();
164 
166  = gammaSnGradCorr(SfGammaCorr, vf);
167 
168  if (this->tsnGradScheme_().corrected())
169  {
170  tfaceFluxCorrection() +=
171  SfGammaSn*this->tsnGradScheme_().correction(vf);
172  }
173 
174  fvm.source() -= mesh.V()*fvc::div(tfaceFluxCorrection())().internalField();
175 
176  if (mesh.fluxRequired(vf.name()))
177  {
178  fvm.faceFluxCorrectionPtr() = tfaceFluxCorrection.ptr();
179  }
180 
181  return tfvm;
182 }
183 
184 
185 template<class Type, class GType>
188 (
191 )
192 {
193  const fvMesh& mesh = this->mesh();
194 
195  surfaceVectorField Sn = mesh.Sf()/mesh.magSf();
196 
197  surfaceVectorField SfGamma = mesh.Sf() & gamma;
199  surfaceVectorField SfGammaCorr = SfGamma - SfGammaSn*Sn;
200 
202  (
203  fvc::div
204  (
205  SfGammaSn*this->tsnGradScheme_().snGrad(vf)
206  + gammaSnGradCorr(SfGammaCorr, vf)
207  )
208  );
209 
210  tLaplacian().rename("laplacian(" + gamma.name() + ',' + vf.name() + ')');
211 
212  return tLaplacian;
213 }
214 
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 } // End namespace fv
219 
220 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 
222 } // End namespace Foam
223 
224 // ************************ vim: set sw=4 sts=4 et: ************************ //