MFEM  v3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
tesla_solver.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_TESLA_SOLVER
13 #define MFEM_TESLA_SOLVER
14 
15 #include "../../config/config.hpp"
16 
17 #ifdef MFEM_USE_MPI
18 
19 #include "../common/pfem_extras.hpp"
20 #include <string>
21 #include <map>
22 
23 namespace mfem
24 {
25 
26 using miniapps::H1_ParFESpace;
27 using miniapps::ND_ParFESpace;
28 using miniapps::RT_ParFESpace;
29 using miniapps::ParDiscreteGradOperator;
30 using miniapps::ParDiscreteCurlOperator;
31 using miniapps::DivergenceFreeProjector;
32 
33 namespace electromagnetics
34 {
35 
36 // Physical Constants
37 // Permeability of Free Space (units H/m)
38 static double mu0_ = 4.0e-7*M_PI;
39 
40 class SurfaceCurrent;
42 {
43 public:
44  TeslaSolver(ParMesh & pmesh, int order, Array<int> & kbcs,
45  Array<int> & vbcs, Vector & vbcv,
46  Coefficient & muInvCoef,
47  void (*a_bc )(const Vector&, Vector&),
48  void (*j_src)(const Vector&, Vector&),
49  void (*m_src)(const Vector&, Vector&));
50  ~TeslaSolver();
51 
52  HYPRE_Int GetProblemSize();
53 
54  void PrintSizes();
55 
56  void Assemble();
57 
58  void Update();
59 
60  void Solve();
61 
62  void GetErrorEstimates(Vector & errors);
63 
65 
66  void WriteVisItFields(int it = 0);
67 
68  void InitializeGLVis();
69 
70  void DisplayToGLVis();
71 
72  const ParGridFunction & GetVectorPotential() { return *a_; }
73 
74 private:
75 
76  int myid_;
77  int num_procs_;
78  int order_;
79 
80  ParMesh * pmesh_;
81 
82  VisItDataCollection * visit_dc_;
83 
84  H1_ParFESpace * H1FESpace_;
85  ND_ParFESpace * HCurlFESpace_;
86  RT_ParFESpace * HDivFESpace_;
87 
88  ParBilinearForm * curlMuInvCurl_;
89  ParBilinearForm * hCurlMass_;
90  ParBilinearForm * hDivMassMuInv_;
91  ParMixedBilinearForm * hDivHCurlMuInv_;
92 
95 
96  ParGridFunction * a_;
97  ParGridFunction * b_;
98  ParGridFunction * h_;
99  ParGridFunction * j_;
100  ParGridFunction * k_;
101  ParGridFunction * m_;
102 
103  DivergenceFreeProjector * DivFreeProj_;
104  SurfaceCurrent * SurfCur_;
105 
106  Coefficient * muInvCoef_; // Dia/Paramagnetic Material Coefficient
107  VectorCoefficient * aBCCoef_; // Vector Potential BC Function
108  VectorCoefficient * jCoef_; // Volume Current Density Function
109  VectorCoefficient * mCoef_; // Magnetization Vector Function
110 
111  void (*a_bc_ )(const Vector&, Vector&);
112  void (*j_src_)(const Vector&, Vector&);
113  void (*m_src_)(const Vector&, Vector&);
114 
115  Array<int> ess_bdr_;
116  Array<int> non_k_bdr_;
117 
118  std::map<std::string,socketstream*> socks_;
119 };
120 
122 {
123 public:
125  ParFiniteElementSpace & HCurlFESpace,
127  Array<int> & kbcs, Array<int> & vbcs, Vector & vbcv);
128  ~SurfaceCurrent();
129 
131 
132  void Update();
133 
134  ParGridFunction * GetPsi() { return psi_; }
135 
136 private:
137  int myid_;
138 
139  ParFiniteElementSpace * H1FESpace_;
140  ParFiniteElementSpace * HCurlFESpace_;
141  ParDiscreteGradOperator * Grad_;
142  Array<int> * kbcs_;
143  Array<int> * vbcs_;
144  Vector * vbcv_;
145 
146  ParBilinearForm * s0_;
147  ParGridFunction * psi_;
148 
149  HypreBoomerAMG * amg_;
150  HyprePCG * pcg_;
151  HypreParMatrix * S0_;
152  HypreParVector * PSI_;
153  HypreParVector * RHS_;
154  HypreParVector * K_;
155 
156  Array<int> ess_bdr_;
157  Array<int> non_k_bdr_;
158 };
159 
160 } // namespace electromagnetics
161 
162 } // namespace mfem
163 
164 #endif // MFEM_USE_MPI
165 
166 #endif // MFEM_TESLA_SOLVER
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void GetErrorEstimates(Vector &errors)
The BoomerAMG solver in hypre.
Definition: hypre.hpp:705
Data collection with VisIt I/O routines.
void RegisterVisItFields(VisItDataCollection &visit_dc)
SurfaceCurrent(ParFiniteElementSpace &H1FESpace, ParFiniteElementSpace &HCurlFESpace, ParDiscreteGradOperator &Grad, Array< int > &kbcs, Array< int > &vbcs, Vector &vbcv)
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:58
PCG solver in hypre.
Definition: hypre.hpp:565
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
Class for parallel bilinear form.
const ParGridFunction & GetVectorPotential()
Class for parallel bilinear form.
void ComputeSurfaceCurrent(ParGridFunction &k)
Vector data type.
Definition: vector.hpp:33
Class for parallel grid function.
Definition: pgridfunc.hpp:31
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:150
Class for parallel meshes.
Definition: pmesh.hpp:28
TeslaSolver(ParMesh &pmesh, int order, Array< int > &kbcs, Array< int > &vbcs, Vector &vbcv, Coefficient &muInvCoef, void(*a_bc)(const Vector &, Vector &), void(*j_src)(const Vector &, Vector &), void(*m_src)(const Vector &, Vector &))