MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
pfem_extras.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_PFEM_EXTRAS
13 #define MFEM_PFEM_EXTRAS
14 
15 #include "mfem.hpp"
16 
17 #ifdef MFEM_USE_MPI
18 
19 #include <cstddef>
20 
21 namespace mfem
22 {
23 
24 namespace common
25 {
26 
27 /** The H1_ParFESpace class is a ParFiniteElementSpace which automatically
28  allocates and destroys its own FiniteElementCollection, in this
29  case an H1_FECollection object.
30 */
32 {
33 public:
35  const int p, const int space_dim = 3,
36  const int type = BasisType::GaussLobatto,
37  int vdim = 1, int order = Ordering::byNODES);
39 private:
40  const FiniteElementCollection *FEC_;
41 };
42 
43 /** The ND_ParFESpace class is a ParFiniteElementSpace which automatically
44  allocates and destroys its own FiniteElementCollection, in this
45  case an ND_FECollection object.
46 */
48 {
49 public:
50  ND_ParFESpace(ParMesh *m, const int p, const int space_dim,
51  int vdim = 1, int order = Ordering::byNODES);
53 private:
54  const FiniteElementCollection *FEC_;
55 };
56 
57 /** The RT_ParFESpace class is a ParFiniteElementSpace which automatically
58  allocates and destroys its own FiniteElementCollection, in this
59  case an RT_FECollection object.
60 */
62 {
63 public:
64  RT_ParFESpace(ParMesh *m, const int p, const int space_dim,
65  int vdim = 1, int order = Ordering::byNODES);
67 private:
68  const FiniteElementCollection *FEC_;
69 };
70 
71 /** The L2_ParFESpace class is a ParFiniteElementSpace which automatically
72  allocates and destroys its own FiniteElementCollection, in this
73  case an L2_FECollection object.
74 */
76 {
77 public:
78  L2_ParFESpace(ParMesh *m, const int p, const int space_dim,
79  int vdim = 1, int order = Ordering::byNODES);
81 private:
82  const FiniteElementCollection *FEC_;
83 };
84 
86 {
87 public:
90  : ParDiscreteLinearOperator(dfes, rfes) {}
92 };
93 
95 {
96 public:
98  ParFiniteElementSpace *rfes);
99 };
100 
102 {
103 public:
105  ParFiniteElementSpace *rfes);
106 };
107 
109 {
110 public:
112  ParFiniteElementSpace *rfes);
113 };
114 
115 /// This class computes the irrotational portion of a vector field.
116 /// This vector field must be discretized using Nedelec basis
117 /// functions.
119 {
120 public:
122  ParFiniteElementSpace & HCurlFESpace,
123  const int & irOrder,
124  ParBilinearForm * s0 = NULL,
125  ParMixedBilinearForm * weakDiv = NULL,
126  ParDiscreteGradOperator * grad = NULL);
127  virtual ~IrrotationalProjector();
128 
129  // Given a GridFunction 'x' of Nedelec DoFs for an arbitrary vector field,
130  // compute the Nedelec DoFs of the irrotational portion, 'y', of
131  // this vector field. The resulting GridFunction will satisfy Curl y = 0
132  // to machine precision.
133  virtual void Mult(const Vector &x, Vector &y) const;
134 
135  void Update();
136 
137 private:
138  void InitSolver() const;
139 
140  ParFiniteElementSpace * H1FESpace_;
141  ParFiniteElementSpace * HCurlFESpace_;
142 
143  ParBilinearForm * s0_;
144  ParMixedBilinearForm * weakDiv_;
145  ParDiscreteGradOperator * grad_;
146 
147  ParGridFunction * psi_;
148  ParGridFunction * xDiv_;
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 
159  bool ownsS0_;
160  bool ownsWeakDiv_;
161  bool ownsGrad_;
162 };
163 
164 /// This class computes the divergence free portion of a vector field.
165 /// This vector field must be discretized using Nedelec basis
166 /// functions.
168 {
169 public:
171  ParFiniteElementSpace & HCurlFESpace,
172  const int & irOrder,
173  ParBilinearForm * s0 = NULL,
174  ParMixedBilinearForm * weakDiv = NULL,
175  ParDiscreteGradOperator * grad = NULL);
176  virtual ~DivergenceFreeProjector();
177 
178  // Given a vector 'x' of Nedelec DoFs for an arbitrary vector field,
179  // compute the Nedelec DoFs of the divergence free portion, 'y', of
180  // this vector field. The resulting vector will satisfy Div y = 0
181  // in a weak sense.
182  virtual void Mult(const Vector &x, Vector &y) const;
183 
184  void Update();
185 };
186 
187 
188 /// Visualize the given parallel mesh object, using a GLVis server on the
189 /// specified host and port. Set the visualization window title, and optionally,
190 /// its geometry.
191 void VisualizeMesh(socketstream &sock, const char *vishost, int visport,
192  ParMesh &pmesh, const char *title,
193  int x = 0, int y = 0, int w = 400, int h = 400,
194  const char *keys = NULL);
195 
196 /// Visualize the given parallel grid function, using a GLVis server on the
197 /// specified host and port. Set the visualization window title, and optionally,
198 /// its geometry.
199 void VisualizeField(socketstream &sock, const char *vishost, int visport,
200  ParGridFunction &gf, const char *title,
201  int x = 0, int y = 0, int w = 400, int h = 400,
202  const char *keys = NULL, bool vec = false);
203 
204 } // namespace common
205 
206 } // namespace mfem
207 
208 #endif // MFEM_USE_MPI
209 #endif
IrrotationalProjector(ParFiniteElementSpace &H1FESpace, ParFiniteElementSpace &HCurlFESpace, const int &irOrder, ParBilinearForm *s0=NULL, ParMixedBilinearForm *weakDiv=NULL, ParDiscreteGradOperator *grad=NULL)
Definition: pfem_extras.cpp:98
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
int vdim
Vector dimension (number of unknowns per degree of freedom).
Definition: fespace.hpp:99
Abstract parallel finite element space.
Definition: pfespace.hpp:28
The BoomerAMG solver in hypre.
Definition: hypre.hpp:951
ParDiscreteCurlOperator(ParFiniteElementSpace *dfes, ParFiniteElementSpace *rfes)
Definition: pfem_extras.cpp:83
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
L2_ParFESpace(ParMesh *m, const int p, const int space_dim, int vdim=1, int order=Ordering::byNODES)
Definition: pfem_extras.cpp:61
PCG solver in hypre.
Definition: hypre.hpp:743
ParDiscreteGradOperator(ParFiniteElementSpace *dfes, ParFiniteElementSpace *rfes)
Definition: pfem_extras.cpp:76
ND_ParFESpace(ParMesh *m, const int p, const int space_dim, int vdim=1, int order=Ordering::byNODES)
Definition: pfem_extras.cpp:37
void VisualizeMesh(socketstream &sock, const char *vishost, int visport, Mesh &mesh, const char *title, int x, int y, int w, int h, const char *keys, bool vec)
Definition: fem_extras.cpp:59
RT_ParFESpace(ParMesh *m, const int p, const int space_dim, int vdim=1, int order=Ordering::byNODES)
Definition: pfem_extras.cpp:49
Class for parallel bilinear form using different test and trial FE spaces.
ParDiscreteDivOperator(ParFiniteElementSpace *dfes, ParFiniteElementSpace *rfes)
Definition: pfem_extras.cpp:90
Class for parallel bilinear form.
Vector data type.
Definition: vector.hpp:48
void VisualizeField(socketstream &sock, const char *vishost, int visport, GridFunction &gf, const char *title, int x, int y, int w, int h, const char *keys, bool vec)
Definition: fem_extras.cpp:94
Class for parallel grid function.
Definition: pgridfunc.hpp:32
H1_ParFESpace(ParMesh *m, const int p, const int space_dim=3, const int type=BasisType::GaussLobatto, int vdim=1, int order=Ordering::byNODES)
Definition: pfem_extras.cpp:24
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:187
Class for parallel meshes.
Definition: pmesh.hpp:32
ParDiscreteInterpolationOperator(ParFiniteElementSpace *dfes, ParFiniteElementSpace *rfes)
Definition: pfem_extras.hpp:88
DivergenceFreeProjector(ParFiniteElementSpace &H1FESpace, ParFiniteElementSpace &HCurlFESpace, const int &irOrder, ParBilinearForm *s0=NULL, ParMixedBilinearForm *weakDiv=NULL, ParDiscreteGradOperator *grad=NULL)