Go to the documentation of this file.
61 int main(
int argc,
char *argv[])
65 "Solver for magnetohydrodynamics (MHD):"
66 " incompressible, laminar flow of a conducting fluid"
67 " under the influence of a magnetic field."
76 #include "createControl.H"
77 #include "createFields.H"
78 #include "initContinuityErrs.H"
88 #include "CourantNo.H"
100 if (
piso.momentumPredictor())
107 while (
piso.correct())
122 while (
piso.correctNonOrthogonal())
130 pEqn.solve(
mesh.solver(
p.select(
piso.finalInnerIter())));
132 if (
piso.finalNonOrthogonalIter())
138 #include "continuityErrs.H"
141 U.correctBoundaryConditions();
146 while (
bpiso.correct())
163 while (
bpiso.correctNonOrthogonal())
170 pBEqn.solve(
mesh.solver(pB.select(
bpiso.finalInnerIter())));
172 if (
bpiso.finalNonOrthogonalIter())
174 phiB -= pBEqn.flux();
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
tmp< volScalarField > rAU
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Ostream & endl(Ostream &os)
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
fvMatrix< vector > fvVectorMatrix
fvMatrix< scalar > fvScalarMatrix
GeometricField< scalar, fvPatchField, volMesh > volScalarField
surfaceScalarField & phiB
Execute application functionObjects to post-process existing results.
pisoControl bpiso(mesh, "BPISO")
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
tmp< surfaceScalarField > flux(const volVectorField &vvf)
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))