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) 2011-2013 OpenFOAM Foundation
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"
27 #include "surfaceInterpolate.H"
28 #include "fvcDiv.H"
29 #include "fvcGrad.H"
30 #include "fvMatrices.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> >
47 (
48  const surfaceScalarField& gammaMagSf,
49  const surfaceScalarField& deltaCoeffs,
51 )
52 {
53  tmp<fvMatrix<Type> > tfvm
54  (
55  new fvMatrix<Type>
56  (
57  vf,
58  deltaCoeffs.dimensions()*gammaMagSf.dimensions()*vf.dimensions()
59  )
60  );
61  fvMatrix<Type>& fvm = tfvm();
62 
63  fvm.upper() = deltaCoeffs.internalField()*gammaMagSf.internalField();
64  fvm.negSumDiag();
65 
67  {
68  const fvPatchField<Type>& pvf = vf.boundaryField()[patchi];
69  const fvsPatchScalarField& pGamma = gammaMagSf.boundaryField()[patchi];
70  const fvsPatchScalarField& pDeltaCoeffs =
71  deltaCoeffs.boundaryField()[patchi];
72 
73  if (pvf.coupled())
74  {
75  fvm.internalCoeffs()[patchi] =
76  pGamma*pvf.gradientInternalCoeffs(pDeltaCoeffs);
77  fvm.boundaryCoeffs()[patchi] =
78  -pGamma*pvf.gradientBoundaryCoeffs(pDeltaCoeffs);
79  }
80  else
81  {
82  fvm.internalCoeffs()[patchi] = pGamma*pvf.gradientInternalCoeffs();
83  fvm.boundaryCoeffs()[patchi] = -pGamma*pvf.gradientBoundaryCoeffs();
84  }
85  }
86 
87  return tfvm;
88 }
89 
90 
91 template<class Type, class GType>
94 (
95  const surfaceVectorField& SfGammaCorr,
97 )
98 {
99  const fvMesh& mesh = this->mesh();
100 
102  (
104  (
105  IOobject
106  (
107  "gammaSnGradCorr("+vf.name()+')',
108  vf.instance(),
109  mesh,
112  ),
113  mesh,
114  SfGammaCorr.dimensions()
115  *vf.dimensions()*mesh.deltaCoeffs().dimensions()
116  )
117  );
118 
119  for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; cmpt++)
120  {
121  tgammaSnGradCorr().replace
122  (
123  cmpt,
124  SfGammaCorr & fvc::interpolate(fvc::grad(vf.component(cmpt)))
125  );
126  }
127 
128  return tgammaSnGradCorr;
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133 
134 template<class Type, class GType>
137 (
139 )
140 {
141  const fvMesh& mesh = this->mesh();
142 
144  (
145  fvc::div(this->tsnGradScheme_().snGrad(vf)*mesh.magSf())
146  );
147 
148  tLaplacian().rename("laplacian(" + vf.name() + ')');
149 
150  return tLaplacian;
151 }
152 
153 
154 template<class Type, class GType>
157 (
160 )
161 {
162  const fvMesh& mesh = this->mesh();
163 
164  const surfaceVectorField Sn(mesh.Sf()/mesh.magSf());
165 
166  const surfaceVectorField SfGamma(mesh.Sf() & gamma);
168  (
169  SfGamma & Sn
170  );
171  const surfaceVectorField SfGammaCorr(SfGamma - SfGammaSn*Sn);
172 
173  tmp<fvMatrix<Type> > tfvm = fvmLaplacianUncorrected
174  (
175  SfGammaSn,
176  this->tsnGradScheme_().deltaCoeffs(vf),
177  vf
178  );
179  fvMatrix<Type>& fvm = tfvm();
180 
182  = gammaSnGradCorr(SfGammaCorr, vf);
183 
184  if (this->tsnGradScheme_().corrected())
185  {
186  tfaceFluxCorrection() +=
187  SfGammaSn*this->tsnGradScheme_().correction(vf);
188  }
189 
190  fvm.source() -= mesh.V()*fvc::div(tfaceFluxCorrection())().internalField();
191 
192  if (mesh.fluxRequired(vf.name()))
193  {
194  fvm.faceFluxCorrectionPtr() = tfaceFluxCorrection.ptr();
195  }
196 
197  return tfvm;
198 }
199 
200 
201 template<class Type, class GType>
204 (
207 )
208 {
209  const fvMesh& mesh = this->mesh();
210 
211  const surfaceVectorField Sn(mesh.Sf()/mesh.magSf());
212  const surfaceVectorField SfGamma(mesh.Sf() & gamma);
214  (
215  SfGamma & Sn
216  );
217  const surfaceVectorField SfGammaCorr(SfGamma - SfGammaSn*Sn);
218 
220  (
221  fvc::div
222  (
223  SfGammaSn*this->tsnGradScheme_().snGrad(vf)
224  + gammaSnGradCorr(SfGammaCorr, vf)
225  )
226  );
227 
228  tLaplacian().rename("laplacian(" + gamma.name() + ',' + vf.name() + ')');
229 
230  return tLaplacian;
231 }
232 
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 } // End namespace fv
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 } // End namespace Foam
241 
242 // ************************************************************************* //
Foam::fvPatchField< Type >
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::GeometricField::component
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::fvc::interpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Definition: surfaceInterpolate.C:110
Foam::fv::gaussLaplacianScheme::fvcLaplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcLaplacian(const GeometricField< Type, fvPatchField, volMesh > &)
Definition: gaussLaplacianScheme.C:137
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
fvcDiv.H
Calculate the divergence of the given field.
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::fvMatrix::internalCoeffs
FieldField< Field, Type > & internalCoeffs()
fvBoundary scalar field containing pseudo-matrix coeffs
Definition: fvMatrix.H:303
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::fvMatrix::faceFluxCorrectionPtr
surfaceTypeFieldPtr & faceFluxCorrectionPtr()
Return pointer to face-flux non-orthogonal correction field.
Definition: fvMatrix.H:321
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::fvPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:342
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::fvPatchField::gradientInternalCoeffs
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
Definition: fvPatchField.H:452
Foam::fv::gaussLaplacianScheme::fvmLaplacian
tmp< fvMatrix< Type > > fvmLaplacian(const GeometricField< GType, fvsPatchField, surfaceMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Definition: gaussLaplacianScheme.C:157
Foam::fv::gaussLaplacianScheme::fvmLaplacianUncorrected
static tmp< fvMatrix< Type > > fvmLaplacianUncorrected(const surfaceScalarField &gammaMagSf, const surfaceScalarField &deltaCoeffs, const GeometricField< Type, fvPatchField, volMesh > &)
Definition: gaussLaplacianScheme.C:47
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
surfaceInterpolate.H
Surface Interpolation.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::tmp::ptr
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:146
Foam::fvMatrix::boundaryCoeffs
FieldField< Field, Type > & boundaryCoeffs()
fvBoundary scalar field containing pseudo-matrix coeffs
Definition: fvMatrix.H:310
internalField
conserve internalField()+
fv
labelList fv(nPoints)
Foam::fv::gaussLaplacianScheme::gammaSnGradCorr
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > gammaSnGradCorr(const surfaceVectorField &SfGammaCorr, const GeometricField< Type, fvPatchField, volMesh > &)
Definition: gaussLaplacianScheme.C:94
Foam::fvPatchField::gradientBoundaryCoeffs
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
Definition: fvPatchField.H:472
fvcGrad.H
Calculate the gradient of the given field.
Foam::fvMatrix::source
Field< Type > & source()
Definition: fvMatrix.H:291
patchi
label patchi
Definition: getPatchFieldScalar.H:1
gaussLaplacianScheme.H
Foam::fvMatrix< Type >
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52