Solute equation solver¶
The multicomponent alloy class for solving the solute transport equation.
Installation¶
Pre-requisites:
A working installation of OpenFOAM 9.
The Solute model library.
In the directChillFoam/applications/solver/heatTransfer/directChillFoam/multicomponentAlloy directory, run:
$ wmake libso
Usage¶
Create the alloy in the CreateFields.H file:
multicomponentAlloy alloy("alloy", U, rho, rho2, phi, fvModels, fvConstraints);
Solve the solute transport equations in the energy corrector loop:
alloy.solve(Us, fvModels, fvConstraints);
C++ Classes¶
-
class multicomponentAlloy : public IOdictionary¶
Public Functions
-
multicomponentAlloy(const word &alloyName, const volVectorField &U, const volScalarField &rho, const dimensionedScalar &rho2, const surfaceScalarField &phi, const fvModels &fvModels, const fvConstraints &fvConstraints)¶
Construct from components.
-
inline virtual ~multicomponentAlloy()¶
Destructor.
-
inline const PtrDictionary<soluteModel> &solutes() const¶
Return the solutes.
-
inline PtrDictionary<soluteModel> &solutes()¶
Return the solutes.
-
void solve(const volVectorField &Us_, const fvModels &fvModels_, const fvConstraints &fvConstraints_)¶
Solve each solute equation.
-
bool read()¶
Read base transportProperties dictionary.
-
multicomponentAlloy(const word &alloyName, const volVectorField &U, const volScalarField &rho, const dimensionedScalar &rho2, const surfaceScalarField &phi, const fvModels &fvModels, const fvConstraints &fvConstraints)¶
Testing¶
The unit tests for the multicomponentAlloy libraries require the Boost libraries (>= 1.69.0). The tests are run in the tests/multicomponentAlloy/case folder:
$ wmake .. # Build the test executable
$ blockMesh # Create the mesh for the test case
$ ./test_multicomponentAlloy --log_level=all # Run the tests