1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156
|
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2015 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"
simpleControl simple(mesh);
#include "createFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "Calculating the magnetic field potential" << endl;
runTime++;
while (simple.correctNonOrthogonal())
{
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;
}
// ************************************************************************* //
|