advectionDiffusionPatchDistMethod.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) 2015-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "surfaceInterpolate.H"
31#include "fvcGrad.H"
32#include "fvcDiv.H"
33#include "fvmDiv.H"
34#include "fvmLaplacian.H"
35#include "fvmSup.H"
37
40
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
43namespace Foam
44{
45namespace patchDistMethods
46{
49}
50}
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
55(
56 const dictionary& dict,
57 const fvMesh& mesh,
58 const labelHashSet& patchIDs
59)
60:
61 patchDistMethod(mesh, patchIDs),
62 //- We do not want to recurse into 'advectionDiffusion' again so
63 // make sure we pick up 'method' always from the subdictionary.
64 //coeffs_(dict.optionalSubDict(type() + "Coeffs")),
65 coeffs_(dict.subDict(type() + "Coeffs")),
66 pdmPredictor_
67 (
69 (
70 coeffs_,
71 mesh,
72 patchIDs
73 )
74 ),
75 epsilon_(coeffs_.getOrDefault<scalar>("epsilon", 0.1)),
76 tolerance_(coeffs_.getOrDefault<scalar>("tolerance", 1e-3)),
77 maxIter_(coeffs_.getOrDefault<int>("maxIter", 10)),
78 predicted_(false)
79{}
80
81
82// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
83
85{
86 return correct(y, const_cast<volVectorField&>(volVectorField::null()));
87}
88
89
91(
94)
95{
96 if (!predicted_)
97 {
98 pdmPredictor_->correct(y);
99 predicted_ = true;
100 }
101
103 (
105 (
106 "ny",
107 mesh_.time().timeName(),
108 mesh_,
111 false
112 ),
113 mesh_,
115 patchTypes<vector>(mesh_, patchIDs_)
116 );
117
118 const fvPatchList& patches = mesh_.boundary();
120
121 for (const label patchi : patchIDs_)
122 {
123 nybf[patchi] == -patches[patchi].nf();
124 }
125
126 int iter = 0;
127 scalar initialResidual = 0;
128
129 do
130 {
131 ny = fvc::grad(y);
132 ny /= (mag(ny) + SMALL);
133
135 nf /= (mag(nf) + SMALL);
136
137 surfaceScalarField yPhi("yPhi", mesh_.Sf() & nf);
138
139 fvScalarMatrix yEqn
140 (
141 fvm::div(yPhi, y)
142 - fvm::Sp(fvc::div(yPhi), y)
143 - epsilon_*y*fvm::laplacian(y)
144 ==
145 dimensionedScalar("1", dimless, 1.0)
146 );
147
148 yEqn.relax();
149 initialResidual = yEqn.solve().initialResidual();
150
151 } while (initialResidual > tolerance_ && ++iter < maxIter_);
152
153 // Need to stabilise the y for overset meshes since the holed cells
154 // keep the initial value (0.0) so the gradient of that will be
155 // zero as well. Turbulence models do not like zero wall distance.
156 y.max(SMALL);
157
158 // Only calculate n if the field is defined
159 if (notNull(n))
160 {
161 n = -ny;
162 }
163
164 return true;
165}
166
167
168// ************************************************************************* //
scalar y
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1092
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Specialisation of patchDist for wall distance calculation.
Calculation of approximate distance to nearest patch for all cells and boundary by solving the Eikona...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
thermo correct()
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
Calculate the divergence of the given field.
Calculate the gradient of the given field.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
Calculate the finiteVolume matrix for implicit and explicit sources.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
bool notNull(const T *ptr)
True if ptr is not a pointer (of type T) to the nullObject.
Definition: nullObject.H:207
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
dictionary dict
volScalarField & e
Definition: createFields.H:11