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 | Copyright (C) 2011-2015 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 Group
28  grpElectroMagneticsSolvers
29 
30 Description
31  Solver for the magnetic field generated by permanent magnets.
32 
33  A Poisson's equation for the magnetic scalar potential psi is solved
34  from which the magnetic field intensity H and magnetic flux density B
35  are obtained. The paramagnetic particle force field (H dot grad(H))
36  is optionally available.
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #include "fvCFD.H"
41 #include "OSspecific.H"
42 #include "magnet.H"
44 #include "simpleControl.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 int main(int argc, char *argv[])
49 {
50  argList::addBoolOption
51  (
52  "noH",
53  "do not write the magnetic field intensity field"
54  );
55 
56  argList::addBoolOption
57  (
58  "noB",
59  "do not write the magnetic flux density field"
60  );
61 
62  argList::addBoolOption
63  (
64  "HdotGradH",
65  "write the paramagnetic particle force field"
66  );
67 
68  #include "setRootCase.H"
69  #include "createTime.H"
70  #include "createMesh.H"
71 
72  simpleControl simple(mesh);
73 
74  #include "createFields.H"
75 
76  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
77 
78  Info<< "Calculating the magnetic field potential" << endl;
79 
80  runTime++;
81 
82  while (simple.correctNonOrthogonal())
83  {
85  }
86 
87  psi.write();
88 
89  if (!args.optionFound("noH") || args.optionFound("HdotGradH"))
90  {
92  (
93  IOobject
94  (
95  "H",
96  runTime.timeName(),
97  mesh
98  ),
100  );
101 
102  if (!args.optionFound("noH"))
103  {
104  Info<< nl
105  << "Creating field H for time "
106  << runTime.timeName() << endl;
107 
108  H.write();
109  }
110 
111  if (args.optionFound("HdotGradH"))
112  {
113  Info<< nl
114  << "Creating field HdotGradH for time "
115  << runTime.timeName() << endl;
116 
117  volVectorField HdotGradH
118  (
119  IOobject
120  (
121  "HdotGradH",
122  runTime.timeName(),
123  mesh
124  ),
125  H & fvc::grad(H)
126  );
127 
128  HdotGradH.write();
129  }
130  }
131 
132  if (!args.optionFound("noB"))
133  {
134  Info<< nl
135  << "Creating field B for time "
136  << runTime.timeName() << endl;
137 
139  (
140  IOobject
141  (
142  "B",
143  runTime.timeName(),
144  mesh
145  ),
148  );
149 
150  B.write();
151  }
152 
153  Info<< "\nEnd\n" << endl;
154 
155  return 0;
156 }
157 
158 
159 // ************************************************************************* //
Foam::fvc::reconstruct
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcReconstruct.C:54
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
OSspecific.H
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
simple
Simple relative velocity model.
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
solve
rhoEqn solve()
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
simpleControl.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
setRootCase.H
electromagneticConstants.H
psi
const volScalarField & psi
Definition: setRegionFluidFields.H:13
simple
const dictionary & simple
Definition: readFluidMultiRegionSIMPLEControls.H:1
createMesh.H
murf
surfaceScalarField murf(IOobject("murf", runTime.timeName(), mesh), mesh, 1)
createTime.H
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
Foam::constant::electromagnetic::mu0
const dimensionedScalar mu0
Magnetic constant/permeability of free space: default SI units: [H/m].
fvCFD.H
Mrf
surfaceScalarField Mrf(IOobject("Mrf", runTime.timeName(), mesh), mesh, dimensionedScalar("Mr", dimensionSet(0, 1, 0, 0, 0, 1, 0), 0))
args
Foam::argList args(argc, argv)
magnet.H