fvMesh: fvSchemes and fvSolution are now demand-driven

fvMesh is no longer derived from fvSchemes and fvSolution, these are now
demand-driven and accessed by the member functions schemes() and solution()
respectively.  This means that the system/fvSchemes and system/fvSolution files
are no longer required during fvMesh constructions simplifying the mesh
generation and manipulation phase; theses files are read on the first call of
their access functions.

The fvSchemes member function names have also been simplified taking advantage
of the context in which they are called, for example

    mesh.ddtScheme(fieldName) -> mesh.schemes().ddt(fieldName)
This commit is contained in:
Henry Weller
2022-03-23 16:23:55 +00:00
parent acd5528557
commit ddbf2d7853
122 changed files with 418 additions and 296 deletions

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2018-2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2018-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -156,9 +156,9 @@ bool Foam::functionObjects::age::execute()
// Set under-relaxation coeff
scalar relaxCoeff = 0.0;
if (mesh_.relaxEquation(schemesField_))
if (mesh_.solution().relaxEquation(schemesField_))
{
relaxCoeff = mesh_.equationRelaxationFactor(schemesField_);
relaxCoeff = mesh_.solution().equationRelaxationFactor(schemesField_);
}
const Foam::fvModels& fvModels(Foam::fvModels::New(mesh_));

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2013-2018 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2013-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -42,7 +42,7 @@ bool Foam::functionObjects::blendingFactor::calcBF()
const FieldType& field = lookupObject<FieldType>(fieldName_);
const word divScheme("div(" + phiName_ + ',' + fieldName_ + ')');
ITstream& its = mesh_.divScheme(divScheme);
ITstream& its = mesh_.schemes().div(divScheme);
const surfaceScalarField& phi = lookupObject<surfaceScalarField>(phiName_);

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2012-2018 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2012-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -39,7 +39,7 @@ bool Foam::functionObjects::grad::calcGrad()
(
resultName_,
fvc::grad(lookupObject<VolFieldType>(fieldName_)),
mesh_.changing() && mesh_.cache(resultName_)
mesh_.changing() && mesh_.solution().cache(resultName_)
);
}
else if (foundObject<SurfaceFieldType>(fieldName_))
@ -48,7 +48,7 @@ bool Foam::functionObjects::grad::calcGrad()
(
resultName_,
fvc::grad(lookupObject<SurfaceFieldType>(fieldName_)),
mesh_.changing() && mesh_.cache(resultName_)
mesh_.changing() && mesh_.solution().cache(resultName_)
);
}
else

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2019-2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2019-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -102,7 +102,7 @@ Foam::volScalarField& Foam::functionObjects::phaseScalarTransport::Phi()
)
);
mesh_.setFluxRequired(PhiPtr_->name());
mesh_.schemes().setFluxRequired(PhiPtr_->name());
}
return PhiPtr_();
@ -392,8 +392,8 @@ bool Foam::functionObjects::phaseScalarTransport::execute()
// Get the relaxation coefficient
const scalar relaxCoeff =
mesh_.relaxEquation(schemesField_)
? mesh_.equationRelaxationFactor(schemesField_)
mesh_.solution().relaxEquation(schemesField_)
? mesh_.solution().equationRelaxationFactor(schemesField_)
: 0;
const Foam::fvModels& fvModels(Foam::fvModels::New(mesh_));

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2012-2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2012-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -147,7 +147,7 @@ Foam::functionObjects::scalarTransport::scalarTransport
{
read(dict);
const dictionary& controls = mesh_.solverDict(s_.name());
const dictionary& controls = mesh_.solution().solverDict(s_.name());
if (controls.found("nSubCycles"))
{
@ -238,9 +238,9 @@ bool Foam::functionObjects::scalarTransport::execute()
// Set under-relaxation coeff
scalar relaxCoeff = 0.0;
if (mesh_.relaxEquation(schemesField_))
if (mesh_.solution().relaxEquation(schemesField_))
{
relaxCoeff = mesh_.equationRelaxationFactor(schemesField_);
relaxCoeff = mesh_.solution().equationRelaxationFactor(schemesField_);
}
const Foam::fvModels& fvModels(Foam::fvModels::New(mesh_));
@ -329,7 +329,7 @@ bool Foam::functionObjects::scalarTransport::execute()
void Foam::functionObjects::scalarTransport::subCycleMULES()
{
const dictionary& controls = mesh_.solverDict(s_.name());
const dictionary& controls = mesh_.solution().solverDict(s_.name());
const label nSubCycles(controls.lookup<label>("nSubCycles"));
const bool LTS = fv::localEulerDdt::enabled(mesh_);
@ -381,7 +381,7 @@ void Foam::functionObjects::scalarTransport::subCycleMULES()
void Foam::functionObjects::scalarTransport::solveMULES()
{
const dictionary& controls = mesh_.solverDict(s_.name());
const dictionary& controls = mesh_.solution().solverDict(s_.name());
const label nCorr(controls.lookup<label>("nCorr"));
const label nSubCycles(controls.lookup<label>("nSubCycles"));
const bool MULESCorr(controls.lookupOrDefault<Switch>("MULESCorr", false));
@ -403,7 +403,7 @@ void Foam::functionObjects::scalarTransport::solveMULES()
surfaceScalarField& sPhi_ = tsPhi_.ref();
const word sScheme(mesh_.divScheme(divScheme)[1].wordToken());
const word sScheme(mesh_.schemes().div(divScheme)[1].wordToken());
// If a compressive convection scheme is used
// the interface normal must be cached
@ -433,7 +433,7 @@ void Foam::functionObjects::scalarTransport::solveMULES()
fv::ddtScheme<scalar>::New
(
mesh_,
mesh_.ddtScheme("ddt(s)")
mesh_.schemes().ddt("ddt(s)")
)
);
const fv::ddtScheme<scalar>& ddtS = tddtS();
@ -547,7 +547,7 @@ void Foam::functionObjects::scalarTransport::solveMULES()
(
phiCN(),
(cnCoeff*s_ + (1.0 - cnCoeff)*s_.oldTime())(),
mesh_.divScheme(divScheme)
mesh_.schemes().div(divScheme)
)
);
@ -605,7 +605,7 @@ void Foam::functionObjects::scalarTransport::solveMULES()
if
(
word(mesh_.ddtScheme("ddt(s)"))
word(mesh_.schemes().ddt("ddt(s)"))
== fv::CrankNicolsonDdtScheme<scalar>::typeName
)
{