block_solve.F90       coverage:  75.00 %func     51.92 %block


     1) module Block_Solve_module
     2) 
     3) #include "petsc/finclude/petscsys.h"
     4)   
     5)   private
     6)   
     7)   public :: bl3dfac, &
     8)             bl3dsolf, &
     9)             bl3dsolb
    10)       
    11) contains
    12) 
    13) ! ************************************************************************** !
    14) 
    15) subroutine bl3dfac(n, k, E, D, F, pivot)
    16) 
    17) !	This version of bl3dfac dated March 15, 2000.
    18) 
    19) ! Modified from bl3dfac.f, bl3dsol.f - Fortran codes of
    20) ! P. Keast, http://www.mscs.dal.ca/ keast/research/pubs.html, 
    21) ! Mathematics and Statistics Department, Dalhousie University, Canada.
    22) 
    23) !************************************************************************
    24) 
    25) !	bl3dfac performs an LU factorization of the block tridiagonal
    26) !	matrix given by:
    27) 
    28) !	D(1) E(1)  O   O    O ............  O
    29) !	F(1) D(2) E(2) O    O ............  O
    30) !	 O   F(2) D(3) E(3) O ............  O
    31) !	 ....................................
    32) !
    33) !        O    O   ........F(n-1) D(n-1) E(n-1)
    34) !        O    O   .............. F(n-1)   D(n)
    35) 
    36) !	where each block, E(j), D(j), F(j) is k by k.
    37) !
    38) !	Pivoting is done only within the diagonal blocks.
    39) 
    40) !	The L matrix is NOT unit lower triangular, but the U matrix
    41) !	is unit UPPER triangular, with the identity matrix on the
    42) !	diagonal. 
    43) 
    44) !************************************************************************
    45) 
    46) !	Uses: Lapack routines dgetrf and dgetrs for Gauss LU factorization
    47) !	and solution, and the BLAS routine dgemm for matrix-matrix 
    48) !       multiplication.
    49) 
    50) !************************************************************************
    51) 
    52) !	Input variables:
    53) 
    54)   PetscInt :: n, k
    55)   PetscInt :: pivot(k,n)
    56)   PetscReal :: E(k,k,n-1), D(k,k,n), F(k,k,n-1)
    57) 
    58) ! Input:
    59) 
    60) ! k           Integer, the size of the blocks.
    61) ! n           Integer, the number of blocks.
    62) ! E(k,k,n-1)  Double precision, the sub-diagonal blocks.
    63) ! D(k,k,n)    Double precision, the diagonal blocks.
    64) ! F(k,k,n-1)  Double precision, the super-diagonal blocks.
    65) 
    66) ! Output:
    67) 
    68) ! E(k,k,n-1)  The lower triangular blocks of L in the LU factorization.
    69) ! F(k,k,n-1)  The upper triangular blocks of U in the LU factorization.
    70) ! D(k,k,n)    The diagonal blocks of L in the LU factorization.
    71) !                   The diagonal blocks of U are unit upper triangular.
    72) ! pivot(k,n)  Integer pivot vector used in factoring the diagonal blocks.
    73) ! 	    pivot(1:k,p) records the pivoting done in diagonal block p.
    74) 
    75) ! Local variables:
    76) 
    77)   PetscInt :: j, i, l, ll, info, info1, iwork(k)
    78)   character(len=1) :: trans, norm = '1'
    79)   PetscReal, parameter :: one = 1.d0
    80)   PetscReal anorm, sum, rcond, work(4*k)
    81) 
    82) !************************************************************************
    83) 
    84) ! Executable statements.
    85) 
    86) !************************************************************************
    87) 
    88) ! First, in main blocks 1 to n-1:
    89) 
    90)   trans = 'N'
    91)         
    92)   do j = 1, n-1
    93) 
    94) !  compute 1-norm: only needed for computing condition nr.
    95) #ifdef CONDNR
    96)     anorm = 0.d0
    97)     do l = 1, k
    98)       sum = 0.d0
    99)       do i = 1, k
   100)         sum = sum + abs(d(i,l,j))
   101)       enddo
   102)       anorm = max(anorm,sum)
   103)     enddo
   104) #endif
   105) 
   106) !  First, factor D(j).
   107)     call dgetrf(k, k, D(1,1,j), k, pivot(1,j), info)
   108) 
   109)     if ( info .ne. 0 ) then
   110) !        write(*,'("node ",3i3,1p15e12.4)') j,n,k,(d(i,ll,j),ll=1,k)
   111) 
   112)       print *,'At block ',j,',',' 1-norm = ',anorm
   113)       print *,'problem, info = ', info   ! make this better later.
   114)       return
   115)     endif
   116) 
   117) !  Estimate condition number
   118) #ifdef CONDNR
   119)     CALL dgecon(norm,k,D(1,1,j),k,anorm,rcond,work,iwork,info1)
   120)     rcond = 1.e0/rcond
   121) !   if (rcond > 1.e10)
   122)         WRITE (*,999) 'Estimate of condition number =', &
   123)         rcond,info1
   124)     999 FORMAT (1X,A,1P,e11.4,i3)
   125) #endif
   126) 
   127) !  Now, compute E(j) from D(j) * E(j) = E(j).
   128)     call dgetrs(trans, k, k, D(1,1,j), k, pivot(1,j), &
   129)                        E(1,1,j), k, info)
   130)     if ( info .lt. 0 ) then
   131)       print *,'Illegal parameter number ',-info
   132)       return
   133)     endif
   134) 
   135) !  Finally, compute D(j+1) = D(j+1) - F(j) * E(j).
   136)     call dgemm(trans, trans, k, k, k, -one, F(1,1,j), k, &
   137)                       E(1,1,j), k, one, D(1,1,j+1), k)
   138)   enddo
   139) 
   140) ! Finally, obtain the LU factorization of D(n).
   141)   call dgetrf(k, k, D(1,1,n), k, pivot(1,n), info)
   142)   if ( info .ne. 0 ) then
   143)      print *,'At block ',j,','
   144)      print *,'problem, info = ', info   ! make this better later.
   145)     return
   146)   endif
   147) 
   148)   return
   149)   
   150) end subroutine bl3dfac
   151) 
   152) ! ************************************************************************** !
   153) 
   154) subroutine bl3dsol(n, k, E, D, F, pivot, nrhs, rhs)
   155) 
   156) !       This version of bl3dsol dated April 20, 2000.
   157) 
   158) !************************************************************************
   159) 
   160) !       bl3dsol solves the system of linear equations whose coefficient
   161) !       matrix is the block tridiagonal matrix given by:
   162) 
   163) !       D(1) E(1)  O   O    O ............  O
   164) !       F(1) D(2) E(2) O    O ............  O
   165) !        O   F(2) D(3) E(3) O ............  O
   166) !        ....................................
   167) !
   168) !        O    O   ........F(n-1) D(n-1) E(n-1)
   169) !        O    O   .............. F(n-1)   D(n)
   170) 
   171) !       where each block, E(j), D(j), F(j) is k by k, and which has already
   172) !       been factored using bl3dfac.
   173) 
   174) !       The right hand sides are stored by blocks (see below).
   175) 
   176) !       The L matrix is NOT unit lower triangular, but the U matrix
   177) !       is unit UPPER triangular.
   178) 
   179) !************************************************************************
   180) 
   181) !       Uses: 
   182) 
   183) !************************************************************************
   184) 
   185) !       Input variables:
   186) 
   187)         character(len=1) :: trans
   188)         PetscInt :: n, k
   189)         PetscReal :: E(k,k,n-1), D(k,k,n), F(k,k,n-1)
   190)         PetscInt :: pivot(k,n), nrhs
   191)         PetscReal :: rhs(k,nrhs,n)
   192) 
   193) !       Input:
   194) 
   195) !       trans         character*1, 'n' or 'N' is A.x = b is to be solved,
   196) !                     't' or 'T' is A^T.x = b is to be solved.
   197) !       n             Integer, the number of blocks.
   198) !       k             Integer, the size of the blocks.
   199) !       E(k,k,n-1)    Double precision, the sub-diagonal blocks.
   200) !       D(k,k,n)      Double precision, the diagonal blocks, after factoring 
   201) !                     by bl3dfac.
   202) !       F(k,k,n-1)    Double precision, the super-diagonal blocks, after 
   203) !                     modification in bl3dfac.
   204) !       pivot(k,n)    Integer pivot vector returned by bl3dfac, used in 
   205) !                     factoring the diagonal blocks.
   206) !                     pivot(1:k,p) records the pivoting done in diagonal 
   207) !                     block p.
   208) !       nrhs          integer, the number of right hand sides.
   209) !       rhs(k,nrhs,n) Double precision, right hand sides. These are stored
   210) !                     by blocks - that is, the nrhs right hand sides 
   211) !                     associated with each k by k block are stored 
   212) !                     consecutively. For example, for n = 3, k = 2 and
   213) !                     nrhs = 3, the elements of the right hand sides
   214) !                     are stored in the following order:
   215) 
   216) !                     [1 ] [3 ] [5 ]
   217) !                     [2 ] [4 ] [6 ]  First block.
   218) 
   219) !                     [7 ] [9 ] [11]
   220) !                     [8 ] [10] [12]  Second block.
   221) 
   222) !                     [13] [15] [17]  Third block.
   223) !                     [14] [16] [18]
   224) 
   225) !       Output:
   226) 
   227) !       rhs(k,nrhs,n) double precision, the solutions.
   228) !       
   229) !       Local variables:
   230) 
   231)         PetscInt :: j, info
   232)         PetscReal, parameter :: one = 1.d0
   233) 
   234) !************************************************************************
   235) 
   236) !       Executable statements.
   237) 
   238) !************************************************************************
   239)         
   240) !       Modification of the right hand sides, that is, solve the
   241) !       block bi-diagonal lower triangular system L.x = b, overwriting
   242) !       b with the solution.
   243) 
   244)         trans = 'N'
   245) 
   246) 
   247)         if ( trans .eq. 'n' .or. trans .eq. 'N' ) then 
   248) 
   249) !          First, solve D(1).x(1..k) = b(1..k).
   250)            call dgetrs(trans, k, nrhs, D(1,1,1), k, pivot(1,1), &
   251)                        rhs(1,1,1), k, info)
   252) 
   253)            if ( info .lt. 0 ) then
   254)               print *,'Illegal parameter number ',-info
   255)               return
   256)            endif
   257) 
   258) !          Then, for j = 2 to n solve 
   259) !          D(j)x((j-1)*k+1..j*k) = b((j-1)*k+1..j*k) 
   260) !                                - F(j-1)x((j-2)*k+1..(j-1)*k)
   261) 
   262)            do j = 2, n
   263)            
   264)               call dgemm(trans, trans, k, nrhs, k, -one, F(1,1,j-1), k, &
   265)                          rhs(1,1,j-1), k, one, rhs(1,1,j), k)
   266) 
   267)               call dgetrs(trans, k, nrhs, D(1,1,j), k, pivot(1,j), &
   268)                           rhs(1,1,j), k, info)
   269)               if ( info .lt. 0 ) then
   270)                  print *,'Illegal parameter number ',-info
   271)                  return
   272)               endif
   273) 
   274)            enddo
   275) 
   276) !        Now, the back substitution.
   277) 
   278)            do j = n-1,1,-1
   279)    
   280) !             Form rhs(:,:,j) = rhs(:,:,j) - E(j)*rhs(:,:,j+1)
   281) 
   282)               call dgemm(trans, trans, k, nrhs, k, -one, E(1,1,j), k, &
   283)                          rhs(1,1,j+1), k, one, rhs(1,1,j), k)
   284) 
   285)            enddo
   286) 
   287)      else
   288) 
   289) !    Solve the transposed system.
   290) 
   291) !    First, the right hand side modification.
   292) 
   293)          do j = 2,n
   294) 
   295) !        Form rhs(:,:,j) = rhs(:,:,j) - E(j-1)^T*rhs(:,:,j-1)
   296)             call dgemm(trans, 'N', k, nrhs, k, -one, E(1,1,j-1), k, &
   297)                            rhs(1,1,j-1), k, one, rhs(1,1,j), k)
   298)          enddo
   299) 
   300) !      Now, the back substitution.
   301) 
   302) !      First, solve D(n)^T.rhs(:,:,n) = rhs(:,:,n)
   303) 
   304)          call dgetrs(trans, k, nrhs, D(1,1,n), k, pivot(1,n), &
   305)                      rhs(1,1,n), k, info)
   306)          if ( info .lt. 0 ) then
   307)             print *,'Illegal parameter number ',-info
   308)             return
   309)          endif
   310) 
   311) !          For j = n-1 down to 1, solve 
   312) !          D(j)^T.rhs(:,:,j) = rhs(:,:,j) - F(j)^T.rhs(:,:,j+1)
   313) 
   314)          do j = n-1,1,-1
   315) 
   316) !        Form rhs(:,:,j) = rhs(:,:,j) - F(j)^T*rhs(:,:,j+1)
   317)             call dgemm(trans, 'N', k, nrhs, k, -one, F(1,1,j), k, &
   318)                            rhs(1,1,j+1), k, one, rhs(1,1,j), k)
   319) 
   320)             call dgetrs(trans, k, nrhs, D(1,1,j), k, pivot(1,j), &
   321)                         rhs(1,1,j), k, info)
   322)             if ( info .lt. 0 ) then
   323)                print *,'Illegal parameter number ',-info
   324)                return
   325)             endif
   326)          enddo
   327)       endif
   328) 
   329)       return
   330)       
   331) end subroutine bl3dsol
   332) 
   333) ! ************************************************************************** !
   334) 
   335) subroutine bl3dsolf(n, k, E, D, F, pivot, nrhs, rhs)
   336) 
   337) !       This version of bl3dsol dated April 20, 2000.
   338) 
   339) !************************************************************************
   340) 
   341) !       bl3dsol solves the system of linear equations whose coefficient
   342) !       matrix is the block tridiagonal matrix given by:
   343) 
   344) !       D(1) E(1)  O   O    O ............  O
   345) !       F(1) D(2) E(2) O    O ............  O
   346) !        O   F(2) D(3) E(3) O ............  O
   347) !        ....................................
   348) !
   349) !        O    O   ........F(n-1) D(n-1) E(n-1)
   350) !        O    O   .............. F(n-1)   D(n)
   351) 
   352) !       where each block, E(j), D(j), F(j) is k by k, and which has already
   353) !       been factored using bl3dfac.
   354) 
   355) !       The right hand sides are stored by blocks (see below).
   356) 
   357) !       The L matrix is NOT unit lower triangular, but the U matrix
   358) !       is unit UPPER triangular.
   359) 
   360) !************************************************************************
   361) 
   362) !       Uses: 
   363) 
   364) !************************************************************************
   365) 
   366) !       Input variables:
   367) 
   368)         character(len=1) :: trans
   369)         PetscInt :: n, k
   370)         PetscReal :: E(k,k,n-1), D(k,k,n), F(k,k,n-1)
   371)         PetscInt :: pivot(k,n), nrhs
   372)         PetscReal :: rhs(k,nrhs,n)
   373) 
   374) !       Input:
   375) 
   376) !       trans         character*1, 'n' or 'N' is A.x = b is to be solved,
   377) !                     't' or 'T' is A^T.x = b is to be solved.
   378) !       n             Integer, the number of blocks.
   379) !       k             Integer, the size of the blocks.
   380) !       E(k,k,n-1)    Double precision, the sub-diagonal blocks.
   381) !       D(k,k,n)      Double precision, the diagonal blocks, after factoring 
   382) !                     by bl3dfac.
   383) !       F(k,k,n-1)    Double precision, the super-diagonal blocks, after 
   384) !                     modification in bl3dfac.
   385) !       pivot(k,n)    Integer pivot vector returned by bl3dfac, used in 
   386) !                     factoring the diagonal blocks.
   387) !                     pivot(1:k,p) records the pivoting done in diagonal 
   388) !                     block p.
   389) !       nrhs          integer, the number of right hand sides.
   390) !       rhs(k,nrhs,n) Double precision, right hand sides. These are stored
   391) !                     by blocks - that is, the nrhs right hand sides 
   392) !                     associated with each k by k block are stored 
   393) !                     consecutively. For example, for n = 3, k = 2 and
   394) !                     nrhs = 3, the elements of the right hand sides
   395) !                     are stored in the following order:
   396) 
   397) !                     [1 ] [3 ] [5 ]
   398) !                     [2 ] [4 ] [6 ]  First block.
   399) 
   400) !                     [7 ] [9 ] [11]
   401) !                     [8 ] [10] [12]  Second block.
   402) 
   403) !                     [13] [15] [17]  Third block.
   404) !                     [14] [16] [18]
   405) 
   406) !       Output:
   407) 
   408) !       rhs(k,nrhs,n) double precision, the solutions.
   409) !       
   410) !       Local variables:
   411) 
   412)         PetscInt :: j, info
   413)         PetscReal, parameter :: one = 1.d0
   414) 
   415) !************************************************************************
   416) 
   417) !       Executable statements.
   418) 
   419) !************************************************************************
   420)         
   421) !       Modification of the right hand sides, that is, solve the
   422) !       block bi-diagonal lower triangular system L.x = b, overwriting
   423) !       b with the solution.
   424) 
   425)           trans = 'N'
   426) 
   427) 
   428) !       if ( trans .eq. 'n' .or. trans .eq. 'N' ) then 
   429) 
   430) !          First, solve D(1).x(1..k) = b(1..k).
   431)            call dgetrs(trans, k, nrhs, D(1,1,1), k, pivot(1,1), &
   432)                        rhs(1,1,1), k, info)
   433) 
   434)            if ( info .lt. 0 ) then
   435)               print *,'Illegal parameter number ',-info
   436)               return
   437)            endif
   438) 
   439) !          Then, for j = 2 to n solve 
   440) !          D(j)x((j-1)*k+1..j*k) = b((j-1)*k+1..j*k) 
   441) !                                - F(j-1)x((j-2)*k+1..(j-1)*k)
   442) 
   443)            do j = 2, n
   444)            
   445)               call dgemm(trans, trans, k, nrhs, k, -one, F(1,1,j-1), k, &
   446)                          rhs(1,1,j-1), k, one, rhs(1,1,j), k)
   447) 
   448)               call dgetrs(trans, k, nrhs, D(1,1,j), k, pivot(1,j), &
   449)                           rhs(1,1,j), k, info)
   450)               if ( info .lt. 0 ) then
   451)                  print *,'Illegal parameter number ',-info
   452)                  return
   453)               endif
   454) 
   455)            enddo
   456) 
   457)       return
   458) 
   459) end subroutine bl3dsolf
   460) 
   461) ! ************************************************************************** !
   462) 
   463) subroutine bl3dsolb(n, k, E, D, F, pivot, nrhs, rhs)
   464) 
   465) !       This version of bl3dsol dated April 20, 2000.
   466) 
   467) !************************************************************************
   468) 
   469) !       bl3dsol solves the system of linear equations whose coefficient
   470) !       matrix is the block tridiagonal matrix given by:
   471) 
   472) !       D(1) E(1)  O   O    O ............  O
   473) !       F(1) D(2) E(2) O    O ............  O
   474) !        O   F(2) D(3) E(3) O ............  O
   475) !        ....................................
   476) !
   477) !        O    O   ........F(n-1) D(n-1) E(n-1)
   478) !        O    O   .............. F(n-1)   D(n)
   479) 
   480) !       where each block, E(j), D(j), F(j) is k by k, and which has already
   481) !       been factored using bl3dfac.
   482) 
   483) !       The right hand sides are stored by blocks (see below).
   484) 
   485) !       The L matrix is NOT unit lower triangular, but the U matrix
   486) !       is unit UPPER triangular.
   487) 
   488) !************************************************************************
   489) 
   490) !       Uses: 
   491) 
   492) !************************************************************************
   493) 
   494) !       Input variables:
   495) 
   496)         character(len=1) :: trans
   497)         PetscInt :: n, k
   498)         PetscReal :: E(k,k,n-1), D(k,k,n), F(k,k,n-1)
   499)         PetscInt :: pivot(k,n), nrhs
   500)         PetscReal :: rhs(k,nrhs,n)
   501) 
   502) !       Input:
   503) 
   504) !       trans         character*1, 'n' or 'N' is A.x = b is to be solved,
   505) !                     't' or 'T' is A^T.x = b is to be solved.
   506) !       n             Integer, the number of blocks.
   507) !       k             Integer, the size of the blocks.
   508) !       E(k,k,n-1)    Double precision, the sub-diagonal blocks.
   509) !       D(k,k,n)      Double precision, the diagonal blocks, after factoring 
   510) !                     by bl3dfac.
   511) !       F(k,k,n-1)    Double precision, the super-diagonal blocks, after 
   512) !                     modification in bl3dfac.
   513) !       pivot(k,n)    Integer pivot vector returned by bl3dfac, used in 
   514) !                     factoring the diagonal blocks.
   515) !                     pivot(1:k,p) records the pivoting done in diagonal 
   516) !                     block p.
   517) !       nrhs          integer, the number of right hand sides.
   518) !       rhs(k,nrhs,n) Double precision, right hand sides. These are stored
   519) !                     by blocks - that is, the nrhs right hand sides 
   520) !                     associated with each k by k block are stored 
   521) !                     consecutively. For example, for n = 3, k = 2 and
   522) !                     nrhs = 3, the elements of the right hand sides
   523) !                     are stored in the following order:
   524) 
   525) !                     [1 ] [3 ] [5 ]
   526) !                     [2 ] [4 ] [6 ]  First block.
   527) 
   528) !                     [7 ] [9 ] [11]
   529) !                     [8 ] [10] [12]  Second block.
   530) 
   531) !                     [13] [15] [17]  Third block.
   532) !                     [14] [16] [18]
   533) 
   534) !       Output:
   535) 
   536) !       rhs(k,nrhs,n) double precision, the solutions.
   537) !       
   538) !       Local variables:
   539) 
   540)         PetscInt :: j !, info
   541)         PetscReal, parameter :: one = 1.d0
   542) 
   543)         trans = 'N'
   544) 
   545) !        Now, the back substitution.
   546) 
   547)            do j = n-1,1,-1
   548)    
   549) !             Form rhs(:,:,j) = rhs(:,:,j) - E(j)*rhs(:,:,j+1)
   550) 
   551)               call dgemm(trans, trans, k, nrhs, k, -one, E(1,1,j), k, &
   552)                          rhs(1,1,j+1), k, one, rhs(1,1,j), k)
   553) 
   554)            enddo
   555) 
   556) 
   557)       return
   558) 
   559) end subroutine bl3dsolb
   560) 
   561) !************************ END OF bl3dsol ********************************
   562) end module Block_Solve_module
   563)   

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