Dune-Functions 2.11
Loading...
Searching...
No Matches
poisson-pq2.py

This is the raw source code of the poisson-pq2.py example. There is also a commented version of this example in the Examples section.

1# SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file AUTHORS.md
2# SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception OR LGPL-3.0-or-later
3
4#! [imports]
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#! [imports]
13
14# Compute element stiffness matrix and element load vector
15#
16# TODO: This assembler loop is very inefficient in terms of run time and should be improved using Python vectorization.
17# See discussion at https://gitlab.dune-project.org/staging/dune-functions/-/merge_requests/295 for hints and code pointers.
18#! [localAssembler]
19def localAssembler(localView, volumeTerm):
20
21 # Number of degrees of freedom on this element
22 n = len(localView)
23
24 # The grid element
25 element = localView.element()
26
27 # The set of shape functions
28 localBasis = localView.tree().finiteElement().localBasis
29
30 # Initialize element matrix and load vector
31 localA = np.zeros((n,n))
32 localB = np.zeros(n)
33
34 # Create a quadrature rule and integrate
35 quadRule = dune.geometry.quadratureRule(element.type, order=4)
36 for pt in quadRule:
37
38 # The determinant that appears in the integral transformation formula
39 integrationElement = element.geometry.integrationElement(pt.position)
40
41 # Evaluate all shape functions (The return value is an array!)
42 values = localBasis.evaluateFunction(pt.position)
43
44 # Evaluate the shape function Jacobians on the reference element (array of arrays)
45 referenceJacobians = localBasis.evaluateJacobian(pt.position)
46
47 # Transform the reference Jacobians to the actual element
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#! [localAssembler]
61
62
63# The assembler for the global stiffness matrix
64#! [assembleLaplaceMatrix]
65def assembleLaplaceMatrix(basis, volumeTerm):
66
67 # Total number of degrees of freedom
68 n = len(basis)
69
70 # Make empty sparse matrix
71 A = lil_matrix((n,n))
72
73 # Make empty vector
74 b = np.zeros(n)
75
76 # View on the finite element basis on a single element
77 localView = basis.localView()
78
79 # Loop over all grid elements
80 grid = basis.gridView
81 for element in grid.elements:
82
83 # Bind the localView to the current element
84 localView.bind(element)
85
86 # Number of degrees of freedom on the current element
87 localN = len(localView)
88
89 # Assemble the local stiffness matrix and load vector
90 localA, localb = localAssembler(localView, volumeTerm)
91
92 # Copy the local entries into the global matrix and vector
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 # Convert matrix to CSR format
104 return A.tocsr(), b
105#! [assembleLaplaceMatrix]
106
107
108
109
110#! [createGrid]
111# Number of grid elements in one direction
112gridSize = 32
113
114# Create a grid of the unit square
115grid = grid.structuredGrid([0,0],[1,1],[gridSize,gridSize])
116#! [createGrid]
117
118#! [createBasis]
119# Create a second-order Lagrange FE basis
120basis = functions.defaultGlobalBasis(grid, functions.Lagrange(order=2))
121#! [createBasis]
122
123#! [assembly]
124# Source term
125f = lambda x : 10
126
127# Compute stiffness matrix and load vector
128A,b = assembleLaplaceMatrix(basis, f)
129#! [assembly]
130
131#! [dirichletDOFs]
132# Determine all coefficients that are on the boundary
133isDirichlet = np.zeros(len(basis))
134
135def markDOF(boundaryDOFNumber):
136 isDirichlet[boundaryDOFNumber] = True
137
138functions.boundarydofs.forEachBoundaryDOF(basis,markDOF)
139#! [dirichletDOFs]
140
141#! [dirichletValues]
142# The function that implements the Dirichlet values
143dirichletValueFunction = lambda x : np.sin(2*np.pi*x[0])
144
145# Get coefficients of a Lagrange-FE approximation of the Dirichlet values
146dirichletValues = np.zeros(len(basis))
147basis.interpolate(dirichletValues, dirichletValueFunction)
148#! [dirichletValues]
149
150#! [dirichletIntegration]
151# Integrate Dirichlet conditions into the matrix and load vector
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#! [dirichletIntegration]
162
163#! [solving]
164# Solve linear system!
165x = scipy.sparse.linalg.spsolve(A, b)
166#! [solving]
167
168#! [vtkWriting]
169# Write result as vtu file
170uh = basis.asFunction(x)
171
172vtkWriter = grid.vtkWriter(subsampling=2)
173uh.addToVTKWriter("u", vtkWriter, dune.grid.DataType.PointData)
174vtkWriter.write("poisson-pq2-result")
175#! [vtkWriting]