regularisationRadiusIsotropic.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) 2021 PCOpt/NTUA
9  Copyright (C) 2021 FOSS GP
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "fvm.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
39 }
40 
41 
42 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
43 
45 {
46  scalar averageVol(gAverage(mesh_.V().field()));
47  const Vector<label>& geometricD = mesh_.geometricD();
48  const boundBox& bounds = mesh_.bounds();
49  forAll(geometricD, iDir)
50  {
51  if (geometricD[iDir] == -1)
52  {
53  averageVol /= bounds.span()[iDir];
54  }
55  }
56  scalar radius = pow(averageVol, scalar(1)/scalar(mesh_.nGeometricD()));
57 
58  scalar multMeanRadius = dict.getOrDefault<scalar>("meanRadiusMult", 10);
59  Info<< "Computed a mean radius of " << radius
60  << " and multiplying with " << multMeanRadius << endl;
61  return multMeanRadius*radius;
62 }
63 
64 
65 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
66 
67 Foam::isotropic::isotropic
68 (
69  const fvMesh& mesh,
70  const dictionary& dict,
71  bool adjustWallThickness
72 )
73 :
74  regularisationRadius(mesh, dict, adjustWallThickness),
75  radius_
76  (
77  "radius",
78  dimLength,
79  scalar
80  (
81  dict_.getOrDefault<scalar>("radius", computeRadius(dict))
82  /(2.*::sqrt(3.))
83  )
84  )
85 {
86  if (adjustWallThickness)
87  {
88  const scalar wallThicknessMult =
89  dict.getOrDefault<scalar>("wallThicknessMult", 0.75);
90  DebugInfo<<
91  "Adjusting wall thickness by " << wallThicknessMult << endl;
92  radius_ *= wallThicknessMult;
93  }
94  DebugInfo
95  << "Using radius " << radius_ << endl;
96 }
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
102 (
103  fvScalarMatrix& matrix,
104  bool isTopoField
105 ) const
106 {
107  const volScalarField& field = matrix.psi();
108  matrix -= fvm::laplacian(sqr(radius_), field);
109 }
110 
111 
112 // ************************************************************************* //
Base class for selecting the regulatisation radius.
dictionary dict
rDeltaTY field()
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
An isotropic regularisationRadius (same in all spatial directions)
virtual void addRegularisationTerm(fvScalarMatrix &matrix, bool isTopoField) const
Add a Laplacian term with an isotropic diffusivity.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
dimensionedScalar radius_
Smoothing radius of the first regulatisation.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:876
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:865
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:485
dynamicFvMesh & mesh
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
#define DebugInfo
Report an information message using Foam::Info.
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Type gAverage(const FieldField< Field, Type > &f)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:192
messageStream Info
Information stream (stdout output on master, null elsewhere)
const boundBox & bounds() const noexcept
Return mesh bounding box.
Definition: polyMesh.H:616
scalar computeRadius(const dictionary &dict)
Compute smoothing radius, if not directly given.
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)