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 | Copyright (C) 2015 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.subDict(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().timeName(),
102  mesh_,
105  false
106  ),
107  mesh_,
109  patchTypes<vector>(mesh_, patchIDs_)
110  );
111 
112  const fvPatchList& patches = mesh_.boundary();
113 
114  forAllConstIter(labelHashSet, patchIDs_, iter)
115  {
116  label patchi = iter.key();
117  ny.boundaryField()[patchi] == -patches[patchi].nf();
118  }
119 
120  int iter = 0;
121  scalar initialResidual = 0;
122 
123  do
124  {
125  ny = fvc::grad(y);
126  ny /= (mag(ny) + SMALL);
127 
129  nf /= (mag(nf) + SMALL);
130 
131  surfaceScalarField yPhi("yPhi", mesh_.Sf() & nf);
132 
133  fvScalarMatrix yEqn
134  (
135  fvm::div(yPhi, y)
136  - fvm::Sp(fvc::div(yPhi), y)
137  - epsilon_*y*fvm::laplacian(y)
138  ==
139  dimensionedScalar("1", dimless, 1.0)
140  );
141 
142  yEqn.relax();
143  initialResidual = yEqn.solve().initialResidual();
144 
145  } while (initialResidual > tolerance_ && ++iter < maxIter_);
146 
147  // Only calculate n if the field is defined
148  if (notNull(n))
149  {
150  n = -ny;
151  }
152 
153  return true;
154 }
155 
156 
157 // ************************************************************************* //
Foam::patchDistMethods::advectionDiffusion::advectionDiffusion
advectionDiffusion(const advectionDiffusion &)
Disallow default bitwise copy construct.
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::patchDistMethods::advectionDiffusion
Calculation of approximate distance to nearest patch for all cells and boundary by solving the Eikona...
Definition: advectionDiffusionPatchDistMethod.H:180
Foam::notNull
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:46
Foam::fvc::interpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Definition: surfaceInterpolate.C:110
fvcDiv.H
Calculate the divergence of the given field.
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::HashSet< label, Hash< label > >
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::patchDistMethods::defineTypeNameAndDebug
defineTypeNameAndDebug(advectionDiffusion, 0)
fvmDiv.H
Calculate the matrix for the divergence of the given field and flux.
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::patchDistMethods::advectionDiffusion::correct
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
Definition: advectionDiffusionPatchDistMethod.C:78
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::label
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
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:519
fvmLaplacian.H
Calculate the matrix for the laplacian of the field.
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
dict
dictionary dict
Definition: searchingEngine.H:14
advectionDiffusionPatchDistMethod.H
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::patchDistMethod
Specialisation of patchDist for wall distance calculation.
Definition: patchDistMethod.H:56
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
surfaceInterpolate.H
Surface Interpolation.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::patchDistMethods::addToRunTimeSelectionTable
addToRunTimeSelectionTable(patchDistMethod, advectionDiffusion, dictionary)
fvcGrad.H
Calculate the gradient of the given field.
fixedValueFvPatchFields.H
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::GeometricField::null
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
Definition: GeometricFieldI.H:30
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
zeroGradientFvPatchFields.H
y
scalar y
Definition: LISASMDCalcMethod1.H:14