168 lines
4.9 KiB
C
168 lines
4.9 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | foam-extend: Open Source CFD
|
|
\\ / O peration | Version: 4.0
|
|
\\ / A nd | Web: http://www.foam-extend.org
|
|
\\/ M anipulation | For copyright notice see file Copyright
|
|
-------------------------------------------------------------------------------
|
|
License
|
|
This file is part of foam-extend.
|
|
|
|
foam-extend is free software: you can redistribute it and/or modify it
|
|
under the terms of the GNU General Public License as published by the
|
|
Free Software Foundation, either version 3 of the License, or (at your
|
|
option) any later version.
|
|
|
|
foam-extend is distributed in the hope that it will be useful, but
|
|
WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
|
|
General Public License for more details.
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
Application
|
|
mhdFoam
|
|
|
|
Description
|
|
Solver for magnetohydrodynamics (MHD): incompressible, laminar flow of a
|
|
conducting fluid under the influence of a magnetic field.
|
|
|
|
An applied magnetic field H acts as a driving force,
|
|
at present boundary conditions cannot be set via the
|
|
electric field E or current density J. The fluid viscosity nu,
|
|
conductivity sigma and permeability mu are read in as uniform
|
|
constants.
|
|
|
|
A fictitous magnetic flux pressure pH is introduced in order to
|
|
compensate for discretisation errors and create a magnetic face flux
|
|
field which is divergence free as required by Maxwell's equations.
|
|
|
|
However, in this formulation discretisation error prevents the normal
|
|
stresses in UB from cancelling with those from BU, but it is unknown
|
|
whether this is a serious error. A correction could be introduced
|
|
whereby the normal stresses in the discretised BU term are replaced
|
|
by those from the UB term, but this would violate the boundedness
|
|
constraint presently observed in the present numerics which
|
|
guarantees div(U) and div(H) are zero.
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "fvCFD.H"
|
|
#include "OSspecific.H"
|
|
#include "pisoControl.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
int main(int argc, char *argv[])
|
|
{
|
|
# include "setRootCase.H"
|
|
|
|
# include "createTime.H"
|
|
# include "createMesh.H"
|
|
|
|
pisoControl piso(mesh);
|
|
|
|
# include "createFields.H"
|
|
# include "initContinuityErrs.H"
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
Info<< nl << "Starting time loop" << endl;
|
|
|
|
while (runTime.loop())
|
|
{
|
|
# include "readBPISOControls.H"
|
|
|
|
Info<< "Time = " << runTime.timeName() << nl << endl;
|
|
|
|
# include "CourantNo.H"
|
|
|
|
{
|
|
fvVectorMatrix UEqn
|
|
(
|
|
fvm::ddt(U)
|
|
+ fvm::div(phi, U)
|
|
- fvc::div(phiB, 2.0*DBU*B)
|
|
- fvm::laplacian(nu, U)
|
|
+ fvc::grad(DBU*magSqr(B))
|
|
);
|
|
|
|
solve(UEqn == -fvc::grad(p));
|
|
|
|
|
|
// --- PISO loop
|
|
|
|
while (piso.correct())
|
|
{
|
|
volScalarField rUA = 1.0/UEqn.A();
|
|
|
|
U = rUA*UEqn.H();
|
|
|
|
phi = (fvc::interpolate(U) & mesh.Sf())
|
|
+ fvc::ddtPhiCorr(rUA, U, phi);
|
|
|
|
while (piso.correctNonOrthogonal())
|
|
{
|
|
fvScalarMatrix pEqn
|
|
(
|
|
fvm::laplacian(rUA, p) == fvc::div(phi)
|
|
);
|
|
|
|
pEqn.setReference(pRefCell, pRefValue);
|
|
pEqn.solve();
|
|
|
|
if (piso.finalNonOrthogonalIter())
|
|
{
|
|
phi -= pEqn.flux();
|
|
}
|
|
}
|
|
|
|
# include "continuityErrs.H"
|
|
|
|
U -= rUA*fvc::grad(p);
|
|
U.correctBoundaryConditions();
|
|
}
|
|
}
|
|
|
|
// --- B-PISO loop
|
|
|
|
for (int Bcorr=0; Bcorr<nBcorr; Bcorr++)
|
|
{
|
|
fvVectorMatrix BEqn
|
|
(
|
|
fvm::ddt(B)
|
|
+ fvm::div(phi, B)
|
|
- fvc::div(phiB, U)
|
|
- fvm::laplacian(DB, B)
|
|
);
|
|
|
|
BEqn.solve();
|
|
|
|
volScalarField rBA = 1.0/BEqn.A();
|
|
|
|
phiB = (fvc::interpolate(B) & mesh.Sf())
|
|
+ fvc::ddtPhiCorr(rBA, B, phiB);
|
|
|
|
fvScalarMatrix pBEqn
|
|
(
|
|
fvm::laplacian(rBA, pB) == fvc::div(phiB)
|
|
);
|
|
pBEqn.solve();
|
|
|
|
phiB -= pBEqn.flux();
|
|
|
|
# include "magneticFieldErr.H"
|
|
}
|
|
|
|
runTime.write();
|
|
}
|
|
|
|
Info<< "End\n" << endl;
|
|
|
|
return 0;
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|