averaging in isoSurfaces

This commit is contained in:
mattijs
2008-11-07 18:14:06 +00:00
parent 856c44c264
commit 29d67c3240
3 changed files with 65 additions and 51 deletions

View File

@ -40,7 +40,7 @@ surfaceFormat vtk;
// 1] vertex values determined from neighbouring cell-centre values
// 2] face values determined using the current face interpolation scheme
// for the field (linear, gamma, etc.)
interpolationScheme cellPointFace;
interpolationScheme cellPoint;
// Fields to sample.
fields
@ -155,21 +155,21 @@ surfaces
// Optional: whether to leave as faces (=default) or triangulate
}
/* not yet (re)implemented --
constantIso
{
name iso;
field rho;
value 0.5;
name isoSurface;
isoField rho;
isoValue 0.5;
interpolate false
}
someIso
interpolatedIso
{
type iso;
field rho;
value 0.5;
name isoSurface;
isoField rho;
isoValue 0.5;
interpolate true;
//regularise false; //optional: do not simplify
}
*/
);

View File

@ -105,6 +105,44 @@ void Foam::sampledIsoSurface::createGeometry() const
volPointInterpolation::New(fvm).interpolate(cellFld)
);
if (average_)
{
//- From point field and interpolated cell.
scalarField cellAvg(fvm.nCells(), scalar(0.0));
labelField nPointCells(fvm.nCells(), 0);
{
for (label pointI = 0; pointI < fvm.nPoints(); pointI++)
{
const labelList& pCells = fvm.pointCells(pointI);
forAll(pCells, i)
{
label cellI = pCells[i];
cellAvg[cellI] += pointFld().internalField()[pointI];
nPointCells[cellI]++;
}
}
}
forAll(cellAvg, cellI)
{
cellAvg[cellI] /= nPointCells[cellI];
}
const isoSurface iso
(
fvm,
cellAvg,
pointFld().internalField(),
isoVal_,
regularise_
);
const_cast<sampledIsoSurface&>(*this).triSurface::operator=(iso);
meshCells_ = iso.meshCells();
}
else
{
//- Direct from cell field and point field. Gives bad continuity.
const isoSurface iso
(
@ -115,45 +153,17 @@ void Foam::sampledIsoSurface::createGeometry() const
regularise_
);
////- From point field and interpolated cell.
//scalarField cellAvg(fvm.nCells(), scalar(0.0));
//labelField nPointCells(fvm.nCells(), 0);
//{
// for (label pointI = 0; pointI < fvm.nPoints(); pointI++)
// {
// const labelList& pCells = fvm.pointCells(pointI);
//
// forAll(pCells, i)
// {
// label cellI = pCells[i];
//
// cellAvg[cellI] += pointFld().internalField()[pointI];
// nPointCells[cellI]++;
// }
// }
//}
//forAll(cellAvg, cellI)
//{
// cellAvg[cellI] /= nPointCells[cellI];
//}
//
//const isoSurface iso
//(
// fvm,
// cellAvg,
// pointFld().internalField(),
// isoVal_,
// regularise_
//);
const_cast<sampledIsoSurface&>(*this).triSurface::operator=(iso);
meshCells_ = iso.meshCells();
}
if (debug)
{
Pout<< "sampledIsoSurface::createGeometry() : constructed iso:"
<< nl
<< " regularise : " << regularise_ << nl
<< " average : " << average_ << nl
<< " isoField : " << isoField_ << nl
<< " isoValue : " << isoVal_ << nl
<< " points : " << points().size() << nl
@ -177,6 +187,7 @@ Foam::sampledIsoSurface::sampledIsoSurface
isoField_(dict.lookup("isoField")),
isoVal_(readScalar(dict.lookup("isoValue"))),
regularise_(dict.lookupOrDefault("regularise", true)),
average_(dict.lookupOrDefault("average", true)),
zoneName_(word::null),
facesPtr_(NULL),
storedTimeIndex_(-1),

View File

@ -66,6 +66,9 @@ class sampledIsoSurface
//- Whether to coarse
const Switch regularise_;
//- Whether to recalculate cell values as average of point values
const Switch average_;
//- zone name (if restricted to zones)
word zoneName_;