MFEM  v3.1
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  double (*muInv)(const Vector&),
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 Update();
57 
58  void Solve();
59 
60  void GetErrorEstimates(Vector & errors);
61 
63 
64  void WriteVisItFields(int it = 0);
65 
66  void InitializeGLVis();
67 
68  void DisplayToGLVis();
69 
70  const ParGridFunction & GetVectorPotential() { return *a_; }
71 
72 private:
73 
74  int myid_;
75  int num_procs_;
76  int order_;
77 
78  ParMesh * pmesh_;
79 
80  VisItDataCollection * visit_dc_;
81 
82  H1_ParFESpace * H1FESpace_;
83  ND_ParFESpace * HCurlFESpace_;
84  RT_ParFESpace * HDivFESpace_;
85 
86  ParBilinearForm * curlMuInvCurl_;
87  ParBilinearForm * hCurlMass_;
88  ParBilinearForm * hDivMassMuInv_;
89  ParMixedBilinearForm * hDivHCurlMuInv_;
90 
93 
94  ParGridFunction * a_;
95  ParGridFunction * b_;
96  ParGridFunction * h_;
97  ParGridFunction * j_;
98  ParGridFunction * k_;
99  ParGridFunction * m_;
100 
101  DivergenceFreeProjector * DivFreeProj_;
102  SurfaceCurrent * SurfCur_;
103 
104  Coefficient * muInvCoef_; // Dia/Paramagnetic Material Coefficient
105  VectorCoefficient * aBCCoef_; // Vector Potential BC Function
106  VectorCoefficient * jCoef_; // Volume Current Density Function
107  VectorCoefficient * mCoef_; // Magnetization Vector Function
108 
109  double (*muInv_)(const Vector&);
110  void (*a_bc_ )(const Vector&, Vector&);
111  void (*j_src_)(const Vector&, Vector&);
112  void (*m_src_)(const Vector&, Vector&);
113 
114  Array<int> ess_bdr_;
115  Array<int> non_k_bdr_;
116 
117  std::map<std::string,socketstream*> socks_;
118 };
119 
121 {
122 public:
124  ParFiniteElementSpace & HCurlFESpace,
126  Array<int> & kbcs, Array<int> & vbcs, Vector & vbcv);
127  ~SurfaceCurrent();
128 
130 
131  void Update();
132 
133  ParGridFunction * GetPsi() { return psi_; }
134 
135 private:
136  int myid_;
137 
138  ParFiniteElementSpace * H1FESpace_;
139  ParFiniteElementSpace * HCurlFESpace_;
140  ParDiscreteGradOperator * Grad_;
141  Array<int> * kbcs_;
142  Array<int> * vbcs_;
143  Vector * vbcv_;
144 
145  ParBilinearForm * s0_;
146  ParGridFunction * psi_;
147 
148  HypreBoomerAMG * amg_;
149  HyprePCG * pcg_;
150  HypreParMatrix * S0_;
151  HypreParVector * PSI_;
152  HypreParVector * RHS_;
153  HypreParVector * K_;
154 
155  Array<int> ess_bdr_;
156  Array<int> non_k_bdr_;
157 };
158 
159 } // namespace electromagnetics
160 
161 } // namespace mfem
162 
163 #endif // MFEM_USE_MPI
164 
165 #endif // MFEM_TESLA_SOLVER
double muInv(const Vector &x)
Definition: tesla.cpp:70
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void GetErrorEstimates(Vector &errors)
TeslaSolver(ParMesh &pmesh, int order, Array< int > &kbcs, Array< int > &vbcs, Vector &vbcv, double(*muInv)(const Vector &), void(*a_bc)(const Vector &, Vector &), void(*j_src)(const Vector &, Vector &), void(*m_src)(const Vector &, Vector &))
The BoomerAMG solver in hypre.
Definition: hypre.hpp:698
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:558
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:143
Class for parallel meshes.
Definition: pmesh.hpp:28