relaxedNonOrthoGaussLaplacianScheme.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 "surfaceInterpolate.H"
30 #include "fvcDiv.H"
31 #include "fvcGrad.H"
32 #include "fvMatrices.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace fv
42 {
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 template<class Type, class GType>
49 (
50  const surfaceScalarField& gammaMagSf,
51  const surfaceScalarField& deltaCoeffs,
53 )
54 {
55  tmp<fvMatrix<Type>> tfvm
56  (
57  new fvMatrix<Type>
58  (
59  vf,
60  deltaCoeffs.dimensions()*gammaMagSf.dimensions()*vf.dimensions()
61  )
62  );
63  fvMatrix<Type>& fvm = tfvm.ref();
64 
65  fvm.upper() = deltaCoeffs.primitiveField()*gammaMagSf.primitiveField();
66  fvm.negSumDiag();
67 
68  forAll(vf.boundaryField(), patchi)
69  {
70  const fvPatchField<Type>& pvf = vf.boundaryField()[patchi];
71  const fvsPatchScalarField& pGamma = gammaMagSf.boundaryField()[patchi];
72  const fvsPatchScalarField& pDeltaCoeffs =
73  deltaCoeffs.boundaryField()[patchi];
74 
75  if (pvf.coupled())
76  {
77  fvm.internalCoeffs()[patchi] =
78  pGamma*pvf.gradientInternalCoeffs(pDeltaCoeffs);
79  fvm.boundaryCoeffs()[patchi] =
80  -pGamma*pvf.gradientBoundaryCoeffs(pDeltaCoeffs);
81  }
82  else
83  {
84  fvm.internalCoeffs()[patchi] = pGamma*pvf.gradientInternalCoeffs();
85  fvm.boundaryCoeffs()[patchi] = -pGamma*pvf.gradientBoundaryCoeffs();
86  }
87  }
88 
89  return tfvm;
90 }
91 
92 
93 template<class Type, class GType>
96 (
97  const surfaceVectorField& SfGammaCorr,
99 )
100 {
101  const fvMesh& mesh = this->mesh();
102 
104  (
106  (
107  IOobject
108  (
109  "gammaSnGradCorr("+vf.name()+')',
110  vf.instance(),
111  mesh,
114  ),
115  mesh,
116  SfGammaCorr.dimensions()
117  *vf.dimensions()*mesh.deltaCoeffs().dimensions()
118  )
119  );
120 
121  for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; cmpt++)
122  {
123  tgammaSnGradCorr.ref().replace
124  (
125  cmpt,
126  fvc::dotInterpolate(SfGammaCorr, fvc::grad(vf.component(cmpt)))
127  );
128  }
129 
130  return tgammaSnGradCorr;
131 }
133 
134 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
135 
136 template<class Type, class GType>
139 (
141 )
142 {
143  const fvMesh& mesh = this->mesh();
144 
146  (
147  fvc::div(this->tsnGradScheme_().snGrad(vf)*mesh.magSf())
148  );
149 
150  tLaplacian.ref().rename("laplacian(" + vf.name() + ')');
151 
152  return tLaplacian;
153 }
154 
155 
156 template<class Type, class GType>
159 (
162 )
163 {
164  const fvMesh& mesh = this->mesh();
165 
167 
168  const surfaceVectorField Sn(mesh.Sf()/mesh.magSf());
169 
170  const surfaceVectorField SfGamma(mesh.Sf() & gamma);
172  (
173  SfGamma & Sn
174  );
175  const surfaceVectorField SfGammaCorr(SfGamma - SfGammaSn*Sn);
176 
177  tmp<fvMatrix<Type>> tfvm = fvmLaplacianUncorrected
178  (
179  SfGammaSn,
180  this->tsnGradScheme_().deltaCoeffs(vf),
181  vf
182  );
183  fvMatrix<Type>& fvm = tfvm.ref();
184 
185  tmp<SType> tfaceFluxCorrection = gammaSnGradCorr(SfGammaCorr, vf);
186 
187  if (this->tsnGradScheme_().corrected())
188  {
189  tfaceFluxCorrection.ref() +=
190  SfGammaSn*this->tsnGradScheme_().correction(vf);
191  }
192 
193  const word corrName(tfaceFluxCorrection().name());
194 
195  tmp<SType> trelaxedCorrection(new SType(tfaceFluxCorrection()));
196 
197  const word oldName(corrName + "_0");
198  const scalar relax(vf.mesh().equationRelaxationFactor(oldName));
199 
200  const objectRegistry& obr = vf.db();
201  if (obr.foundObject<SType>(oldName))
202  {
203  SType& oldCorrection = obr.lookupObjectRef<SType>(oldName);
204 
205  trelaxedCorrection.ref() *= relax;
206  trelaxedCorrection.ref() += (1.0-relax)*oldCorrection;
207 
208  oldCorrection = tfaceFluxCorrection;
209  }
210  else
211  {
212  SType* s = new SType(oldName, tfaceFluxCorrection);
213  s->store();
214  }
215 
216  fvm.source() -=
217  mesh.V()
218  *fvc::div
219  (
220  trelaxedCorrection()
221  )().primitiveField();
222 
223  if (mesh.fluxRequired(vf.name()))
224  {
225  fvm.faceFluxCorrectionPtr() = trelaxedCorrection.ptr();
226  }
227 
228  return tfvm;
229 }
230 
231 
232 template<class Type, class GType>
235 (
238 )
239 {
240  const fvMesh& mesh = this->mesh();
241 
242  const surfaceVectorField Sn(mesh.Sf()/mesh.magSf());
243  const surfaceVectorField SfGamma(mesh.Sf() & gamma);
245  (
246  SfGamma & Sn
247  );
248  const surfaceVectorField SfGammaCorr(SfGamma - SfGammaSn*Sn);
249 
251  (
252  fvc::div
253  (
254  SfGammaSn*this->tsnGradScheme_().snGrad(vf)
255  + gammaSnGradCorr(SfGammaCorr, vf)
256  )
257  );
258 
259  tLaplacian.ref().rename
260  (
261  "laplacian(" + gamma.name() + ',' + vf.name() + ')'
262  );
263 
264  return tLaplacian;
265 }
266 
267 
268 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
269 
270 } // End namespace fv
271 
272 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
273 
274 } // End namespace Foam
275 
276 // ************************************************************************* //
virtual bool coupled() const
True if the patch field is coupled.
Definition: fvPatchField.H:253
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
uint8_t direction
Definition: direction.H:46
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
UEqn relax()
static scalar Sn(const scalar a, const scalar x)
Definition: invIncGamma.C:64
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
static tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &tvf)
Interpolate field onto faces.
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
Ignore writing from objectRegistry::writeObject()
scalarField & upper()
Definition: lduMatrix.C:208
Basic second-order laplacian using face-gradients and Gauss&#39; theorem.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Calculate the gradient of the given field.
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
Definition: fvPatchField.H:828
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:450
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcLaplacian(const GeometricField< Type, fvPatchField, volMesh > &)
Calculate the divergence of the given field.
const Mesh & mesh() const noexcept
Return mesh.
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:266
faceFluxFieldPtrType & faceFluxCorrectionPtr()
Return pointer to face-flux non-orthogonal correction field.
Definition: fvMatrix.H:588
Field< Type > & source() noexcept
Definition: fvMatrix.H:533
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
Definition: fvPatchField.H:854
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
const FieldField< Field, Type > & internalCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for internal cells
Definition: fvMatrix.H:547
const FieldField< Field, Type > & boundaryCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for boundary cells
Definition: fvMatrix.H:565
const scalar gamma
Definition: EEqn.H:9
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static tmp< fvMatrix< Type > > fvmLaplacianUncorrected(const surfaceScalarField &gammaMagSf, const surfaceScalarField &deltaCoeffs, const GeometricField< Type, fvPatchField, volMesh > &)
Registry of regIOobjects.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
tmp< fvMatrix< Type > > fvmLaplacian(const GeometricField< GType, fvsPatchField, surfaceMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
Namespace for OpenFOAM.
const dimensionSet & dimensions() const noexcept
Return dimensions.