PoissonPatchDistMethod.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-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2017 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 "fvcGrad.H"
31#include "fvmLaplacian.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace patchDistMethods
39{
42}
43}
44
45// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46
48(
49 const dictionary& dict,
50 const fvMesh& mesh,
51 const labelHashSet& patchIDs
52)
53:
54 patchDistMethod(mesh, patchIDs)
55{}
56
57
59(
60 const fvMesh& mesh,
61 const labelHashSet& patchIDs
62)
63:
64 patchDistMethod(mesh, patchIDs)
65{}
66
67
68// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69
71{
72 return correct(y, const_cast<volVectorField&>(volVectorField::null()));
73}
74
75
77(
80)
81{
82 if (!tyPsi_)
83 {
84 tyPsi_ = tmp<volScalarField>
85 (
87 (
89 (
90 "yPsi",
91 mesh_.time().timeName(),
92 mesh_
93 ),
94 mesh_,
96 y.boundaryFieldRef().types()
97 )
98 );
99 }
100 volScalarField& yPsi = tyPsi_.ref();
101
102 solve(fvm::laplacian(yPsi) == dimensionedScalar("1", dimless, -1.0));
103
104 volVectorField gradyPsi(fvc::grad(yPsi));
105 volScalarField magGradyPsi(mag(gradyPsi));
106
107 // Need to stabilise the y for overset meshes since the holed cells
108 // keep the initial value (0.0) so the gradient of that will be
109 // zero as well. Turbulence models do not like zero wall distance.
110 y = max
111 (
112 sqrt(magSqr(gradyPsi) + 2*yPsi) - magGradyPsi,
113 dimensionedScalar("smallY", dimLength, SMALL)
114 );
115
116 // For overset: enforce smooth y field (yPsi is smooth, magGradyPsi is
117 // not)
118 mesh_.interpolate(y);
119
120 // Need to stabilise the y for overset meshes since the holed cells
121 // keep the initial value (0.0) so the gradient of that will be
122 // zero as well. Turbulence models do not like zero wall distance.
123 y.max(SMALL);
124
125 // For overset: enforce smooth y field (yPsi is smooth, magGradyPsi is
126 // not)
127 mesh_.interpolate(y);
128
129 // Cache yPsi if the mesh is moving otherwise delete
130 if (!mesh_.changing())
131 {
132 tyPsi_.clear();
133 }
134
135 // Only calculate n if the field is defined
136 if (notNull(n))
137 {
138 n =
139 -gradyPsi
140 /max
141 (
142 magGradyPsi,
143 dimensionedScalar("smallMagGradyPsi", dimLength, SMALL)
144
145 );
146
147 // For overset: enforce smooth field
148 mesh_.interpolate(n);
149 }
150
151 return true;
152}
153
154
155// ************************************************************************* //
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.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
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
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 Poisson's ...
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
thermo correct()
dynamicFvMesh & mesh
Calculate the gradient of the given field.
Calculate the matrix for the laplacian of the field.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar sqrt(const dimensionedScalar &ds)
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
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dictionary dict
CEqn solve()