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