All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
mhdFoam.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-2021 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  mhdFoam
26 
27 Description
28  Solver for magnetohydrodynamics (MHD): incompressible, laminar flow of a
29  conducting fluid under the influence of a magnetic field.
30 
31  An applied magnetic field H acts as a driving force,
32  at present boundary conditions cannot be set via the
33  electric field E or current density J. The fluid viscosity nu,
34  conductivity sigma and permeability mu are read in as uniform
35  constants.
36 
37  A fictitious magnetic flux pressure pB is introduced in order to
38  compensate for discretisation errors and create a magnetic face flux
39  field which is divergence free as required by Maxwell's equations.
40 
41  However, in this formulation discretisation error prevents the normal
42  stresses in UB from cancelling with those from BU, but it is unknown
43  whether this is a serious error. A correction could be introduced
44  whereby the normal stresses in the discretised BU term are replaced
45  by those from the UB term, but this would violate the boundedness
46  constraint presently observed in the present numerics which
47  guarantees div(U) and div(H) are zero.
48 
49 \*---------------------------------------------------------------------------*/
50 
51 #include "fvCFD.H"
52 #include "pisoControl.H"
53 #include "pressureReference.H"
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 int main(int argc, char *argv[])
58 {
59  #include "postProcess.H"
60 
61  #include "setRootCaseLists.H"
62  #include "createTime.H"
63  #include "createMesh.H"
64  #include "createControl.H"
65  #include "createFields.H"
66  #include "initContinuityErrs.H"
67 
68  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70  Info<< nl << "Starting time loop" << endl;
71 
72  while (runTime.loop())
73  {
74  Info<< "Time = " << runTime.userTimeName() << nl << endl;
75 
76  #include "CourantNo.H"
77 
78  {
80  (
81  fvm::ddt(U)
82  + fvm::div(phi, U)
83  - fvc::div(phiB, 2.0*DBU*B)
84  - fvm::laplacian(nu, U)
85  + fvc::grad(DBU*magSqr(B))
86  );
87 
88  if (piso.momentumPredictor())
89  {
90  solve(UEqn == -fvc::grad(p));
91  }
92 
93 
94  // --- PISO loop
95  while (piso.correct())
96  {
97  volScalarField rAU(1.0/UEqn.A());
101  (
102  "phiHbyA",
103  fvc::flux(HbyA)
104  + rAUf*fvc::ddtCorr(U, phi)
105  );
106 
107  // Update the pressure BCs to ensure flux consistency
109 
110  while (piso.correctNonOrthogonal())
111  {
112  fvScalarMatrix pEqn
113  (
115  );
116 
117  pEqn.setReference
118  (
119  pressureReference.refCell(),
120  pressureReference.refValue()
121  );
122 
123  pEqn.solve();
124 
125  if (piso.finalNonOrthogonalIter())
126  {
127  phi = phiHbyA - pEqn.flux();
128  }
129  }
130 
131  #include "continuityErrs.H"
132 
133  U = HbyA - rAU*fvc::grad(p);
134  U.correctBoundaryConditions();
135  }
136  }
137 
138  // --- B-PISO loop
139  while (bpiso.correct())
140  {
141  fvVectorMatrix BEqn
142  (
143  fvm::ddt(B)
144  + fvm::div(phi, B)
145  - fvc::div(phiB, U)
146  - fvm::laplacian(DB, B)
147  );
148 
149  BEqn.solve();
150 
151  volScalarField rAB(1.0/BEqn.A());
152  surfaceScalarField rABf("rABf", fvc::interpolate(rAB));
153 
154  phiB = fvc::flux(B);
155 
156  while (bpiso.correctNonOrthogonal())
157  {
158  fvScalarMatrix pBEqn
159  (
160  fvm::laplacian(rABf, pB) == fvc::div(phiB)
161  );
162 
163  pBEqn.solve();
164 
165  if (bpiso.finalNonOrthogonalIter())
166  {
167  phiB -= pBEqn.flux();
168  }
169  }
170 
171  #include "magneticFieldErr.H"
172  }
173 
174  runTime.write();
175  }
176 
177  Info<< "End\n" << endl;
178 
179  return 0;
180 }
181 
182 
183 // ************************************************************************* //
pressureReference & pressureReference
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
U
Definition: pEqn.H:72
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:170
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
pisoControl piso(mesh)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
dimensionedScalar rAUf("rAUf", dimTime, 1.0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
phiHbyA
Definition: pEqn.H:30
phi
Definition: correctPhi.H:3
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
pisoControl bpiso(mesh, "BPISO")
dimensioned< scalar > magSqr(const dimensioned< Type > &)
static const char nl
Definition: Ostream.H:260
surfaceScalarField & phiB
Definition: createPhiB.H:47
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
rhoEqn solve()
tmp< volScalarField > rAU
volVectorField & HbyA
Definition: pEqn.H:13
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
messageStream Info
volScalarField & p
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Execute application functionObjects to post-process existing results.
fvVectorMatrix & UEqn
Definition: UEqn.H:13