ENH: finiteArea - faMesh now derived from faSolution, faSchemes and data classes

This commit is contained in:
Andrew Heather
2017-12-18 10:50:37 +00:00
parent f385e3b984
commit 2728a96b9c
24 changed files with 79 additions and 122 deletions

View File

@ -155,4 +155,4 @@
dimensionedScalar("one", dimless, 0.01) dimensionedScalar("one", dimless, 0.01)
); );
aMesh.schemesDict().setFluxRequired("h"); aMesh.setFluxRequired("h");

View File

@ -532,9 +532,9 @@ void Foam::faMatrix<Type>::relax(const scalar alpha)
template<class Type> template<class Type>
void Foam::faMatrix<Type>::relax() void Foam::faMatrix<Type>::relax()
{ {
if (psi_.mesh().solutionDict().relaxEquation(psi_.name())) if (psi_.mesh().relaxEquation(psi_.name()))
{ {
relax(psi_.mesh().solutionDict().equationRelaxationFactor(psi_.name())); relax(psi_.mesh().equationRelaxationFactor(psi_.name()));
} }
else else
{ {
@ -628,7 +628,7 @@ template<class Type>
Foam::tmp<Foam::GeometricField<Type, Foam::faePatchField, Foam::edgeMesh>> Foam::tmp<Foam::GeometricField<Type, Foam::faePatchField, Foam::edgeMesh>>
Foam::faMatrix<Type>::flux() const Foam::faMatrix<Type>::flux() const
{ {
if (!psi_.mesh().schemesDict().fluxRequired(psi_.name())) if (!psi_.mesh().fluxRequired(psi_.name()))
{ {
FatalErrorInFunction FatalErrorInFunction
<< "flux requested but " << psi_.name() << "flux requested but " << psi_.name()

View File

@ -143,6 +143,8 @@ Foam::SolverPerformance<Type> Foam::faMatrix<Type>::solve
psi.correctBoundaryConditions(); psi.correctBoundaryConditions();
psi.mesh().setSolverPerformance(psi.name(), solverPerfVec);
return solverPerfVec; return solverPerfVec;
} }
@ -150,26 +152,14 @@ Foam::SolverPerformance<Type> Foam::faMatrix<Type>::solve
template<class Type> template<class Type>
Foam::SolverPerformance<Type> Foam::faMatrix<Type>::faSolver::solve() Foam::SolverPerformance<Type> Foam::faMatrix<Type>::faSolver::solve()
{ {
return solvei return solve(faMat_.psi().mesh().solverDict(faMat_.psi().name()));
(
faMat_.psi().mesh().solutionDict().solverDict
(
faMat_.psi().name()
)
);
} }
template<class Type> template<class Type>
Foam::SolverPerformance<Type> Foam::faMatrix<Type>::solve() Foam::SolverPerformance<Type> Foam::faMatrix<Type>::solve()
{ {
return solve return solve(this->psi().mesh().solverDict(this->psi().name()));
(
this->psi().mesh().solutionDict().solverDict
(
this->psi().name()
)
);
} }

View File

@ -89,6 +89,8 @@ Foam::solverPerformance Foam::faMatrix<Foam::scalar>::solve
psi.correctBoundaryConditions(); psi.correctBoundaryConditions();
psi.mesh().setSolverPerformance(psi.name(), solverPerf);
return solverPerf; return solverPerf;
} }

View File

@ -55,8 +55,7 @@ const int Foam::faMesh::quadricsFit_ = 0;
void Foam::faMesh::setPrimitiveMeshData() void Foam::faMesh::setPrimitiveMeshData()
{ {
DebugInFunction DebugInFunction << "Setting primitive data" << endl;
<< "Setting primitive data" << endl;
const indirectPrimitivePatch& bp = patch(); const indirectPrimitivePatch& bp = patch();
@ -121,11 +120,7 @@ void Foam::faMesh::setPrimitiveMeshData()
void Foam::faMesh::clearGeomNotAreas() const void Foam::faMesh::clearGeomNotAreas() const
{ {
if (debug) DebugInFunction << "Clearing geometry" << endl;
{
Info<< "void faMesh::clearGeomNotAreas() const : "
<< "Clearing geometry" << endl;
}
deleteDemandDrivenData(SPtr_); deleteDemandDrivenData(SPtr_);
deleteDemandDrivenData(patchPtr_); deleteDemandDrivenData(patchPtr_);
@ -144,11 +139,7 @@ void Foam::faMesh::clearGeomNotAreas() const
void Foam::faMesh::clearGeom() const void Foam::faMesh::clearGeom() const
{ {
if (debug) DebugInFunction << "Clearing geometry" << endl;
{
Info<< "void faMesh::clearGeom() const : "
<< "Clearing geometry" << endl;
}
clearGeomNotAreas(); clearGeomNotAreas();
deleteDemandDrivenData(S0Ptr_); deleteDemandDrivenData(S0Ptr_);
@ -159,11 +150,7 @@ void Foam::faMesh::clearGeom() const
void Foam::faMesh::clearAddressing() const void Foam::faMesh::clearAddressing() const
{ {
if (debug) DebugInFunction << "Clearing addressing" << endl;
{
Info<< "void faMesh::clearAddressing() const : "
<< "Clearing addressing" << endl;
}
deleteDemandDrivenData(lduPtr_); deleteDemandDrivenData(lduPtr_);
} }
@ -184,6 +171,9 @@ Foam::faMesh::faMesh(const polyMesh& pMesh)
GeoMesh<polyMesh>(pMesh), GeoMesh<polyMesh>(pMesh),
MeshObject<polyMesh, Foam::UpdateableMeshObject, faMesh>(pMesh), MeshObject<polyMesh, Foam::UpdateableMeshObject, faMesh>(pMesh),
edgeInterpolation(*this), edgeInterpolation(*this),
faSchemes(mesh()),
faSolution(mesh()),
data(mesh()),
faceLabels_ faceLabels_
( (
IOobject IOobject
@ -229,11 +219,7 @@ Foam::faMesh::faMesh(const polyMesh& pMesh)
correctPatchPointNormalsPtr_(nullptr), correctPatchPointNormalsPtr_(nullptr),
globalMeshDataPtr_(nullptr) globalMeshDataPtr_(nullptr)
{ {
if (debug) DebugInFunction << "Creating faMesh from IOobject" << endl;
{
Info<< "faMesh::faMesh(...) : "
<< "Creating faMesh from IOobject" << endl;
}
setPrimitiveMeshData(); setPrimitiveMeshData();
@ -277,6 +263,9 @@ Foam::faMesh::faMesh
GeoMesh<polyMesh>(pMesh), GeoMesh<polyMesh>(pMesh),
MeshObject<polyMesh, Foam::UpdateableMeshObject, faMesh>(pMesh), MeshObject<polyMesh, Foam::UpdateableMeshObject, faMesh>(pMesh),
edgeInterpolation(*this), edgeInterpolation(*this),
faSchemes(mesh()),
faSolution(mesh()),
data(mesh()),
faceLabels_ faceLabels_
( (
IOobject IOobject
@ -324,11 +313,7 @@ Foam::faMesh::faMesh
correctPatchPointNormalsPtr_(nullptr), correctPatchPointNormalsPtr_(nullptr),
globalMeshDataPtr_(nullptr) globalMeshDataPtr_(nullptr)
{ {
if (debug) DebugInFunction << "Creating faMesh from components" << endl;
{
Info<< "faMesh::faMesh(...) : "
<< "Creating faMesh from components" << endl;
}
} }
@ -341,6 +326,9 @@ Foam::faMesh::faMesh
GeoMesh<polyMesh>(pMesh), GeoMesh<polyMesh>(pMesh),
MeshObject<polyMesh, Foam::UpdateableMeshObject, faMesh>(pMesh), MeshObject<polyMesh, Foam::UpdateableMeshObject, faMesh>(pMesh),
edgeInterpolation(*this), edgeInterpolation(*this),
faSchemes(mesh()),
faSolution(mesh()),
data(mesh()),
faceLabels_ faceLabels_
( (
IOobject IOobject
@ -388,8 +376,7 @@ Foam::faMesh::faMesh
correctPatchPointNormalsPtr_(nullptr), correctPatchPointNormalsPtr_(nullptr),
globalMeshDataPtr_(nullptr) globalMeshDataPtr_(nullptr)
{ {
DebugInFunction DebugInFunction << "Creating faMesh from definition file" << endl;
<< "Creating faMesh from definition file" << endl;
// Reading faMeshDefinition dictionary // Reading faMeshDefinition dictionary
IOdictionary faMeshDefinition IOdictionary faMeshDefinition
@ -791,6 +778,9 @@ Foam::faMesh::faMesh
GeoMesh<polyMesh>(pMesh), GeoMesh<polyMesh>(pMesh),
MeshObject<polyMesh, Foam::UpdateableMeshObject, faMesh>(pMesh), MeshObject<polyMesh, Foam::UpdateableMeshObject, faMesh>(pMesh),
edgeInterpolation(*this), edgeInterpolation(*this),
faSchemes(mesh()),
faSolution(mesh()),
data(mesh()),
faceLabels_ faceLabels_
( (
IOobject IOobject
@ -838,8 +828,7 @@ Foam::faMesh::faMesh
correctPatchPointNormalsPtr_(nullptr), correctPatchPointNormalsPtr_(nullptr),
globalMeshDataPtr_(nullptr) globalMeshDataPtr_(nullptr)
{ {
DebugInFunction DebugInFunction << "Creating faMesh from polyPatch" << endl;
<< "Creating faMesh from polyPatch" << endl;
const polyBoundaryMesh& pbm = pMesh.boundaryMesh(); const polyBoundaryMesh& pbm = pMesh.boundaryMesh();
@ -977,8 +966,7 @@ const Foam::faceList& Foam::faMesh::faces() const
void Foam::faMesh::addFaPatches(const List<faPatch*>& p) void Foam::faMesh::addFaPatches(const List<faPatch*>& p)
{ {
DebugInFunction DebugInFunction << "Adding patches to faMesh" << endl;
<< "Adding patches to faMesh" << endl;
if (boundary().size() > 0) if (boundary().size() > 0)
{ {
@ -1304,7 +1292,7 @@ Foam::boolList& Foam::faMesh::correctPatchPointNormals() const
} }
bool Foam::faMesh::write() const bool Foam::faMesh::write(const bool valid) const
{ {
faceLabels_.write(); faceLabels_.write();
boundary_.write(); boundary_.write();

View File

@ -58,6 +58,9 @@ Author
#include "labelIOList.H" #include "labelIOList.H"
#include "FieldFields.H" #include "FieldFields.H"
#include "faGlobalMeshData.H" #include "faGlobalMeshData.H"
#include "faSchemes.H"
#include "faSolution.H"
#include "data.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -77,7 +80,10 @@ class faMesh
public GeoMesh<polyMesh>, public GeoMesh<polyMesh>,
public MeshObject<polyMesh, Foam::UpdateableMeshObject, faMesh>, public MeshObject<polyMesh, Foam::UpdateableMeshObject, faMesh>,
public lduMesh, public lduMesh,
public edgeInterpolation public edgeInterpolation,
public faSchemes,
public faSolution,
public data
{ {
// Private data // Private data
@ -402,6 +408,13 @@ public:
//- Return reference to the mesh database //- Return reference to the mesh database
virtual const objectRegistry& thisDb() const; virtual const objectRegistry& thisDb() const;
//- Name function is needed to disambiguate those inherited
// from base classes
const word& name() const
{
return thisDb().name();
}
//- Return constant reference to boundary mesh //- Return constant reference to boundary mesh
const faBoundaryMesh& boundary() const; const faBoundaryMesh& boundary() const;
@ -531,7 +544,7 @@ public:
boolList& correctPatchPointNormals() const; boolList& correctPatchPointNormals() const;
//- Write mesh //- Write mesh
virtual bool write() const; virtual bool write(const bool valid = true) const;
// Member Operators // Member Operators

View File

@ -52,7 +52,7 @@ ddt
return fa::faDdtScheme<Type>::New return fa::faDdtScheme<Type>::New
( (
mesh, mesh,
mesh.schemesDict().ddtScheme("ddt(" + dt.name() + ')') mesh.ddtScheme("ddt(" + dt.name() + ')')
).ref().facDdt(dt); ).ref().facDdt(dt);
} }
@ -67,7 +67,7 @@ ddt
return fa::faDdtScheme<Type>::New return fa::faDdtScheme<Type>::New
( (
vf.mesh(), vf.mesh(),
vf.mesh().schemesDict().ddtScheme("ddt(" + vf.name() + ')') vf.mesh().ddtScheme("ddt(" + vf.name() + ')')
).ref().facDdt(vf); ).ref().facDdt(vf);
} }
@ -83,7 +83,7 @@ ddt
return fa::faDdtScheme<Type>::New return fa::faDdtScheme<Type>::New
( (
vf.mesh(), vf.mesh(),
vf.mesh().schemesDict().ddtScheme vf.mesh().ddtScheme
( (
"ddt(" + rho.name() + ',' + vf.name() + ')' "ddt(" + rho.name() + ',' + vf.name() + ')'
) )
@ -102,7 +102,7 @@ ddt
return fa::faDdtScheme<Type>::New return fa::faDdtScheme<Type>::New
( (
vf.mesh(), vf.mesh(),
vf.mesh().schemesDict().ddtScheme vf.mesh().ddtScheme
( (
"ddt(" + rho.name() + ',' + vf.name() + ')' "ddt(" + rho.name() + ',' + vf.name() + ')'
) )

View File

@ -107,7 +107,7 @@ div
( (
fa::divScheme<Type>::New fa::divScheme<Type>::New
( (
vf.mesh(), vf.mesh().schemesDict().divScheme(name) vf.mesh(), vf.mesh().divScheme(name)
).ref().facDiv(vf) ).ref().facDiv(vf)
); );
GeometricField GeometricField
@ -207,7 +207,7 @@ div
( (
vf.mesh(), vf.mesh(),
flux, flux,
vf.mesh().schemesDict().divScheme(name) vf.mesh().divScheme(name)
).ref().facDiv(flux, vf) ).ref().facDiv(flux, vf)
); );
GeometricField<Type, faPatchField, areaMesh>& Div = tDiv.ref(); GeometricField<Type, faPatchField, areaMesh>& Div = tDiv.ref();

View File

@ -115,7 +115,7 @@ grad
fa::gradScheme<Type>::New fa::gradScheme<Type>::New
( (
vf.mesh(), vf.mesh(),
vf.mesh().schemesDict().gradScheme(name) vf.mesh().gradScheme(name)
).ref().grad(vf); ).ref().grad(vf);
GeometricField<GradType, faPatchField, areaMesh>& gGrad = tgGrad.ref(); GeometricField<GradType, faPatchField, areaMesh>& gGrad = tgGrad.ref();

View File

@ -52,7 +52,7 @@ laplacian
return fa::laplacianScheme<Type>::New return fa::laplacianScheme<Type>::New
( (
vf.mesh(), vf.mesh(),
vf.mesh().schemesDict().laplacianScheme(name) vf.mesh().laplacianScheme(name)
).ref().facLaplacian(vf); ).ref().facLaplacian(vf);
} }
@ -172,7 +172,7 @@ laplacian
return fa::laplacianScheme<Type>::New return fa::laplacianScheme<Type>::New
( (
vf.mesh(), vf.mesh(),
vf.mesh().schemesDict().laplacianScheme(name) vf.mesh().laplacianScheme(name)
).ref().facLaplacian(gamma, vf); ).ref().facLaplacian(gamma, vf);
} }
@ -306,7 +306,7 @@ laplacian
return fa::laplacianScheme<Type>::New return fa::laplacianScheme<Type>::New
( (
vf.mesh(), vf.mesh(),
vf.mesh().schemesDict().laplacianScheme(name) vf.mesh().laplacianScheme(name)
).ref().facLaplacian(gamma, vf); ).ref().facLaplacian(gamma, vf);
} }

View File

@ -52,7 +52,7 @@ lnGrad
return fa::lnGradScheme<Type>::New return fa::lnGradScheme<Type>::New
( (
vf.mesh(), vf.mesh(),
vf.mesh().schemesDict().lnGradScheme(name) vf.mesh().lnGradScheme(name)
).ref().lnGrad(vf); ).ref().lnGrad(vf);
} }

View File

@ -98,7 +98,7 @@ ndiv
( (
fa::divScheme<Type>::New fa::divScheme<Type>::New
( (
vf.mesh(), vf.mesh().schemesDict().divScheme(name) vf.mesh(), vf.mesh().divScheme(name)
).ref().facDiv(vf) ).ref().facDiv(vf)
); );
GeometricField<Type, faPatchField, areaMesh>& Div = tDiv.ref(); GeometricField<Type, faPatchField, areaMesh>& Div = tDiv.ref();
@ -192,7 +192,7 @@ ndiv
( (
vf.mesh(), vf.mesh(),
flux, flux,
vf.mesh().schemesDict().divScheme(name) vf.mesh().divScheme(name)
).ref().facDiv(flux, vf) ).ref().facDiv(flux, vf)
); );

View File

@ -118,7 +118,7 @@ ngrad
fa::gradScheme<Type>::New fa::gradScheme<Type>::New
( (
vf.mesh(), vf.mesh(),
vf.mesh().schemesDict().gradScheme(name) vf.mesh().gradScheme(name)
).ref().grad(vf); ).ref().grad(vf);
GeometricField<GradType, faPatchField, areaMesh>& gGrad = tgGrad.ref(); GeometricField<GradType, faPatchField, areaMesh>& gGrad = tgGrad.ref();

View File

@ -52,7 +52,7 @@ ddt
return fa::faDdtScheme<Type>::New return fa::faDdtScheme<Type>::New
( (
vf.mesh(), vf.mesh(),
vf.mesh().schemesDict().ddtScheme("ddt(" + vf.name() + ')') vf.mesh().ddtScheme("ddt(" + vf.name() + ')')
).ref().famDdt(vf); ).ref().famDdt(vf);
} }
@ -68,7 +68,7 @@ ddt
return fa::faDdtScheme<Type>::New return fa::faDdtScheme<Type>::New
( (
vf.mesh(), vf.mesh(),
vf.mesh().schemesDict().ddtScheme vf.mesh().ddtScheme
( (
"ddt(" + rho.name() + ',' + vf.name() + ')' "ddt(" + rho.name() + ',' + vf.name() + ')'
) )
@ -87,7 +87,7 @@ ddt
return fa::faDdtScheme<Type>::New return fa::faDdtScheme<Type>::New
( (
vf.mesh(), vf.mesh(),
vf.mesh().schemesDict().ddtScheme vf.mesh().ddtScheme
( (
"ddt(" + rho.name() + ',' + vf.name() + ')' "ddt(" + rho.name() + ',' + vf.name() + ')'
) )

View File

@ -59,7 +59,7 @@ div
( (
vf.mesh(), vf.mesh(),
flux, flux,
vf.mesh().schemesDict().divScheme(name) vf.mesh().divScheme(name)
).ref().famDiv(flux, vf) ).ref().famDiv(flux, vf)
); );
faMatrix<Type>& M = tM.ref(); faMatrix<Type>& M = tM.ref();
@ -70,7 +70,7 @@ div
( (
vf.mesh(), vf.mesh(),
flux, flux,
vf.mesh().schemesDict().divScheme(name) vf.mesh().divScheme(name)
).ref().facDiv(flux, vf) ).ref().facDiv(flux, vf)
); );

View File

@ -172,7 +172,7 @@ laplacian
return fa::laplacianScheme<Type>::New return fa::laplacianScheme<Type>::New
( (
vf.mesh(), vf.mesh(),
vf.mesh().schemesDict().laplacianScheme(name) vf.mesh().laplacianScheme(name)
).ref().famLaplacian(gamma, vf); ).ref().famLaplacian(gamma, vf);
} }
@ -218,7 +218,7 @@ laplacian
return fa::laplacianScheme<Type>::New return fa::laplacianScheme<Type>::New
( (
vf.mesh(), vf.mesh(),
vf.mesh().schemesDict().laplacianScheme(name) vf.mesh().laplacianScheme(name)
).ref().famLaplacian(gamma, vf); ).ref().famLaplacian(gamma, vf);
} }

View File

@ -55,7 +55,7 @@ ndiv
( (
vf.mesh(), vf.mesh(),
flux, flux,
vf.mesh().schemesDict().divScheme(name) vf.mesh().divScheme(name)
).ref().famDiv(flux, vf);//TODO calculate normal ).ref().famDiv(flux, vf);//TODO calculate normal
} }

View File

@ -55,7 +55,7 @@ div
( (
vf.mesh(), vf.mesh(),
flux, flux,
vf.mesh().schemesDict().divScheme(name) vf.mesh().divScheme(name)
).ref().famDiv(flux, vf); ).ref().famDiv(flux, vf);
} }

View File

@ -50,8 +50,6 @@ SourceFiles
namespace Foam namespace Foam
{ {
class mapPolyMesh;
/*---------------------------------------------------------------------------*\ /*---------------------------------------------------------------------------*\
Class leastSquaresFaVectors Declaration Class leastSquaresFaVectors Declaration
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/

View File

@ -80,7 +80,7 @@ gaussLaplacianScheme<Type>::famLaplacian
if (this->tlnGradScheme_().corrected()) if (this->tlnGradScheme_().corrected())
{ {
if (this->mesh().schemesDict().fluxRequired(vf.name())) if (this->mesh().fluxRequired(vf.name()))
{ {
fam.faceFluxCorrectionPtr() = new fam.faceFluxCorrectionPtr() = new
GeometricField<Type, faePatchField, edgeMesh> GeometricField<Type, faePatchField, edgeMesh>
@ -89,8 +89,8 @@ gaussLaplacianScheme<Type>::famLaplacian
); );
fam.source() -= fam.source() -=
this->mesh().S()* this->mesh().S()
fac::div *fac::div
( (
*fam.faceFluxCorrectionPtr() *fam.faceFluxCorrectionPtr()
)().internalField(); )().internalField();
@ -98,8 +98,8 @@ gaussLaplacianScheme<Type>::famLaplacian
else else
{ {
fam.source() -= fam.source() -=
this->mesh().S()* this->mesh().S()
fac::div *fac::div
( (
gammaMagSf*this->tlnGradScheme_().correction(vf) gammaMagSf*this->tlnGradScheme_().correction(vf)
)().internalField(); )().internalField();

View File

@ -90,7 +90,7 @@ correctedLnGrad<Type>::correction
gradScheme<typename pTraits<Type>::cmptType>::New gradScheme<typename pTraits<Type>::cmptType>::New
( (
mesh, mesh,
mesh.schemesDict().gradScheme(ssf.name()) mesh.gradScheme(ssf.name())
)() )()
.grad(vf.component(cmpt)) .grad(vf.component(cmpt))
) )

View File

@ -56,7 +56,7 @@ Foam::tmp<Foam::edgeInterpolationScheme<Type>> Foam::fac::scheme
( (
faceFlux.mesh(), faceFlux.mesh(),
faceFlux, faceFlux,
faceFlux.mesh().schemesDict().interpolationScheme(name) faceFlux.mesh().interpolationScheme(name)
); );
} }
@ -86,7 +86,7 @@ Foam::tmp<Foam::edgeInterpolationScheme<Type>> Foam::fac::scheme
return edgeInterpolationScheme<Type>::New return edgeInterpolationScheme<Type>::New
( (
mesh, mesh,
mesh.schemesDict().interpolationScheme(name) mesh.interpolationScheme(name)
); );
} }

View File

@ -57,8 +57,6 @@ void Foam::edgeInterpolation::clearOut()
Foam::edgeInterpolation::edgeInterpolation(const faMesh& fam) Foam::edgeInterpolation::edgeInterpolation(const faMesh& fam)
: :
faMesh_(fam), faMesh_(fam),
schemesDict_(fam()),
solutionDict_(fam()),
lPN_(nullptr), lPN_(nullptr),
weightingFactors_(nullptr), weightingFactors_(nullptr),
differenceFactors_(nullptr), differenceFactors_(nullptr),
@ -219,7 +217,7 @@ void Foam::edgeInterpolation::makeLPN() const
( (
"lPN", "lPN",
faMesh_.time().constant(), faMesh_.time().constant(),
faMesh_.db(), faMesh_(),
IOobject::NO_READ, IOobject::NO_READ,
IOobject::NO_WRITE, IOobject::NO_WRITE,
false false

View File

@ -43,8 +43,6 @@ SourceFiles
#include "tmp.H" #include "tmp.H"
#include "scalar.H" #include "scalar.H"
#include "faSchemes.H"
#include "faSolution.H"
#include "areaFieldsFwd.H" #include "areaFieldsFwd.H"
#include "edgeFieldsFwd.H" #include "edgeFieldsFwd.H"
#include "className.H" #include "className.H"
@ -67,12 +65,6 @@ class edgeInterpolation
// Reference to faMesh // Reference to faMesh
const faMesh& faMesh_; const faMesh& faMesh_;
//- Discretisation schemes
faSchemes schemesDict_;
//- Solver settings
faSolution solutionDict_;
// Demand-driven data // Demand-driven data
@ -153,30 +145,6 @@ public:
return faMesh_; return faMesh_;
} }
//- Return schemes
const faSchemes& schemesDict() const
{
return schemesDict_;
}
//- Return access to schemes
faSchemes& schemesDict()
{
return schemesDict_;
}
//- Return solver settings
const faSolution& solutionDict() const
{
return solutionDict_;
}
//- Return access to solver settings
faSolution& solutionDict()
{
return solutionDict_;
}
//- Return reference to PN geodesic distance //- Return reference to PN geodesic distance
const edgeScalarField& lPN() const; const edgeScalarField& lPN() const;