MFEM  v3.4
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 "../common/pfem_extras.hpp"
16 #include "electromagnetics.hpp"
17 
18 #ifdef MFEM_USE_MPI
19 
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 class SurfaceCurrent;
38 {
39 public:
40  TeslaSolver(ParMesh & pmesh, int order, Array<int> & kbcs,
41  Array<int> & vbcs, Vector & vbcv,
42  Coefficient & muInvCoef,
43  void (*a_bc )(const Vector&, Vector&),
44  void (*j_src)(const Vector&, Vector&),
45  void (*m_src)(const Vector&, Vector&));
46  ~TeslaSolver();
47 
48  HYPRE_Int GetProblemSize();
49 
50  void PrintSizes();
51 
52  void Assemble();
53 
54  void Update();
55 
56  void Solve();
57 
58  void GetErrorEstimates(Vector & errors);
59 
61 
62  void WriteVisItFields(int it = 0);
63 
64  void InitializeGLVis();
65 
66  void DisplayToGLVis();
67 
68  const ParGridFunction & GetVectorPotential() { return *a_; }
69 
70 private:
71 
72  int myid_;
73  int num_procs_;
74  int order_;
75 
76  ParMesh * pmesh_;
77 
78  VisItDataCollection * visit_dc_;
79 
80  H1_ParFESpace * H1FESpace_;
81  ND_ParFESpace * HCurlFESpace_;
82  RT_ParFESpace * HDivFESpace_;
83 
84  ParBilinearForm * curlMuInvCurl_;
85  ParBilinearForm * hCurlMass_;
86  ParMixedBilinearForm * hDivHCurlMuInv_;
87  ParMixedBilinearForm * weakCurlMuInv_;
88 
91 
92  ParGridFunction * a_; // Vector Potential (HCurl)
93  ParGridFunction * b_; // Magnetic Flux (HDiv)
94  ParGridFunction * h_; // Magnetic Field (HCurl)
95  ParGridFunction * jr_; // Raw Volumetric Current Density (HCurl)
96  ParGridFunction * j_; // Volumetric Current Density (HCurl)
97  ParGridFunction * k_; // Surface Current Density (HCurl)
98  ParGridFunction * m_; // Magnetization (HDiv)
99  ParGridFunction * bd_; // Dual of B (HCurl)
100  ParGridFunction * jd_; // Dual of J, the rhs vector (HCurl)
101 
102  DivergenceFreeProjector * DivFreeProj_;
103  SurfaceCurrent * SurfCur_;
104 
105  Coefficient * muInvCoef_; // Dia/Paramagnetic Material Coefficient
106  VectorCoefficient * aBCCoef_; // Vector Potential BC Function
107  VectorCoefficient * jCoef_; // Volume Current Density Function
108  VectorCoefficient * mCoef_; // Magnetization Vector Function
109 
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> ess_bdr_tdofs_;
116  Array<int> non_k_bdr_;
117 
118  std::map<std::string,socketstream*> socks_;
119 };
120 
122 {
123 public:
126  Array<int> & kbcs, Array<int> & vbcs, Vector & vbcv);
127  ~SurfaceCurrent();
128 
129  void InitSolver() const;
130 
132 
133  void Update();
134 
135  ParGridFunction * GetPsi() { return psi_; }
136 
137 private:
138  int myid_;
139 
140  ParFiniteElementSpace * H1FESpace_;
141  ParDiscreteGradOperator * grad_;
142  Array<int> * kbcs_;
143  Array<int> * vbcs_;
144  Vector * vbcv_;
145 
146  ParBilinearForm * s0_;
147  ParGridFunction * psi_;
148  ParGridFunction * rhs_;
149 
150  HypreParMatrix * S0_;
151  mutable Vector Psi_;
152  mutable Vector RHS_;
153 
154  mutable HypreBoomerAMG * amg_;
155  mutable HyprePCG * pcg_;
156 
157  Array<int> ess_bdr_, ess_bdr_tdofs_;
158  Array<int> non_k_bdr_;
159 };
160 
161 } // namespace electromagnetics
162 
163 } // namespace mfem
164 
165 #endif // MFEM_USE_MPI
166 
167 #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:796
SurfaceCurrent(ParFiniteElementSpace &H1FESpace, ParDiscreteGradOperator &Grad, Array< int > &kbcs, Array< int > &vbcs, Vector &vbcv)
Data collection with VisIt I/O routines.
void RegisterVisItFields(VisItDataCollection &visit_dc)
void j_src(const Vector &x, double t, Vector &j)
Definition: maxwell.cpp:108
PCG solver in hypre.
Definition: hypre.hpp:656
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
Class for parallel bilinear form using different test and trial FE spaces.
const ParGridFunction & GetVectorPotential()
Class for parallel bilinear form.
void ComputeSurfaceCurrent(ParGridFunction &k)
Vector data type.
Definition: vector.hpp:48
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:175
Class for parallel meshes.
Definition: pmesh.hpp:32
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 &))