Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
ElmerCSC
GitHub Repository: ElmerCSC/elmerfem
Path: blob/devel/fem/src/LinearAlgebra.F90
3203 views
1
!/*****************************************************************************/
2
! *
3
! * Elmer, A Finite Element Software for Multiphysical Problems
4
! *
5
! * Copyright 1st April 1995 - , CSC - IT Center for Science Ltd., Finland
6
! *
7
! * This library is free software; you can redistribute it and/or
8
! * modify it under the terms of the GNU Lesser General Public
9
! * License as published by the Free Software Foundation; either
10
! * version 2.1 of the License, or (at your option) any later version.
11
! *
12
! * This library is distributed in the hope that it will be useful,
13
! * but WITHOUT ANY WARRANTY; without even the implied warranty of
14
! * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
15
! * Lesser General Public License for more details.
16
! *
17
! * You should have received a copy of the GNU Lesser General Public
18
! * License along with this library (in file ../LGPL-2.1); if not, write
19
! * to the Free Software Foundation, Inc., 51 Franklin Street,
20
! * Fifth Floor, Boston, MA 02110-1301 USA
21
! *
22
! *****************************************************************************/
23
!
24
!/******************************************************************************
25
! *
26
! * Authors: Juha Ruokolainen
27
! * Email: [email protected]
28
! * Web: http://www.csc.fi/elmer
29
! * Address: CSC - IT Center for Science Ltd.
30
! * Keilaranta 14
31
! * 02101 Espoo, Finland
32
! *
33
! * Original Date: 01 Oct 1996
34
! *
35
! *****************************************************************************/
36
37
!> \ingroup ElmerLib
38
!> \{
39
40
!------------------------------------------------------------------------------
41
!> Linear Algebra: LU-decomposition & matrix inverse & nonsymmetric full
42
!> matrix eigenvalues (don't use this for anything big, use for example
43
!> LAPACK routines instead...)
44
!------------------------------------------------------------------------------
45
MODULE LinearAlgebra
46
47
USE Types
48
IMPLICIT NONE
49
50
CONTAINS
51
52
SUBROUTINE InvertMatrix( A,n )
53
54
INTEGER :: n
55
REAL(KIND=dp) :: A(:,:)
56
57
REAL(KIND=dp) :: s
58
INTEGER :: i,j,k
59
INTEGER :: pivot(n)
60
LOGICAL :: erroneous
61
62
! /*
63
! * AP = LU
64
! */
65
CALL LUDecomp( a,n,pivot,erroneous )
66
67
IF (erroneous) CALL Fatal('InvertMatrix', 'inversion needs successfull LU-decomposition')
68
69
! /*
70
! * INV(U)
71
! */
72
DO i=N-1,1,-1
73
DO j=N,i+1,-1
74
s = -A(i,j)
75
DO k=i+1,j-1
76
s = s - A(i,k)*A(k,j)
77
END DO
78
A(i,j) = s
79
END DO
80
END DO
81
82
! /*
83
! * INV(L)
84
! */
85
DO i=n-1,1,-1
86
DO j=n,i+1,-1
87
s = 0.D00
88
DO k=i+1,j
89
s = s - A(j,k)*A(k,i)
90
END DO
91
A(j,i) = A(i,i)*s
92
END DO
93
END DO
94
95
! /*
96
! * A = INV(AP)
97
! */
98
DO i=1,n
99
DO j=1,n
100
s = 0.0D0
101
DO k=MAX(i,j),n
102
IF ( k /= i ) THEN
103
s = s + A(i,k)*A(k,j)
104
ELSE
105
s = s + A(k,j)
106
END IF
107
END DO
108
A(i,j) = s
109
END DO
110
END DO
111
112
! /*
113
! * A = INV(A) (at last)
114
! */
115
DO i=n,1,-1
116
IF ( pivot(i) /= i ) THEN
117
DO j = 1,n
118
s = A(i,j)
119
A(i,j) = A(pivot(i),j)
120
A(pivot(i),j) = s
121
END DO
122
END IF
123
END DO
124
125
END SUBROUTINE InvertMatrix
126
127
128
SUBROUTINE LUSolve( n,A,x,pivot_in )
129
REAL(KIND=dp) :: A(:,:)
130
REAL(KIND=dp) :: x(n)
131
INTEGER :: n
132
INTEGER, OPTIONAL :: pivot_in(:)
133
134
REAL(KIND=dp) :: s
135
INTEGER :: i,j
136
INTEGER :: pivot(n)
137
LOGICAL :: erroneous
138
139
! /*
140
! * AP = LU
141
! */
142
IF(PRESENT(pivot_in)) THEN
143
Pivot(1:n) = Pivot_in(1:n)
144
ELSE
145
CALL LUDecomp( A,n,pivot,erroneous )
146
IF (erroneous) CALL Fatal('LUSolve', 'LU-decomposition fails')
147
END IF
148
149
! Forward substitute
150
DO i=1,n
151
s = x(i)
152
DO j=1,i-1
153
s = s - A(i,j) * x(j)
154
END DO
155
x(i) = A(i,i) * s
156
END DO
157
158
!
159
! Backward substitute (solve x from Ux = z)
160
DO i=n,1,-1
161
s = x(i)
162
DO j=i+1,n
163
s = s - A(i,j) * x(j)
164
END DO
165
x(i) = s
166
END DO
167
168
DO i=n,1,-1
169
IF ( pivot(i) /= i ) THEN
170
s = x(i)
171
x(i) = x(pivot(i))
172
x(pivot(i)) = s
173
END IF
174
END DO
175
176
END SUBROUTINE LUSolve
177
178
179
!----------------------------------------------------------------------
180
!> LU- decomposition by gaussian elimination. Row pivoting is used.
181
!
182
!> result : AP = L'U ; L' = LD; pivot[i] is the swapped column number
183
!> for column i.
184
!
185
!> Result is stored in place of original matrix.
186
!----------------------------------------------------------------------
187
SUBROUTINE LUDecomp( a,n,pivot,erroneous )
188
189
REAL(KIND=dp), DIMENSION (:,:) :: a
190
INTEGER :: n
191
INTEGER, DIMENSION (:) :: pivot
192
LOGICAL, OPTIONAL :: erroneous
193
194
INTEGER :: i,j,k,l
195
REAL(KIND=dp) :: swap
196
197
IF (PRESENT(erroneous)) erroneous = .FALSE.
198
199
DO i=1,n
200
j = i
201
DO k=i+1,n
202
IF ( ABS(A(i,k)) > ABS(A(i,j)) ) j = k
203
END DO
204
205
IF ( ABS(A(i,j)) == 0.0d0) THEN
206
CALL Error( 'LUDecomp', 'Matrix is singular.' )
207
IF (PRESENT(erroneous)) erroneous = .TRUE.
208
RETURN
209
END IF
210
211
pivot(i) = j
212
213
IF ( j /= i ) THEN
214
DO k=1,i
215
swap = A(k,j)
216
A(k,j) = A(k,i)
217
A(k,i) = swap
218
END DO
219
END IF
220
221
DO k=i+1,n
222
A(i,k) = A(i,k) / A(i,i)
223
END DO
224
225
DO k=i+1,n
226
IF ( j /= i ) THEN
227
swap = A(k,i)
228
A(k,i) = A(k,j)
229
A(k,j) = swap
230
END IF
231
232
DO l=i+1,n
233
A(k,l) = A(k,l) - A(k,i) * A(i,l)
234
END DO
235
END DO
236
END DO
237
238
pivot(n) = n
239
IF ( ABS(A(n,n)) == 0.0d0 ) THEN
240
CALL Error( 'LUDecomp', 'Matrix is (at least almost) singular.' )
241
END IF
242
243
DO i=1,n
244
IF ( ABS(A(i,i)) == 0.0d0 ) THEN
245
CALL Error( 'LUSolve', 'Matrix is singular.' )
246
IF (PRESENT(erroneous)) erroneous = .TRUE.
247
RETURN
248
END IF
249
A(i,i) = 1.0_dp / A(i,i)
250
END DO
251
252
END SUBROUTINE LUDecomp
253
254
255
256
SUBROUTINE ComplexInvertMatrix( A,n )
257
258
COMPLEX(KIND=dp), DIMENSION(:,:) :: A
259
INTEGER :: n
260
261
COMPLEX(KIND=dp) :: s
262
INTEGER :: i,j,k
263
INTEGER :: pivot(n)
264
LOGICAL :: erroneous
265
266
! /*
267
! * AP = LU
268
! */
269
CALL ComplexLUDecomp( a,n,pivot,erroneous )
270
271
IF (erroneous) CALL Fatal('ComplexInvertMatrix', 'inversion needs successfull LU-decomposition')
272
273
DO i=1,n
274
IF ( ABS(A(i,i))==0.0d0 ) THEN
275
CALL Error( 'ComplexInvertMatrix', 'Matrix is singular.' )
276
RETURN
277
END IF
278
A(i,i) = 1.0D0/A(i,i)
279
END DO
280
281
! /*
282
! * INV(U)
283
! */
284
DO i=N-1,1,-1
285
DO j=N,i+1,-1
286
s = -A(i,j)
287
DO k=i+1,j-1
288
s = s - A(i,k)*A(k,j)
289
END DO
290
A(i,j) = s
291
END DO
292
END DO
293
294
! /*
295
! * INV(L)
296
! */
297
DO i=n-1,1,-1
298
DO j=n,i+1,-1
299
s = 0.D00
300
DO k=i+1,j
301
s = s - A(j,k)*A(k,i)
302
END DO
303
A(j,i) = A(i,i)*s
304
END DO
305
END DO
306
307
! /*
308
! * A = INV(AP)
309
! */
310
DO i=1,n
311
DO j=1,n
312
s = 0.0D0
313
DO k=MAX(i,j),n
314
IF ( k /= i ) THEN
315
s = s + A(i,k)*A(k,j)
316
ELSE
317
s = s + A(k,j)
318
END IF
319
END DO
320
A(i,j) = s
321
END DO
322
END DO
323
324
! /*
325
! * A = INV(A) (at last)
326
! */
327
DO i=n,1,-1
328
IF ( pivot(i) /= i ) THEN
329
DO j = 1,n
330
s = A(i,j)
331
A(i,j) = A(pivot(i),j)
332
A(pivot(i),j) = s
333
END DO
334
END IF
335
END DO
336
337
END SUBROUTINE ComplexInvertMatrix
338
339
!-------------------------------------------------------------------------
340
!> LU- decomposition by gaussian elimination for complex valued linear system.
341
!> Row pivoting is used.
342
!
343
!> result : AP = L'U ; L' = LD; pivot[i] is the swapped column number
344
!> for column i.
345
!
346
!> Result is stored in place of original matrix.
347
!-------------------------------------------------------------------------
348
349
SUBROUTINE ComplexLUDecomp( a,n,pivot,erroneous )
350
351
COMPLEX(KIND=dp), DIMENSION (:,:) :: a
352
INTEGER :: n
353
INTEGER, DIMENSION (:) :: pivot
354
LOGICAL, OPTIONAL :: erroneous
355
356
INTEGER :: i,j,k,l
357
COMPLEX(KIND=dp) :: swap
358
359
IF (PRESENT(erroneous)) erroneous = .FALSE.
360
361
DO i=1,n
362
j = i
363
DO k=i+1,n
364
IF ( ABS(A(i,k)) > ABS(A(i,j)) ) j = k
365
END DO
366
367
IF ( ABS(A(i,j))==0.0d0 ) THEN
368
CALL Error( 'ComplexLUDecomp', 'Matrix is singular.' )
369
IF (PRESENT(erroneous)) erroneous = .TRUE.
370
RETURN
371
END IF
372
373
pivot(i) = j
374
375
IF ( j /= i ) THEN
376
DO k=1,i
377
swap = A(k,j)
378
A(k,j) = A(k,i)
379
A(k,i) = swap
380
END DO
381
END IF
382
383
DO k=i+1,n
384
A(i,k) = A(i,k) / A(i,i)
385
END DO
386
387
DO k=i+1,n
388
IF ( j /= i ) THEN
389
swap = A(k,i)
390
A(k,i) = A(k,j)
391
A(k,j) = swap
392
END IF
393
394
DO l=i+1,n
395
A(k,l) = A(k,l) - A(k,i) * A(i,l)
396
END DO
397
END DO
398
END DO
399
400
pivot(n) = n
401
IF ( ABS(A(n,n))==0.0d0 ) THEN
402
CALL Error( 'ComplexLUDecomp', 'Matrix is (at least almost) singular.' )
403
END IF
404
END SUBROUTINE ComplexLUDecomp
405
406
407
!------------------------------------------------------------------------------
408
!> Solves a 2 x 2 linear system.
409
!------------------------------------------------------------------------------
410
SUBROUTINE SolveLinSys2x2( A, x, b )
411
!------------------------------------------------------------------------------
412
REAL(KIND=dp), INTENT(out) :: x(:)
413
REAL(KIND=dp), INTENT(in) :: A(:,:),b(:)
414
!------------------------------------------------------------------------------
415
REAL(KIND=dp) :: detA
416
!------------------------------------------------------------------------------
417
detA = A(1,1) * A(2,2) - A(1,2) * A(2,1)
418
419
IF ( detA == 0.0d0 ) THEN
420
WRITE( Message, * ) 'Singular matrix, sorry!'
421
CALL Error( 'SolveLinSys2x2', Message )
422
RETURN
423
END IF
424
425
detA = 1.0d0 / detA
426
x(1) = detA * (A(2,2) * b(1) - A(1,2) * b(2))
427
x(2) = detA * (A(1,1) * b(2) - A(2,1) * b(1))
428
!------------------------------------------------------------------------------
429
END SUBROUTINE SolveLinSys2x2
430
!------------------------------------------------------------------------------
431
432
433
!------------------------------------------------------------------------------
434
!> Solves a 3 x 3 linear system.
435
!------------------------------------------------------------------------------
436
SUBROUTINE SolveLinSys3x3( A, x, b )
437
!------------------------------------------------------------------------------
438
REAL(KIND=dp), INTENT(out) :: x(:)
439
REAL(KIND=dp), INTENT(in) :: A(:,:),b(:)
440
!------------------------------------------------------------------------------
441
REAL(KIND=dp) :: C(2,2),y(2),g(2),s,t,q
442
!------------------------------------------------------------------------------
443
444
IF ( ABS(A(1,1))>ABS(A(1,2)) .AND. ABS(A(1,1))>ABS(A(1,3)) ) THEN
445
q = 1.0d0 / A(1,1)
446
s = q * A(2,1)
447
t = q * A(3,1)
448
C(1,1) = A(2,2) - s * A(1,2)
449
C(1,2) = A(2,3) - s * A(1,3)
450
C(2,1) = A(3,2) - t * A(1,2)
451
C(2,2) = A(3,3) - t * A(1,3)
452
453
g(1) = b(2) - s * b(1)
454
g(2) = b(3) - t * b(1)
455
CALL SolveLinSys2x2( C,y,g )
456
457
x(2) = y(1)
458
x(3) = y(2)
459
x(1) = q * ( b(1) - A(1,2) * x(2) - A(1,3) * x(3) )
460
ELSE IF ( ABS(A(1,2)) > ABS(A(1,3)) ) THEN
461
q = 1.0d0 / A(1,2)
462
s = q * A(2,2)
463
t = q * A(3,2)
464
C(1,1) = A(2,1) - s * A(1,1)
465
C(1,2) = A(2,3) - s * A(1,3)
466
C(2,1) = A(3,1) - t * A(1,1)
467
C(2,2) = A(3,3) - t * A(1,3)
468
469
g(1) = b(2) - s * b(1)
470
g(2) = b(3) - t * b(1)
471
CALL SolveLinSys2x2( C,y,g )
472
473
x(1) = y(1)
474
x(3) = y(2)
475
x(2) = q * ( b(1) - A(1,1) * x(1) - A(1,3) * x(3) )
476
ELSE
477
q = 1.0d0 / A(1,3)
478
s = q * A(2,3)
479
t = q * A(3,3)
480
C(1,1) = A(2,1) - s * A(1,1)
481
C(1,2) = A(2,2) - s * A(1,2)
482
C(2,1) = A(3,1) - t * A(1,1)
483
C(2,2) = A(3,2) - t * A(1,2)
484
485
g(1) = b(2) - s * b(1)
486
g(2) = b(3) - t * b(1)
487
CALL SolveLinSys2x2( C,y,g )
488
489
x(1) = y(1)
490
x(2) = y(2)
491
x(3) = q * ( b(1) - A(1,1) * x(1) - A(1,2) * x(2) )
492
END IF
493
!------------------------------------------------------------------------------
494
END SUBROUTINE SolveLinSys3x3
495
!------------------------------------------------------------------------------
496
497
498
!> Solves a small dense linear system using Lapack routines
499
!------------------------------------------------------------------------------
500
SUBROUTINE SolveLinSys( A, x, n )
501
!------------------------------------------------------------------------------
502
INTEGER :: n
503
REAL(KIND=dp) :: A(n,n), x(n), b(n)
504
505
INTERFACE
506
SUBROUTINE SolveLapack( N,A,x )
507
INTEGER N
508
DOUBLE PRECISION A(n*n),x(n)
509
END SUBROUTINE
510
END INTERFACE
511
512
!------------------------------------------------------------------------------
513
SELECT CASE(n)
514
CASE(1)
515
x(1) = x(1) / A(1,1)
516
CASE(2)
517
b = x
518
CALL SolveLinSys2x2(A,x,b)
519
CASE(3)
520
b = x
521
CALL SolveLinSys3x3(A,x,b)
522
CASE DEFAULT
523
CALL SolveLapack(n,A,x)
524
END SELECT
525
!------------------------------------------------------------------------------
526
END SUBROUTINE SolveLinSys
527
!------------------------------------------------------------------------------
528
529
530
!------------------------------------------------------------------------------
531
SUBROUTINE InvertMatrix3x3( G,GI,detG )
532
!------------------------------------------------------------------------------
533
REAL(KIND=dp) :: G(3,3),GI(3,3)
534
REAL(KIND=dp) :: detG, s
535
!------------------------------------------------------------------------------
536
s = 1.0 / DetG
537
538
GI(1,1) = s * (G(2,2)*G(3,3) - G(3,2)*G(2,3));
539
GI(2,1) = -s * (G(2,1)*G(3,3) - G(3,1)*G(2,3));
540
GI(3,1) = s * (G(2,1)*G(3,2) - G(3,1)*G(2,2));
541
542
GI(1,2) = -s * (G(1,2)*G(3,3) - G(3,2)*G(1,3));
543
GI(2,2) = s * (G(1,1)*G(3,3) - G(3,1)*G(1,3));
544
GI(3,2) = -s * (G(1,1)*G(3,2) - G(3,1)*G(1,2));
545
546
GI(1,3) = s * (G(1,2)*G(2,3) - G(2,2)*G(1,3));
547
GI(2,3) = -s * (G(1,1)*G(2,3) - G(2,1)*G(1,3));
548
GI(3,3) = s * (G(1,1)*G(2,2) - G(2,1)*G(1,2));
549
!------------------------------------------------------------------------------
550
END SUBROUTINE InvertMatrix3x3
551
!------------------------------------------------------------------------------
552
553
554
!------------------------------------------------------------------------------
555
! Quadratic precision version of the previous routine!
556
!------------------------------------------------------------------------------
557
SUBROUTINE InvertMatrix3x3QP( G,GI,detG )
558
!------------------------------------------------------------------------------
559
INTEGER, PARAMETER :: qp = SELECTED_REAL_KIND(24)
560
REAL(KIND=qp) :: G(3,3),GI(3,3)
561
REAL(KIND=qp) :: detG, s
562
!------------------------------------------------------------------------------
563
s = 1.0 / DetG
564
565
GI(1,1) = s * (G(2,2)*G(3,3) - G(3,2)*G(2,3));
566
GI(2,1) = -s * (G(2,1)*G(3,3) - G(3,1)*G(2,3));
567
GI(3,1) = s * (G(2,1)*G(3,2) - G(3,1)*G(2,2));
568
569
GI(1,2) = -s * (G(1,2)*G(3,3) - G(3,2)*G(1,3));
570
GI(2,2) = s * (G(1,1)*G(3,3) - G(3,1)*G(1,3));
571
GI(3,2) = -s * (G(1,1)*G(3,2) - G(3,1)*G(1,2));
572
573
GI(1,3) = s * (G(1,2)*G(2,3) - G(2,2)*G(1,3));
574
GI(2,3) = -s * (G(1,1)*G(2,3) - G(2,1)*G(1,3));
575
GI(3,3) = s * (G(1,1)*G(2,2) - G(2,1)*G(1,2));
576
!------------------------------------------------------------------------------
577
END SUBROUTINE InvertMatrix3x3QP
578
!------------------------------------------------------------------------------
579
580
!------------------------------------------------------------------------------
581
!> Compute the value of 3x3 determinant
582
!------------------------------------------------------------------------------
583
FUNCTION Det3x3( A ) RESULT ( val )
584
!------------------------------------------------------------------------------
585
REAL(KIND=dp) :: A(3,3)
586
REAL(KIND=dp) :: val
587
588
val = A(1,1) * ( A(2,2) * A(3,3) - A(2,3) * A(3,2) ) &
589
- A(1,2) * ( A(2,1) * A(3,3) - A(2,3) * A(3,1) ) &
590
+ A(1,3) * ( A(2,1) * A(3,2) - A(2,2) * A(3,1) )
591
!------------------------------------------------------------------------------
592
END FUNCTION Det3x3
593
!------------------------------------------------------------------------------
594
595
! --------------------------------------------------
596
!> Solve eigenvalues of a nonsymmetric matrix A(n,n)
597
!> The matrix is modified in the process.
598
! --------------------------------------------------
599
SUBROUTINE EigenValues( A, n, Vals )
600
IMPLICIT NONE
601
REAL(KIND=dp) :: A(:,:)
602
INTEGER :: n
603
COMPLEX(KIND=dp) :: Vals(:)
604
605
INTEGER :: iter, i, j, k
606
REAL(KIND=dp) :: b, s, t
607
608
IF (n==1 ) THEN
609
Vals(1) = A(1,1)
610
RETURN
611
END IF
612
613
CALL Hesse(A, n)
614
615
DO iter=1,1000
616
DO i=1,n-1
617
s = AEPS * ( ABS(A(i,i))+ABS(A(i+1,i+1)) )
618
IF ( ABS(A(i+1,i)) < s ) A(i+1,i) = 0.0
619
END DO
620
621
i = 1
622
DO WHILE( .TRUE. )
623
DO j=i,n-1
624
IF (A(j+1,j) /= 0 ) EXIT
625
END DO
626
627
DO k=j,n-1
628
IF (A(k+1,k) == 0 ) EXIT
629
END DO
630
i = k;
631
IF ( i >= n .OR. k-j+1 >= 3 ) EXIT
632
END DO
633
634
IF (k-j+1 < 3) EXIT
635
CALL Francis(A(j:,j:), k-j+1);
636
END DO
637
638
j = 0
639
i = 1
640
DO WHILE( i<n )
641
IF (A(i+1,i) == 0) THEN
642
j = j + 1
643
Vals(j) = A(i,i)
644
ELSE
645
b = ( A(i,i)+A(i+1,i+1) ); s=b*b
646
t = A(i,i)*A(i+1,i+1) - A(i,i+1)*A(i+1,i)
647
s = s - 4*t
648
IF (s >= 0) THEN
649
s = SQRT(s)
650
j = j + 1
651
IF ( b>0 ) THEN
652
Vals(j) = (b+s)/2
653
ELSE
654
Vals(j) = 2*t/(b-s)
655
END IF
656
j = j + 1
657
IF ( b > 0 ) THEN
658
Vals(j) = 2*t/(b+s)
659
ELSE
660
Vals(j) = (b-s)/2
661
END IF
662
ELSE
663
j = j + 1
664
Vals(j) = CMPLX(b, SQRT(-s),KIND=dp)/2
665
j = j + 1
666
Vals(j) = CMPLX(b,-SQRT(-s),KIND=dp)/2
667
END IF
668
i = i + 1
669
END IF
670
i = i + 1
671
END DO
672
673
IF (A(n, n-1) == 0) THEN
674
j = j + 1
675
Vals(j) = A(n,n)
676
END IF
677
678
CONTAINS
679
680
SUBROUTINE vbcalc( x,v,b,beg,end )
681
IMPLICIT NONE
682
REAL(KIND=dp) :: x(:),v(:),b
683
INTEGER :: beg, end
684
685
REAL(KIND=dp) :: alpha,m
686
INTEGER :: i
687
688
m = MAXVAL( ABS(x(beg:end)) )
689
690
IF ( m == 0 ) THEN
691
v(beg:end) = 0
692
ELSE
693
alpha = 0
694
m = 1 / m
695
DO i=beg,end
696
v(i) = m*x(i)
697
alpha = alpha+v(i)*v(i)
698
END DO
699
alpha = SQRT(alpha)
700
b = 1/(alpha*(alpha+ABS(v(beg))))
701
IF ( v(beg) < 0 ) THEN
702
v(beg) = v(beg) - alpha
703
ELSE
704
v(beg) = v(beg) + alpha
705
END IF
706
END IF
707
END SUBROUTINE vbcalc
708
709
710
SUBROUTINE Hesse(H, dim)
711
IMPLICIT NONE
712
INTEGER :: dim
713
REAL(KIND=dp) :: H(:,:)
714
715
REAL(KIND=dp) :: b, s
716
INTEGER :: i, j, k;
717
REAL(KIND=dp), ALLOCATABLE :: v(:), x(:)
718
719
ALLOCATE( x(dim), v(dim) )
720
721
DO i=1,dim-1
722
DO j=i+1,dim
723
x(j) = H(j,i)
724
END DO
725
726
CALL vbcalc(x, v, b, i+1, dim )
727
728
IF (v(i+1) == 0) EXIT
729
730
DO j=i+2, dim
731
x(j) = v(j)/v(i+1)
732
v(j) = b*v(i+1)*v(j)
733
END DO
734
v(i+1) = b*v(i+1)*v(i+1)
735
736
DO j=1,dim
737
s = 0.0
738
DO k=i+1,dim
739
s = s + H(j,k)*v(k)
740
END DO
741
H(j,i+1) = H(j,i+1) - s
742
DO k=i+2,dim
743
H(j,k) = H(j,k) - s*x(k)
744
END DO
745
END DO
746
747
DO j=1,dim
748
s = H(i+1,j)
749
DO k=i+2,dim
750
s = s + H(k,j)*x(k)
751
END DO
752
DO k=i+1,dim
753
H(k,j) = H(k,j) - s*v(k)
754
END DO
755
END DO
756
H(i+2:dim,i) = 0
757
END DO
758
759
DEALLOCATE( x, v )
760
END SUBROUTINE Hesse
761
762
763
SUBROUTINE Francis( H, dim )
764
IMPLICIT NONE
765
INTEGER :: dim
766
REAL(KIND=dp) :: H(:,:)
767
768
INTEGER :: i, i1, j, k, m, n, end
769
REAL(KIND=dp) :: x(3), v(3), b, s, t, v0i
770
771
n = dim; m = n-1
772
773
t = H(m,m)*H(n,n) - H(m,n)*H(n,m)
774
s = H(m,m)+H(n,n)
775
776
x(1) = H(1,1)*H(1,1)+H(1,2)*H(2,1)-s*H(1,1)+t
777
x(2) = H(2,1)*(H(1,1)+H(2,2)-s)
778
x(3) = H(2,1)*H(3,2);
779
780
CALL vbcalc(x, v, b, 1, 3)
781
782
IF (v(1) == 0) RETURN
783
784
! DO i=2,3
785
! x(i) = v(i)/v(1);
786
! v(i) = b * v(1) * v(i);
787
! END DO
788
x(2) = v(2)/v(1)
789
v(2) = b*v(1)*v(2)
790
x(3) = v(3)/v(1)
791
v(3) = b*v(1)*v(3)
792
x(1) = 1; v(1) = b*v(1)*v(1)
793
794
DO i=1,dim
795
! DO j=1,3
796
! s = s + H(i,j) * v(j)
797
! END DO
798
s = H(i,1)*v(1) + H(i,2)*v(2) + H(i,3)*v(3)
799
800
! DO j=1,3
801
! H(i,j) = H(i,j) - s * x(j)
802
! END DO
803
H(i,1) = H(i,1) - s
804
H(i,2) = H(i,2) - s*x(2)
805
H(i,3) = H(i,3) - s*x(3)
806
END DO
807
808
DO i=1,dim
809
! s = 0.0d0
810
! DO j=1,3
811
! s = s + H(j,i) * x(j)
812
! END DO
813
s = h(1,i) + H(2,i)*x(2) + H(3,i)*x(3)
814
815
! DO j=1,3
816
! H(j,i) = H(j,i) - s * v(j)
817
! END DO
818
H(1,i) = H(1,i) - s*v(1)
819
H(2,i) = H(2,i) - s*v(2)
820
H(3,i) = H(3,i) - s*v(3)
821
END DO
822
823
DO i=1,dim-1
824
end = MIN(2, dim-i-1)
825
DO j=1,end+1
826
x(j) = H(i+j,i)
827
END DO
828
829
CALL vbcalc(x, v, b, 1, end+1)
830
831
IF (v(1) == 0) EXIT
832
833
DO j=2,end+1
834
x(j) = v(j)/v(1);
835
v(j) = b*v(1)*v(j)
836
END DO
837
x(1) = 1; v(1) = b*v(1)*v(1)
838
839
DO j=1,dim
840
s = 0.0
841
DO k=1,end+1
842
s = s + H(j,i+k)*v(k)
843
END DO
844
DO k=1,end+1
845
H(j,i+k) = H(j,i+k) - s*x(k)
846
END DO
847
END DO
848
849
DO j=1,dim
850
s = H(i+1,j)
851
DO k=2,end+1
852
s = s + H(i+k,j)*x(k)
853
END DO
854
DO k=1,end+1
855
H(i+k,j) = H(i+k,j) - s*v(k)
856
END DO
857
END DO
858
H(i+2:dim,i) = 0
859
END DO
860
END SUBROUTINE Francis
861
END SUBROUTINE EigenValues
862
863
END MODULE LinearAlgebra
864
865
!> \} ElmerLib
866
867