mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
Merge branch 'master' of /home/noisy3/OpenFOAM/OpenFOAM-dev
This commit is contained in:
@ -60,7 +60,7 @@ class magnet
|
||||
// Private data
|
||||
|
||||
word name_;
|
||||
scalar relativPermeability_;
|
||||
scalar relativePermeability_;
|
||||
dimensionedScalar remanence_;
|
||||
vector orientation_;
|
||||
|
||||
@ -85,7 +85,7 @@ public:
|
||||
)
|
||||
:
|
||||
name_(name),
|
||||
relativPermeability_(mur),
|
||||
relativePermeability_(mur),
|
||||
remanence_("Mr", dimensionSet(0, -1, 0, 0, 0, 1, 0), Mr),
|
||||
orientation_(orientation)
|
||||
{}
|
||||
@ -111,7 +111,7 @@ public:
|
||||
//- Return relative permeability
|
||||
inline scalar mur() const
|
||||
{
|
||||
return relativPermeability_;
|
||||
return relativePermeability_;
|
||||
}
|
||||
|
||||
//- Return remenance
|
||||
@ -133,7 +133,7 @@ public:
|
||||
{
|
||||
is.readBegin("magnet");
|
||||
is >> m.name_
|
||||
>> m.relativPermeability_
|
||||
>> m.relativePermeability_
|
||||
>> m.remanence_.value()
|
||||
>> m.orientation_;
|
||||
is.readEnd("magnet");
|
||||
@ -148,7 +148,7 @@ public:
|
||||
{
|
||||
os << token::BEGIN_LIST
|
||||
<< m.name_ << token::SPACE
|
||||
<< m.relativPermeability_ << token::SPACE
|
||||
<< m.relativePermeability_ << token::SPACE
|
||||
<< m.remanence_.value()
|
||||
<< m.orientation_
|
||||
<< token::END_LIST;
|
||||
|
||||
@ -29,7 +29,8 @@ Description
|
||||
|
||||
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.
|
||||
are obtained. The paramagnetic particle force field (H dot grad(H))
|
||||
is optionally available.
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
@ -42,8 +43,23 @@ Description
|
||||
|
||||
int main(int argc, char *argv[])
|
||||
{
|
||||
argList::addBoolOption("noH", "do not write the magnetic field");
|
||||
argList::addBoolOption("noB", "do not write the magnetic field");
|
||||
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"
|
||||
@ -64,10 +80,8 @@ int main(int argc, char *argv[])
|
||||
|
||||
psi.write();
|
||||
|
||||
if (!args.optionFound("noH"))
|
||||
if (!args.optionFound("noH") || args.optionFound("HdotGradH"))
|
||||
{
|
||||
Info<< nl
|
||||
<< "Creating field H for time " << runTime.timeName() << endl;
|
||||
volVectorField H
|
||||
(
|
||||
IOobject
|
||||
@ -79,13 +93,42 @@ int main(int argc, char *argv[])
|
||||
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;
|
||||
<< "Creating field B for time "
|
||||
<< runTime.timeName() << endl;
|
||||
|
||||
volVectorField B
|
||||
(
|
||||
IOobject
|
||||
|
||||
Reference in New Issue
Block a user