relaxedNonOrthoGaussLaplacianSchemes.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) 2019 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "fvMesh.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 makeFvLaplacianScheme(relaxedNonOrthoGaussLaplacianScheme)
34 
35 #define declareFvmLaplacianScalarGamma(Type) \
36  \
37 template<> \
38 Foam::tmp<Foam::fvMatrix<Foam::Type>> \
39 Foam::fv::relaxedNonOrthoGaussLaplacianScheme<Foam::Type, Foam::scalar>:: \
40 fvmLaplacian \
41 ( \
42  const GeometricField<scalar, fvsPatchField, surfaceMesh>& gamma, \
43  const GeometricField<Type, fvPatchField, volMesh>& vf \
44 ) \
45 { \
46  const fvMesh& mesh = this->mesh(); \
47  \
48  typedef GeometricField<Type, fvsPatchField, surfaceMesh> SType; \
49  \
50  GeometricField<scalar, fvsPatchField, surfaceMesh> gammaMagSf \
51  ( \
52  gamma*mesh.magSf() \
53  ); \
54  \
55  tmp<fvMatrix<Type>> tfvm = fvmLaplacianUncorrected \
56  ( \
57  gammaMagSf, \
58  this->tsnGradScheme_().deltaCoeffs(vf), \
59  vf \
60  ); \
61  fvMatrix<Type>& fvm = tfvm.ref(); \
62  \
63  if (this->tsnGradScheme_().corrected()) \
64  { \
65  tmp<SType> tCorr(this->tsnGradScheme_().correction(vf)); \
66  const word corrName(tCorr().name()); \
67  tmp<SType> tfaceFluxCorrection(gammaMagSf*tCorr); \
68  \
69  tmp<SType> trelaxedCorrection(new SType(tfaceFluxCorrection())); \
70  \
71  const word oldName(corrName + "_0"); \
72  const scalar relax(vf.mesh().equationRelaxationFactor(corrName)); \
73  const objectRegistry& obr = vf.db(); \
74  if (obr.foundObject<SType>(oldName)) \
75  { \
76  SType& oldCorrection = obr.lookupObjectRef<SType>(oldName); \
77  trelaxedCorrection.ref() *= relax; \
78  trelaxedCorrection.ref() += (1.0-relax)*oldCorrection; \
79  \
80  oldCorrection = trelaxedCorrection(); \
81  } \
82  else \
83  { \
84  SType* s = new SType(oldName, tfaceFluxCorrection); \
85  s->store(); \
86  } \
87  \
88  tmp<Field<Type>> tcorr \
89  ( \
90  mesh.V() \
91  *fvc::div \
92  ( \
93  trelaxedCorrection() \
94  )().primitiveField() \
95  ); \
96  \
97  fvm.source() -= tcorr(); \
98  \
99  if (mesh.fluxRequired(vf.name())) \
100  { \
101  fvm.faceFluxCorrectionPtr() = trelaxedCorrection.ptr(); \
102  } \
103  } \
104  \
105  return tfvm; \
106 } \
107  \
108  \
109 template<> \
110 Foam::tmp<Foam::GeometricField<Foam::Type, Foam::fvPatchField, Foam::volMesh>> \
111 Foam::fv::relaxedNonOrthoGaussLaplacianScheme<Foam::Type, Foam::scalar>::fvcLaplacian \
112 ( \
113  const GeometricField<scalar, fvsPatchField, surfaceMesh>& gamma, \
114  const GeometricField<Type, fvPatchField, volMesh>& vf \
115 ) \
116 { \
117  const fvMesh& mesh = this->mesh(); \
118  \
119  tmp<GeometricField<Type, fvPatchField, volMesh>> tLaplacian \
120  ( \
121  fvc::div(gamma*this->tsnGradScheme_().snGrad(vf)*mesh.magSf()) \
122  ); \
123  \
124  tLaplacian.ref().rename \
125  ( \
126  "laplacian(" + gamma.name() + ',' + vf.name() + ')' \
127  ); \
128  \
129  return tLaplacian; \
130 }
131 
132 
138 
139 
140 // ************************************************************************* //
Foam::Tensor< scalar >
Foam::SymmTensor< scalar >
relaxedNonOrthoGaussLaplacianScheme.H
fvMesh.H
Foam::SphericalTensor< scalar >
declareFvmLaplacianScalarGamma
#define declareFvmLaplacianScalarGamma(Type)
Foam::Vector< scalar >
makeFvLaplacianScheme
makeFvLaplacianScheme(relaxedNonOrthoGaussLaplacianScheme)
Definition: relaxedNonOrthoGaussLaplacianSchemes.C:33