advectionDiffusionPatchDistMethod.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "surfaceInterpolate.H"
28 #include "fvcGrad.H"
29 #include "fvcDiv.H"
30 #include "fvmDiv.H"
31 #include "fvmLaplacian.H"
32 #include "fvmSup.H"
34 
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace patchDistMethods
43 {
46 }
47 }
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const dictionary& dict,
54  const fvMesh& mesh,
55  const labelHashSet& patchIDs
56 )
57 :
58  patchDistMethod(mesh, patchIDs),
59  coeffs_(dict.optionalSubDict(type() + "Coeffs")),
60  pdmPredictor_
61  (
63  (
64  coeffs_,
65  mesh,
66  patchIDs
67  )
68  ),
69  epsilon_(coeffs_.lookupOrDefault<scalar>("epsilon", 0.1)),
70  tolerance_(coeffs_.lookupOrDefault<scalar>("tolerance", 1e-3)),
71  maxIter_(coeffs_.lookupOrDefault<int>("maxIter", 10)),
72  predicted_(false)
73 {}
74 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
79 {
80  return correct(y, const_cast<volVectorField&>(volVectorField::null()));
81 }
82 
83 
85 (
88 )
89 {
90  if (!predicted_)
91  {
92  pdmPredictor_->correct(y);
93  predicted_ = true;
94  }
95 
97  (
98  IOobject
99  (
100  "ny",
101  mesh_.time().name(),
102  mesh_,
105  false
106  ),
107  mesh_,
109  patchTypes<vector>(mesh_, patchIDs_)
110  );
111 
112  const fvPatchList& patches = mesh_.boundary();
114 
115  forAllConstIter(labelHashSet, patchIDs_, iter)
116  {
117  label patchi = iter.key();
118  nybf[patchi] == -patches[patchi].nf();
119  }
120 
121  int iter = 0;
122  scalar initialResidual = 0;
123 
124  do
125  {
126  ny = fvc::grad(y);
127  ny /= (mag(ny) + small);
128 
130  nf /= (mag(nf) + small);
131 
132  surfaceScalarField yPhi("yPhi", mesh_.Sf() & nf);
133 
134  fvScalarMatrix yEqn
135  (
136  fvm::div(yPhi, y)
137  - fvm::Sp(fvc::div(yPhi), y)
138  - epsilon_*y*fvm::laplacian(y)
139  ==
141  );
142 
143  yEqn.relax();
144  initialResidual = yEqn.solve().initialResidual();
145 
146  } while (initialResidual > tolerance_ && ++iter < maxIter_);
147 
148  // Only calculate n if the field is defined
149  if (notNull(n))
150  {
151  n = -ny;
152  }
153 
154  return true;
155 }
156 
157 
158 // ************************************************************************* //
scalar y
label n
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
Generic GeometricBoundaryField class.
Generic GeometricField class.
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:604
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Specialisation of patchDist for wall distance calculation.
Calculation of approximate distance to nearest patch for all cells and boundary by solving the Eikona...
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
advectionDiffusion(const dictionary &dict, const fvMesh &mesh, const labelHashSet &patchIDs)
Construct from coefficients dictionary, 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 matrix for implicit and explicit sources.
label patchi
const fvPatchList & patches
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar e
Elementary charge.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
defineTypeNameAndDebug(advectionDiffusion, 0)
addToRunTimeSelectionTable(patchDistMethod, advectionDiffusion, dictionary)
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static const zero Zero
Definition: zero.H:97
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimless
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:58
dimensioned< scalar > mag(const dimensioned< Type > &)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict