STYLE: scalarTransport/energyTransport: modernise the code

- Remove redundant copy ctor and assignment operator (already deleted in base class)
- Remove unused header files
- Use default destructor
- Reorder member variables
This commit is contained in:
Kutalmis Bercin
2024-07-25 11:45:17 +01:00
committed by Andrew Heather
parent 559f13d450
commit 0ff5eb5687
4 changed files with 94 additions and 121 deletions

View File

@ -26,11 +26,6 @@ License
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
#include "energyTransport.H" #include "energyTransport.H"
#include "surfaceFields.H"
#include "fvmDdt.H"
#include "fvmDiv.H"
#include "fvmLaplacian.H"
#include "fvmSup.H"
#include "turbulentTransportModel.H" #include "turbulentTransportModel.H"
#include "turbulentFluidThermoModel.H" #include "turbulentFluidThermoModel.H"
#include "addToRunTimeSelectionTable.H" #include "addToRunTimeSelectionTable.H"
@ -188,24 +183,6 @@ Foam::functionObjects::energyTransport::energyTransport
) )
: :
fvMeshFunctionObject(name, runTime, dict), fvMeshFunctionObject(name, runTime, dict),
fieldName_(dict.getOrDefault<word>("field", "T")),
phiName_(dict.getOrDefault<word>("phi", "phi")),
rhoName_(dict.getOrDefault<word>("rho", "rho")),
nCorr_(0),
tol_(1),
schemesField_("unknown-schemesField"),
fvOptions_(mesh_),
multiphaseThermo_(dict.subOrEmptyDict("phaseThermos")),
Cp_("Cp", dimEnergy/dimMass/dimTemperature, 0, dict),
kappa_
(
"kappa",
dimEnergy/dimTime/dimLength/dimTemperature,
0,
dict
),
rho_("rhoInf", dimDensity, 0, dict),
Prt_("Prt", dimless, 1, dict),
rhoCp_ rhoCp_
( (
IOobject IOobject
@ -219,7 +196,25 @@ Foam::functionObjects::energyTransport::energyTransport
), ),
mesh_, mesh_,
dimensionedScalar(dimEnergy/dimTemperature/dimVolume, Zero) dimensionedScalar(dimEnergy/dimTemperature/dimVolume, Zero)
) ),
fvOptions_(mesh_),
multiphaseThermo_(dict.subOrEmptyDict("phaseThermos")),
Cp_("Cp", dimEnergy/dimMass/dimTemperature, 0, dict),
kappa_
(
"kappa",
dimEnergy/dimTime/dimLength/dimTemperature,
0,
dict
),
rho_("rhoInf", dimDensity, 0, dict),
Prt_("Prt", dimless, 1, dict),
fieldName_(dict.getOrDefault<word>("field", "T")),
schemesField_("unknown-schemesField"),
phiName_(dict.getOrDefault<word>("phi", "phi")),
rhoName_(dict.getOrDefault<word>("rho", "rho")),
tol_(1),
nCorr_(0)
{ {
read(dict); read(dict);
@ -306,17 +301,14 @@ Foam::functionObjects::energyTransport::energyTransport
} }
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::functionObjects::energyTransport::~energyTransport()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
bool Foam::functionObjects::energyTransport::read(const dictionary& dict) bool Foam::functionObjects::energyTransport::read(const dictionary& dict)
{ {
fvMeshFunctionObject::read(dict); if (!fvMeshFunctionObject::read(dict))
{
return false;
}
dict.readIfPresent("phi", phiName_); dict.readIfPresent("phi", phiName_);
dict.readIfPresent("rho", rhoName_); dict.readIfPresent("rho", rhoName_);
@ -359,14 +351,14 @@ bool Foam::functionObjects::energyTransport::execute()
// Convergence monitor parameters // Convergence monitor parameters
bool converged = false; bool converged = false;
label iter = 0; int iter = 0;
if (phi.dimensions() == dimMass/dimTime) if (phi.dimensions() == dimMass/dimTime)
{ {
rhoCp_ = rho()*Cp(); rhoCp_ = rho()*Cp();
const surfaceScalarField rhoCpPhi(fvc::interpolate(Cp())*phi); const surfaceScalarField rhoCpPhi(fvc::interpolate(Cp())*phi);
for (label i = 0; i <= nCorr_; i++) for (int i = 0; i <= nCorr_; ++i)
{ {
fvScalarMatrix sEqn fvScalarMatrix sEqn
( (
@ -401,7 +393,7 @@ bool Foam::functionObjects::energyTransport::execute()
rhoCp rhoCp
); );
for (label i = 0; i <= nCorr_; i++) for (int i = 0; i <= nCorr_; ++i)
{ {
fvScalarMatrix sEqn fvScalarMatrix sEqn
( (

View File

@ -70,7 +70,7 @@ Usage
Prt <scalar>; Prt <scalar>;
schemesField <word>; schemesField <word>;
tolerance <scalar>; tolerance <scalar>;
nCorr <label>; nCorr <int>;
fvOptions <dict>; fvOptions <dict>;
phaseThermos <dict>; phaseThermos <dict>;
@ -93,7 +93,7 @@ Usage
Prt | Turbulent Prandtl number | scalar | no | 1 Prt | Turbulent Prandtl number | scalar | no | 1
schemesField | Name of field to specify schemes | word | no | field schemesField | Name of field to specify schemes | word | no | field
tolerance | Outer-loop initial-residual tolerance | scalar | no | 1 tolerance | Outer-loop initial-residual tolerance | scalar | no | 1
nCorr | Number of outer-loop correctors | label | no | 0 nCorr | Number of outer-loop correctors | int | no | 0
fvOptions | List of finite-volume options | dict | no | - fvOptions | List of finite-volume options | dict | no | -
phaseThermos | Dictionary for multi-phase thermo | dict | no | null phaseThermos | Dictionary for multi-phase thermo | dict | no | null
\endtable \endtable
@ -232,25 +232,10 @@ class energyTransport
: :
public fvMeshFunctionObject public fvMeshFunctionObject
{ {
// Private data // Private Data
//- Name of the transport field. //- Volumetric heat capacity field [J/m^3/K]
word fieldName_; volScalarField rhoCp_;
//- Name of flux field
word phiName_;
//- Name of density field
word rhoName_;
//- Number of corrector iterations (optional)
label nCorr_;
//- Outer-loop initial-residual tolerance
scalar tol_;
//- Name of field whose schemes are used (optional)
word schemesField_;
//- Run-time selectable finite volume options, e.g. sources, constraints //- Run-time selectable finite volume options, e.g. sources, constraints
fv::optionList fvOptions_; fv::optionList fvOptions_;
@ -261,7 +246,7 @@ class energyTransport
//- List of phase names //- List of phase names
wordList phaseNames_; wordList phaseNames_;
//- List of phase heat capacities //- List of phase specific heat capacities at constant pressure
PtrList<dimensionedScalar> Cps_; PtrList<dimensionedScalar> Cps_;
//- List of phase thermal diffusivity for temperature [J/m/s/K] //- List of phase thermal diffusivity for temperature [J/m/s/K]
@ -270,7 +255,7 @@ class energyTransport
//- Unallocated phase list //- Unallocated phase list
UPtrList<volScalarField> phases_; UPtrList<volScalarField> phases_;
//- Heat capacity for single phase flows //- Specific heat capacity at constant pressure for single phase flows
dimensionedScalar Cp_; dimensionedScalar Cp_;
//- Thermal diffusivity for temperature for single phase flows //- Thermal diffusivity for temperature for single phase flows
@ -282,8 +267,23 @@ class energyTransport
//- Turbulent Prandt number //- Turbulent Prandt number
dimensionedScalar Prt_; dimensionedScalar Prt_;
//- rhoCp //- Name of the transport field
volScalarField rhoCp_; word fieldName_;
//- Name of field whose schemes are used
word schemesField_;
//- Name of flux field
word phiName_;
//- Name of density field
word rhoName_;
//- Outer-loop initial-residual tolerance
scalar tol_;
//- Number of corrector iterations
int nCorr_;
// Private Member Functions // Private Member Functions
@ -294,21 +294,15 @@ class energyTransport
//- Return the diffusivity field //- Return the diffusivity field
tmp<volScalarField> kappaEff() const; tmp<volScalarField> kappaEff() const;
//- Return rho field //- Return the density field, rho
tmp<volScalarField> rho() const; tmp<volScalarField> rho() const;
//- Return Cp //- Return the specific heat capacity at constant pressure field, Cp
tmp<volScalarField> Cp() const; tmp<volScalarField> Cp() const;
//- Return kappa //- Return the thermal diffusivity field
tmp<volScalarField> kappa() const; tmp<volScalarField> kappa() const;
//- No copy construct
energyTransport(const energyTransport&) = delete;
//- No copy assignment
void operator=(const energyTransport&) = delete;
public: public:
@ -328,7 +322,7 @@ public:
//- Destructor //- Destructor
virtual ~energyTransport(); virtual ~energyTransport() = default;
// Member Functions // Member Functions

View File

@ -27,11 +27,6 @@ License
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
#include "scalarTransport.H" #include "scalarTransport.H"
#include "surfaceFields.H"
#include "fvmDdt.H"
#include "fvmDiv.H"
#include "fvmLaplacian.H"
#include "fvmSup.H"
#include "CMULES.H" #include "CMULES.H"
#include "turbulentTransportModel.H" #include "turbulentTransportModel.H"
#include "turbulentFluidThermoModel.H" #include "turbulentFluidThermoModel.H"
@ -172,7 +167,9 @@ Foam::functionObjects::scalarTransport::scalarTransport
) )
: :
fvMeshFunctionObject(name, runTime, dict), fvMeshFunctionObject(name, runTime, dict),
fvOptions_(mesh_),
fieldName_(dict.getOrDefault<word>("field", "s")), fieldName_(dict.getOrDefault<word>("field", "s")),
schemesField_("unknown-schemesField"),
phiName_(dict.getOrDefault<word>("phi", "phi")), phiName_(dict.getOrDefault<word>("phi", "phi")),
rhoName_(dict.getOrDefault<word>("rho", "rho")), rhoName_(dict.getOrDefault<word>("rho", "rho")),
nutName_(dict.getOrDefault<word>("nut", "none")), nutName_(dict.getOrDefault<word>("nut", "none")),
@ -182,12 +179,12 @@ Foam::functionObjects::scalarTransport::scalarTransport
dict.getOrDefault<word>("phasePhiCompressed", "alphaPhiUn") dict.getOrDefault<word>("phasePhiCompressed", "alphaPhiUn")
), ),
D_(0), D_(0),
constantD_(false), alphaD_(1),
alphaDt_(1),
tol_(1), tol_(1),
nCorr_(0), nCorr_(0),
resetOnStartUp_(false), resetOnStartUp_(false),
schemesField_("unknown-schemesField"), constantD_(false),
fvOptions_(mesh_),
bounded01_(dict.getOrDefault("bounded01", true)) bounded01_(dict.getOrDefault("bounded01", true))
{ {
read(dict); read(dict);
@ -203,32 +200,30 @@ Foam::functionObjects::scalarTransport::scalarTransport
} }
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::functionObjects::scalarTransport::~scalarTransport()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
bool Foam::functionObjects::scalarTransport::read(const dictionary& dict) bool Foam::functionObjects::scalarTransport::read(const dictionary& dict)
{ {
fvMeshFunctionObject::read(dict); if (!fvMeshFunctionObject::read(dict))
{
return false;
}
dict.readIfPresent("phi", phiName_); dict.readIfPresent("phi", phiName_);
dict.readIfPresent("rho", rhoName_); dict.readIfPresent("rho", rhoName_);
dict.readIfPresent("nut", nutName_); dict.readIfPresent("nut", nutName_);
dict.readIfPresent("phase", phaseName_); dict.readIfPresent("phase", phaseName_);
dict.readIfPresent("bounded01", bounded01_); dict.readIfPresent("phasePhiCompressed", phasePhiCompressedName_);
schemesField_ = dict.getOrDefault("schemesField", fieldName_); schemesField_ = dict.getOrDefault("schemesField", fieldName_);
constantD_ = dict.readIfPresent("D", D_);
alphaD_ = dict.getOrDefault<scalar>("alphaD", 1);
alphaDt_ = dict.getOrDefault<scalar>("alphaDt", 1);
dict.readIfPresent("alphaD", alphaD_);
dict.readIfPresent("alphaDt", alphaDt_);
dict.readIfPresent("tolerance", tol_); dict.readIfPresent("tolerance", tol_);
dict.readIfPresent("nCorr", nCorr_); dict.readIfPresent("nCorr", nCorr_);
dict.readIfPresent("resetOnStartUp", resetOnStartUp_); dict.readIfPresent("resetOnStartUp", resetOnStartUp_);
constantD_ = dict.readIfPresent("D", D_);
dict.readIfPresent("bounded01", bounded01_);
if (dict.found("fvOptions")) if (dict.found("fvOptions"))
{ {
@ -260,7 +255,7 @@ bool Foam::functionObjects::scalarTransport::execute()
// Convergence monitor parameters // Convergence monitor parameters
bool converged = false; bool converged = false;
label iter = 0; int iter = 0;
// Two phase scalar transport // Two phase scalar transport
if (phaseName_ != "none") if (phaseName_ != "none")
@ -278,7 +273,7 @@ bool Foam::functionObjects::scalarTransport::execute()
// Solve // Solve
tmp<surfaceScalarField> tTPhiUD; tmp<surfaceScalarField> tTPhiUD;
for (label i = 0; i <= nCorr_; i++) for (int i = 0; i <= nCorr_; ++i)
{ {
fvScalarMatrix sEqn fvScalarMatrix sEqn
( (
@ -317,9 +312,8 @@ bool Foam::functionObjects::scalarTransport::execute()
{ {
const volScalarField& rho = lookupObject<volScalarField>(rhoName_); const volScalarField& rho = lookupObject<volScalarField>(rhoName_);
for (label i = 0; i <= nCorr_; i++) for (int i = 0; i <= nCorr_; ++i)
{ {
fvScalarMatrix sEqn fvScalarMatrix sEqn
( (
fvm::ddt(rho, s) fvm::ddt(rho, s)
@ -340,7 +334,7 @@ bool Foam::functionObjects::scalarTransport::execute()
} }
else if (phi.dimensions() == dimVolume/dimTime) else if (phi.dimensions() == dimVolume/dimTime)
{ {
for (label i = 0; i <= nCorr_; i++) for (int i = 0; i <= nCorr_; ++i)
{ {
fvScalarMatrix sEqn fvScalarMatrix sEqn
( (

View File

@ -73,7 +73,7 @@ Usage
alphaD <scalar>; alphaD <scalar>;
alphaDt <scalar>; alphaDt <scalar>;
tolerance <scalar>; tolerance <scalar>;
nCorr <label>; nCorr <int>;
resetOnStartUp <bool>; resetOnStartUp <bool>;
fvOptions <dict>; fvOptions <dict>;
@ -98,7 +98,7 @@ Usage
alphaD | Laminar diffusivity coefficient | scalar | no | 1 alphaD | Laminar diffusivity coefficient | scalar | no | 1
alphaDt | Turbulent diffusivity coefficient | scalar | no | 1 alphaDt | Turbulent diffusivity coefficient | scalar | no | 1
tolerance | Outer-loop initial-residual tolerance | scalar | no | 1 tolerance | Outer-loop initial-residual tolerance | scalar | no | 1
nCorr | Number of outer-loop correctors | label | no | 0 nCorr | Number of outer-loop correctors | int | no | 0
resetOnStartUp | Flag to reset field to zero on start-up | bool | no | no resetOnStartUp | Flag to reset field to zero on start-up | bool | no | no
fvOptions | List of finite-volume options | dict | no | - fvOptions | List of finite-volume options | dict | no | -
\endtable \endtable
@ -187,52 +187,52 @@ class scalarTransport
: :
public fvMeshFunctionObject public fvMeshFunctionObject
{ {
// Private data // Private Data
//- Run-time selectable finite volume options, e.g. sources, constraints
fv::optionList fvOptions_;
//- Name of the transport field. //- Name of the transport field.
word fieldName_; word fieldName_;
//- Name of flux field (optional) //- Name of field whose schemes are used
word schemesField_;
//- Name of flux field
word phiName_; word phiName_;
//- Name of density field (optional) //- Name of density field
word rhoName_; word rhoName_;
//- Name of turbulent viscosity field (optional) //- Name of turbulent viscosity field
word nutName_; word nutName_;
//- Name of phase field (optional) //- Name of phase field
word phaseName_; word phaseName_;
//- Name of phase field compressed flux (optional) //- Name of phase field compressed flux
word phasePhiCompressedName_; word phasePhiCompressedName_;
//- Diffusion coefficient (optional) //- Diffusion coefficient
scalar D_; scalar D_;
//- Flag to indicate whether a constant, uniform D_ is specified //- Laminar diffusion coefficient
bool constantD_;
//- Laminar diffusion coefficient (optional)
scalar alphaD_; scalar alphaD_;
//- Turbulent diffusion coefficient (optional) //- Turbulent diffusion coefficient
scalar alphaDt_; scalar alphaDt_;
//- Outer-loop initial-residual tolerance //- Outer-loop initial-residual tolerance
scalar tol_; scalar tol_;
//- Number of corrector iterations (optional) //- Number of corrector iterations
label nCorr_; int nCorr_;
//- Flag to reset the scalar to zero on start-up //- Flag to reset the scalar to zero on start-up
bool resetOnStartUp_; bool resetOnStartUp_;
//- Name of field whose schemes are used (optional) //- Flag to indicate whether a constant, uniform D_ is specified
word schemesField_; bool constantD_;
//- Run-time selectable finite volume options, e.g. sources, constraints
fv::optionList fvOptions_;
//- Bound scalar between 0-1 using MULES for multiphase case //- Bound scalar between 0-1 using MULES for multiphase case
bool bounded01_; bool bounded01_;
@ -251,13 +251,6 @@ class scalarTransport
) const; ) const;
//- No copy construct
scalarTransport(const scalarTransport&) = delete;
//- No copy assignment
void operator=(const scalarTransport&) = delete;
public: public:
//- Runtime type information //- Runtime type information
@ -276,7 +269,7 @@ public:
//- Destructor //- Destructor
virtual ~scalarTransport(); virtual ~scalarTransport() = default;
// Member Functions // Member Functions