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 -------------------------------------------------------------------------------
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 
29 #include "PoissonPatchDistMethod.H"
30 #include "fvcGrad.H"
31 #include "fvmLaplacian.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace patchDistMethods
39 {
40  defineTypeNameAndDebug(Poisson, 0);
42 }
43 }
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
47 Foam::patchDistMethods::Poisson::Poisson
48 (
49  const dictionary& dict,
50  const fvMesh& mesh,
51  const labelHashSet& patchIDs
52 )
53 :
55 {}
56 
57 
58 Foam::patchDistMethods::Poisson::Poisson
59 (
60  const fvMesh& mesh,
61  const labelHashSet& patchIDs
62 )
63 :
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_ = volScalarField::New
85  (
86  "yPsi",
88  mesh_,
90  y.boundaryFieldRef().types()
91  );
92  }
93  auto& yPsi = tyPsi_.ref();
94 
95  solve(fvm::laplacian(yPsi) == dimensionedScalar("1", dimless, -1.0));
96 
97  volVectorField gradyPsi(fvc::grad(yPsi));
98  volScalarField magGradyPsi(mag(gradyPsi));
99 
100  // Need to stabilise the y for overset meshes since the holed cells
101  // keep the initial value (0.0) so the gradient of that will be
102  // zero as well. Turbulence models do not like zero wall distance.
103  y = max
104  (
105  sqrt(magSqr(gradyPsi) + 2*yPsi) - magGradyPsi,
106  dimensionedScalar("smallY", dimLength, SMALL)
107  );
108 
109  // For overset: enforce smooth y field (yPsi is smooth, magGradyPsi is
110  // not)
111  mesh_.interpolate(y);
112 
113  // Need to stabilise the y for overset meshes since the holed cells
114  // keep the initial value (0.0) so the gradient of that will be
115  // zero as well. Turbulence models do not like zero wall distance.
116  y.clamp_min(SMALL);
117 
118  // For overset: enforce smooth y field (yPsi is smooth, magGradyPsi is
119  // not)
120  mesh_.interpolate(y);
121 
122  // Cache yPsi if the mesh is moving otherwise delete
123  if (!mesh_.changing())
124  {
125  tyPsi_.clear();
126  }
127 
128  // Only calculate n if the field is defined
129  if (notNull(n))
130  {
131  n =
132  -gradyPsi
133  /max
134  (
135  magGradyPsi,
136  dimensionedScalar("smallMagGradyPsi", dimLength, SMALL)
137 
138  );
139 
140  // For overset: enforce smooth field
141  mesh_.interpolate(n);
142  }
143 
144  return true;
145 }
146 
147 
148 // ************************************************************************* //
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
const labelList patchIDs(pbm.indices(polyPatchNames, true))
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
addToRunTimeSelectionTable(patchDistMethod, advectionDiffusion, dictionary)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Calculate the matrix for the laplacian of the field.
const dimensionSet dimless
Dimensionless.
Macros for easy insertion into run-time selection tables.
scalar y
dynamicFvMesh & mesh
Calculate the gradient of the given field.
static const GeometricField< vector, fvPatchField, volMesh > & null() noexcept
Return a null GeometricField (reference to a nullObject).
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2406/OpenFOAM-v2406/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
Calculation of approximate distance to nearest patch for all cells and boundary by solving Poisson&#39;s ...
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
bool notNull(const T *ptr) noexcept
True if ptr is not a pointer (of type T) to the nullObject.
Definition: nullObject.H:267
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
defineTypeNameAndDebug(advectionDiffusion, 0)
label n
Specialisation of patchDist for wall distance calculation.
Do not request registration (bool: false)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127