156 lines
4.1 KiB
C
156 lines
4.1 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
|
|
\\/ M anipulation |
|
|
-------------------------------------------------------------------------------
|
|
License
|
|
This file is part of OpenFOAM.
|
|
|
|
OpenFOAM 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.
|
|
|
|
OpenFOAM 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 OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
Application
|
|
magneticFoam
|
|
|
|
Description
|
|
Solver for the magnetic field generated by permanent magnets.
|
|
|
|
A Poisson's equation for the magnetic scalar potential psi is solved
|
|
from which the magnetic field intensity H and magnetic flux density B
|
|
are obtained. The paramagnetic particle force field (H dot grad(H))
|
|
is optionally available.
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "fvCFD.H"
|
|
#include "OSspecific.H"
|
|
#include "magnet.H"
|
|
#include "electromagneticConstants.H"
|
|
#include "simpleControl.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
int main(int argc, char *argv[])
|
|
{
|
|
argList::addBoolOption
|
|
(
|
|
"noH",
|
|
"do not write the magnetic field intensity field"
|
|
);
|
|
|
|
argList::addBoolOption
|
|
(
|
|
"noB",
|
|
"do not write the magnetic flux density field"
|
|
);
|
|
|
|
argList::addBoolOption
|
|
(
|
|
"HdotGradH",
|
|
"write the paramagnetic particle force field"
|
|
);
|
|
|
|
#include "setRootCase.H"
|
|
#include "createTime.H"
|
|
#include "createMesh.H"
|
|
#include "createFields.H"
|
|
|
|
simpleControl simple(mesh);
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
Info<< "Calculating the magnetic field potential" << endl;
|
|
|
|
runTime++;
|
|
|
|
for (int nonOrth=0; nonOrth<=simple.nNonOrthCorr(); nonOrth++)
|
|
{
|
|
solve(fvm::laplacian(murf, psi) + fvc::div(murf*Mrf));
|
|
}
|
|
|
|
psi.write();
|
|
|
|
if (!args.optionFound("noH") || args.optionFound("HdotGradH"))
|
|
{
|
|
volVectorField H
|
|
(
|
|
IOobject
|
|
(
|
|
"H",
|
|
runTime.timeName(),
|
|
mesh
|
|
),
|
|
fvc::reconstruct(fvc::snGrad(psi)*mesh.magSf())
|
|
);
|
|
|
|
if (!args.optionFound("noH"))
|
|
{
|
|
Info<< nl
|
|
<< "Creating field H for time "
|
|
<< runTime.timeName() << endl;
|
|
|
|
H.write();
|
|
}
|
|
|
|
if (args.optionFound("HdotGradH"))
|
|
{
|
|
Info<< nl
|
|
<< "Creating field HdotGradH for time "
|
|
<< runTime.timeName() << endl;
|
|
|
|
volVectorField HdotGradH
|
|
(
|
|
IOobject
|
|
(
|
|
"HdotGradH",
|
|
runTime.timeName(),
|
|
mesh
|
|
),
|
|
H & fvc::grad(H)
|
|
);
|
|
|
|
HdotGradH.write();
|
|
}
|
|
}
|
|
|
|
if (!args.optionFound("noB"))
|
|
{
|
|
Info<< nl
|
|
<< "Creating field B for time "
|
|
<< runTime.timeName() << endl;
|
|
|
|
volVectorField B
|
|
(
|
|
IOobject
|
|
(
|
|
"B",
|
|
runTime.timeName(),
|
|
mesh
|
|
),
|
|
constant::electromagnetic::mu0
|
|
*fvc::reconstruct(murf*fvc::snGrad(psi)*mesh.magSf() + murf*Mrf)
|
|
);
|
|
|
|
B.write();
|
|
}
|
|
|
|
Info<< "\nEnd\n" << endl;
|
|
|
|
return 0;
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|