nearWallDist, wallDist: Parallel consistentency and finer control
The near wall distances calculated for use in wall functions and the
corrections applied to near wall cells as part of the meshWave
wall/patch distance method have been made consistent across processor
and cyclic boundaries.
The extent to which these corrections are performed in the meshWave
method is now controllable by an nCorrectors entry. This defaults to 2,
which produces a result rougly equivalent to the previous correction
procedure. A higher level of correction could be specified as follows,
in system/fvSchemes:
wallDist
{
method meshWave;
nCorrectors 3;
}
Corrections replace basic cell-centre-face-centre distances with more
accurate cell-centre-face-polygon calculations in which all the points
and edges of the wall face are taken into account. The number of
correctors represents the number of layers of cells that these
corrections propagate into the mesh from the wall faces in question.
Note that correctors are expensive, and returns diminish as the number
of corrections increase, as the error in the basic calculation reduces
with distance from the wall faces. It is unlikely that more than 2 or
3 correctors would ever be warranted. Indeed, this control is
potentially more useful in reducing the number of corrections to 1 or 0
in order to reduce the computational expense in dynamic mesh cases where
distances are being constantly recomputed.
This commit is contained in:
@ -26,7 +26,6 @@ License
|
||||
#include "inverseDistanceDiffusivity.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
#include "fvPatchDistWave.H"
|
||||
#include "fvWallPoint.H"
|
||||
#include "HashSet.H"
|
||||
#include "surfaceInterpolate.H"
|
||||
#include "zeroGradientFvPatchFields.H"
|
||||
@ -85,12 +84,11 @@ Foam::inverseDistanceDiffusivity::operator()() const
|
||||
|
||||
if (patchNames_.size())
|
||||
{
|
||||
fvPatchDistWave::wave<fvWallPoint>
|
||||
fvPatchDistWave::calculate
|
||||
(
|
||||
mesh(),
|
||||
mesh().boundaryMesh().patchSet(patchNames_),
|
||||
y,
|
||||
false
|
||||
y
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
@ -25,8 +25,7 @@ License
|
||||
|
||||
#include "inverseFaceDistanceDiffusivity.H"
|
||||
#include "surfaceFields.H"
|
||||
#include "FvFaceCellWave.H"
|
||||
#include "fvWallPoint.H"
|
||||
#include "fvPatchDistWave.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
|
||||
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||
@ -68,112 +67,51 @@ Foam::inverseFaceDistanceDiffusivity::~inverseFaceDistanceDiffusivity()
|
||||
Foam::tmp<Foam::surfaceScalarField>
|
||||
Foam::inverseFaceDistanceDiffusivity::operator()() const
|
||||
{
|
||||
surfaceScalarField y
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"y",
|
||||
mesh().time().timeName(),
|
||||
mesh()
|
||||
),
|
||||
mesh(),
|
||||
dimensionedScalar(dimless, 1)
|
||||
);
|
||||
|
||||
const labelHashSet patchIDs(mesh().boundaryMesh().patchSet(patchNames_));
|
||||
|
||||
const surfaceVectorField::Boundary& CfBf = mesh().Cf().boundaryField();
|
||||
|
||||
// Create changed face information
|
||||
List<labelPair> changedPatchAndFaces;
|
||||
List<fvWallPoint> changedFacesInfo;
|
||||
if (patchNames_.size())
|
||||
{
|
||||
label changedFacei = 0;
|
||||
forAllConstIter(labelHashSet, patchIDs, iter)
|
||||
{
|
||||
const label patchi = iter.key();
|
||||
changedFacei += mesh().boundary()[patchi].size();
|
||||
}
|
||||
fvPatchDistWave::calculate
|
||||
(
|
||||
mesh(),
|
||||
patchIDs,
|
||||
y
|
||||
);
|
||||
|
||||
changedPatchAndFaces.resize(changedFacei);
|
||||
changedFacesInfo.resize(changedFacei);
|
||||
|
||||
changedFacei = 0;
|
||||
// Use cell distance on faces that are part of the patch set. This
|
||||
// avoids divide-by-zero issues.
|
||||
forAllConstIter(labelHashSet, patchIDs, iter)
|
||||
{
|
||||
const label patchi = iter.key();
|
||||
|
||||
forAll(mesh().boundary()[patchi], patchFacei)
|
||||
const labelUList& patchCells =
|
||||
mesh().boundary()[patchi].faceCells();
|
||||
|
||||
forAll(patchCells, patchFacei)
|
||||
{
|
||||
changedPatchAndFaces[changedFacei] =
|
||||
labelPair(patchi, patchFacei);
|
||||
changedFacesInfo[changedFacei] =
|
||||
fvWallPoint(CfBf[patchi][patchFacei], 0);
|
||||
|
||||
changedFacei ++;
|
||||
y.boundaryFieldRef()[patchi][patchFacei] =
|
||||
mag
|
||||
(
|
||||
mesh().Cf().boundaryField()[patchi][patchFacei]
|
||||
- mesh().C()[patchCells[patchFacei]]
|
||||
);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
// Initialise wave storage
|
||||
List<fvWallPoint> internalFaceInfo(mesh().nInternalFaces());
|
||||
List<List<fvWallPoint>> patchFaceInfo
|
||||
(
|
||||
FvFaceCellWave<fvWallPoint>::template
|
||||
sizesListList<List<List<fvWallPoint>>>
|
||||
(
|
||||
FvFaceCellWave<fvWallPoint>::template
|
||||
listListSizes(mesh().boundary()),
|
||||
fvWallPoint()
|
||||
)
|
||||
);
|
||||
List<fvWallPoint> cellInfo(mesh().nCells());
|
||||
|
||||
// Wave through the mesh
|
||||
FvFaceCellWave<fvWallPoint> wave
|
||||
(
|
||||
mesh(),
|
||||
changedPatchAndFaces,
|
||||
changedFacesInfo,
|
||||
internalFaceInfo,
|
||||
patchFaceInfo,
|
||||
cellInfo,
|
||||
mesh().globalData().nTotalCells() + 1 // max iterations
|
||||
);
|
||||
|
||||
// Create the diffusivity field
|
||||
tmp<surfaceScalarField> tfaceDiffusivity
|
||||
(
|
||||
new surfaceScalarField
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"faceDiffusivity",
|
||||
mesh().time().timeName(),
|
||||
mesh(),
|
||||
IOobject::NO_READ,
|
||||
IOobject::NO_WRITE
|
||||
),
|
||||
mesh(),
|
||||
dimensionedScalar(dimless, 1.0)
|
||||
)
|
||||
);
|
||||
surfaceScalarField& faceDiffusivity = tfaceDiffusivity.ref();
|
||||
|
||||
// Convert waved distance data into diffusivities
|
||||
forAll(internalFaceInfo, facei)
|
||||
{
|
||||
const scalar dist = internalFaceInfo[facei].dist(wave.data());
|
||||
faceDiffusivity[facei] = 1/dist;
|
||||
}
|
||||
forAll(patchFaceInfo, patchi)
|
||||
{
|
||||
// Use cell distance on faces that are part of the patch set. This
|
||||
// avoids divide-by-zero issues.
|
||||
const bool useCellDist = patchIDs.found(patchi);
|
||||
|
||||
const labelUList& patchCells = mesh().boundary()[patchi].faceCells();
|
||||
|
||||
forAll(patchFaceInfo[patchi], patchFacei)
|
||||
{
|
||||
const scalar dist =
|
||||
useCellDist
|
||||
? cellInfo[patchCells[patchFacei]].dist(wave.data())
|
||||
: patchFaceInfo[patchi][patchFacei].dist(wave.data());
|
||||
|
||||
faceDiffusivity.boundaryFieldRef()[patchi][patchFacei] = 1/dist;
|
||||
}
|
||||
}
|
||||
|
||||
return tfaceDiffusivity;
|
||||
return surfaceScalarField::New("faceDiffusivity", 1/y);
|
||||
}
|
||||
|
||||
|
||||
|
||||
Reference in New Issue
Block a user