mirror of
https://github.com/OpenFOAM/OpenFOAM-6.git
synced 2025-12-08 06:57:46 +00:00
wallDist/patchDistMethods/Poisson: New method for fast calculation of an approximate wall-distance field
by solving Poisson's equation
This commit is contained in:
@ -44,6 +44,7 @@ $(wallDist)/nearWallDist/nearWallDist.C
|
|||||||
$(wallDist)/wallDist/wallDist.C
|
$(wallDist)/wallDist/wallDist.C
|
||||||
$(wallDist)/patchDistMethods/patchDistMethod/patchDistMethod.C
|
$(wallDist)/patchDistMethods/patchDistMethod/patchDistMethod.C
|
||||||
$(wallDist)/patchDistMethods/meshWave/meshWavePatchDistMethod.C
|
$(wallDist)/patchDistMethods/meshWave/meshWavePatchDistMethod.C
|
||||||
|
$(wallDist)/patchDistMethods/Poisson/PoissonPatchDistMethod.C
|
||||||
|
|
||||||
|
|
||||||
fvMeshMapper = fvMesh/fvMeshMapper
|
fvMeshMapper = fvMesh/fvMeshMapper
|
||||||
|
|||||||
@ -0,0 +1,129 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / 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 "PoissonPatchDistMethod.H"
|
||||||
|
#include "fvcGrad.H"
|
||||||
|
#include "fvmLaplacian.H"
|
||||||
|
#include "addToRunTimeSelectionTable.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
namespace patchDistMethods
|
||||||
|
{
|
||||||
|
defineTypeNameAndDebug(Poisson, 0);
|
||||||
|
addToRunTimeSelectionTable(patchDistMethod, Poisson, dictionary);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
Foam::patchDistMethods::Poisson::Poisson
|
||||||
|
(
|
||||||
|
const dictionary& dict,
|
||||||
|
const fvMesh& mesh,
|
||||||
|
const labelHashSet& patchIDs
|
||||||
|
)
|
||||||
|
:
|
||||||
|
patchDistMethod(mesh, patchIDs)
|
||||||
|
{}
|
||||||
|
|
||||||
|
|
||||||
|
Foam::patchDistMethods::Poisson::Poisson
|
||||||
|
(
|
||||||
|
const fvMesh& mesh,
|
||||||
|
const labelHashSet& patchIDs
|
||||||
|
)
|
||||||
|
:
|
||||||
|
patchDistMethod(mesh, patchIDs)
|
||||||
|
{}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
bool Foam::patchDistMethods::Poisson::correct(volScalarField& y)
|
||||||
|
{
|
||||||
|
return correct(y, const_cast<volVectorField&>(volVectorField::null()));
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
bool Foam::patchDistMethods::Poisson::correct
|
||||||
|
(
|
||||||
|
volScalarField& y,
|
||||||
|
volVectorField& n
|
||||||
|
)
|
||||||
|
{
|
||||||
|
if (!tyPsi_.valid())
|
||||||
|
{
|
||||||
|
tyPsi_ = tmp<volScalarField>
|
||||||
|
(
|
||||||
|
new volScalarField
|
||||||
|
(
|
||||||
|
IOobject
|
||||||
|
(
|
||||||
|
"yPsi",
|
||||||
|
mesh_.time().timeName(),
|
||||||
|
mesh_
|
||||||
|
),
|
||||||
|
mesh_,
|
||||||
|
dimensionedScalar("yPsi", sqr(dimLength), 0.0),
|
||||||
|
y.boundaryField().types()
|
||||||
|
)
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
volScalarField& yPsi = tyPsi_();
|
||||||
|
|
||||||
|
solve(fvm::laplacian(yPsi) == dimensionedScalar("1", dimless, -1.0));
|
||||||
|
|
||||||
|
volVectorField gradyPsi(fvc::grad(yPsi));
|
||||||
|
volScalarField magGradyPsi(mag(gradyPsi));
|
||||||
|
|
||||||
|
y = sqrt(magSqr(gradyPsi) + 2*yPsi) - magGradyPsi;
|
||||||
|
|
||||||
|
// Cache yPsi if the mesh is moving otherwise delete
|
||||||
|
if (!mesh_.changing())
|
||||||
|
{
|
||||||
|
tyPsi_.clear();
|
||||||
|
}
|
||||||
|
|
||||||
|
// Only calculate n if the field is defined
|
||||||
|
if (notNull(n))
|
||||||
|
{
|
||||||
|
n =
|
||||||
|
-gradyPsi
|
||||||
|
/max
|
||||||
|
(
|
||||||
|
magGradyPsi,
|
||||||
|
dimensionedScalar("smallMagGradyPsi", dimLength, SMALL)
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -0,0 +1,177 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / 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::Poisson
|
||||||
|
|
||||||
|
Description
|
||||||
|
Calculation of approximate distance to nearest patch for all cells and
|
||||||
|
boundary using a Poisson equation.
|
||||||
|
|
||||||
|
References:
|
||||||
|
\verbatim
|
||||||
|
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,
|
||||||
|
"Differential Equation for Approximate Wall Distance",
|
||||||
|
Int.J.Numer.Meth., 39:743-762, (2002).
|
||||||
|
|
||||||
|
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)
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
Example of the wallDist specification in fvSchemes:
|
||||||
|
\verbatim
|
||||||
|
laplacianSchemes
|
||||||
|
{
|
||||||
|
.
|
||||||
|
.
|
||||||
|
laplacian(yPsi) Gauss linear corrected;
|
||||||
|
.
|
||||||
|
.
|
||||||
|
}
|
||||||
|
|
||||||
|
wallDist
|
||||||
|
{
|
||||||
|
method Poisson;
|
||||||
|
}
|
||||||
|
\endverbatim
|
||||||
|
Also the solver specification for yPsi is required in fvSolution, e.g.
|
||||||
|
for simple cases:
|
||||||
|
\verbatim
|
||||||
|
yPsi
|
||||||
|
{
|
||||||
|
solver PCG;
|
||||||
|
preconditioner DIC;
|
||||||
|
tolerance 1e-5;
|
||||||
|
relTol 0;
|
||||||
|
}
|
||||||
|
|
||||||
|
or for more complex cases:
|
||||||
|
|
||||||
|
yPsi
|
||||||
|
{
|
||||||
|
solver GAMG;
|
||||||
|
smoother GaussSeidel;
|
||||||
|
cacheAgglomeration true;
|
||||||
|
nCellsInCoarsestLevel 10;
|
||||||
|
agglomerator faceAreaPair;
|
||||||
|
mergeLevels 1;
|
||||||
|
tolerance 1e-5;
|
||||||
|
relTol 0;
|
||||||
|
}
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
SeeAlso
|
||||||
|
Foam::patchDistMethod::meshWave
|
||||||
|
Foam::wallDist
|
||||||
|
|
||||||
|
SourceFiles
|
||||||
|
PoissonPatchDistMethod.C
|
||||||
|
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
#ifndef PoissonPatchDistMethod_H
|
||||||
|
#define PoissonPatchDistMethod_H
|
||||||
|
|
||||||
|
#include "patchDistMethod.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
namespace patchDistMethods
|
||||||
|
{
|
||||||
|
|
||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
Class Poisson Declaration
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
class Poisson
|
||||||
|
:
|
||||||
|
public patchDistMethod
|
||||||
|
{
|
||||||
|
// Private Member Data
|
||||||
|
|
||||||
|
//- Cache yPsi for moving meshes
|
||||||
|
tmp<volScalarField> tyPsi_;
|
||||||
|
|
||||||
|
|
||||||
|
// Private Member Functions
|
||||||
|
|
||||||
|
//- Disallow default bitwise copy construct
|
||||||
|
Poisson(const Poisson&);
|
||||||
|
|
||||||
|
//- Disallow default bitwise assignment
|
||||||
|
void operatorconst Poisson&);
|
||||||
|
|
||||||
|
|
||||||
|
public:
|
||||||
|
|
||||||
|
//- Runtime type information
|
||||||
|
TypeName("Poisson");
|
||||||
|
|
||||||
|
|
||||||
|
// Constructors
|
||||||
|
|
||||||
|
//- Construct from coefficients dictionary, mesh
|
||||||
|
// and fixed-value patch set
|
||||||
|
Poisson
|
||||||
|
(
|
||||||
|
const dictionary& dict,
|
||||||
|
const fvMesh& mesh,
|
||||||
|
const labelHashSet& patchIDs
|
||||||
|
);
|
||||||
|
|
||||||
|
//- Construct from mesh and fixed-value patch set
|
||||||
|
Poisson
|
||||||
|
(
|
||||||
|
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
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -75,6 +75,8 @@ Foam::patchDistMethods::meshWave::meshWave
|
|||||||
|
|
||||||
bool Foam::patchDistMethods::meshWave::correct(volScalarField& y)
|
bool Foam::patchDistMethods::meshWave::correct(volScalarField& y)
|
||||||
{
|
{
|
||||||
|
y = dimensionedScalar("yWall", dimLength, GREAT);
|
||||||
|
|
||||||
// Calculate distance starting from patch faces
|
// Calculate distance starting from patch faces
|
||||||
patchWave wave(mesh_, patchIDs_, correctWalls_);
|
patchWave wave(mesh_, patchIDs_, correctWalls_);
|
||||||
|
|
||||||
@ -105,6 +107,8 @@ bool Foam::patchDistMethods::meshWave::correct
|
|||||||
volVectorField& n
|
volVectorField& n
|
||||||
)
|
)
|
||||||
{
|
{
|
||||||
|
y = dimensionedScalar("yWall", dimLength, GREAT);
|
||||||
|
|
||||||
// Collect pointers to data on patches
|
// Collect pointers to data on patches
|
||||||
UPtrList<vectorField> patchData(mesh_.boundaryMesh().size());
|
UPtrList<vectorField> patchData(mesh_.boundaryMesh().size());
|
||||||
|
|
||||||
|
|||||||
@ -22,21 +22,29 @@ License
|
|||||||
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
|
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
|
||||||
|
|
||||||
Class
|
Class
|
||||||
Foam::patchDist
|
Foam::patchDistMethods::meshWave
|
||||||
|
|
||||||
Description
|
Description
|
||||||
Calculation of distance to nearest patch for all cells and boundary
|
Fast topological mesh-wave method for calculating the distance to nearest
|
||||||
using meshWave.
|
patch for all cells and boundary.
|
||||||
|
|
||||||
Distance correction (correctWalls = true):
|
For regular/un-distorted meshes this method is accurate but for skewed,
|
||||||
For each cell with face on wall calculate the true nearest point (by
|
non-orthogonal meshes it is approximate with the error increasing with the
|
||||||
triangle decomposition) on that face and do the same for that face's
|
degree of mesh distortion. The distance from the near-wall cells to the
|
||||||
pointNeighbours. This will find the true nearest distance in almost all
|
boundary may optionally be corrected for mesh distortion by setting
|
||||||
cases. Only very skewed cells or cells close to another wall might be
|
correctWalls = true.
|
||||||
missed.
|
|
||||||
|
|
||||||
For each cell with only one point on wall the same is done except now it
|
Example of the wallDist specification in fvSchemes:
|
||||||
takes the pointFaces() of the wall point to look for the nearest point.
|
\verbatim
|
||||||
|
wallDist
|
||||||
|
{
|
||||||
|
method meshWave;
|
||||||
|
}
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
SeeAlso
|
||||||
|
Foam::patchDistMethod::Poisson
|
||||||
|
Foam::wallDist
|
||||||
|
|
||||||
SourceFiles
|
SourceFiles
|
||||||
meshWavePatchDistMethod.C
|
meshWavePatchDistMethod.C
|
||||||
@ -89,8 +97,8 @@ public:
|
|||||||
|
|
||||||
// Constructors
|
// Constructors
|
||||||
|
|
||||||
//- Construct from mesh and flag whether or not to correct wall.
|
//- Construct from coefficients dictionary, mesh
|
||||||
// Calculate for all cells.
|
// and fixed-value patch set
|
||||||
meshWave
|
meshWave
|
||||||
(
|
(
|
||||||
const dictionary& dict,
|
const dictionary& dict,
|
||||||
@ -98,7 +106,8 @@ public:
|
|||||||
const labelHashSet& patchIDs
|
const labelHashSet& patchIDs
|
||||||
);
|
);
|
||||||
|
|
||||||
//- Construct from mesh and flag whether or not to correct wall.
|
//- Construct from mesh, fixed-value patch set and flag specifying
|
||||||
|
// whether or not to correct wall.
|
||||||
// Calculate for all cells. correctWalls : correct wall (face&point)
|
// Calculate for all cells. correctWalls : correct wall (face&point)
|
||||||
// cells for correct distance, searching neighbours.
|
// cells for correct distance, searching neighbours.
|
||||||
meshWave
|
meshWave
|
||||||
|
|||||||
@ -145,7 +145,7 @@ public:
|
|||||||
//- Correct the given distance-to-patch field
|
//- Correct the given distance-to-patch field
|
||||||
virtual bool correct(volScalarField& y) = 0;
|
virtual bool correct(volScalarField& y) = 0;
|
||||||
|
|
||||||
//- Correct the given distance-to-patch and reflection vector fields
|
//- Correct the given distance-to-patch and normal-to-patch fields
|
||||||
virtual bool correct(volScalarField& y, volVectorField& n) = 0;
|
virtual bool correct(volScalarField& y, volVectorField& n) = 0;
|
||||||
};
|
};
|
||||||
|
|
||||||
|
|||||||
@ -79,7 +79,7 @@ Foam::wallDist::wallDist(const fvMesh& mesh)
|
|||||||
mesh
|
mesh
|
||||||
),
|
),
|
||||||
mesh,
|
mesh,
|
||||||
dimensionedScalar("yWall", dimLength, GREAT),
|
dimensionedScalar("yWall", dimLength, SMALL),
|
||||||
patchTypes<scalar>(pdm_->patchIDs())
|
patchTypes<scalar>(pdm_->patchIDs())
|
||||||
),
|
),
|
||||||
n_(NULL)
|
n_(NULL)
|
||||||
|
|||||||
Reference in New Issue
Block a user