]> gitweb.michael.orlitzky.com - numerical-analysis.git/blob - src/FEM/R1.hs
098dcf34530ed1ccb26c5f79532f1a9f4cfd9156
[numerical-analysis.git] / src / FEM / R1.hs
1 {-# LANGUAGE GADTs #-}
2 {-# LANGUAGE ScopedTypeVariables #-}
3
4 -- | Finite Element Method (FEM) to solve the problem,
5 --
6 -- -(A(x)*u'(x))' + c(x)*u*(x) = f(x) on (a,b)
7 --
8 -- with one of the two types of boundary conditions,
9 --
10 -- * Dirichlet: u(a) = u(b) = 0
11 --
12 -- * A(a)u'(a) = alpha, A(b)u'(b) = beta
13 --
14 -- if c(x) = 0, then it is assumed that the boundary conditions are
15 -- Dirichlet.
16 --
17 -- The code creates a linear system,
18 --
19 -- (K + M)x = F
20 --
21 -- to be solved for @x@, where @K@ (\"big_K\") is the stiffness
22 -- matrix, @M@ (\"big_M\") is the mass matrix, and @F@ (\"big_F\")
23 -- is the load vector.
24 --
25 -- Warning: until PLU factorization is implemented, we can only
26 -- solve the resulting system if it's positive definite!
27 --
28 module FEM.R1
29 where
30
31 import qualified Algebra.Algebraic as Algebraic ( C )
32 import qualified Algebra.Field as Field ( C )
33 import qualified Algebra.RealField as RealField ( C )
34 import qualified Algebra.ToRational as ToRational ( C )
35 import Data.Vector.Fixed ( Arity, S )
36 import NumericPrelude
37 import qualified Prelude as P
38
39 import Integration.Gaussian ( gaussian )
40 import Linear.Matrix (
41 Col,
42 Mat(..),
43 Row,
44 (!!!),
45 construct,
46 fromList,
47 ifoldl2,
48 nrows )
49 import Polynomials.Orthogonal ( legendre )
50
51 -- | Dirichlet boundary conditions. Since u(a)=u(b)=0 are fixed,
52 -- there's no additional information conveyed by this type.
53 data Dirichlet a = Dirichlet { domain_dirichlet :: Interval a }
54
55 -- | Neumann boundary conditions. @alpha@ specifies A(a)u'(b) and
56 -- @beta@ specifies A(b)u'(b).
57 data Neumann a =
58 Neumann { domain_neumann :: Interval a,
59 alpha :: a,
60 beta :: a }
61
62 -- | Boundary conditions can be either Dirichlet or Neumann.
63 type BoundaryConditions a = Either (Dirichlet a) (Neumann a)
64
65 type Interval a = (a,a)
66
67 -- | The data for the PDE that we are attempting to solve.
68 data PDE a =
69 PDE {
70 -- | A(x)
71 big_A :: (a -> a),
72 -- | c(x)
73 c :: (a -> a),
74 -- | f(x)
75 f :: (a -> a),
76
77 -- | The boundary conditions. The domain also specifies the
78 -- boundary in R^1.
79 bdy :: BoundaryConditions a }
80
81
82
83 -- | Non-PDE parameters for the finite element method. The additional
84 -- type parameter @n@ should be a type-level representation of the
85 -- largest element in @max_degrees@. It needs to be known statically
86 -- for the dimensions of the pointer matrix. The parameter @l@ is
87 -- the number of global basis functions. It's equal to the number of
88 -- /internal/ mesh nodes (i.e. m-1), plus the sum of (p_i - 1) for
89 -- each p_i in max_degrees.
90 data Params m n l a =
91 Params {
92 -- | A partition of the domain.
93 mesh :: Col m (Interval a),
94
95 -- | A list of the highest-degree polynomial that we will use over
96 -- each interval in our mesh.
97 max_degrees :: Col m Int }
98
99
100 -- | The pointer matrix mapping local basis elements to global
101 -- ones. The (i,j)th entry of the matrix contains the index of the
102 -- global basis function corresponding to the jth local basis
103 -- function over element i.
104 --
105 -- TODO: This works for Dirichlet boundary conditions only.
106 --
107 -- Note that the number of columns in the matrix depends on the
108 --
109 -- Examples:
110 --
111 -- >>> import Data.Vector.Fixed ( N5, N6 )
112 -- >>> import Linear.Matrix ( Col5, fromList )
113 -- >>> import Naturals ( N19 )
114 --
115 -- >>> let p = fromList [[3],[3],[5],[4],[5]] :: Col5 Int
116 -- >>> let mesh = undefined :: Col5 (Int,Int)
117 -- >>> let params = Params mesh p :: Params N5 N5 N19 Int
118 -- >>> let row1 = [0,1,5,6,0,0] :: [Int]
119 -- >>> let row2 = [1,2,7,8,0,0] :: [Int]
120 -- >>> let row3 = [2,3,9,10,11,12] :: [Int]
121 -- >>> let row4 = [3,4,13,14,15,0] :: [Int]
122 -- >>> let row5 = [4,0,16,17,18,19] :: [Int]
123 -- >>> let expected = fromList [row1,row2,row3,row4,row5] :: Mat N5 N6 Int
124 -- >>> pointer params == expected
125 -- True
126 --
127 pointer :: (Arity m, Arity n, Arity l) => Params m n l a -> Mat m (S n) Int
128 pointer params =
129 construct lambda
130 where
131 -- | The number of polynomials in the local basis for element i.
132 numpolys :: Int -> Int
133 numpolys i = ((max_degrees params) !!! (i,0)) + 1
134
135 lambda i j
136
137 -- The first two (linear) basis functions are easy to do.
138 | j < 2 = if i + j >= (nrows $ mesh params) then 0 else i + j
139
140 -- No local basis function exists for this j.
141 | j >= (numpolys i) = 0
142
143 -- The first row we handle as a special case.
144 | i == 0 = (nrows $ mesh params) + j - 2
145
146 -- The first entry for this row not corresponding to one of the
147 -- linear polynomials (which we did in the first guard). This
148 -- grabs the biggest number in the previous row and begins counting
149 -- from there.
150 | j == 2 = (lambda (i-1) (numpolys (i-1) - 1)) + 1
151
152 -- If none of the other cases apply, we can compute our value by
153 -- looking to the left in the matrix.
154 | otherwise = (lambda i (j-1)) + 1
155
156
157
158 -- | An affine transform taking the interval [x1,x2] into [-1,1].
159 --
160 -- Examples:
161 --
162 -- >>> let phi = affine (-6,7)
163 -- >>> phi (-6)
164 -- -1.0
165 -- >>> phi 7
166 -- 1.0
167 --
168 affine :: Field.C a => (a,a) -> (a -> a)
169 affine (x1,x2) x = (fromInteger 2)*(x - x1)/(x2 - x1) - (fromInteger 1)
170
171 -- | The inverse of 'affine'. It should send [-1,1] into [x1,x2].
172 --
173 -- Examples:
174 --
175 -- >>> let phi = affine_inv (-6,7)
176 -- >>> phi (-1)
177 -- -6.0
178 -- >>> phi 1
179 -- 7.0
180 --
181 affine_inv :: Field.C a => (a,a) -> (a -> a)
182 affine_inv (x1,x2) x =
183 x*(x2 - x1)/two + (x1 + x2)/two
184 where
185 two = fromInteger 2
186
187
188 -- | Normalized integrals of orthogonal basis functions over
189 -- n[-1,1]. The test case below comes from Sage where the
190 -- orthogonality of the polynomials' derivatives can easily be
191 -- tested.
192 --
193 -- Examples:
194 --
195 -- >>> import qualified Algebra.Absolute as Absolute ( abs )
196 --
197 -- >>> let expected = 2.99624907925257
198 -- >>> let actual = big_N 4 1.5 :: Double
199 -- >>> Absolute.abs (actual - expected) < 1e-12
200 -- True
201 --
202 big_N :: forall a. (Algebraic.C a, RealField.C a) => Integer -> a -> a
203 big_N k x
204 | k < 0 = error "requested a negative basis function"
205 | k == 0 = (one - x) / (fromInteger 2)
206 | k == 1 = (one + x) / (fromInteger 2)
207 | otherwise =
208 coeff * ( legendre k x - legendre (k-2) x )
209 where
210 two = fromInteger 2
211 four = fromInteger 4
212 coeff = one / (sqrt (four*(fromInteger k) - two)) :: a
213
214
215 -- | A matrix containing 'big_N' functions indexed by their
216 -- element/number. Each row in the matrix represents a finite element
217 -- (i.e. an interval in the mesh). Within row @i@, column @j@ contains
218 -- the @j@th 'big_N' basis function.
219 --
220 -- Any given 'big_N' will probably wind up in this matrix multiple
221 -- times; the basis functions don't change depending on the
222 -- interval. Only the /number/ of basis functions does. Computing
223 -- them this way allows us to easily construct a lookup \"table\" of
224 -- the proper dimensions.
225 --
226 -- The second example below relies on the fact that @big_N 3@ and
227 -- @big_N 6@ expand to Legendre polynomials (2,4) and (5,7)
228 -- respectively and so should be orthogonal over [-1,1].
229 --
230 -- Examples:
231 --
232 -- >>> import Data.Vector.Fixed ( N5, N6 )
233 -- >>> import Integration.Gaussian ( gaussian )
234 -- >>> import Linear.Matrix ( Col5, fromList )
235 -- >>> import Naturals ( N19 )
236 --
237 -- >>> let p = fromList [[3],[3],[5],[4],[5]] :: Col5 Int
238 -- >>> let mesh = undefined :: Col5 (Double,Double)
239 -- >>> let params = Params mesh p :: Params N5 N5 N19 Double
240 -- >>> let big_ns = big_N_matrix :: Mat N5 N6 (Double -> Double)
241 -- >>> let n1 = big_ns !!! (1,0)
242 -- >>> let n4 = big_ns !!! (4,0)
243 -- >>> n1 1.5 == n4 1.5
244 -- True
245 -- >>> let n1 = big_ns !!! (1,3)
246 -- >>> let n2 = big_ns !!! (2,4)
247 -- >>> gaussian (\x -> (n1 x) * (n2 x)) < 1e-12
248 -- True
249 --
250 big_N_matrix :: (Arity m, Arity n, Algebraic.C a, RealField.C a)
251 => Mat m n (a -> a)
252 big_N_matrix =
253 construct lambda
254 where
255 lambda _ j x = big_N (toInteger j) x
256
257
258
259
260 -- | Derivatives of the 'big_N's, that is, orthogonal basis functions
261 -- over [-1,1]. The test case below comes from Sage where the
262 -- orthogonality of the polynomials' derivatives can easily be
263 -- tested. The indices are shifted by one so that k=0 is the first
264 -- basis function.
265 --
266 -- Examples:
267 --
268 -- >>> import qualified Algebra.Absolute as Absolute ( abs )
269 --
270 -- >>> let expected = 11.5757525403319
271 -- >>> let actual = big_N' 3 1.5 :: Double
272 -- >>> Absolute.abs (actual - expected) < 1e-10
273 -- True
274 --
275 big_N' :: forall a. (Algebraic.C a, RealField.C a) => Integer -> a -> a
276 big_N' k x
277 | k < 0 = error "requested a negative basis function"
278 | k == 0 = negate ( one / (fromInteger 2))
279 | k == 1 = one / (fromInteger 2)
280 | otherwise = coeff * ( legendre k x )
281 where
282 two = fromInteger 2
283 coeff = sqrt ((two*(fromInteger k) + one) / two) :: a
284
285
286 -- | The matrix of (N_i' * N_j') functions used in the integrand of
287 -- the stiffness/mass matrices.
288 big_N's_matrix :: (Arity m, Arity n, Algebraic.C a, RealField.C a)
289 => Mat m n (a -> a)
290 big_N's_matrix =
291 construct lambda
292 where
293 lambda i j x = (big_N' (toInteger i) x) * (big_N' (toInteger j) x)
294
295
296 -- | Compute the global load vector F.
297 --
298 -- Examples:
299 --
300 -- >>> import Data.Vector.Fixed ( N3, N4 )
301 -- >>> import Linear.Matrix ( Col4, frobenius_norm, fromList )
302 -- >>> import Naturals ( N7 )
303 --
304 -- >>> let big_A = const (1::Double)
305 -- >>> let c x = sin x
306 -- >>> let f x = x*(sin x)
307 -- >>> let bdy = Left (Dirichlet (0,1::Double))
308 -- >>> let pde = PDE big_A c f bdy
309 --
310 -- >>> let i1 = (0.0,1/3)
311 -- >>> let i2 = (1/3,2/3)
312 -- >>> let i3 = (2/3,4/5)
313 -- >>> let i4 = (4/5,1.0)
314 -- >>> let mesh = fromList [[i1], [i2], [i3], [i4]] :: Col4 (Double,Double)
315 -- >>> let pvec = fromList [[2],[3],[2],[1]] :: Col4 Int
316 -- >>> let params = Params mesh pvec :: Params N4 N3 N7 Double
317 --
318 -- >>> let f1 = [0.0418]
319 -- >>> let f2 = [0.0805]
320 -- >>> let f3 = [0.1007]
321 -- >>> let f4 = [-0.0045]
322 -- >>> let f5 = [-0.0332]
323 -- >>> let f6 = [-0.0054]
324 -- >>> let f7 = [-0.0267]
325 -- >>> let expected = fromList [f1,f2,f3,f4,f5,f6,f7] :: Col N7 Double
326 -- >>> let actual = big_F pde params
327 -- >>> frobenius_norm (actual - expected) < 1e-4
328 -- True
329 --
330 big_F :: forall m n l a.
331 (Arity l, Arity m, Arity n,
332 Algebraic.C a, RealField.C a, ToRational.C a)
333 => PDE a
334 -> Params m n l a
335 -> Col l a
336 big_F pde params =
337 ifoldl2 accum zero (big_N_matrix :: Mat m (S n) (a -> a))
338 where
339 accum :: Int -> Int -> Col l a -> (a -> a) -> Col l a
340 accum i j prev_F this_N =
341 prev_F + this_F
342 where
343 two = fromInteger 2
344 (x1,x2) = (mesh params) !!! (i,0)
345 q = affine_inv (x1,x2)
346 integrand x = ((f pde) (q x)) * (this_N x)
347
348 -- The pointer matrix numbers from 1 so subtract one here to
349 -- get the right index.
350 global_idx = ((pointer params) !!! (i,j)) - 1
351 this_F = construct lambda
352 lambda k _ = if k == global_idx
353 then (gaussian integrand)*(x2 - x1) / two
354 else fromInteger 0
355
356
357
358 big_K_elem :: forall m n l a b.
359 (Arity l, Arity m, Arity n,
360 Algebraic.C a, RealField.C a, ToRational.C a)
361 => PDE a
362 -> Params m n l a
363 -> Int
364 -> Int
365 -> Mat l l a
366 -> b
367 -> Mat l l a
368 big_K_elem pde params _ k cur_K _ =
369 ifoldl2 accum cur_K (big_N's_matrix :: Mat m (S n) (a -> a))
370 where
371 accum :: Int -> Int -> Mat l l a -> (a -> a) -> Mat l l a
372 accum i j prev_K these_N's =
373 prev_K + this_K
374 where
375 two = fromInteger 2
376 (x1,x2) = (mesh params) !!! (k,0)
377 q = affine_inv (x1,x2)
378 integrand x = ((big_A pde) (q x)) * (these_N's x)
379 -- The pointer matrix numbers from 1 so subtract one here to
380 -- get the right index.
381 global_row_idx = ((pointer params) !!! (k,i)) - 1
382 global_col_idx = ((pointer params) !!! (k,j)) - 1
383 this_K = construct lambda
384 lambda v w = if v == global_row_idx && w == global_col_idx
385 then (two/(x2 - x1))* (gaussian integrand)
386 else fromInteger 0
387
388
389
390 -- | Compute the \"big K\" stiffness matrix. There are three
391 -- parameters needed for K, namely i,j,k so a fold over a matrix will
392 -- not do. This little gimmick simulates a three-index fold by doing a
393 -- two-index fold over a row of the proper dimensions.
394 --
395 -- Examples:
396 --
397 -- >>> import Data.Vector.Fixed ( N3, N4 )
398 -- >>> import Linear.Matrix ( Col4, frobenius_norm, fromList )
399 -- >>> import Naturals ( N7 )
400 --
401 -- >>> let big_A = const (1::Double)
402 -- >>> let c x = sin x
403 -- >>> let f x = x*(sin x)
404 -- >>> let bdy = Left (Dirichlet (0,1::Double))
405 -- >>> let pde = PDE big_A c f bdy
406 --
407 -- >>> let i1 = (0.0,1/3)
408 -- >>> let i2 = (1/3,2/3)
409 -- >>> let i3 = (2/3,4/5)
410 -- >>> let i4 = (4/5,1.0)
411 -- >>> let mesh = fromList [[i1], [i2], [i3], [i4]] :: Col4 (Double,Double)
412 -- >>> let pvec = fromList [[2],[3],[2],[1]] :: Col4 Int
413 -- >>> let params = Params mesh pvec :: Params N4 N3 N7 Double
414 --
415 -- >>> let k1 = [6, -3, 0, 0, 0, 0, 0] :: [Double]
416 -- >>> let k2 = [-3, 10.5, -7.5, 0, 0, 0, 0] :: [Double]
417 -- >>> let k3 = [0, -7.5, 12.5, 0, 0, 0, 0] :: [Double]
418 -- >>> let k4 = [0, 0, 0, 6, 0, 0, 0] :: [Double]
419 -- >>> let k5 = [0, 0, 0, 0, 6, 0, 0] :: [Double]
420 -- >>> let k6 = [0, 0, 0, 0, 0, 6, 0] :: [Double]
421 -- >>> let k7 = [0, 0, 0, 0, 0, 0, 15] :: [Double]
422 -- >>> let expected = fromList [k1,k2,k3,k4,k5,k6,k7] :: Mat N7 N7 Double
423 -- >>> let actual = big_K pde params
424 -- >>> frobenius_norm (actual - expected) < 1e-10
425 -- True
426 --
427 big_K :: forall m n l a.
428 (Arity l, Arity m, Arity n,
429 Algebraic.C a, RealField.C a, ToRational.C a)
430 => PDE a
431 -> Params m n l a
432 -> Mat l l a
433 big_K pde params =
434 ifoldl2 (big_K_elem pde params) zero col_idxs
435 where
436 col_idxs = fromList [map fromInteger [0..]] :: Row m a
437