mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
New spatialLinear cellSizeFunction, allowing a linear grade in space in a given
direction. Not a function of surface distance.
This commit is contained in:
@ -15,6 +15,7 @@ $(cellSiseFunctions)/uniform/uniform.C
|
|||||||
$(cellSiseFunctions)/uniformDistance/uniformDistance.C
|
$(cellSiseFunctions)/uniformDistance/uniformDistance.C
|
||||||
$(cellSiseFunctions)/linearDistance/linearDistance.C
|
$(cellSiseFunctions)/linearDistance/linearDistance.C
|
||||||
$(cellSiseFunctions)/surfaceOffsetLinearDistance/surfaceOffsetLinearDistance.C
|
$(cellSiseFunctions)/surfaceOffsetLinearDistance/surfaceOffsetLinearDistance.C
|
||||||
|
$(cellSiseFunctions)/linearSpatial/linearSpatial.C
|
||||||
|
|
||||||
initialPointsMethod/initialPointsMethod/initialPointsMethod.C
|
initialPointsMethod/initialPointsMethod/initialPointsMethod.C
|
||||||
initialPointsMethod/uniformGrid/uniformGrid.C
|
initialPointsMethod/uniformGrid/uniformGrid.C
|
||||||
|
|||||||
@ -0,0 +1,145 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | Copyright (C) 2009-2009 OpenCFD Ltd.
|
||||||
|
\\/ 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 2 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, write to the Free Software Foundation,
|
||||||
|
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
||||||
|
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
#include "linearSpatial.H"
|
||||||
|
#include "addToRunTimeSelectionTable.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
defineTypeNameAndDebug(linearSpatial, 0);
|
||||||
|
addToRunTimeSelectionTable(cellSizeFunction, linearSpatial, dictionary);
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
linearSpatial::linearSpatial
|
||||||
|
(
|
||||||
|
const dictionary& initialPointsDict,
|
||||||
|
const conformalVoronoiMesh& cvMesh,
|
||||||
|
const searchableSurface& surface
|
||||||
|
)
|
||||||
|
:
|
||||||
|
cellSizeFunction(typeName, initialPointsDict, cvMesh, surface),
|
||||||
|
referencePoint_(coeffsDict().lookup("referencePoint")),
|
||||||
|
referenceCellSize_(readScalar(coeffsDict().lookup("referenceCellSize"))),
|
||||||
|
direction_(coeffsDict().lookup("direction")),
|
||||||
|
cellSizeGradient_(readScalar(coeffsDict().lookup("cellSizeGradient")))
|
||||||
|
{
|
||||||
|
direction_ /= mag(direction_);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
scalar linearSpatial::sizeFunction(const point& pt) const
|
||||||
|
{
|
||||||
|
return
|
||||||
|
referenceCellSize_
|
||||||
|
+ ((pt - referencePoint_) & direction_)*cellSizeGradient_;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
bool linearSpatial::cellSize
|
||||||
|
(
|
||||||
|
const point& pt,
|
||||||
|
scalar& size,
|
||||||
|
bool isSurfacePoint
|
||||||
|
) const
|
||||||
|
{
|
||||||
|
if (sideMode_ == BOTHSIDES || isSurfacePoint)
|
||||||
|
{
|
||||||
|
size = sizeFunction(pt);
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
size = 0;
|
||||||
|
|
||||||
|
List<pointIndexHit> hits;
|
||||||
|
|
||||||
|
surface_.findNearest
|
||||||
|
(
|
||||||
|
pointField(1, pt),
|
||||||
|
scalarField(1, sqr(snapToSurfaceTol_)),
|
||||||
|
hits
|
||||||
|
);
|
||||||
|
|
||||||
|
const pointIndexHit& hitInfo = hits[0];
|
||||||
|
|
||||||
|
// If the nearest point is essentially on the surface, do not do a
|
||||||
|
// getVolumeType calculation, as it will be prone to error.
|
||||||
|
if (hitInfo.hit())
|
||||||
|
{
|
||||||
|
size = sizeFunction(pt);
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
pointField ptF(1, pt);
|
||||||
|
List<searchableSurface::volumeType> vTL;
|
||||||
|
|
||||||
|
surface_.getVolumeType(ptF, vTL);
|
||||||
|
|
||||||
|
bool functionApplied = false;
|
||||||
|
|
||||||
|
if
|
||||||
|
(
|
||||||
|
sideMode_ == INSIDE
|
||||||
|
&& vTL[0] == searchableSurface::INSIDE
|
||||||
|
)
|
||||||
|
{
|
||||||
|
size = sizeFunction(pt);
|
||||||
|
|
||||||
|
functionApplied = true;
|
||||||
|
}
|
||||||
|
else if
|
||||||
|
(
|
||||||
|
sideMode_ == OUTSIDE
|
||||||
|
&& vTL[0] == searchableSurface::OUTSIDE
|
||||||
|
)
|
||||||
|
{
|
||||||
|
size = sizeFunction(pt);
|
||||||
|
|
||||||
|
functionApplied = true;
|
||||||
|
}
|
||||||
|
|
||||||
|
return functionApplied;
|
||||||
|
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
} // End namespace Foam
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -0,0 +1,122 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | Copyright (C) 2009-2009 OpenCFD Ltd.
|
||||||
|
\\/ 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 2 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, write to the Free Software Foundation,
|
||||||
|
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
||||||
|
|
||||||
|
Class
|
||||||
|
Foam::linearSpatial
|
||||||
|
|
||||||
|
Description
|
||||||
|
|
||||||
|
SourceFiles
|
||||||
|
linearSpatial.C
|
||||||
|
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
#ifndef linearSpatial_H
|
||||||
|
#define linearSpatial_H
|
||||||
|
|
||||||
|
#include "cellSizeFunction.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
|
||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
Class linearSpatial Declaration
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
class linearSpatial
|
||||||
|
:
|
||||||
|
public cellSizeFunction
|
||||||
|
{
|
||||||
|
|
||||||
|
private:
|
||||||
|
|
||||||
|
// Private data
|
||||||
|
|
||||||
|
//- Reference point for spatial size grading
|
||||||
|
point referencePoint_;
|
||||||
|
|
||||||
|
//- Cell size at reference point
|
||||||
|
scalar referenceCellSize_;
|
||||||
|
|
||||||
|
//- Direction of cell size grading, stored as unit vector, may be
|
||||||
|
// supplied with any magnitude
|
||||||
|
vector direction_;
|
||||||
|
|
||||||
|
//- Gradient of cell size change in direction of direction_
|
||||||
|
scalar cellSizeGradient_;
|
||||||
|
|
||||||
|
|
||||||
|
// Private Member Functions
|
||||||
|
|
||||||
|
//- Calculate the cell size as a function of the given position
|
||||||
|
scalar sizeFunction(const point& pt) const;
|
||||||
|
|
||||||
|
|
||||||
|
public:
|
||||||
|
|
||||||
|
//- Runtime type information
|
||||||
|
TypeName("linearSpatial");
|
||||||
|
|
||||||
|
// Constructors
|
||||||
|
|
||||||
|
//- Construct from components
|
||||||
|
linearSpatial
|
||||||
|
(
|
||||||
|
const dictionary& initialPointsDict,
|
||||||
|
const conformalVoronoiMesh& cvMesh,
|
||||||
|
const searchableSurface& surface
|
||||||
|
);
|
||||||
|
|
||||||
|
|
||||||
|
//- Destructor
|
||||||
|
virtual ~linearSpatial()
|
||||||
|
{}
|
||||||
|
|
||||||
|
|
||||||
|
// Member Functions
|
||||||
|
|
||||||
|
//- Modify scalar argument to the cell size specified by function.
|
||||||
|
// Return a boolean specifying if the function was used, i.e. false if
|
||||||
|
// the point was not in range of the surface for a spatially varying
|
||||||
|
// size.
|
||||||
|
virtual bool cellSize
|
||||||
|
(
|
||||||
|
const point& pt,
|
||||||
|
scalar& size,
|
||||||
|
bool isSurfacePoint = false
|
||||||
|
) const;
|
||||||
|
};
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
} // End namespace Foam
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
#endif
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -60,14 +60,7 @@ bool uniform::cellSize
|
|||||||
bool isSurfacePoint
|
bool isSurfacePoint
|
||||||
) const
|
) const
|
||||||
{
|
{
|
||||||
if (isSurfacePoint)
|
if (sideMode_ == BOTHSIDES || isSurfacePoint)
|
||||||
{
|
|
||||||
size = cellSize_;
|
|
||||||
|
|
||||||
return true;
|
|
||||||
}
|
|
||||||
|
|
||||||
if (sideMode_ == BOTHSIDES)
|
|
||||||
{
|
{
|
||||||
size = cellSize_;
|
size = cellSize_;
|
||||||
|
|
||||||
|
|||||||
Reference in New Issue
Block a user