mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
magneticFoam: New magnetic-field solver
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.
This commit is contained in:
@ -0,0 +1,3 @@
|
||||
magneticFoam.C
|
||||
|
||||
EXE = $(FOAM_APPBIN)/magneticFoam
|
||||
@ -0,0 +1,4 @@
|
||||
EXE_INC = \
|
||||
-I$(LIB_SRC)/finiteVolume/lnInclude
|
||||
|
||||
EXE_LIBS = -lfiniteVolume
|
||||
@ -0,0 +1,81 @@
|
||||
Info<< "Reading field psi\n" << endl;
|
||||
volScalarField psi
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"psi",
|
||||
runTime.timeName(),
|
||||
mesh,
|
||||
IOobject::MUST_READ,
|
||||
IOobject::AUTO_WRITE
|
||||
),
|
||||
mesh
|
||||
);
|
||||
|
||||
Info<< "Reading transportProperties\n" << endl;
|
||||
|
||||
IOdictionary transportProperties
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"transportProperties",
|
||||
runTime.constant(),
|
||||
mesh,
|
||||
IOobject::MUST_READ,
|
||||
IOobject::NO_WRITE
|
||||
)
|
||||
);
|
||||
|
||||
List<magnet> magnets(transportProperties.lookup("magnets"));
|
||||
|
||||
surfaceScalarField murf
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"murf",
|
||||
runTime.timeName(),
|
||||
mesh
|
||||
),
|
||||
mesh,
|
||||
1
|
||||
);
|
||||
|
||||
surfaceScalarField Mrf
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"Mrf",
|
||||
runTime.timeName(),
|
||||
mesh
|
||||
),
|
||||
mesh,
|
||||
dimensionedScalar("Mr", dimensionSet(0, 1, 0, 0, 0, 1, 0), 0)
|
||||
);
|
||||
|
||||
forAll(magnets, i)
|
||||
{
|
||||
label magnetZonei = mesh.faceZones().findZoneID(magnets[i].name());
|
||||
|
||||
if (magnetZonei == -1)
|
||||
{
|
||||
FatalIOErrorIn(args.executable().c_str(), transportProperties)
|
||||
<< "Cannot find faceZone for magnet " << magnets[i].name()
|
||||
<< exit(FatalIOError);
|
||||
}
|
||||
|
||||
const labelList& faces =
|
||||
mesh.faceZones()[magnetZonei];
|
||||
|
||||
const scalar muri = magnets[i].mur();
|
||||
const scalar Mri = magnets[i].Mr().value();
|
||||
const vector& orientationi = magnets[i].orientation();
|
||||
|
||||
const surfaceVectorField& Sf = mesh.Sf();
|
||||
|
||||
forAll(faces, i)
|
||||
{
|
||||
label facei = faces[i];
|
||||
murf[facei] = muri;
|
||||
Mrf[facei] = Mri*(orientationi & Sf[facei]);
|
||||
}
|
||||
}
|
||||
169
applications/solvers/electromagnetics/magneticFoam/magnet.H
Normal file
169
applications/solvers/electromagnetics/magneticFoam/magnet.H
Normal file
@ -0,0 +1,169 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2010-2010 OpenCFD Ltd.
|
||||
\\/ 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/>.
|
||||
|
||||
Class
|
||||
Foam::magnet
|
||||
|
||||
Description
|
||||
Class to hold the defining data for a permanent magnet, in particular
|
||||
the name, relative permeability and remanence.
|
||||
|
||||
SourceFiles
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#ifndef magnet_H
|
||||
#define magnet_H
|
||||
|
||||
#include "dimensionedVector.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
|
||||
// Forward declaration of classes
|
||||
class Istream;
|
||||
class Ostream;
|
||||
|
||||
// Forward declaration of friend functions and operators
|
||||
class magnet;
|
||||
Istream& operator>>(Istream&, magnet&);
|
||||
Ostream& operator<<(Ostream&, const magnet&);
|
||||
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class magnet Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
class magnet
|
||||
{
|
||||
// Private data
|
||||
|
||||
word name_;
|
||||
scalar relativPermeability_;
|
||||
dimensionedScalar remanence_;
|
||||
vector orientation_;
|
||||
|
||||
public:
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Null constructor for lists
|
||||
inline magnet()
|
||||
:
|
||||
remanence_("Mr", dimensionSet(0, -1, 0, 0, 0, 1, 0), 0),
|
||||
orientation_(vector::zero)
|
||||
{}
|
||||
|
||||
//- Construct from components
|
||||
inline magnet
|
||||
(
|
||||
const word& name,
|
||||
const scalar mur,
|
||||
const scalar Mr,
|
||||
const vector& orientation
|
||||
)
|
||||
:
|
||||
name_(name),
|
||||
relativPermeability_(mur),
|
||||
remanence_("Mr", dimensionSet(0, -1, 0, 0, 0, 1, 0), Mr),
|
||||
orientation_(orientation)
|
||||
{}
|
||||
|
||||
//- Construct from Istream
|
||||
inline magnet(Istream& is)
|
||||
:
|
||||
remanence_("Mr", dimensionSet(0, -1, 0, 0, 0, 1, 0), 0),
|
||||
orientation_(vector::zero)
|
||||
{
|
||||
is >> *this;
|
||||
}
|
||||
|
||||
|
||||
// Member Functions
|
||||
|
||||
//- Return name
|
||||
inline const word& name() const
|
||||
{
|
||||
return name_;
|
||||
}
|
||||
|
||||
//- Return relative permeability
|
||||
inline scalar mur() const
|
||||
{
|
||||
return relativPermeability_;
|
||||
}
|
||||
|
||||
//- Return remenance
|
||||
inline const dimensionedScalar& Mr() const
|
||||
{
|
||||
return remanence_;
|
||||
}
|
||||
|
||||
//- Return orientation
|
||||
inline const vector& orientation() const
|
||||
{
|
||||
return orientation_;
|
||||
}
|
||||
|
||||
|
||||
// IOstream operators
|
||||
|
||||
inline friend Istream& operator>>(Istream& is, magnet& m)
|
||||
{
|
||||
is.readBegin("magnet");
|
||||
is >> m.name_
|
||||
>> m.relativPermeability_
|
||||
>> m.remanence_.value()
|
||||
>> m.orientation_;
|
||||
is.readEnd("magnet");
|
||||
|
||||
// Check state of Istream
|
||||
is.check("operator>>(Istream&, magnet&)");
|
||||
|
||||
return is;
|
||||
}
|
||||
|
||||
inline friend Ostream& operator<<(Ostream& os, const magnet& m)
|
||||
{
|
||||
os << token::BEGIN_LIST
|
||||
<< m.name_ << token::SPACE
|
||||
<< m.relativPermeability_ << token::SPACE
|
||||
<< m.remanence_.value()
|
||||
<< m.orientation_
|
||||
<< token::END_LIST;
|
||||
|
||||
return os;
|
||||
}
|
||||
};
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#endif
|
||||
|
||||
// ************************************************************************* //
|
||||
@ -0,0 +1,111 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2010-2010 OpenCFD Ltd.
|
||||
\\/ 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 2 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, write to the Free Software Foundation,
|
||||
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
||||
|
||||
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.
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "fvCFD.H"
|
||||
#include "OSspecific.H"
|
||||
#include "magnet.H"
|
||||
#include "electromagneticConstants.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
int main(int argc, char *argv[])
|
||||
{
|
||||
argList::addBoolOption("noH", "do not write the magnetic field");
|
||||
argList::addBoolOption("noB", "do not write the magnetic field");
|
||||
|
||||
#include "setRootCase.H"
|
||||
#include "createTime.H"
|
||||
#include "createMesh.H"
|
||||
#include "createFields.H"
|
||||
#include "readSIMPLEControls.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
Info<< "Calculating the magnetic field potential" << endl;
|
||||
|
||||
runTime++;
|
||||
|
||||
for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
|
||||
{
|
||||
solve(fvm::laplacian(murf, psi) + fvc::div(murf*Mrf));
|
||||
}
|
||||
|
||||
psi.write();
|
||||
|
||||
if (!args.optionFound("noH"))
|
||||
{
|
||||
Info<< nl
|
||||
<< "Creating field H for time " << runTime.timeName() << endl;
|
||||
volVectorField H
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"H",
|
||||
runTime.timeName(),
|
||||
mesh
|
||||
),
|
||||
fvc::reconstruct(fvc::snGrad(psi)*mesh.magSf())
|
||||
);
|
||||
|
||||
H.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;
|
||||
}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
||||
Reference in New Issue
Block a user