magneticFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2015 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Application
27  magneticFoam
28 
29 Group
30  grpElectroMagneticsSolvers
31 
32 Description
33  Solver for the magnetic field generated by permanent magnets.
34 
35  A Poisson's equation for the magnetic scalar potential psi is solved
36  from which the magnetic field intensity H and magnetic flux density B
37  are obtained. The paramagnetic particle force field (H dot grad(H))
38  is optionally available.
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #include "fvCFD.H"
43 #include "OSspecific.H"
44 #include "magnet.H"
46 #include "simpleControl.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 int main(int argc, char *argv[])
51 {
52  argList::addNote
53  (
54  "Solver for the magnetic field generated by permanent magnets."
55  );
56 
57  argList::addBoolOption
58  (
59  "noH",
60  "Do not write the magnetic field intensity field"
61  );
62 
63  argList::addBoolOption
64  (
65  "noB",
66  "Do not write the magnetic flux density field"
67  );
68 
69  argList::addBoolOption
70  (
71  "HdotGradH",
72  "Write the paramagnetic particle force field"
73  );
74 
75  #include "addCheckCaseOptions.H"
76  #include "setRootCaseLists.H"
77  #include "createTime.H"
78  #include "createMesh.H"
79 
80  simpleControl simple(mesh);
81 
82  #include "createFields.H"
83 
84  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
85 
86  Info<< "Calculating the magnetic field potential" << endl;
87 
88  ++runTime;
89 
90  while (simple.correctNonOrthogonal())
91  {
93  }
94 
95  psi.write();
96 
97  if (!args.found("noH") || args.found("HdotGradH"))
98  {
100  (
101  IOobject
102  (
103  "H",
104  runTime.timeName(),
105  mesh
106  ),
108  );
109 
110  if (!args.found("noH"))
111  {
112  Info<< nl
113  << "Creating field H for time "
114  << runTime.timeName() << endl;
115 
116  H.write();
117  }
118 
119  if (args.found("HdotGradH"))
120  {
121  Info<< nl
122  << "Creating field HdotGradH for time "
123  << runTime.timeName() << endl;
124 
125  volVectorField HdotGradH
126  (
127  IOobject
128  (
129  "HdotGradH",
130  runTime.timeName(),
131  mesh
132  ),
133  H & fvc::grad(H)
134  );
135 
136  HdotGradH.write();
137  }
138  }
139 
140  if (!args.found("noB"))
141  {
142  Info<< nl
143  << "Creating field B for time "
144  << runTime.timeName() << endl;
145 
147  (
148  IOobject
149  (
150  "B",
151  runTime.timeName(),
152  mesh
153  ),
156  );
157 
158  B.write();
159  }
160 
161  Info<< "\nEnd\n" << endl;
162 
163  return 0;
164 }
165 
166 
167 // ************************************************************************* //
const dimensionedScalar mu0
Magnetic constant/permeability of free space: default SI units: [H/m].
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:43
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:85
surfaceScalarField murf(IOobject("murf", runTime.timeName(), mesh), mesh, dimensionedScalar("one", dimless, 1.0))
CEqn solve()
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
dynamicFvMesh & mesh
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:40
surfaceScalarField Mrf(IOobject("Mrf", runTime.timeName(), mesh), mesh, dimensionedScalar(dimensionSet(0, 1, 0, 0, 0, 1, 0), Zero))
Required Variables.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
const dictionary & simple
messageStream Info
Information stream (stdout output on master, null elsewhere)
const volScalarField & psi
virtual bool write(const bool valid=true) const
Write using setting from DB.
Foam::argList args(argc, argv)
Required Classes.
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:49