LCOV - code coverage report
Current view: top level - src/fm - cp_cfm_dlaf_api.F (source / functions) Coverage Total Hit
Test: CP2K Regtests (git:936074a) Lines: 0.0 % 46 0
Test Date: 2025-12-04 06:27:48 Functions: 0.0 % 6 0

            Line data    Source code
       1              : !--------------------------------------------------------------------------------------------------!
       2              : !   CP2K: A general program to perform molecular dynamics simulations                              !
       3              : !   Copyright 2000-2025 CP2K developers group <https://cp2k.org>                                   !
       4              : !                                                                                                  !
       5              : !   SPDX-License-Identifier: GPL-2.0-or-later                                                      !
       6              : !--------------------------------------------------------------------------------------------------!
       7              : 
       8              : MODULE cp_cfm_dlaf_api
       9              : 
      10              :    USE cp_cfm_basic_linalg, ONLY: cp_cfm_uplo_to_full
      11              :    USE cp_cfm_types, ONLY: cp_cfm_type
      12              : #if defined(__DLAF)
      13              :    USE cp_dlaf_utils_api, ONLY: cp_dlaf_create_grid
      14              :    USE dlaf_fortran, ONLY: dlaf_pzheevd, &
      15              :                            dlaf_pzhegvd, &
      16              :                            dlaf_pzpotrf, &
      17              :                            dlaf_pzpotri
      18              : #endif
      19              :    USE kinds, ONLY: dp
      20              : #include "../base/base_uses.f90"
      21              : 
      22              :    IMPLICIT NONE
      23              : 
      24              :    PRIVATE
      25              : 
      26              :    CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'cp_cfm_dlaf_api'
      27              : 
      28              :    PUBLIC :: cp_cfm_pzpotrf_dlaf, cp_cfm_pzpotri_dlaf
      29              :    PUBLIC :: cp_cfm_diag_dlaf, cp_cfm_diag_gen_dlaf
      30              : 
      31              : CONTAINS
      32              : 
      33              : !***************************************************************************************************
      34              : !> \brief Cholesky factorization using DLA-Future
      35              : !> \param uplo ...
      36              : !> \param n Matrix size
      37              : !> \param a Local matrix
      38              : !> \param ia Row index of first row (has to be 1)
      39              : !> \param ja Col index of first column (has to be 1)
      40              : !> \param desca ScaLAPACK matrix descriptor
      41              : !> \param info 0 if factorization completed normally
      42              : !> \author Rocco Meli
      43              : ! **************************************************************************************************
      44            0 :    SUBROUTINE cp_cfm_pzpotrf_dlaf(uplo, n, a, ia, ja, desca, info)
      45              :       CHARACTER, INTENT(IN)                              :: uplo
      46              :       INTEGER, INTENT(IN)                                :: n
      47              :       COMPLEX(KIND=dp), DIMENSION(:, :), TARGET          :: a
      48              :       INTEGER, INTENT(IN)                                :: ia, ja
      49              :       INTEGER, DIMENSION(9)                              :: desca
      50              :       INTEGER, TARGET                                    :: info
      51              : 
      52              :       CHARACTER(len=*), PARAMETER :: routineN = 'cp_cfm_pzpotrf_dlaf'
      53              : 
      54              :       INTEGER                                            :: handle
      55              : 
      56            0 :       CALL timeset(routineN, handle)
      57              : #if defined(__DLAF)
      58              :       CALL dlaf_pzpotrf(uplo, n, a, ia, ja, desca, info)
      59              : #else
      60              :       MARK_USED(uplo)
      61              :       MARK_USED(n)
      62              :       MARK_USED(a)
      63              :       MARK_USED(ia)
      64              :       MARK_USED(ja)
      65              :       MARK_USED(desca)
      66              :       MARK_USED(info)
      67            0 :       CPABORT("CP2K compiled without the DLA-Future library.")
      68              : #endif
      69            0 :       CALL timestop(handle)
      70            0 :    END SUBROUTINE cp_cfm_pzpotrf_dlaf
      71              : 
      72              : !***************************************************************************************************
      73              : !> \brief Inverse from Cholesky factorization using DLA-Future
      74              : !> \param uplo ...
      75              : !> \param n Matrix size
      76              : !> \param a Local matrix
      77              : !> \param ia Row index of first row (has to be 1)
      78              : !> \param ja Col index of first column (has to be 1)
      79              : !> \param desca ScaLAPACK matrix descriptor
      80              : !> \param info 0 if factorization completed normally
      81              : !> \author Rocco Meli
      82              : ! **************************************************************************************************
      83            0 :    SUBROUTINE cp_cfm_pzpotri_dlaf(uplo, n, a, ia, ja, desca, info)
      84              :       CHARACTER, INTENT(IN)                              :: uplo
      85              :       INTEGER, INTENT(IN)                                :: n
      86              :       COMPLEX(KIND=dp), DIMENSION(:, :), TARGET          :: a
      87              :       INTEGER, INTENT(IN)                                :: ia, ja
      88              :       INTEGER, DIMENSION(9)                              :: desca
      89              :       INTEGER, TARGET                                    :: info
      90              : 
      91              :       CHARACTER(len=*), PARAMETER :: routineN = 'cp_cfm_pzpotri_dlaf'
      92              : 
      93              :       INTEGER                                            :: handle
      94              : 
      95            0 :       CALL timeset(routineN, handle)
      96              : #if defined(__DLAF)
      97              :       CALL dlaf_pzpotri(uplo, n, a, ia, ja, desca, info)
      98              : #else
      99              :       MARK_USED(uplo)
     100              :       MARK_USED(n)
     101              :       MARK_USED(a)
     102              :       MARK_USED(ia)
     103              :       MARK_USED(ja)
     104              :       MARK_USED(desca)
     105              :       MARK_USED(info)
     106            0 :       CPABORT("CP2K compiled without the DLA-Future library.")
     107              : #endif
     108            0 :       CALL timestop(handle)
     109            0 :    END SUBROUTINE cp_cfm_pzpotri_dlaf
     110              : 
     111              : ! **************************************************************************************************
     112              : !> \brief DLA-Future eigensolver for complex Hermitian matrices
     113              : !> \param matrix ...
     114              : !> \param eigenvectors ...
     115              : !> \param eigenvalues ...
     116              : !> \author Rocco Meli
     117              : ! **************************************************************************************************
     118            0 :    SUBROUTINE cp_cfm_diag_dlaf(matrix, eigenvectors, eigenvalues)
     119              : 
     120              :       TYPE(cp_cfm_type), INTENT(IN)                      :: matrix, eigenvectors
     121              :       REAL(KIND=dp), DIMENSION(:), INTENT(OUT)           :: eigenvalues
     122              : 
     123              :       CHARACTER(LEN=*), PARAMETER                        :: routineN = 'cp_cfm_diag_dlaf'
     124              : 
     125              :       INTEGER                                            :: handle, n, nmo
     126              :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:), TARGET   :: eig
     127              : 
     128            0 :       CALL timeset(routineN, handle)
     129              : 
     130            0 :       n = matrix%matrix_struct%nrow_global
     131            0 :       ALLOCATE (eig(n))
     132              : 
     133            0 :       CALL cp_cfm_diag_dlaf_base(matrix, eigenvectors, eig)
     134              : 
     135            0 :       nmo = SIZE(eigenvalues, 1)
     136            0 :       IF (nmo > n) THEN
     137            0 :          eigenvalues(1:n) = eig(1:n)
     138              :       ELSE
     139            0 :          eigenvalues(1:nmo) = eig(1:nmo)
     140              :       END IF
     141              : 
     142            0 :       DEALLOCATE (eig)
     143              : 
     144            0 :       CALL timestop(handle)
     145              : 
     146            0 :    END SUBROUTINE cp_cfm_diag_dlaf
     147              : 
     148              : ! **************************************************************************************************
     149              : !> \brief DLA-Future generalized eigensolver for complex Hermitian matrices
     150              : !> \param amatrix ...
     151              : !> \param bmatrix ...
     152              : !> \param eigenvectors ...
     153              : !> \param eigenvalues ...
     154              : !> \author Rocco Meli
     155              : ! **************************************************************************************************
     156            0 :    SUBROUTINE cp_cfm_diag_gen_dlaf(amatrix, bmatrix, eigenvectors, eigenvalues)
     157              : 
     158              :       TYPE(cp_cfm_type), INTENT(IN)                      :: amatrix, bmatrix, eigenvectors
     159              :       REAL(KIND=dp), DIMENSION(:), INTENT(OUT)           :: eigenvalues
     160              : 
     161              :       CHARACTER(LEN=*), PARAMETER :: routineN = 'cp_cfm_diag_gen_dlaf'
     162              : 
     163              :       INTEGER                                            :: handle, n, nmo
     164              :       REAL(KIND=dp), ALLOCATABLE, DIMENSION(:), TARGET   :: eig
     165              : 
     166            0 :       CALL timeset(routineN, handle)
     167              : 
     168            0 :       n = amatrix%matrix_struct%nrow_global
     169            0 :       ALLOCATE (eig(n))
     170              : 
     171            0 :       CALL cp_cfm_diag_gen_dlaf_base(amatrix, bmatrix, eigenvectors, eig)
     172              : 
     173            0 :       nmo = SIZE(eigenvalues, 1)
     174            0 :       IF (nmo > n) THEN
     175            0 :          eigenvalues(1:n) = eig(1:n)
     176              :       ELSE
     177            0 :          eigenvalues(1:nmo) = eig(1:nmo)
     178              :       END IF
     179              : 
     180            0 :       DEALLOCATE (eig)
     181              : 
     182            0 :       CALL timestop(handle)
     183              : 
     184            0 :    END SUBROUTINE cp_cfm_diag_gen_dlaf
     185              : 
     186              :    !***************************************************************************************************
     187              : !> \brief DLA-Future standard eigensolver for complex Hermitian matrices
     188              : !> \param matrix ...
     189              : !> \param eigenvectors ...
     190              : !> \param eigenvalues ...
     191              : !> \author Rocco Meli
     192              : ! **************************************************************************************************
     193            0 :    SUBROUTINE cp_cfm_diag_dlaf_base(matrix, eigenvectors, eigenvalues)
     194              :       TYPE(cp_cfm_type), INTENT(IN)                      :: matrix, eigenvectors
     195              :       REAL(kind=dp), DIMENSION(:), INTENT(OUT), TARGET   :: eigenvalues
     196              : 
     197              :       CHARACTER(len=*), PARAMETER :: dlaf_name = 'pzheevd_dlaf', routineN = 'cp_cfm_diag_dlaf_base'
     198              :       CHARACTER, PARAMETER                               :: uplo = 'L'
     199              : 
     200              :       CHARACTER(LEN=100)                                 :: message
     201            0 :       COMPLEX(KIND=dp), DIMENSION(:, :), POINTER         :: a, z
     202              :       INTEGER                                            :: blacs_context, dlaf_handle, handle, n
     203              :       INTEGER, DIMENSION(9)                              :: desca, descz
     204              :       INTEGER, TARGET                                    :: info
     205              : 
     206            0 :       CALL timeset(routineN, handle)
     207              : 
     208              : #if defined(__DLAF)
     209              :       ! DLAF needs the lower triangular part
     210              :       ! Use eigenvectors matrix as workspace
     211              :       CALL cp_cfm_uplo_to_full(matrix, eigenvectors)
     212              : 
     213              :       ! Create DLAF grid from BLACS context (if already present, does nothing)
     214              :       blacs_context = matrix%matrix_struct%context%get_handle()
     215              :       CALL cp_dlaf_create_grid(blacs_context)
     216              : 
     217              :       n = matrix%matrix_struct%nrow_global
     218              : 
     219              :       a => matrix%local_data
     220              :       z => eigenvectors%local_data
     221              : 
     222              :       desca(:) = matrix%matrix_struct%descriptor(:)
     223              :       descz(:) = eigenvectors%matrix_struct%descriptor(:)
     224              : 
     225              :       info = -1
     226              :       CALL timeset(dlaf_name, dlaf_handle)
     227              :       CALL dlaf_pzheevd(uplo, n, a, 1, 1, desca, eigenvalues, z, 1, 1, descz, info)
     228              :       CALL timestop(dlaf_handle)
     229              : 
     230              :       IF (info /= 0) THEN
     231              :          WRITE (message, "(A,I0,A)") "ERROR in DLAF_PZHEEVD: Eigensolver failed (INFO = ", info, ")"
     232              :          CPABORT(TRIM(message))
     233              :       END IF
     234              : #else
     235              :       MARK_USED(a)
     236              :       MARK_USED(z)
     237              :       MARK_USED(desca)
     238              :       MARK_USED(descz)
     239              :       MARK_USED(matrix)
     240              :       MARK_USED(eigenvectors)
     241              :       MARK_USED(eigenvalues)
     242              :       MARK_USED(uplo)
     243              :       MARK_USED(n)
     244              :       MARK_USED(info)
     245              :       MARK_USED(dlaf_handle)
     246              :       MARK_USED(dlaf_name)
     247              :       MARK_USED(message)
     248              :       MARK_USED(blacs_context)
     249            0 :       CPABORT("CP2K compiled without DLAF library.")
     250              : #endif
     251              : 
     252            0 :       CALL timestop(handle)
     253              : 
     254            0 :    END SUBROUTINE cp_cfm_diag_dlaf_base
     255              : 
     256              : !***************************************************************************************************
     257              : !> \brief DLA-Future generalized eigensolver for complex Hermitian matrices
     258              : !> \param amatrix ...
     259              : !> \param bmatrix ...
     260              : !> \param eigenvectors ...
     261              : !> \param eigenvalues ...
     262              : !> \author Rocco Meli
     263              : ! **************************************************************************************************
     264            0 :    SUBROUTINE cp_cfm_diag_gen_dlaf_base(amatrix, bmatrix, eigenvectors, eigenvalues)
     265              :       TYPE(cp_cfm_type), INTENT(IN)                      :: amatrix, bmatrix, eigenvectors
     266              :       REAL(kind=dp), DIMENSION(:), INTENT(OUT), TARGET   :: eigenvalues
     267              : 
     268              :       CHARACTER(len=*), PARAMETER :: dlaf_name = 'pzhegvd_dlaf', &
     269              :          routineN = 'cp_cfm_diag_gen_dlaf_base'
     270              :       CHARACTER, PARAMETER                               :: uplo = 'L'
     271              : 
     272              :       CHARACTER(LEN=100)                                 :: message
     273            0 :       COMPLEX(KIND=dp), DIMENSION(:, :), POINTER         :: a, b, z
     274              :       INTEGER                                            :: blacs_context, dlaf_handle, handle, n
     275              :       INTEGER, DIMENSION(9)                              :: desca, descb, descz
     276              :       INTEGER, TARGET                                    :: info
     277              : 
     278            0 :       CALL timeset(routineN, handle)
     279              : 
     280              : #if defined(__DLAF)
     281              :       ! DLAF needs the lower triangular part
     282              :       ! Use eigenvectors matrix as workspace
     283              :       CALL cp_cfm_uplo_to_full(amatrix, eigenvectors)
     284              :       CALL cp_cfm_uplo_to_full(bmatrix, eigenvectors)
     285              : 
     286              :       ! Create DLAF grid from BLACS context (if already present, does nothing)
     287              :       blacs_context = amatrix%matrix_struct%context%get_handle()
     288              :       CALL cp_dlaf_create_grid(blacs_context)
     289              : 
     290              :       n = amatrix%matrix_struct%nrow_global
     291              : 
     292              :       a => amatrix%local_data
     293              :       b => bmatrix%local_data
     294              :       z => eigenvectors%local_data
     295              : 
     296              :       desca(:) = amatrix%matrix_struct%descriptor(:)
     297              :       descb(:) = bmatrix%matrix_struct%descriptor(:)
     298              :       descz(:) = eigenvectors%matrix_struct%descriptor(:)
     299              : 
     300              :       info = -1
     301              :       CALL timeset(dlaf_name, dlaf_handle)
     302              :       CALL dlaf_pzhegvd(uplo, n, a, 1, 1, desca, b, 1, 1, descb, eigenvalues, z, 1, 1, descz, info)
     303              :       CALL timestop(dlaf_handle)
     304              : 
     305              :       IF (info /= 0) THEN
     306              :          WRITE (message, "(A,I0,A)") "ERROR in DLAF_PZHEGVD: Eigensolver failed (INFO = ", info, ")"
     307              :          CPABORT(TRIM(message))
     308              :       END IF
     309              : #else
     310              :       MARK_USED(a)
     311              :       MARK_USED(b)
     312              :       MARK_USED(z)
     313              :       MARK_USED(desca)
     314              :       MARK_USED(descb)
     315              :       MARK_USED(descz)
     316              :       MARK_USED(amatrix)
     317              :       MARK_USED(bmatrix)
     318              :       MARK_USED(eigenvectors)
     319              :       MARK_USED(eigenvalues)
     320              :       MARK_USED(uplo)
     321              :       MARK_USED(n)
     322              :       MARK_USED(info)
     323              :       MARK_USED(dlaf_handle)
     324              :       MARK_USED(dlaf_name)
     325              :       MARK_USED(message)
     326              :       MARK_USED(blacs_context)
     327            0 :       CPABORT("CP2K compiled without DLAF library.")
     328              : #endif
     329              : 
     330            0 :       CALL timestop(handle)
     331              : 
     332            0 :    END SUBROUTINE cp_cfm_diag_gen_dlaf_base
     333              : 
     334              : END MODULE cp_cfm_dlaf_api
        

Generated by: LCOV version 2.0-1