advectionDiffusionPatchDistMethod: Calculation of approximate distance to nearest patch for all cells

and boundary by solving the Eikonal equation in advection form with diffusion smoothing.
This commit is contained in:
Henry
2015-01-09 16:12:58 +00:00
parent 113f7a9d7b
commit cb3cd0afb0
4 changed files with 381 additions and 3 deletions

View File

@ -45,6 +45,7 @@ $(wallDist)/wallDist/wallDist.C
$(wallDist)/patchDistMethods/patchDistMethod/patchDistMethod.C
$(wallDist)/patchDistMethods/meshWave/meshWavePatchDistMethod.C
$(wallDist)/patchDistMethods/Poisson/PoissonPatchDistMethod.C
$(wallDist)/patchDistMethods/advectionDiffusion/advectionDiffusionPatchDistMethod.C
fvMeshMapper = fvMesh/fvMeshMapper

View File

@ -30,15 +30,15 @@ Description
References:
\verbatim
D.B.Spalding,
D.B. Spalding,
"Calculation of turbulent heat transfer in cluttered spaces",
Proc. 10th Int. Heat Transfer Conference, Brighton, UK, (1994).
E.Fares and W.Schroder,
E. Fares and W. Schroder,
"Differential Equation for Approximate Wall Distance",
Int.J.Numer.Meth., 39:743-762, (2002).
P. G. Tucker,
P.G. Tucker,
"Differential equation based wall distance computation for DES and
RANS",
J.Comp.Phys., Vol. 190, Issue 1, 1 st September, pp. 229-248 (2003)

View File

@ -0,0 +1,171 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "advectionDiffusionPatchDistMethod.H"
#include "surfaceInterpolate.H"
#include "fvcGrad.H"
#include "fvcDiv.H"
#include "fvmDiv.H"
#include "fvmLaplacian.H"
#include "fvmSup.H"
#include "addToRunTimeSelectionTable.H"
#include "fixedValueFvPatchFields.H"
#include "zeroGradientFvPatchFields.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace patchDistMethods
{
defineTypeNameAndDebug(advectionDiffusion, 0);
addToRunTimeSelectionTable(patchDistMethod, advectionDiffusion, dictionary);
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::patchDistMethods::advectionDiffusion::advectionDiffusion
(
const dictionary& dict,
const fvMesh& mesh,
const labelHashSet& patchIDs
)
:
patchDistMethod(mesh, patchIDs),
coeffs_(dict.subDict(type() + "Coeffs")),
pdmPredictor_
(
patchDistMethod::New
(
coeffs_,
mesh,
patchIDs
)
),
epsilon_(coeffs_.lookupOrDefault<scalar>("epsilon", 0.1)),
maxIter_(coeffs_.lookupOrDefault<int>("maxIter", 10)),
tolerance_(coeffs_.lookupOrDefault<scalar>("tolerance", 1e-3))
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
bool Foam::patchDistMethods::advectionDiffusion::correct(volScalarField& y)
{
return correct(y, const_cast<volVectorField&>(volVectorField::null()));
}
namespace Foam
{
template<class Type>
wordList patchTypes
(
const fvMesh& mesh,
const labelHashSet& patchIDs
)
{
wordList yTypes
(
mesh.boundary().size(),
zeroGradientFvPatchField<Type>::typeName
);
forAllConstIter(labelHashSet, patchIDs, iter)
{
yTypes[iter.key()] = fixedValueFvPatchField<Type>::typeName;
}
return yTypes;
}
}
bool Foam::patchDistMethods::advectionDiffusion::correct
(
volScalarField& y,
volVectorField& n
)
{
pdmPredictor_->correct(y);
int iter = 0;
scalar initialResidual = 0;
do
{
volVectorField ny
(
IOobject
(
"ny",
mesh_.time().timeName(),
mesh_
),
mesh_,
dimensionedVector("nWall", dimless, vector::zero),
patchTypes<vector>(mesh_, patchIDs_)
);
const fvPatchList& patches = mesh_.boundary();
forAllConstIter(labelHashSet, patchIDs_, iter)
{
label patchi = iter.key();
ny.boundaryField()[patchi] == -patches[patchi].nf();
}
ny = fvc::grad(y);
ny /= (mag(ny) + SMALL);
surfaceVectorField nf(fvc::interpolate(ny));
nf /= (mag(nf) + SMALL);
surfaceScalarField yPhi("yPhi", mesh_.Sf() & nf);
fvScalarMatrix yEqn
(
fvm::div(yPhi, y)
- fvm::Sp(fvc::div(yPhi), y)
- epsilon_*y*fvm::laplacian(y)
==
dimensionedScalar("1", dimless, 1.0)
);
yEqn.relax();
initialResidual = yEqn.solve().initialResidual();
// Only calculate n if the field is defined
if (notNull(n))
{
n = -ny;
}
} while (initialResidual > tolerance_ && ++iter < maxIter_);
return true;
}
// ************************************************************************* //

View File

@ -0,0 +1,206 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::patchDistMethods::advectionDiffusion
Description
Calculation of approximate distance to nearest patch for all cells and
boundary by solving the Eikonal equation in advection form with diffusion
smoothing.
References:
\verbatim
P.G. Tucker, C.L. Rumsey, R.E. Bartels, R.T. Biedron,
"Transport equation based wall distance computations aimed at flows
with time-dependent geometry",
NASA/TM-2003-212680, December 2003.
\endverbatim
Example of the wallDist specification in fvSchemes:
\verbatim
laplacianSchemes
{
.
.
laplacian(yPsi) Gauss linear corrected;
laplacian(yWall) Gauss linear corrected;
.
.
}
wallDist
{
method advectionDiffusion;
advectionDiffusionCoeffs
{
method Poisson;
epsilon 0.1;
tolerance 1e-3;
maxIter 10;
}
}
\endverbatim
Also the solver specification for yWall is required in fvSolution, e.g.
for simple cases:
\verbatim
yPsi
{
solver PCG;
preconditioner DIC;
tolerance 1e-4;
relTol 0;
}
yWall
{
solver PBiCG;
preconditioner DILU;
tolerance 1e-4;
relTol 0;
}
or for more complex cases:
yPsi
{
solver GAMG;
smoother GaussSeidel;
cacheAgglomeration true;
nCellsInCoarsestLevel 10;
agglomerator faceAreaPair;
mergeLevels 1;
tolerance 1e-4;
relTol 0;
}
yWall
{
solver GAMG;
smoother symGaussSeidel;
cacheAgglomeration true;
nCellsInCoarsestLevel 10;
agglomerator faceAreaPair;
mergeLevels 1;
tolerance 1e-4;
relTol 0;
}
\endverbatim
SeeAlso
Foam::patchDistMethod::meshWave
Foam::wallDist
SourceFiles
advectionDiffusionPatchDistMethod.C
\*---------------------------------------------------------------------------*/
#ifndef advectionDiffusionPatchDistMethod_H
#define advectionDiffusionPatchDistMethod_H
#include "patchDistMethod.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace patchDistMethods
{
/*---------------------------------------------------------------------------*\
Class advectionDiffusion Declaration
\*---------------------------------------------------------------------------*/
class advectionDiffusion
:
public patchDistMethod
{
// Private Member Data
//- Sub-dictionary of coefficients
dictionary coeffs_;
//- Run-time selected method to predict the distance-to-wall field
autoPtr<patchDistMethod> pdmPredictor_;
//- Diffusion coefficient multiplying y*laplacian(y)
scalar epsilon_;
//- Convergence tolerance for the iterations of the advection-diffusion
// equation to correct the distance-to-patch and normal-to-patch fields
scalar tolerance_;
//- Maximum number of iterations of the advection-diffusion equation
// to correct the distance-to-patch and normal-to-patch fields
int maxIter_;
// Private Member Functions
//- Disallow default bitwise copy construct
advectionDiffusion(const advectionDiffusion&);
//- Disallow default bitwise assignment
void operator=(const advectionDiffusion&);
public:
//- Runtime type information
TypeName("advectionDiffusion");
// Constructors
//- Construct from coefficients dictionary, mesh
// and fixed-value patch set
advectionDiffusion
(
const dictionary& dict,
const fvMesh& mesh,
const labelHashSet& patchIDs
);
// Member Functions
//- Correct the given distance-to-patch field
virtual bool correct(volScalarField& y);
//- Correct the given distance-to-patch and normal-to-patch fields
virtual bool correct(volScalarField& y, volVectorField& n);
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace patchDistMethods
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //