block_tridiag.F90       coverage:  0.00 %func     0.00 %block


     1) module Block_Tridiag_module
     2) 
     3) #include "petsc/finclude/petscsys.h"
     4) 
     5)   private
     6) 
     7)   public :: decbt, &
     8)             solbtf, &
     9)             solbtb
    10)      
    11) contains
    12) 
    13) ! ************************************************************************** !
    14) 
    15) subroutine decbt (m, n, ndim, a, b, c, ip, ier)
    16)       
    17)   use PFLOTRAN_Constants_module
    18) 
    19)   implicit none
    20) 
    21)   PetscInt :: m,n,ndim,ip(ndim,n),ier
    22)   PetscInt :: nm1,nm2,km1,i,j,k,l
    23)   PetscReal :: dp
    24)   PetscReal :: a(ndim,ndim,n), b(ndim,ndim,n), c(ndim,ndim,n)
    25)      
    26) !--------------------------------------
    27) ! block-tridiagonal matrix decomposition routine.
    28) ! written by a. c. hindmarsh.
    29) ! latest revision january 26, 1977  (ag)
    30) ! reference]  ucid-30150
    31) !             solution of block-tridiagonal systems of linear
    32) !             algebraic equations
    33) !             a.c. hindmarsh
    34) !             february 1977
    35) ! the input matrix contains three blocks of elements in each block-row,
    36) ! including blocks in the (1,3) and (n,n-2) block positions.
    37) ! decbt uses block gauss elimination and subroutines dec and sol
    38) ! for solution of blocks.  partial pivoting is done within
    39) ! block-rows only.
    40) ! input..
    41) !     m = order of each block.
    42) !     n = number of blocks in each direction of the matrix.
    43) !         n must be 4 or more.  the complete matrix has order m*n.
    44) !     a = m by m by n array containing diagonal blocks.
    45) !         a(i,j,k) contains the (i,j) element of the k-th block.
    46) !     b = m by m by n array containing the super-diagonal blocks
    47) !         (in b(*,*,k) for k = 1,...,n-1) and the block in the (n,n-2)
    48) !         block position (in b(*,*,n)).
    49) !     c = m by m by n array containing the subdiagonal blocks
    50) !         (in c(*,*,k) for k = 2,3,...,n) and the block in the
    51) !         (1,3) block position (in c(*,*,1)).
    52) !    ip = integer array of length m*n for working storage.
    53) ! output..
    54) ! a,b,c = m by m by n arrays containing the block lu decomposition
    55) !         of the input matrix.
    56) !    ip = m by n array of pivot information.  ip(*,k) contains
    57) !         information for the k-th digonal block.
    58) !   ier = 0  if no trouble occurred, or
    59) !       = -1 if the input value of m or n was illegal, or
    60) !       = k  if a singular matrix was found in the k-th diagonal block.
    61) ! use solbt to solve the associated linear system.
    62) ! decbt calls subroutines  dec(m,m0,a,ip,ier)  and  sol(m,m0,a,y,ip)
    63) ! for solution of m by m linear systems.
    64) !------------------------------------------------------------------
    65)       
    66)   if (m .lt. 1 .or. n .lt. 4) goto 210
    67)     nm1 = n - 1
    68)     nm2 = n - 2
    69) 
    70) ! process the first block-row. ------------------------------
    71)     call dec (m, ndim, a(:,:,1), ip, ier)
    72)     k = 1
    73)     if (ier .ne. 0) goto 200
    74)     do j = 1,m
    75)       call sol (m, ndim, a(:,:,1), b(1,j,1), ip)
    76)       call sol (m, ndim, a(:,:,1), c(1,j,1), ip)
    77)     enddo
    78) ! adjust b(*,*,2). -----------------------------------------------
    79)     do j = 1,m
    80)       do i = 1,m
    81)         dp = 0.d0
    82)         do l = 1,m
    83)           dp = dp + c(i,l,2)*c(l,j,1)
    84)         enddo
    85)         b(i,j,2) = b(i,j,2) - dp
    86)       enddo
    87)     enddo
    88) ! main loop.  process block-rows 2 to n-1. --------------
    89)     do k = 2,nm1
    90)       km1 = k - 1
    91)       do j = 1,m
    92)         do i = 1,m
    93)           dp = 0.d0
    94)           do l = 1,m
    95)             dp = dp + c(i,l,k)*b(l,j,km1)
    96)           enddo
    97)           a(i,j,k) = a(i,j,k) - dp
    98)         enddo
    99)       enddo
   100) 
   101)       call dec (m, ndim, a(1,1,k), ip(1,k), ier)
   102) 
   103)       if (ier .ne. 0) goto 200
   104)       do j = 1,m
   105)         call sol (m, ndim, a(1,1,k), b(1,j,k), ip(1,k))
   106)       enddo
   107)     enddo
   108) 
   109) ! process last block-row and return. ------------------
   110)     do j = 1,m
   111)       do i = 1,m
   112)         dp = 0.d0
   113)         do l = 1,m
   114)           dp = dp + b(i,l,n)*b(l,j,nm2)
   115)         enddo
   116)         c(i,j,n) = c(i,j,n) - dp
   117)       enddo
   118)     enddo
   119)     do j = 1,m
   120)       do i = 1,m
   121)         dp = 0.d0
   122)         do l = 1,m
   123)           dp = dp + c(i,l,n)*b(l,j,nm1)
   124)         enddo
   125)         a(i,j,n) = a(i,j,n) - dp
   126)       enddo
   127)     enddo
   128) 
   129)     call dec (m, ndim, a(1,1,n), ip(1,n), ier)
   130)     k = n
   131)     if (ier .ne. 0) goto 200
   132)     return
   133) 
   134) ! error returns. ------------------------------------
   135) 200  ier = k
   136)      return
   137) 
   138) 210  ier = -1
   139)      return
   140) 
   141) end subroutine decbt
   142) 
   143) ! ************************************************************************** !
   144) 
   145) subroutine solbt (m, n, ndim, a, b, c, ip, y)
   146) 
   147)   implicit none
   148) 
   149)   PetscInt :: m, n, ndim, ip(ndim,n)
   150)   PetscInt :: nm1, nm2, km1, i, j, k, kb, kp1
   151)   PetscReal :: a(ndim,ndim,n),b(ndim,ndim,n),c(ndim,ndim,n)
   152)   PetscReal :: y(ndim,n),dp
   153) 
   154) !---------------------
   155) ! solution of block-tridiagonal linear system.
   156) ! coefficient matrix must have been previously processed by decbt.
   157) ! m, n, a, b, c, and ip  must not have been changed since call to decbt
   158) ! written by a. c. hindmarsh.
   159) ! input..
   160) !     m = order of each block.
   161) !     n = number of blocks in each direction of matrix.
   162) ! a,b,c = m by m by n arrays containing block lu decomposition
   163) !         of coefficient matrix from decbt.
   164) !    ip = m by n integer array of pivot information from decbt.
   165) !     y = array of length m*n containg the right-hand side vector
   166) !         (treated as an m by n array here).
   167) ! output..
   168) !     y = solution vector, of length m*n.
   169) ! solbt makes calls to subroutine sol(m,m0,a,y,ip)
   170) ! for solution of m by m linear systems.
   171) !----------------------------------------------------------------
   172) 
   173)   nm1 = n - 1
   174)   nm2 = n - 2
   175) 
   176) ! forward solution sweep. ---------------------------------
   177)   call sol (m, ndim, a(:,:,1), y, ip)
   178)   do k = 2,nm1
   179)     km1 = k - 1
   180)     do i = 1,m
   181)       dp = 0.d0
   182)       do j = 1,m
   183)         dp = dp + c(i,j,k)*y(j,km1)
   184)       enddo
   185)       y(i,k) = y(i,k) - dp
   186)     enddo
   187)     call sol (m, ndim, a(1,1,k), y(1,k), ip(1,k))
   188)   enddo
   189) 
   190)   do i = 1,m
   191)     dp = 0.d0
   192)     do j = 1,m
   193)       dp = dp + c(i,j,n)*y(j,nm1) + b(i,j,n)*y(j,nm2)
   194)     enddo
   195)     y(i,n) = y(i,n) - dp
   196)   enddo
   197) 
   198)   call sol (m, ndim, a(1,1,n), y(1,n), ip(1,n))
   199) 
   200) ! backward solution sweep. ----------------------------------
   201) 
   202) !     do j=1,n
   203) !       print *,'solbt: ',m,n,ndim,j,y(1,j),ip(1,j),a(1,1,j),b(1,1,j),c(1,1,j)
   204) !     enddo
   205)       
   206)   do kb = 1,nm1
   207)     k = n - kb
   208)     kp1 = k + 1
   209)     do i = 1,m
   210)       dp = 0.d0
   211)       do j = 1,m
   212)         dp = dp + b(i,j,k)*y(j,kp1)
   213)       enddo
   214)       y(i,k) = y(i,k) - dp
   215)     enddo
   216)   enddo
   217) 
   218)   do i = 1,m
   219)     dp = 0.d0
   220)     do j = 1,m
   221)       dp = dp + c(i,j,1)*y(j,3)
   222)     enddo
   223)     y(i,1) = y(i,1) - dp
   224)   enddo
   225) 
   226)   return
   227) end subroutine solbt
   228) 
   229) ! ************************************************************************** !
   230) 
   231) subroutine solbtf (m, n, ndim, a, b, c, ip, y)
   232) 
   233)   implicit none
   234)   PetscInt :: m, n, ndim, ip(ndim,n)
   235)   PetscInt :: nm1, nm2, km1, i, j, k
   236)   PetscReal :: a(ndim,ndim,n),b(ndim,ndim,n),c(ndim,ndim,n)
   237)   PetscReal :: y(ndim,n),dp
   238) 
   239) !---------------------
   240) ! solution of block-tridiagonal linear system.
   241) ! coefficient matrix must have been previously processed by decbt.
   242) ! m, n, a, b, c, and ip  must not have been changed since call to decbt
   243) ! written by a. c. hindmarsh.
   244) ! input..
   245) !     m = order of each block.
   246) !     n = number of blocks in each direction of matrix.
   247) ! a,b,c = m by m by n arrays containing block lu decomposition
   248) !         of coefficient matrix from decbt.
   249) !    ip = m by n integer array of pivot information from decbt.
   250) !     y = array of length m*n containg the right-hand side vector
   251) !         (treated as an m by n array here). 
   252) ! output..
   253) !     y = solution vector, of length m*n.
   254) ! solbt makes calls to subroutine sol(m,m0,a,y,ip)
   255) ! for solution of m by m linear systems.
   256) !----------------------------------------------------------------
   257) 
   258)   nm1 = n - 1
   259)   nm2 = n - 2
   260) 
   261) ! forward solution sweep. ---------------------------------
   262)   call sol (m, ndim, a, y, ip)
   263)   do k = 2,nm1
   264)     km1 = k - 1
   265)     do i = 1,m
   266)       dp = 0.d0
   267)       do j = 1,m
   268)         dp = dp + c(i,j,k)*y(j,km1)
   269)       enddo
   270)       y(i,k) = y(i,k) - dp
   271)     enddo
   272)     call sol (m, ndim, a(1,1,k), y(1,k), ip(1,k))
   273)   enddo
   274) 
   275)   do i = 1,m
   276)     dp = 0.d0
   277)     do j = 1,m
   278)       dp = dp + c(i,j,n)*y(j,nm1) + b(i,j,n)*y(j,nm2)
   279)     enddo
   280)     y(i,n) = y(i,n) - dp
   281)   enddo
   282) 
   283)   call sol (m, ndim, a(1,1,n), y(1,n), ip(1,n))
   284) 
   285)   return
   286) end subroutine solbtf
   287) 
   288) ! ************************************************************************** !
   289) 
   290) subroutine solbtb (m, n, ndim, a, b, c, ip, y)
   291) 
   292)   implicit none
   293) 
   294)   PetscInt :: m, n, ndim, ip(ndim,n)
   295)   PetscInt :: nm1, i, j, k, kb, kp1
   296)   PetscReal :: a(ndim,ndim,n),b(ndim,ndim,n),c(ndim,ndim,n)
   297)   PetscReal :: y(ndim,n),dp
   298) 
   299) !---------------------
   300) ! solution of block-tridiagonal linear system.
   301) ! coefficient matrix must have been previously processed by decbt.
   302) ! m, n, a, b, c, and ip  must not have been changed since call to decbt
   303) ! written by a. c. hindmarsh.
   304) ! input..
   305) !     m = order of each block.
   306) !     n = number of blocks in each direction of matrix.
   307) ! a,b,c = m by m by n arrays containing block lu decomposition
   308) !         of coefficient matrix from decbt.
   309) !    ip = m by n integer array of pivot information from decbt.
   310) !     y = array of length m*n containg the right-hand side vector
   311) !         (treated as an m by n array here).
   312) ! output..
   313) !     y = solution vector, of length m*n.
   314) ! solbt makes calls to subroutine sol(m,m0,a,y,ip)
   315) ! for solution of m by m linear systems.
   316) !----------------------------------------------------------------
   317) 
   318) !     call sol (m, ndim, a(1,1,n), y(1,n), ip(1,n))
   319) 
   320)   nm1 = n - 1
   321) 
   322) !     do j=1,n
   323) !       print *,'solbtb: ',m,n,ndim,j,y(1,j),ip(1,j),a(1,1,j),b(1,1,j),c(1,1,j)
   324) !     enddo
   325) 
   326) ! backward solution sweep. ----------------------------------
   327)   do kb = 1,nm1
   328)     k = n - kb
   329)     kp1 = k + 1
   330)     do i = 1,m
   331)       dp = 0.d0
   332)       do j = 1,m
   333)         dp = dp + b(i,j,k)*y(j,kp1)
   334)       enddo
   335)       y(i,k) = y(i,k) - dp
   336)     enddo
   337)   enddo
   338) 
   339)   do i = 1,m
   340)     dp = 0.d0
   341)     do j = 1,m
   342)       dp = dp + c(i,j,1)*y(j,3)
   343)     enddo
   344)     y(i,1) = y(i,1) - dp
   345)   enddo
   346) 
   347)   return
   348) end subroutine solbtb
   349) 
   350) ! ************************************************************************** !
   351) 
   352) subroutine dec (n, ndim, a, ip, ier)
   353) 
   354)   implicit none
   355) 
   356)   PetscInt :: m, n, ndim, ip(n), ier
   357)   PetscInt :: nm1,i,j,k,kp1
   358)   PetscReal :: a(ndim,n), t
   359) 
   360) !---------------
   361) !  matrix triangularization by gauss elimination with partial pivoting.
   362) !  input..
   363) !     n = order of matrix.
   364) !     ndim = declared first dimension of array  a.
   365) !     a = matrix to be triangularized.
   366) !  output..
   367) !     a(i,j), i.le.j = upper triangular factor, u .
   368) !     a(i,j), i.gt.j = multipliers = lower triangular factor, i - l.
   369) !     ip(k), k.lt.n = index of k-th pivot row.
   370) !     ier = 0 if matrix a is nonsingular, or k if found to be
   371) !           singular at stage k.
   372) !  row interchanges are finished in u, only partly in l.
   373) !  use  sol  to obtain solution of linear system.
   374) !  if ier .ne. 0, a is singular, sol will divide by 0.d0.
   375) !------------------------------------------------------------
   376) !
   377) !  reference:  a. c. hindmarsh, l. j. sloan, k. w. fong, and
   378) !              g. h. rodrigue,
   379) !              dec/sol:  solution of dense systems of linear
   380) !              algebraic equations,
   381) !              lawrence livermore laboratory report ucid-30137,
   382) !              june 1976.
   383) !----------------------------------------------------------------
   384) 
   385)   ier = 0
   386)   if (n .eq. 1) goto 70
   387)   nm1 = n - 1
   388)   do k = 1,nm1
   389)     kp1 = k + 1
   390) 
   391) !  find the pivot in column k.  search rows k to n. -------
   392)     m = k
   393)     do i = kp1,n
   394)       if (abs(a(i,k)) .gt. abs(a(m,k))) m = i
   395)     enddo
   396)     ip(k) = m
   397) 
   398) !  interchange elements in rows k and m. -------------
   399)     t = a(m,k)
   400)     if (m .eq. k) goto 20
   401)     a(m,k) = a(k,k)
   402)     a(k,k) = t
   403)  20     if (t .eq. 0.d0) goto 80
   404) 
   405) !  store multipliers in a(i,k), i = k+1,...,n. -
   406)     t = 1.d0/t
   407)     do i = kp1,n
   408)       a(i,k) = -a(i,k)*t
   409)     enddo
   410) 
   411) !  apply multipliers to other columns of a.
   412)     do j = kp1,n
   413)       t = a(m,j)
   414)       a(m,j) = a(k,j)
   415)       a(k,j) = t
   416)       if (t .ne. 0.d0) then
   417)         do i = kp1,n
   418)           a(i,j) = a(i,j) + a(i,k)*t
   419)         enddo
   420)       endif
   421)     enddo
   422)   enddo
   423) 
   424) 70   k = n
   425)      if (a(n,n) .eq. 0.d0) goto 80
   426)      return
   427) 
   428) 80   ier = k
   429)      return
   430) 
   431) end subroutine dec
   432) 
   433) ! ************************************************************************** !
   434) 
   435) subroutine sol (n, ndim, a, b, ip)
   436) 
   437)   implicit none
   438) 
   439)   PetscInt :: m, n, ndim, ip(n)
   440)   PetscInt :: nm1,i,k,kb,km1,kp1
   441)   PetscReal :: a(ndim,n), b(n), t
   442) 
   443) !--------------------------------------------------
   444) !  solution of linear system a*x = b using output of dec.
   445) !  input..
   446) !     n = order of matrix.
   447) !     ndim = declared first dimension of array  a.
   448) !     a = triangularized matrix obtained from dec.
   449) !     b = right hand side vector.
   450) !     ip = pivot information vector obtained from dec.
   451) !  do not use if dec has set ier .ne. 0.
   452) !  output..
   453) !     b = solution vector, x .
   454) !--------------------------------------------------------------
   455) !
   456) !  reference:  a. c. hindmarsh, l. j. sloan, k. w. fong, and
   457) !              g. h. rodrigue,
   458) !              dec/sol:  solution of dense systems of linear
   459) !              algebraic equations,
   460) !              lawrence livermore laboratory report ucid-30137,
   461) !              june 1976.
   462) !----------------------------------------------------------------
   463) 
   464)   if (n .eq. 1) goto 50
   465)   nm1 = n - 1
   466) 
   467) !  apply row permutations and multipliers to b. --------------
   468)   do k = 1,nm1
   469)     kp1 = k + 1
   470)     m = ip(k)
   471)     t = b(m)
   472)     b(m) = b(k)
   473)     b(k) = t
   474)     do i = kp1,n
   475)       b(i) = b(i) + a(i,k)*t
   476)     enddo
   477)    enddo
   478) 
   479) !  back solve. -------------------------------------
   480)    do kb = 1,nm1
   481)     km1 = n - kb
   482)     k = km1 + 1
   483)     b(k) = b(k)/a(k,k)
   484)     t = -b(k)
   485)     do i = 1,km1
   486)       b(i) = b(i) + a(i,k)*t
   487)     enddo
   488)    enddo
   489) 
   490) 50   b(1) = b(1)/a(1,1)
   491) 
   492)   return
   493) end subroutine sol
   494) 
   495) end module Block_Tridiag_module

generated by
Intel(R) C++/Fortran Compiler code-coverage tool
Web-Page Owner: Nobody