(git:e7e05ae)
qs_condnum.F
Go to the documentation of this file.
1 !--------------------------------------------------------------------------------------------------!
2 ! CP2K: A general program to perform molecular dynamics simulations !
3 ! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4 ! !
5 ! SPDX-License-Identifier: GPL-2.0-or-later !
6 !--------------------------------------------------------------------------------------------------!
7 
8 ! **************************************************************************************************
9 !> \brief Calculation of overlap matrix condition numbers
10 !> \par History
11 !> \author JGH (14.11.2016)
12 ! **************************************************************************************************
13 MODULE qs_condnum
16  USE cp_blacs_env, ONLY: cp_blacs_env_type
19  USE cp_fm_diag, ONLY: cp_fm_power
22  cp_fm_struct_type
23  USE cp_fm_types, ONLY: cp_fm_create,&
24  cp_fm_release,&
25  cp_fm_type
26  USE dbcsr_api, ONLY: &
27  dbcsr_add, dbcsr_copy, dbcsr_create, dbcsr_desymmetrize, dbcsr_gershgorin_norm, &
28  dbcsr_get_block_diag, dbcsr_get_info, dbcsr_get_matrix_type, dbcsr_iterator_blocks_left, &
29  dbcsr_iterator_next_block, dbcsr_iterator_start, dbcsr_iterator_stop, dbcsr_iterator_type, &
30  dbcsr_p_type, dbcsr_release, dbcsr_type, dbcsr_type_no_symmetry, dbcsr_type_symmetric
31  USE kinds, ONLY: dp
32  USE mathlib, ONLY: invmat
33 #include "./base/base_uses.f90"
34 
35  IMPLICIT NONE
36 
37  PRIVATE
38 
39 ! *** Global parameters ***
40 
41  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_condnum'
42 
43 ! *** Public subroutines ***
44 
45  PUBLIC :: overlap_condnum
46 
47 CONTAINS
48 
49 ! **************************************************************************************************
50 !> \brief Calculation of the overlap matrix Condition Number
51 !> \param matrixkp_s The overlap matrices to be calculated (kpoints, optional)
52 !> \param condnum Condition numbers for 1 and 2 norm
53 !> \param iunit output unit
54 !> \param norml1 logical: calculate estimate to 1-norm
55 !> \param norml2 logical: calculate estimate to 1-norm and 2-norm condition number
56 !> \param use_arnoldi logical: use Arnoldi iteration to estimate 2-norm condition number
57 !> \param blacs_env ...
58 !> \date 07.11.2016
59 !> \par History
60 !> \author JHU
61 !> \version 1.0
62 ! **************************************************************************************************
63  SUBROUTINE overlap_condnum(matrixkp_s, condnum, iunit, norml1, norml2, use_arnoldi, blacs_env)
64 
65  TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrixkp_s
66  REAL(kind=dp), DIMENSION(2), INTENT(INOUT) :: condnum
67  INTEGER, INTENT(IN) :: iunit
68  LOGICAL, INTENT(IN) :: norml1, norml2, use_arnoldi
69  TYPE(cp_blacs_env_type), POINTER :: blacs_env
70 
71  CHARACTER(len=*), PARAMETER :: routinen = 'overlap_condnum'
72 
73  INTEGER :: handle, ic, maxiter, nbas, ndep
74  LOGICAL :: converged
75  REAL(kind=dp) :: amnorm, anorm, eps_ev, max_ev, min_ev, &
76  threshold
77  REAL(kind=dp), DIMENSION(2) :: eigvals
78  TYPE(cp_fm_struct_type), POINTER :: matrix_struct
79  TYPE(cp_fm_type) :: fmsmat, fmwork
80  TYPE(dbcsr_type) :: tempmat
81  TYPE(dbcsr_type), POINTER :: smat
82 
83  CALL timeset(routinen, handle)
84 
85  condnum(1:2) = 0.0_dp
86  NULLIFY (smat)
87  IF (SIZE(matrixkp_s, 2) == 1) THEN
88  IF (iunit > 0) WRITE (iunit, '(/,T2,A)') "OVERLAP MATRIX CONDITION NUMBER"
89  smat => matrixkp_s(1, 1)%matrix
90  ELSE
91  IF (iunit > 0) WRITE (iunit, '(/,T2,A)') "OVERLAP MATRIX CONDITION NUMBER AT GAMMA POINT"
92  ALLOCATE (smat)
93  CALL dbcsr_create(smat, template=matrixkp_s(1, 1)%matrix)
94  CALL dbcsr_copy(smat, matrixkp_s(1, 1)%matrix)
95  DO ic = 2, SIZE(matrixkp_s, 2)
96  CALL dbcsr_add(smat, matrixkp_s(1, ic)%matrix, 1.0_dp, 1.0_dp)
97  END DO
98  END IF
99  !
100  IF (ASSOCIATED(smat)) THEN
101  cpassert(dbcsr_get_matrix_type(smat) .EQ. dbcsr_type_symmetric)
102  IF (norml1) THEN
103  ! norm of S
104  anorm = dbcsr_gershgorin_norm(smat)
105  CALL estimate_norm_invmat(smat, amnorm)
106  IF (iunit > 0) THEN
107  WRITE (iunit, '(T2,A)') "1-Norm Condition Number (Estimate)"
108  WRITE (iunit, '(T4,A,ES11.3E3,T32,A,ES11.3E3,A4,ES11.3E3,T63,A,F8.4)') &
109  "CN : |A|*|A^-1|: ", anorm, " * ", amnorm, "=", anorm*amnorm, "Log(1-CN):", log10(anorm*amnorm)
110  END IF
111  condnum(1) = anorm*amnorm
112  END IF
113  IF (norml2) THEN
114  eps_ev = 1.0e-14_dp
115  ! diagonalization
116  CALL dbcsr_get_info(smat, nfullrows_total=nbas)
117  CALL cp_fm_struct_create(fmstruct=matrix_struct, context=blacs_env, &
118  nrow_global=nbas, ncol_global=nbas)
119  CALL cp_fm_create(fmsmat, matrix_struct)
120  CALL cp_fm_create(fmwork, matrix_struct)
121  ! transfer to FM
122  CALL dbcsr_create(tempmat, template=smat, matrix_type=dbcsr_type_no_symmetry)
123  CALL dbcsr_desymmetrize(smat, tempmat)
124  CALL copy_dbcsr_to_fm(tempmat, fmsmat)
125 
126  ! diagonalize
127  anorm = cp_fm_norm(fmsmat, "1")
128  CALL cp_fm_power(fmsmat, fmwork, -1.0_dp, eps_ev, ndep, eigvals=eigvals)
129  min_ev = eigvals(1)
130  max_ev = eigvals(2)
131  amnorm = cp_fm_norm(fmsmat, "1")
132 
133  CALL dbcsr_release(tempmat)
134  CALL cp_fm_release(fmsmat)
135  CALL cp_fm_release(fmwork)
136  CALL cp_fm_struct_release(matrix_struct)
137 
138  IF (iunit > 0) THEN
139  WRITE (iunit, '(T2,A)') "1-Norm and 2-Norm Condition Numbers using Diagonalization"
140  IF (min_ev > 0) THEN
141  WRITE (iunit, '(T4,A,ES11.3E3,T32,A,ES11.3E3,A4,ES11.3E3,T63,A,F8.4)') &
142  "CN : |A|*|A^-1|: ", anorm, " * ", amnorm, "=", anorm*amnorm, "Log(1-CN):", log10(anorm*amnorm)
143  WRITE (iunit, '(T4,A,ES11.3E3,T32,A,ES11.3E3,A4,ES11.3E3,T63,A,F8.4)') &
144  "CN : max/min ev: ", max_ev, " / ", min_ev, "=", max_ev/min_ev, "Log(2-CN):", log10(max_ev/min_ev)
145  ELSE
146  WRITE (iunit, '(T4,A,ES11.3E3,T32,A,ES11.3E3,T63,A)') &
147  "CN : max/min EV: ", max_ev, " / ", min_ev, "Log(CN): infinity"
148  END IF
149  END IF
150  IF (min_ev > 0) THEN
151  condnum(1) = anorm*amnorm
152  condnum(2) = max_ev/min_ev
153  ELSE
154  condnum(1:2) = 0.0_dp
155  END IF
156  END IF
157  IF (use_arnoldi) THEN
158  ! parameters for matrix condition test
159  threshold = 1.0e-6_dp
160  maxiter = 1000
161  eps_ev = 1.0e8_dp
162  CALL arnoldi_extremal(smat, max_ev, min_ev, &
163  threshold=threshold, max_iter=maxiter, converged=converged)
164  IF (iunit > 0) THEN
165  WRITE (iunit, '(T2,A)') "2-Norm Condition Number using Arnoldi iterations"
166  IF (min_ev > 0) THEN
167  WRITE (iunit, '(T4,A,ES11.3E3,T32,A,ES11.3E3,A4,ES11.3E3,T63,A,F8.4)') &
168  "CN : max/min ev: ", max_ev, " / ", min_ev, "=", max_ev/min_ev, "Log(2-CN):", log10(max_ev/min_ev)
169  ELSE
170  WRITE (iunit, '(T4,A,ES11.3E3,T32,A,ES11.3E3,T63,A)') &
171  "CN : max/min ev: ", max_ev, " / ", min_ev, "Log(CN): infinity"
172  END IF
173  END IF
174  IF (min_ev > 0) THEN
175  condnum(2) = max_ev/min_ev
176  ELSE
177  condnum(2) = 0.0_dp
178  END IF
179  IF (converged) THEN
180  IF (min_ev == 0) THEN
181  cpwarn("Ill-conditioned S matrix: basis set is overcomplete.")
182  ELSE IF ((max_ev/min_ev) > eps_ev) THEN
183  cpwarn("Ill-conditioned S matrix: basis set is overcomplete.")
184  END IF
185  ELSE
186  cpwarn("Condition number estimate of overlap matrix is not reliable (not converged).")
187  END IF
188  END IF
189  END IF
190  IF (SIZE(matrixkp_s, 2) == 1) THEN
191  NULLIFY (smat)
192  ELSE
193  CALL dbcsr_release(smat)
194  DEALLOCATE (smat)
195  END IF
196 
197  CALL timestop(handle)
198 
199  END SUBROUTINE overlap_condnum
200 
201 ! **************************************************************************************************
202 !> \brief Calculates an estimate of the 1-norm of the inverse of a matrix
203 !> Uses LAPACK norm estimator algorithm
204 !> NJ Higham, Function of Matrices, Algorithm 3.21, page 66
205 !> \param amat Sparse, symmetric matrix
206 !> \param anorm Estimate of ||INV(A)||
207 !> \date 15.11.2016
208 !> \par History
209 !> \author JHU
210 !> \version 1.0
211 ! **************************************************************************************************
212  SUBROUTINE estimate_norm_invmat(amat, anorm)
213  TYPE(dbcsr_type), POINTER :: amat
214  REAL(kind=dp), INTENT(OUT) :: anorm
215 
216  INTEGER :: i, k, nbas
217  INTEGER, DIMENSION(1) :: r
218  REAL(kind=dp) :: g, gg
219  REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: x, xsi
220  REAL(kind=dp), DIMENSION(2) :: work
221  REAL(kind=dp), EXTERNAL :: dlange
222  TYPE(dbcsr_type) :: pmat
223 
224  ! generate a block diagonal preconditioner
225  CALL dbcsr_create(pmat, name="SMAT Preconditioner", template=amat)
226  ! replicate the diagonal blocks of the overlap matrix
227  CALL dbcsr_get_block_diag(amat, pmat)
228  ! invert preconditioner
229  CALL smat_precon_diag(pmat)
230 
231  anorm = 1.0_dp
232  CALL dbcsr_get_info(amat, nfullrows_total=nbas)
233  ALLOCATE (x(nbas), xsi(nbas))
234  x(1:nbas) = 1._dp/real(nbas, kind=dp)
235  CALL dbcsr_solve(amat, x, pmat)
236  g = dlange("1", nbas, 1, x, nbas, work)
237  xsi(1:nbas) = sign(1._dp, x(1:nbas))
238  x(1:nbas) = xsi(1:nbas)
239  CALL dbcsr_solve(amat, x, pmat)
240  k = 2
241  DO
242  r = maxloc(abs(x))
243  x(1:nbas) = 0._dp
244  x(r) = 1._dp
245  CALL dbcsr_solve(amat, x, pmat)
246  gg = g
247  g = dlange("1", nbas, 1, x, nbas, work)
248  IF (g <= gg) EXIT
249  x(1:nbas) = sign(1._dp, x(1:nbas))
250  IF (sum(abs(x - xsi)) == 0 .OR. sum(abs(x + xsi)) == 0) EXIT
251  xsi(1:nbas) = x(1:nbas)
252  CALL dbcsr_solve(amat, x, pmat)
253  k = k + 1
254  IF (k > 5) EXIT
255  IF (sum(r) == sum(maxloc(abs(x)))) EXIT
256  END DO
257  !
258  IF (nbas > 1) THEN
259  DO i = 1, nbas
260  x(i) = -1._dp**(i + 1)*(1._dp + real(i - 1, dp)/real(nbas - 1, dp))
261  END DO
262  ELSE
263  x(1) = 1.0_dp
264  END IF
265  CALL dbcsr_solve(amat, x, pmat)
266  gg = dlange("1", nbas, 1, x, nbas, work)
267  gg = 2._dp*gg/real(3*nbas, dp)
268  anorm = max(g, gg)
269  DEALLOCATE (x, xsi)
270  CALL dbcsr_release(pmat)
271 
272  END SUBROUTINE estimate_norm_invmat
273 
274 ! **************************************************************************************************
275 !> \brief ...
276 !> \param amat ...
277 !> \param x ...
278 !> \param pmat ...
279 ! **************************************************************************************************
280  SUBROUTINE dbcsr_solve(amat, x, pmat)
281  TYPE(dbcsr_type), POINTER :: amat
282  REAL(kind=dp), DIMENSION(:), INTENT(INOUT) :: x
283  TYPE(dbcsr_type) :: pmat
284 
285  INTEGER :: max_iter, nbas
286  LOGICAL :: converged
287  REAL(kind=dp) :: threshold
288 
289  CALL dbcsr_get_info(amat, nfullrows_total=nbas)
290  max_iter = min(1000, nbas)
291  threshold = 1.e-6_dp
292  CALL arnoldi_conjugate_gradient(amat, x, pmat, converged=converged, threshold=threshold, max_iter=max_iter)
293 
294  END SUBROUTINE dbcsr_solve
295 
296 ! **************************************************************************************************
297 !> \brief ...
298 !> \param pmat ...
299 ! **************************************************************************************************
300  SUBROUTINE smat_precon_diag(pmat)
301  TYPE(dbcsr_type) :: pmat
302 
303  INTEGER :: iatom, info, jatom, n
304  REAL(kind=dp), DIMENSION(:, :), POINTER :: sblock
305  TYPE(dbcsr_iterator_type) :: dbcsr_iter
306 
307  CALL dbcsr_iterator_start(dbcsr_iter, pmat)
308  DO WHILE (dbcsr_iterator_blocks_left(dbcsr_iter))
309  CALL dbcsr_iterator_next_block(dbcsr_iter, iatom, jatom, sblock)
310  cpassert(iatom == jatom)
311  n = SIZE(sblock, 1)
312  CALL invmat(sblock(1:n, 1:n), info)
313  cpassert(info == 0)
314  END DO
315  CALL dbcsr_iterator_stop(dbcsr_iter)
316 
317  END SUBROUTINE smat_precon_diag
318 
319 END MODULE qs_condnum
320 
arnoldi iteration using dbcsr
Definition: arnoldi_api.F:15
subroutine, public arnoldi_extremal(matrix_a, max_ev, min_ev, converged, threshold, max_iter)
simple wrapper to estimate extremal eigenvalues with arnoldi, using the old lanczos interface this hi...
Definition: arnoldi_api.F:254
subroutine, public arnoldi_conjugate_gradient(matrix_a, vec_x, matrix_p, converged, threshold, max_iter)
Wrapper for conjugated gradient algorithm for Ax=b.
Definition: arnoldi_api.F:303
methods related to the blacs parallel environment
Definition: cp_blacs_env.F:15
DBCSR operations in CP2K.
subroutine, public copy_dbcsr_to_fm(matrix, fm)
Copy a DBCSR matrix to a BLACS matrix.
basic linear algebra operations for full matrices
real(kind=dp) function, public cp_fm_norm(matrix, mode)
norm of matrix using (p)dlange
used for collecting some of the diagonalization schemes available for cp_fm_type. cp_fm_power also mo...
Definition: cp_fm_diag.F:17
subroutine, public cp_fm_power(matrix, work, exponent, threshold, n_dependent, verbose, eigvals)
...
Definition: cp_fm_diag.F:896
represent the structure of a full matrix
Definition: cp_fm_struct.F:14
subroutine, public cp_fm_struct_create(fmstruct, para_env, context, nrow_global, ncol_global, nrow_block, ncol_block, descriptor, first_p_pos, local_leading_dimension, template_fmstruct, square_blocks, force_block)
allocates and initializes a full matrix structure
Definition: cp_fm_struct.F:125
subroutine, public cp_fm_struct_release(fmstruct)
releases a full matrix structure
Definition: cp_fm_struct.F:320
represent a full matrix distributed on many processors
Definition: cp_fm_types.F:15
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp)
creates a new full matrix with the given structure
Definition: cp_fm_types.F:167
Defines the basic variable types.
Definition: kinds.F:23
integer, parameter, public dp
Definition: kinds.F:34
Collection of simple mathematical functions and subroutines.
Definition: mathlib.F:15
subroutine, public invmat(a, info)
returns inverse of matrix using the lapack routines DGETRF and DGETRI
Definition: mathlib.F:548
Calculation of overlap matrix condition numbers.
Definition: qs_condnum.F:13
subroutine, public overlap_condnum(matrixkp_s, condnum, iunit, norml1, norml2, use_arnoldi, blacs_env)
Calculation of the overlap matrix Condition Number.
Definition: qs_condnum.F:64