2 {-# LANGUAGE ScopedTypeVariables #-}
4 -- | Finite Element Method (FEM) to solve the problem,
6 -- -(A(x)*u'(x))' + c(x)*u*(x) = f(x) on (a,b)
8 -- with one of the two types of boundary conditions,
10 -- * Dirichlet: u(a) = u(b) = 0
12 -- * A(a)u'(a) = alpha, A(b)u'(b) = beta
14 -- if c(x) = 0, then it is assumed that the boundary conditions are
17 -- The code creates a linear system,
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.
25 -- Warning: until PLU factorization is implemented, we can only
26 -- solve the resulting system if it's positive definite!
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 )
37 import qualified Prelude as P
39 import Integration.Gaussian ( gaussian )
40 import Linear.Matrix (
50 import Linear.System ( solve_positive_definite )
51 import Polynomials.Orthogonal ( legendre )
53 -- | Dirichlet boundary conditions. Since u(a)=u(b)=0 are fixed,
54 -- there's no additional information conveyed by this type.
55 data Dirichlet a = Dirichlet { domain_dirichlet :: Interval a }
57 -- | Neumann boundary conditions. @alpha@ specifies A(a)u'(b) and
58 -- @beta@ specifies A(b)u'(b).
60 Neumann { domain_neumann :: Interval a,
64 -- | Boundary conditions can be either Dirichlet or Neumann.
65 type BoundaryConditions a = Either (Dirichlet a) (Neumann a)
67 type Interval a = (a,a)
69 -- | The data for the PDE that we are attempting to solve.
79 -- | The boundary conditions. The domain also specifies the
81 bdy :: BoundaryConditions a }
85 -- | Non-PDE parameters for the finite element method. The additional
86 -- type parameter @n@ should be a type-level representation of the
87 -- largest element in @max_degrees@. It needs to be known statically
88 -- for the dimensions of the pointer matrix. The parameter @l@ is
89 -- the number of global basis functions. It's equal to the number of
90 -- /internal/ mesh nodes (i.e. m-1), plus the sum of (p_i - 1) for
91 -- each p_i in max_degrees.
94 -- | A partition of the domain.
95 mesh :: Col m (Interval a),
97 -- | A list of the highest-degree polynomial that we will use over
98 -- each interval in our mesh.
99 max_degrees :: Col m Int }
102 -- | The pointer matrix mapping local basis elements to global
103 -- ones. The (i,j)th entry of the matrix contains the index of the
104 -- global basis function corresponding to the jth local basis
105 -- function over element i.
107 -- TODO: This works for Dirichlet boundary conditions only.
109 -- Note that the number of columns in the matrix depends on the
113 -- >>> import Data.Vector.Fixed ( N5, N6 )
114 -- >>> import Linear.Matrix ( Col5, fromList )
115 -- >>> import Naturals ( N19 )
117 -- >>> let p = fromList [[3],[3],[5],[4],[5]] :: Col5 Int
118 -- >>> let mesh = undefined :: Col5 (Int,Int)
119 -- >>> let params = Params mesh p :: Params N5 N5 N19 Int
120 -- >>> let row1 = [0,1,5,6,0,0] :: [Int]
121 -- >>> let row2 = [1,2,7,8,0,0] :: [Int]
122 -- >>> let row3 = [2,3,9,10,11,12] :: [Int]
123 -- >>> let row4 = [3,4,13,14,15,0] :: [Int]
124 -- >>> let row5 = [4,0,16,17,18,19] :: [Int]
125 -- >>> let expected = fromList [row1,row2,row3,row4,row5] :: Mat N5 N6 Int
126 -- >>> pointer params == expected
129 pointer :: (Arity m, Arity n, Arity l) => Params m n l a -> Mat m (S n) Int
133 -- | The number of polynomials in the local basis for element i.
134 numpolys :: Int -> Int
135 numpolys i = ((max_degrees params) !!! (i,0)) + 1
139 -- The first two (linear) basis functions are easy to do.
140 | j < 2 = if i + j >= (nrows $ mesh params) then 0 else i + j
142 -- No local basis function exists for this j.
143 | j >= (numpolys i) = 0
145 -- The first row we handle as a special case.
146 | i == 0 = (nrows $ mesh params) + j - 2
148 -- The first entry for this row not corresponding to one of the
149 -- linear polynomials (which we did in the first guard). This
150 -- grabs the biggest number in the previous row and begins counting
152 | j == 2 = (lambda (i-1) (numpolys (i-1) - 1)) + 1
154 -- If none of the other cases apply, we can compute our value by
155 -- looking to the left in the matrix.
156 | otherwise = (lambda i (j-1)) + 1
160 -- | An affine transform taking the interval [x1,x2] into [-1,1].
164 -- >>> let phi = affine (-6,7)
170 affine :: Field.C a => (a,a) -> (a -> a)
171 affine (x1,x2) x = (fromInteger 2)*(x - x1)/(x2 - x1) - (fromInteger 1)
173 -- | The inverse of 'affine'. It should send [-1,1] into [x1,x2].
177 -- >>> let phi = affine_inv (-6,7)
183 affine_inv :: Field.C a => (a,a) -> (a -> a)
184 affine_inv (x1,x2) x =
185 x*(x2 - x1)/two + (x1 + x2)/two
192 -- | Normalized integrals of orthogonal basis functions over
193 -- n[-1,1]. The test case below comes from Sage where the
194 -- orthogonality of the polynomials' derivatives can easily be
199 -- >>> import qualified Algebra.Absolute as Absolute ( abs )
201 -- >>> let expected = 2.99624907925257
202 -- >>> let actual = big_N 4 1.5 :: Double
203 -- >>> Absolute.abs (actual - expected) < 1e-12
206 big_N :: forall a. (Algebraic.C a, RealField.C a) => Integer -> a -> a
208 | k < 0 = error "requested a negative basis function"
209 | k == 0 = (one - x) / (fromInteger 2)
210 | k == 1 = (one + x) / (fromInteger 2)
212 coeff * ( legendre k x - legendre (k-2) x )
216 coeff = one / (sqrt (four*(fromInteger k) - two)) :: a
219 -- | A matrix containing 'big_N' functions indexed by their
220 -- element/number. Each row in the matrix represents a finite element
221 -- (i.e. an interval in the mesh). Within row @i@, column @j@ contains
222 -- the @j@th 'big_N' basis function.
224 -- Any given 'big_N' will probably wind up in this matrix multiple
225 -- times; the basis functions don't change depending on the
226 -- interval. Only the /number/ of basis functions does. Computing
227 -- them this way allows us to easily construct a lookup \"table\" of
228 -- the proper dimensions.
230 -- The second example below relies on the fact that @big_N 3@ and
231 -- @big_N 6@ expand to Legendre polynomials (2,4) and (5,7)
232 -- respectively and so should be orthogonal over [-1,1].
236 -- >>> import Data.Vector.Fixed ( N5, N6 )
237 -- >>> import Integration.Gaussian ( gaussian )
238 -- >>> import Linear.Matrix ( Col5, fromList )
239 -- >>> import Naturals ( N19 )
241 -- >>> let p = fromList [[3],[3],[5],[4],[5]] :: Col5 Int
242 -- >>> let mesh = undefined :: Col5 (Double,Double)
243 -- >>> let params = Params mesh p :: Params N5 N5 N19 Double
244 -- >>> let big_ns = big_N_matrix :: Mat N5 N6 (Double -> Double)
245 -- >>> let n1 = big_ns !!! (1,0)
246 -- >>> let n4 = big_ns !!! (4,0)
247 -- >>> n1 1.5 == n4 1.5
249 -- >>> let n1 = big_ns !!! (1,3)
250 -- >>> let n2 = big_ns !!! (2,4)
251 -- >>> gaussian (\x -> (n1 x) * (n2 x)) < 1e-12
254 big_N_matrix :: (Arity m, Arity n, Algebraic.C a, RealField.C a)
259 lambda _ j x = big_N (toInteger j) x
262 -- | The matrix of (N_i * N_j) functions used in the integrand of
263 -- the mass matrices.
264 big_Ns_matrix :: (Arity m, Arity n, Algebraic.C a, RealField.C a)
269 lambda i j x = (big_N (toInteger i) x) * (big_N (toInteger j) x)
272 -- | Compute the global load vector F.
276 -- >>> import Linear.Matrix ( Col4, frobenius_norm, fromList )
277 -- >>> import Naturals ( N3, N4, N7 )
279 -- >>> let big_A = const (1::Double)
280 -- >>> let c x = sin x
281 -- >>> let f x = x*(sin x)
282 -- >>> let bdy = Left (Dirichlet (0,1::Double))
283 -- >>> let pde = PDE big_A c f bdy
285 -- >>> let i1 = (0.0,1/3)
286 -- >>> let i2 = (1/3,2/3)
287 -- >>> let i3 = (2/3,4/5)
288 -- >>> let i4 = (4/5,1.0)
289 -- >>> let mesh = fromList [[i1], [i2], [i3], [i4]] :: Col4 (Double,Double)
290 -- >>> let pvec = fromList [[2],[3],[2],[1]] :: Col4 Int
291 -- >>> let params = Params mesh pvec :: Params N4 N3 N7 Double
293 -- >>> let f1 = [0.0418]
294 -- >>> let f2 = [0.0805]
295 -- >>> let f3 = [0.1007]
296 -- >>> let f4 = [-0.0045]
297 -- >>> let f5 = [-0.0332]
298 -- >>> let f6 = [-0.0054]
299 -- >>> let f7 = [-0.0267]
300 -- >>> let expected = fromList [f1,f2,f3,f4,f5,f6,f7] :: Col N7 Double
301 -- >>> let actual = big_F pde params
302 -- >>> frobenius_norm (actual - expected) < 1e-4
305 big_F :: forall m n l a.
306 (Arity l, Arity m, Arity n,
307 Algebraic.C a, RealField.C a, ToRational.C a)
312 ifoldl2 accum zero (big_N_matrix :: Mat m (S n) (a -> a))
314 accum :: Int -> Int -> Col l a -> (a -> a) -> Col l a
315 accum i j prev_F this_N =
319 (x1,x2) = (mesh params) !!! (i,0)
320 q = affine_inv (x1,x2)
321 integrand x = ((f pde) (q x)) * (this_N x)
323 -- The pointer matrix numbers from 1 so subtract one here to
324 -- get the right index.
325 k = ((pointer params) !!! (i,j)) - 1
326 integral = (gaussian integrand)*(x2 - x1) / two
327 this_F = set_idx zero (k,0) integral
330 -- * Stiffness matrix
332 -- | Derivatives of the 'big_N's, that is, orthogonal basis functions
333 -- over [-1,1]. The test case below comes from Sage where the
334 -- orthogonality of the polynomials' derivatives can easily be
335 -- tested. The indices are shifted by one so that k=0 is the first
340 -- >>> import qualified Algebra.Absolute as Absolute ( abs )
342 -- >>> let expected = 11.5757525403319
343 -- >>> let actual = big_N' 3 1.5 :: Double
344 -- >>> Absolute.abs (actual - expected) < 1e-10
347 big_N' :: forall a. (Algebraic.C a, RealField.C a) => Integer -> a -> a
349 | k < 0 = error "requested a negative basis function"
350 | k == 0 = negate ( one / (fromInteger 2))
351 | k == 1 = one / (fromInteger 2)
352 | otherwise = coeff * ( legendre k x )
355 coeff = sqrt ((two*(fromInteger k) + one) / two) :: a
358 -- | The matrix of (N_i' * N_j') functions used in the integrand of
359 -- the stiffness matrix.
360 big_N's_matrix :: (Arity m, Arity n, Algebraic.C a, RealField.C a)
365 lambda i j x = (big_N' (toInteger i) x) * (big_N' (toInteger j) x)
368 big_K_elem :: forall m n l a b.
369 (Arity l, Arity m, Arity n,
370 Algebraic.C a, RealField.C a, ToRational.C a)
378 big_K_elem pde params _ k cur_K _ =
379 ifoldl2 accum cur_K (big_N's_matrix :: Mat m (S n) (a -> a))
381 accum :: Int -> Int -> Mat l l a -> (a -> a) -> Mat l l a
382 accum i j prev_K these_N's =
386 (x1,x2) = (mesh params) !!! (k,0)
387 q = affine_inv (x1,x2)
388 integrand x = ((big_A pde) (q x)) * (these_N's x)
389 -- The pointer matrix numbers from 1 so subtract one here to
390 -- get the right index.
391 row_idx = ((pointer params) !!! (k,i)) - 1
392 col_idx = ((pointer params) !!! (k,j)) - 1
393 integral = (two/(x2 - x1))* (gaussian integrand)
394 this_K = set_idx zero (row_idx, col_idx) integral
398 -- | Compute the \"big K\" stiffness matrix. There are three
399 -- parameters needed for K, namely i,j,k so a fold over a matrix will
400 -- not do. This little gimmick simulates a three-index fold by doing a
401 -- two-index fold over a row of the proper dimensions.
405 -- >>> import Linear.Matrix ( Col4, frobenius_norm, fromList )
406 -- >>> import Naturals ( N3, N4, N7 )
408 -- >>> let big_A = const (1::Double)
409 -- >>> let c x = sin x
410 -- >>> let f x = x*(sin x)
411 -- >>> let bdy = Left (Dirichlet (0,1::Double))
412 -- >>> let pde = PDE big_A c f bdy
414 -- >>> let i1 = (0.0,1/3)
415 -- >>> let i2 = (1/3,2/3)
416 -- >>> let i3 = (2/3,4/5)
417 -- >>> let i4 = (4/5,1.0)
418 -- >>> let mesh = fromList [[i1], [i2], [i3], [i4]] :: Col4 (Double,Double)
419 -- >>> let pvec = fromList [[2],[3],[2],[1]] :: Col4 Int
420 -- >>> let params = Params mesh pvec :: Params N4 N3 N7 Double
422 -- >>> let k1 = [6, -3, 0, 0, 0, 0, 0] :: [Double]
423 -- >>> let k2 = [-3, 10.5, -7.5, 0, 0, 0, 0] :: [Double]
424 -- >>> let k3 = [0, -7.5, 12.5, 0, 0, 0, 0] :: [Double]
425 -- >>> let k4 = [0, 0, 0, 6, 0, 0, 0] :: [Double]
426 -- >>> let k5 = [0, 0, 0, 0, 6, 0, 0] :: [Double]
427 -- >>> let k6 = [0, 0, 0, 0, 0, 6, 0] :: [Double]
428 -- >>> let k7 = [0, 0, 0, 0, 0, 0, 15] :: [Double]
429 -- >>> let expected = fromList [k1,k2,k3,k4,k5,k6,k7] :: Mat N7 N7 Double
430 -- >>> let actual = big_K pde params
431 -- >>> frobenius_norm (actual - expected) < 1e-10
434 big_K :: forall m n l a.
435 (Arity l, Arity m, Arity n,
436 Algebraic.C a, RealField.C a, ToRational.C a)
441 ifoldl2 (big_K_elem pde params) zero col_idxs
443 col_idxs = fromList [map fromInteger [0..]] :: Row m a
448 big_M_elem :: forall m n l a b.
449 (Arity l, Arity m, Arity n,
450 Algebraic.C a, RealField.C a, ToRational.C a)
458 big_M_elem pde params _ k cur_M _ =
459 ifoldl2 accum cur_M (big_Ns_matrix :: Mat m (S n) (a -> a))
461 accum :: Int -> Int -> Mat l l a -> (a -> a) -> Mat l l a
462 accum i j prev_M these_Ns =
466 (x1,x2) = (mesh params) !!! (k,0)
467 q = affine_inv (x1,x2)
468 integrand x = ((c pde) (q x)) * (these_Ns x)
469 -- The pointer matrix numbers from 1 so subtract one here to
470 -- get the right index.
471 row_idx = ((pointer params) !!! (k,i)) - 1
472 col_idx = ((pointer params) !!! (k,j)) - 1
473 integral = (x2 - x1)*(gaussian integrand) / two
474 this_M = set_idx zero (row_idx, col_idx) integral
477 -- | Compute the \"big M\" mass matrix. There are three
478 -- parameters needed for M, namely i,j,k so a fold over a matrix will
479 -- not do. This little gimmick simulates a three-index fold by doing a
480 -- two-index fold over a row of the proper dimensions.
484 -- >>> import Linear.Matrix ( Col4, frobenius_norm, fromList )
485 -- >>> import Naturals ( N3, N4, N7 )
487 -- >>> let big_A = const (1::Double)
488 -- >>> let c x = sin x
489 -- >>> let f x = x*(sin x)
490 -- >>> let bdy = Left (Dirichlet (0,1::Double))
491 -- >>> let pde = PDE big_A c f bdy
493 -- >>> let i1 = (0.0,1/3)
494 -- >>> let i2 = (1/3,2/3)
495 -- >>> let i3 = (2/3,4/5)
496 -- >>> let i4 = (4/5,1.0)
497 -- >>> let mesh = fromList [[i1], [i2], [i3], [i4]] :: Col4 (Double,Double)
498 -- >>> let pvec = fromList [[2],[3],[2],[1]] :: Col4 Int
499 -- >>> let params = Params mesh pvec :: Params N4 N3 N7 Double
501 -- >>> let m1 = [0.0723,0.0266,0,-0.0135,-0.0305,0.0058,0] :: [Double]
502 -- >>> let m2 = [0.0266,0.0897,0.0149,0,-0.0345,-0.0109,-0.0179] :: [Double]
503 -- >>> let m3 = [0,0.0149,0.0809,0,0,0,-0.0185] :: [Double]
504 -- >>> let m4 = [-0.0135,0,0,0.0110,0,0,0] :: [Double]
505 -- >>> let m5 = [-0.0305,-0.0345,0,0,0.0319,0.0018,0] :: [Double]
506 -- >>> let m6 = [0.0058,-0.0109,0,0,0.0018,0.0076,0] :: [Double]
507 -- >>> let m7 = [0,-0.0179,-0.0185,0,0,0,0.0178] :: [Double]
509 -- >>> let expected = fromList [m1,m2,m3,m4,m5,m6,m7] :: Mat N7 N7 Double
510 -- >>> let actual = big_M pde params
511 -- >>> frobenius_norm (actual - expected) < 1e-3
514 big_M :: forall m n l a.
515 (Arity l, Arity m, Arity n,
516 Algebraic.C a, RealField.C a, ToRational.C a)
521 ifoldl2 (big_M_elem pde params) zero col_idxs
523 col_idxs = fromList [map fromInteger [0..]] :: Row m a
527 -- | Determine the coefficient vector @x@ from the system @(K + M)x = F@.
531 -- >>> import Linear.Matrix ( Col4, Col7, frobenius_norm, fromList )
532 -- >>> import Naturals ( N3, N4, N7 )
534 -- >>> let big_A = const (1::Double)
535 -- >>> let c x = sin x
536 -- >>> let f x = x*(sin x)
537 -- >>> let bdy = Left (Dirichlet (0,1::Double))
538 -- >>> let pde = PDE big_A c f bdy
540 -- >>> let i1 = (0.0,1/3)
541 -- >>> let i2 = (1/3,2/3)
542 -- >>> let i3 = (2/3,4/5)
543 -- >>> let i4 = (4/5,1.0)
544 -- >>> let mesh = fromList [[i1], [i2], [i3], [i4]] :: Col4 (Double,Double)
545 -- >>> let pvec = fromList [[2],[3],[2],[1]] :: Col4 Int
546 -- >>> let params = Params mesh pvec :: Params N4 N3 N7 Double
548 -- >>> let c1 = [0.02366220347687] :: [Double]
549 -- >>> let c2 = [0.03431630082636] :: [Double]
550 -- >>> let c3 = [0.02841800893264] :: [Double]
551 -- >>> let c4 = [-0.00069489654996] :: [Double]
552 -- >>> let c5 = [-0.00518637005151] :: [Double]
553 -- >>> let c6 = [-0.00085028505337] :: [Double]
554 -- >>> let c7 = [-0.00170478210110] :: [Double]
555 -- >>> let expected = fromList [c1,c2,c3,c4,c5,c6,c7] :: Col7 Double
556 -- >>> let actual = coefficients pde params
557 -- >>> frobenius_norm (actual - expected) < 1e-8
560 coefficients :: forall m n l a.
561 (Arity m, Arity n, Arity l,
562 Algebraic.C a, Eq a, RealField.C a, ToRational.C a)
564 -> Params m n (S l) a
566 coefficients pde params =
567 solve_positive_definite matrix b
569 matrix = (big_K pde params) + (big_M pde params)