-
Notifications
You must be signed in to change notification settings - Fork 4
/
Constructors.jl
159 lines (126 loc) · 2.97 KB
/
Constructors.jl
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
# ------------------------------------------------------------------------------
# BDDC
# ------------------------------------------------------------------------------
# ------------------------------------------------------------------------------
# Lumped
# ------------------------------------------------------------------------------
"""
```julia
LumpedBDDC(operator)
```
Lumped BDDC preconditioner for finite element operators
# Arguments:
- `operator::Operator`: finite element operator to precondition
# Returns:
- lumped BDDC preconditioner object
# Example:
```jldoctest
# setup
mesh = Mesh2D(1.0, 1.0);
# operator
finediffusion = GalleryOperator("diffusion", 5, 5, mesh);
# preconditioner
bddc = LumpedBDDC(finediffusion);
# verify
println(bddc)
println(bddc.operator)
# output
lumped BDDC preconditioner
finite element operator:
2d mesh:
dx: 1.0
dy: 1.0
2 inputs:
operator field:
tensor product basis:
numbernodes1d: 5
numberquadraturepoints1d: 5
numbercomponents: 1
dimension: 2
evaluation mode:
gradient
operator field:
tensor product basis:
numbernodes1d: 5
numberquadraturepoints1d: 5
numbercomponents: 1
dimension: 2
evaluation mode:
quadratureweights
1 output:
operator field:
tensor product basis:
numbernodes1d: 5
numberquadraturepoints1d: 5
numbercomponents: 1
dimension: 2
evaluation mode:
gradient
```
"""
function LumpedBDDC(operator::Operator)
# common constructor
return BDDC(operator, BDDCInjectionType.scaled)
end
# ------------------------------------------------------------------------------
# Dirichlet
# ------------------------------------------------------------------------------
"""
```julia
DirichletBDDC(operator)
```
Dirichlet BDDC preconditioner for finite element operators
# Arguments:
- `operator::Operator`: finite element operator to precondition
# Returns:
- Dirichlet BDDC preconditioner object
# Example:
```jldoctest
# setup
mesh = Mesh2D(1.0, 1.0);
# operator
finediffusion = GalleryOperator("diffusion", 5, 5, mesh);
# preconditioner
bddc = DirichletBDDC(finediffusion);
# verify
println(bddc)
println(bddc.operator)
# output
Dirichlet BDDC preconditioner
finite element operator:
2d mesh:
dx: 1.0
dy: 1.0
2 inputs:
operator field:
tensor product basis:
numbernodes1d: 5
numberquadraturepoints1d: 5
numbercomponents: 1
dimension: 2
evaluation mode:
gradient
operator field:
tensor product basis:
numbernodes1d: 5
numberquadraturepoints1d: 5
numbercomponents: 1
dimension: 2
evaluation mode:
quadratureweights
1 output:
operator field:
tensor product basis:
numbernodes1d: 5
numberquadraturepoints1d: 5
numbercomponents: 1
dimension: 2
evaluation mode:
gradient
```
"""
function DirichletBDDC(operator::Operator)
# common constructor
return BDDC(operator, BDDCInjectionType.harmonic)
end
# ------------------------------------------------------------------------------