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  Copyright (C) 2020-2023 PCOpt/NTUA
11  Copyright (C) 2020-2023 FOSS GP
12 -------------------------------------------------------------------------------
13 License
14  This file is part of OpenFOAM.
15 
16  OpenFOAM is free software: you can redistribute it and/or modify it
17  under the terms of the GNU General Public License as published by
18  the Free Software Foundation, either version 3 of the License, or
19  (at your option) any later version.
20 
21  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24  for more details.
25 
26  You should have received a copy of the GNU General Public License
27  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28 
29 \*---------------------------------------------------------------------------*/
30 
32 #include "surfaceInterpolate.H"
33 #include "fvcGrad.H"
34 #include "fvcDiv.H"
35 #include "fvmDiv.H"
36 #include "fvmLaplacian.H"
37 #include "fvmSup.H"
38 #include "fvOptions.H"
40 
43 
44 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 namespace patchDistMethods
49 {
52 }
53 }
54 
55 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56 
57 Foam::patchDistMethods::advectionDiffusion::advectionDiffusion
58 (
59  const dictionary& dict,
60  const fvMesh& mesh,
61  const labelHashSet& patchIDs
62 )
63 :
65  //- We do not want to recurse into 'advectionDiffusion' again so
66  // make sure we pick up 'method' always from the subdictionary.
67  //coeffs_(dict.optionalSubDict(type() + "Coeffs")),
68  coeffs_(dict.subDict(type() + "Coeffs")),
69  pdmPredictor_
70  (
72  (
73  coeffs_,
74  mesh,
75  patchIDs
76  )
77  ),
78  epsilon_(coeffs_.getOrDefault<scalar>("epsilon", 0.1)),
79  tolerance_(coeffs_.getOrDefault<scalar>("tolerance", 1e-3)),
80  maxIter_(coeffs_.getOrDefault<int>("maxIter", 10)),
81  predicted_(false),
82  checkAndWriteMesh_
83  (coeffs_.lookupOrDefault<bool>("checkAndWriteMesh", true))
84 {}
85 
86 
87 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88 
90 {
91  return correct(y, const_cast<volVectorField&>(volVectorField::null()));
92 }
93 
94 
96 (
99 )
100 {
101  if (!predicted_)
102  {
103  pdmPredictor_->correct(y);
104  predicted_ = true;
105  }
106 
107  // If the mesh has become invalid, trying to solve the eikonal equation
108  // might lead to crashing and we wont get a chance to see the failed mesh.
109  // Write the mesh points here
110  if (checkAndWriteMesh_)
111  {
112  DebugInfo
113  << "Checking mesh in advectionDiffusion " << endl;
114  if (mesh_.checkMesh(true))
115  {
117  (
118  IOobject
119  (
120  "points",
121  mesh_.pointsInstance(),
122  mesh_.meshSubDir,
123  mesh_,
126  false
127  ),
128  mesh_.points()
129  );
130  points.write();
131  }
132  }
133 
134  volVectorField ny
135  (
136  IOobject
137  (
138  "ny",
139  mesh_.time().timeName(),
140  mesh_,
144  ),
145  mesh_,
147  patchTypes<vector>(mesh_, patchIDs_)
148  );
149 
150  const fvPatchList& patches = mesh_.boundary();
151  volVectorField::Boundary& nybf = ny.boundaryFieldRef();
152 
153  for (const label patchi : patchIDs_)
154  {
155  nybf[patchi] == -patches[patchi].nf();
156  }
157 
158  // Scaling dimensions, to be used with fvOptions
159  volScalarField scaleDims
160  (
161  IOobject
162  (
163  "scaleDims",
164  mesh_.time().timeName(),
165  mesh_,
169  ),
170  mesh_,
171  dimensionedScalar("scaleDims", dimTime/dimLength, scalar(1))
172  );
173 
174  fv::options& fvOptions(fv::options::New(this->mesh_));
175  int iter = 0;
176  scalar initialResidual = 0;
177 
178  do
179  {
180  ny = fvc::grad(y);
181  ny /= (mag(ny) + SMALL);
182 
184  nf /= (mag(nf) + SMALL);
185 
186  surfaceScalarField yPhi("yPhi", mesh_.Sf() & nf);
187 
188  fvScalarMatrix yEqn
189  (
190  fvm::div(yPhi, y)
191  + fvm::SuSp(-fvc::div(yPhi), y)
192  - epsilon_*y*fvm::laplacian(y)
193  ==
194  dimensionedScalar("1", dimless, 1.0)
195  + fvOptions(scaleDims, y)
196  );
197 
198  yEqn.relax();
199  fvOptions.constrain(yEqn);
200  initialResidual = yEqn.solve().initialResidual();
201  fvOptions.correct(y);
202 
203  } while (initialResidual > tolerance_ && ++iter < maxIter_);
204 
205  // Need to stabilise the y for overset meshes since the holed cells
206  // keep the initial value (0.0) so the gradient of that will be
207  // zero as well. Turbulence models do not like zero wall distance.
208  y.clamp_min(SMALL);
209 
210  // Only calculate n if the field is defined
211  if (notNull(n))
212  {
213  n = -ny;
214  }
215 
216  return true;
217 }
218 
219 
220 // ************************************************************************* //
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
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
type
Types of root.
Definition: Roots.H:52
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)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Calculate the matrix for the laplacian of the field.
GeometricBoundaryField< vector, fvPatchField, volMesh > Boundary
Type of boundary fields.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
Macros for easy insertion into run-time selection tables.
fv::options & fvOptions
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
scalar y
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
const pointField & points
Calculate the gradient of the given field.
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
#define DebugInfo
Report an information message using Foam::Info.
Calculate the divergence of the given field.
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
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/v2312/OpenFOAM-v2312/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
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.
static const GeometricField< vector, fvPatchField, volMesh > & null()
Return a null geometric field.
Calculate the matrix for the divergence of the given field and flux.
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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
defineTypeNameAndDebug(advectionDiffusion, 0)
const polyBoundaryMesh & patches
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
label n
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Specialisation of patchDist for wall distance calculation.
A primitive field of type <T> with automated input and output.
Do not request registration (bool: false)
bool notNull(const T *ptr)
True if ptr is not a pointer (of type T) to the nullObject.
Definition: nullObject.H:246
Calculation of approximate distance to nearest patch for all cells and boundary by solving the Eikona...
Calculate the finiteVolume matrix for implicit and explicit sources.
Namespace for OpenFOAM.
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
Definition: fvPatch.H:59
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127