forked from SCOREC/inclusion_solver
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathpfem_extras.hpp
212 lines (174 loc) · 6.75 KB
/
pfem_extras.hpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
#ifndef PFEM_EXTRAS
#define PFEM_EXTRAS
#include <mfem.hpp>
#include <cstddef>
using namespace mfem;
// constants
// Permittivity of Free Space (units F/m)
static const double epsilon0_ = 1.; //8.8541878176e-12;
class ElementMeshStream : public std::stringstream
{
public:
ElementMeshStream(Element::Type e);
};
/** The H1_ParFESpace class is a ParFiniteElementSpace which automatically
allocates and destroys its own FiniteElementCollection, in this
case an H1_FECollection object.
*/
class H1_ParFESpace : public ParFiniteElementSpace
{
public:
H1_ParFESpace(ParMesh *m,
const int p, const int space_dim = 3,
const int type = BasisType::GaussLobatto,
int vdim = 1, int order = Ordering::byNODES);
~H1_ParFESpace();
private:
const FiniteElementCollection *FEC_;
};
/** The ND_ParFESpace class is a ParFiniteElementSpace which automatically
allocates and destroys its own FiniteElementCollection, in this
case an ND_FECollection object.
*/
class ND_ParFESpace : public ParFiniteElementSpace
{
public:
ND_ParFESpace(ParMesh *m, const int p, const int space_dim,
int vdim = 1, int order = Ordering::byNODES);
~ND_ParFESpace();
private:
const FiniteElementCollection *FEC_;
};
/** The RT_ParFESpace class is a ParFiniteElementSpace which automatically
allocates and destroys its own FiniteElementCollection, in this
case an RT_FECollection object.
*/
class RT_ParFESpace : public ParFiniteElementSpace
{
public:
RT_ParFESpace(ParMesh *m, const int p, const int space_dim,
int vdim = 1, int order = Ordering::byNODES);
~RT_ParFESpace();
private:
const FiniteElementCollection *FEC_;
};
/** The L2_ParFESpace class is a ParFiniteElementSpace which automatically
allocates and destroys its own FiniteElementCollection, in this
case an L2_FECollection object.
*/
class L2_ParFESpace : public ParFiniteElementSpace
{
public:
L2_ParFESpace(ParMesh *m, const int p, const int space_dim,
int vdim = 1, int order = Ordering::byNODES);
~L2_ParFESpace();
private:
const FiniteElementCollection *FEC_;
};
class ParDiscreteInterpolationOperator : public ParDiscreteLinearOperator
{
public:
ParDiscreteInterpolationOperator(ParFiniteElementSpace *dfes,
ParFiniteElementSpace *rfes)
: ParDiscreteLinearOperator(dfes, rfes) {}
virtual ~ParDiscreteInterpolationOperator();
};
class ParDiscreteGradOperator : public ParDiscreteInterpolationOperator
{
public:
ParDiscreteGradOperator(ParFiniteElementSpace *dfes,
ParFiniteElementSpace *rfes);
};
class ParDiscreteCurlOperator : public ParDiscreteInterpolationOperator
{
public:
ParDiscreteCurlOperator(ParFiniteElementSpace *dfes,
ParFiniteElementSpace *rfes);
};
class ParDiscreteDivOperator : public ParDiscreteInterpolationOperator
{
public:
ParDiscreteDivOperator(ParFiniteElementSpace *dfes,
ParFiniteElementSpace *rfes);
};
/// This class computes the irrotational portion of a vector field.
/// This vector field must be discretized using Nedelec basis
/// functions.
class IrrotationalProjector : public Operator
{
public:
IrrotationalProjector(ParFiniteElementSpace & H1FESpace,
ParFiniteElementSpace & HCurlFESpace,
const int & irOrder,
ParBilinearForm * s0 = NULL,
ParMixedBilinearForm * weakDiv = NULL,
ParDiscreteGradOperator * grad = NULL);
virtual ~IrrotationalProjector();
// Given a GridFunction 'x' of Nedelec DoFs for an arbitrary vector field,
// compute the Nedelec DoFs of the irrotational portion, 'y', of
// this vector field. The resulting GridFunction will satisfy Curl y = 0
// to machine precision.
virtual void Mult(const Vector &x, Vector &y) const;
void Update();
private:
void InitSolver() const;
ParFiniteElementSpace * H1FESpace_;
ParFiniteElementSpace * HCurlFESpace_;
ParBilinearForm * s0_;
ParMixedBilinearForm * weakDiv_;
ParDiscreteGradOperator * grad_;
ParGridFunction * psi_;
ParGridFunction * xDiv_;
HypreParMatrix * S0_;
mutable Vector Psi_;
mutable Vector RHS_;
mutable HypreBoomerAMG * amg_;
mutable HyprePCG * pcg_;
Array<int> ess_bdr_, ess_bdr_tdofs_;
bool ownsS0_;
bool ownsWeakDiv_;
bool ownsGrad_;
};
/// This class computes the divergence free portion of a vector field.
/// This vector field must be discretized using Nedelec basis
/// functions.
class DivergenceFreeProjector : public IrrotationalProjector
{
public:
DivergenceFreeProjector(ParFiniteElementSpace & H1FESpace,
ParFiniteElementSpace & HCurlFESpace,
const int & irOrder,
ParBilinearForm * s0 = NULL,
ParMixedBilinearForm * weakDiv = NULL,
ParDiscreteGradOperator * grad = NULL);
virtual ~DivergenceFreeProjector();
// Given a vector 'x' of Nedelec DoFs for an arbitrary vector field,
// compute the Nedelec DoFs of the divergence free portion, 'y', of
// this vector field. The resulting vector will satisfy Div y = 0
// in a weak sense.
virtual void Mult(const Vector &x, Vector &y) const;
void Update();
};
/// Visualize the given parallel mesh object, using a GLVis server on the
/// specified host and port. Set the visualization window title, and optionally,
/// its geometry.
void VisualizeMesh(socketstream &sock, const char *vishost, int visport,
ParMesh &pmesh, const char *title,
int x = 0, int y = 0, int w = 400, int h = 400,
const char *keys = NULL);
/// Visualize the given parallel grid function, using a GLVis server on the
/// specified host and port. Set the visualization window title, and optionally,
/// its geometry.
void VisualizeField(socketstream &sock, const char *vishost, int visport,
const ParGridFunction &gf, const char *title,
int x = 0, int y = 0, int w = 400, int h = 400,
const char *keys = NULL, bool vec = false);
/// Merges vertices which lie at the same location
void MergeMeshNodes(Mesh * mesh, int logging);
/// Convert a set of attribute numbers to a marker array
/** The marker array will be of size max_attr and it will contain only zeroes
and ones. Ones indicate which attribute numbers are present in the attrs
array. In the special case when attrs has a single entry equal to -1 the
marker array will contain all ones. */
void AttrToMarker(int max_attr, const Array<int> &attrs, Array<int> &marker);
#endif