1
2
3
4
5import numpy as np
6import scipy.sparse.linalg
7from scipy.sparse import lil_matrix
8
9import dune.geometry
10import dune.grid as grid
11import dune.functions as functions
12
13
14
15
16
17
18
19def localAssembler(localView, volumeTerm):
20
21
22 n = len(localView)
23
24
25 element = localView.element()
26
27
28 localBasis = localView.tree().finiteElement().localBasis
29
30
31 localA = np.zeros((n,n))
32 localB = np.zeros(n)
33
34
35 quadRule = dune.geometry.quadratureRule(element.type, order=4)
36 for pt in quadRule:
37
38
39 integrationElement = element.geometry.integrationElement(pt.position)
40
41
42 values = localBasis.evaluateFunction(pt.position)
43
44
45 referenceJacobians = localBasis.evaluateJacobian(pt.position)
46
47
48 geometryJacobianInverse = element.geometry.jacobianInverse(pt.position)
49 jacobians = [ np.dot(np.array(g)[0], geometryJacobianInverse) for g in referenceJacobians ]
50
51 quadPosGlobal = element.geometry.toGlobal(pt.position)
52
53 for i in range( n ):
54 for j in range( n ):
55 localA[i,j] += pt.weight * integrationElement * np.dot(jacobians[i], jacobians[j])
56
57 localB[i] += pt.weight * integrationElement * values[i] * volumeTerm(quadPosGlobal)
58
59 return localA, localB
60
61
62
63
64
65def assembleLaplaceMatrix(basis, volumeTerm):
66
67
68 n = len(basis)
69
70
71 A = lil_matrix((n,n))
72
73
74 b = np.zeros(n)
75
76
77 localView = basis.localView()
78
79
80 grid = basis.gridView
81 for element in grid.elements:
82
83
84 localView.bind(element)
85
86
87 localN = len(localView)
88
89
90 localA, localb = localAssembler(localView, volumeTerm)
91
92
93 for i in range(localN):
94
95 gi = localView.index(i)[0]
96
97 b[gi] += localb[i]
98
99 for j in range(localN):
100 gj = localView.index(j)[0]
101 A[gi, gj] += localA[i, j]
102
103
104 return A.tocsr(), b
105
106
107
108
109
110
111
112gridSize = 32
113
114
115grid = grid.structuredGrid([0,0],[1,1],[gridSize,gridSize])
116
117
118
119
120basis = functions.defaultGlobalBasis(grid, functions.Lagrange(order=2))
121
122
123
124
125f = lambda x : 10
126
127
128A,b = assembleLaplaceMatrix(basis, f)
129
130
131
132
133isDirichlet = np.zeros(len(basis))
134
135def markDOF(boundaryDOFNumber):
136 isDirichlet[boundaryDOFNumber] = True
137
138functions.boundarydofs.forEachBoundaryDOF(basis,markDOF)
139
140
141
142
143dirichletValueFunction = lambda x : np.sin(2*np.pi*x[0])
144
145
146dirichletValues = np.zeros(len(basis))
147basis.interpolate(dirichletValues, dirichletValueFunction)
148
149
150
151
152rows, cols = A.nonzero()
153
154for i,j in zip(rows, cols):
155 if isDirichlet[i]:
156 if i==j:
157 A[i,j] = 1.0
158 else:
159 A[i,j] = 0
160 b[i] = dirichletValues[i]
161
162
163
164
165x = scipy.sparse.linalg.spsolve(A, b)
166
167
168
169
170uh = basis.asFunction(x)
171
172vtkWriter = grid.vtkWriter(subsampling=2)
173uh.addToVTKWriter("u", vtkWriter, dune.grid.DataType.PointData)
174vtkWriter.write("poisson-pq2-result")
175