(git:9a40378)
Loading...
Searching...
No Matches
ri_environment_methods.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Calculates integral matrices for RIGPW method
10!> \par History
11!> created JGH [08.2012]
12!> Dorothea Golze [02.2014] (1) extended, re-structured, cleaned
13!> (2) heavily debugged
14!> updated for RI JGH [08.2017]
15!> \authors JGH
16!> Dorothea Golze
17! **************************************************************************************************
24 USE cp_dbcsr_api, ONLY: &
28 dbcsr_type_antisymmetric, dbcsr_type_no_symmetry, dbcsr_type_symmetric
29 USE cp_dbcsr_contrib, ONLY: dbcsr_dot,&
35 USE kinds, ONLY: dp
42 USE message_passing, ONLY: mp_comm_type,&
44 USE pw_types, ONLY: pw_c1d_gs_type,&
62 USE qs_rho_types, ONLY: qs_rho_get,&
64 USE util, ONLY: get_limit
65
66!$ USE OMP_LIB, ONLY: omp_get_max_threads, omp_get_thread_num
67#include "./base/base_uses.f90"
68
69 IMPLICIT NONE
70
71 PRIVATE
72
73! **************************************************************************************************
74
75 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'ri_environment_methods'
76
78
79! **************************************************************************************************
80
81CONTAINS
82
83! **************************************************************************************************
84!> \brief creates and initializes an lri_env
85!> \param lri_env the lri_environment you want to create
86!> \param qs_env ...
87!> \param calculate_forces ...
88! **************************************************************************************************
89 SUBROUTINE build_ri_matrices(lri_env, qs_env, calculate_forces)
90
91 TYPE(lri_environment_type), POINTER :: lri_env
92 TYPE(qs_environment_type), POINTER :: qs_env
93 LOGICAL, INTENT(IN) :: calculate_forces
94
95 CALL calculate_ri_integrals(lri_env, qs_env, calculate_forces)
96
97 END SUBROUTINE build_ri_matrices
98
99! **************************************************************************************************
100!> \brief calculates integrals needed for the RI density fitting,
101!> integrals are calculated once, before the SCF starts
102!> \param lri_env ...
103!> \param qs_env ...
104!> \param calculate_forces ...
105! **************************************************************************************************
106 SUBROUTINE calculate_ri_integrals(lri_env, qs_env, calculate_forces)
107
108 TYPE(lri_environment_type), POINTER :: lri_env
109 TYPE(qs_environment_type), POINTER :: qs_env
110 LOGICAL, INTENT(IN) :: calculate_forces
111
112 CHARACTER(LEN=*), PARAMETER :: routinen = 'calculate_ri_integrals'
113 REAL(kind=dp), DIMENSION(2), PARAMETER :: fx = [0.0_dp, 1.0_dp]
114
115 INTEGER :: handle, i, i1, i2, iatom, ispin, izero, &
116 j, j1, j2, jatom, m, n, nbas
117 INTEGER, DIMENSION(:, :), POINTER :: bas_ptr
118 REAL(kind=dp) :: fpre
119 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: eval
120 REAL(kind=dp), DIMENSION(:, :), POINTER :: avec, fblk, fout
121 TYPE(cp_blacs_env_type), POINTER :: blacs_env
122 TYPE(dbcsr_iterator_type) :: iter
123 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_p
124 TYPE(dbcsr_type) :: emat
125 TYPE(dbcsr_type), POINTER :: fmat
126 TYPE(dft_control_type), POINTER :: dft_control
127 TYPE(mp_para_env_type), POINTER :: para_env
128 TYPE(qs_ks_env_type), POINTER :: ks_env
129 TYPE(qs_rho_type), POINTER :: rho
130
131 CALL timeset(routinen, handle)
132
133 cpassert(ASSOCIATED(lri_env))
134 cpassert(ASSOCIATED(qs_env))
135
136 ! overlap matrices
137 CALL get_qs_env(qs_env=qs_env, ks_env=ks_env)
138 CALL get_qs_env(qs_env=qs_env, dft_control=dft_control)
139 NULLIFY (rho, matrix_p)
140 ! orbital overlap matrix (needed for N constraints and forces)
141 ! we replicate this in order to not directly interact qith qs_env
142 IF (calculate_forces) THEN
143 ! charge constraint forces are calculated here
144 CALL get_qs_env(qs_env=qs_env, rho=rho)
145 CALL qs_rho_get(rho, rho_ao=matrix_p)
146 ALLOCATE (fmat)
147 CALL dbcsr_create(fmat, template=matrix_p(1)%matrix)
148 DO ispin = 1, dft_control%nspins
149 fpre = lri_env%ri_fit%ftrm1n(ispin)/lri_env%ri_fit%ntrm1n
150 CALL dbcsr_add(fmat, matrix_p(ispin)%matrix, fx(ispin), -fpre)
151 END DO
152 CALL build_overlap_matrix(ks_env=ks_env, matrix_s=lri_env%ob_smat, nderivative=0, &
153 basis_type_a="ORB", basis_type_b="ORB", sab_nl=lri_env%soo_list, &
154 calculate_forces=.true., matrix_p=fmat)
155 CALL dbcsr_release(fmat)
156 DEALLOCATE (fmat)
157 ELSE
158 CALL build_overlap_matrix(ks_env=ks_env, matrix_s=lri_env%ob_smat, nderivative=0, &
159 basis_type_a="ORB", basis_type_b="ORB", sab_nl=lri_env%soo_list)
160 END IF
161 ! RI overlap matrix
162 CALL build_overlap_matrix(ks_env=ks_env, matrix_s=lri_env%ri_smat, nderivative=0, &
163 basis_type_a="RI_HXC", basis_type_b="RI_HXC", sab_nl=lri_env%saa_list)
164 IF (calculate_forces) THEN
165 ! calculate f*a(T) pseudo density matrix for forces
166 bas_ptr => lri_env%ri_fit%bas_ptr
167 avec => lri_env%ri_fit%avec
168 fout => lri_env%ri_fit%fout
169 ALLOCATE (fmat)
170 CALL dbcsr_create(fmat, template=lri_env%ri_smat(1)%matrix)
171 CALL cp_dbcsr_alloc_block_from_nbl(fmat, lri_env%saa_list)
172 CALL dbcsr_set(fmat, 0.0_dp)
173 CALL dbcsr_iterator_start(iter, fmat)
174 DO WHILE (dbcsr_iterator_blocks_left(iter))
175 CALL dbcsr_iterator_next_block(iter, iatom, jatom, fblk)
176 i1 = bas_ptr(1, iatom)
177 i2 = bas_ptr(2, iatom)
178 j1 = bas_ptr(1, jatom)
179 j2 = bas_ptr(2, jatom)
180 IF (iatom <= jatom) THEN
181 DO ispin = 1, dft_control%nspins
182 DO j = j1, j2
183 m = j - j1 + 1
184 DO i = i1, i2
185 n = i - i1 + 1
186 fblk(n, m) = fblk(n, m) + fout(i, ispin)*avec(j, ispin) + avec(i, ispin)*fout(j, ispin)
187 END DO
188 END DO
189 END DO
190 ELSE
191 DO ispin = 1, dft_control%nspins
192 DO i = i1, i2
193 n = i - i1 + 1
194 DO j = j1, j2
195 m = j - j1 + 1
196 fblk(m, n) = fblk(m, n) + fout(i, ispin)*avec(j, ispin) + avec(i, ispin)*fout(j, ispin)
197 END DO
198 END DO
199 END DO
200 END IF
201 IF (iatom == jatom) THEN
202 fblk(:, :) = 0.25_dp*fblk(:, :)
203 ELSE
204 fblk(:, :) = 0.5_dp*fblk(:, :)
205 END IF
206 END DO
207 CALL dbcsr_iterator_stop(iter)
208 !
209 CALL build_overlap_matrix(ks_env=ks_env, matrix_s=lri_env%ri_smat, nderivative=0, &
210 basis_type_a="RI_HXC", basis_type_b="RI_HXC", sab_nl=lri_env%saa_list, &
211 calculate_forces=.true., matrix_p=fmat)
212 CALL dbcsr_release(fmat)
213 DEALLOCATE (fmat)
214 END IF
215 ! approximation (preconditioner) or exact inverse of RI overlap
216 CALL dbcsr_allocate_matrix_set(lri_env%ri_sinv, 1)
217 ALLOCATE (lri_env%ri_sinv(1)%matrix)
218 SELECT CASE (lri_env%ri_sinv_app)
219 CASE ("NONE")
220 !
221 CASE ("INVS")
222 CALL dbcsr_create(lri_env%ri_sinv(1)%matrix, template=lri_env%ri_smat(1)%matrix)
223 CALL invert_hotelling(lri_env%ri_sinv(1)%matrix, lri_env%ri_smat(1)%matrix, &
224 threshold=1.e-10_dp, use_inv_as_guess=.false., &
225 norm_convergence=1.e-10_dp, filter_eps=1.e-12_dp, silent=.false.)
226 CASE ("INVF")
227 CALL dbcsr_create(emat, matrix_type=dbcsr_type_no_symmetry, template=lri_env%ri_smat(1)%matrix)
228 CALL get_qs_env(qs_env=qs_env, para_env=para_env, blacs_env=blacs_env)
229 CALL dbcsr_get_info(lri_env%ri_smat(1)%matrix, nfullrows_total=nbas)
230 ALLOCATE (eval(nbas))
231 CALL cp_dbcsr_syevd(lri_env%ri_smat(1)%matrix, emat, eval, para_env, blacs_env)
232 izero = 0
233 DO i = 1, nbas
234 IF (eval(i) < 1.0e-10_dp) THEN
235 eval(i) = 0.0_dp
236 izero = izero + 1
237 ELSE
238 eval(i) = sqrt(1.0_dp/eval(i))
239 END IF
240 END DO
241 CALL dbcsr_scale_by_vector(emat, eval, side='right')
242 CALL dbcsr_create(lri_env%ri_sinv(1)%matrix, template=lri_env%ri_smat(1)%matrix)
243 CALL dbcsr_multiply("N", "T", 1.0_dp, emat, emat, 0.0_dp, lri_env%ri_sinv(1)%matrix)
244 DEALLOCATE (eval)
245 CALL dbcsr_release(emat)
246 CASE ("AINV")
247 CALL dbcsr_create(lri_env%ri_sinv(1)%matrix, template=lri_env%ri_smat(1)%matrix)
248 CALL invert_hotelling(lri_env%ri_sinv(1)%matrix, lri_env%ri_smat(1)%matrix, &
249 threshold=1.e-5_dp, use_inv_as_guess=.false., &
250 norm_convergence=1.e-4_dp, filter_eps=1.e-4_dp, silent=.false.)
251 CASE DEFAULT
252 cpabort("Unknown RI_SINV type")
253 END SELECT
254
255 ! solve Rx=n
256 CALL ri_metric_solver(mat=lri_env%ri_smat(1)%matrix, &
257 vecr=lri_env%ri_fit%nvec, &
258 vecx=lri_env%ri_fit%rm1n, &
259 matp=lri_env%ri_sinv(1)%matrix, &
260 solver=lri_env%ri_sinv_app, &
261 ptr=lri_env%ri_fit%bas_ptr)
262
263 ! calculate n(t)x
264 lri_env%ri_fit%ntrm1n = sum(lri_env%ri_fit%nvec(:)*lri_env%ri_fit%rm1n(:))
265
266 ! calculate 3c-overlap integrals
267 IF (ASSOCIATED(lri_env%o3c)) THEN
268 CALL release_o3c_container(lri_env%o3c)
269 ELSE
270 ALLOCATE (lri_env%o3c)
271 END IF
272 CALL init_o3c_container(lri_env%o3c, dft_control%nspins, &
273 lri_env%orb_basis, lri_env%orb_basis, lri_env%ri_basis, &
274 lri_env%soo_list, lri_env%soa_list)
275
276 CALL calculate_o3c_integrals(lri_env%o3c, calculate_forces, matrix_p)
277
278 CALL timestop(handle)
279
280 END SUBROUTINE calculate_ri_integrals
281! **************************************************************************************************
282!> \brief solver for RI systems (R*x=n)
283!> \param mat ...
284!> \param vecr ...
285!> \param vecx ...
286!> \param matp ...
287!> \param solver ...
288!> \param ptr ...
289! **************************************************************************************************
290 SUBROUTINE ri_metric_solver(mat, vecr, vecx, matp, solver, ptr)
291
292 TYPE(dbcsr_type) :: mat
293 REAL(kind=dp), DIMENSION(:), INTENT(IN) :: vecr
294 REAL(kind=dp), DIMENSION(:), INTENT(OUT) :: vecx
295 TYPE(dbcsr_type) :: matp
296 CHARACTER(LEN=*), INTENT(IN) :: solver
297 INTEGER, DIMENSION(:, :), INTENT(IN) :: ptr
298
299 CHARACTER(LEN=*), PARAMETER :: routinen = 'ri_metric_solver'
300
301 INTEGER :: handle, max_iter, n
302 LOGICAL :: converged
303 REAL(kind=dp) :: rerror, threshold
304 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: vect
305
306 CALL timeset(routinen, handle)
307
308 threshold = 1.e-10_dp
309 max_iter = 100
310
311 vecx(:) = vecr(:)
312
313 SELECT CASE (solver)
314 CASE ("NONE")
315 CALL arnoldi_conjugate_gradient(mat, vecx, &
316 converged=converged, threshold=threshold, max_iter=max_iter)
317 CASE ("INVS", "INVF")
318 converged = .false.
319 CALL ri_matvec(matp, vecr, vecx, ptr)
320 CASE ("AINV")
321 CALL arnoldi_conjugate_gradient(mat, vecx, matrix_p=matp, &
322 converged=converged, threshold=threshold, max_iter=max_iter)
323 CASE DEFAULT
324 cpabort("Unknown RI solver")
325 END SELECT
326
327 IF (.NOT. converged) THEN
328 ! get error
329 rerror = 0.0_dp
330 n = SIZE(vecr)
331 ALLOCATE (vect(n))
332 CALL ri_matvec(mat, vecx, vect, ptr)
333 vect(:) = vect(:) - vecr(:)
334 rerror = maxval(abs(vect(:)))
335 DEALLOCATE (vect)
336 cpwarn_if(rerror > threshold, "RI solver: CG did not converge properly")
337 END IF
338
339 CALL timestop(handle)
340
341 END SUBROUTINE ri_metric_solver
342
343! **************************************************************************************************
344!> \brief performs the fitting of the density and distributes the fitted
345!> density on the grid
346!> \param lri_env the lri environment
347!> lri_density the environment for the fitting
348!> pmatrix density matrix
349!> lri_rho_struct where the fitted density is stored
350!> \param qs_env ...
351!> \param pmatrix ...
352!> \param lri_rho_struct ...
353!> \param atomic_kind_set ...
354!> \param para_env ...
355! **************************************************************************************************
356 SUBROUTINE calculate_ri_densities(lri_env, qs_env, pmatrix, &
357 lri_rho_struct, atomic_kind_set, para_env)
358
359 TYPE(lri_environment_type), POINTER :: lri_env
360 TYPE(qs_environment_type), POINTER :: qs_env
361 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: pmatrix
362 TYPE(qs_rho_type), INTENT(IN) :: lri_rho_struct
363 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
364 TYPE(mp_para_env_type), POINTER :: para_env
365
366 CHARACTER(LEN=*), PARAMETER :: routinen = 'calculate_ri_densities'
367
368 INTEGER :: atom_a, handle, i1, i2, iatom, ikind, &
369 ispin, n, natom, nspin
370 INTEGER, ALLOCATABLE, DIMENSION(:) :: atom_of_kind, kind_of
371 INTEGER, DIMENSION(:, :), POINTER :: bas_ptr
372 REAL(kind=dp), DIMENSION(:), POINTER :: tot_rho_r
373 REAL(kind=dp), DIMENSION(:, :), POINTER :: avec
374 TYPE(lri_density_type), POINTER :: lri_density
375 TYPE(lri_kind_type), DIMENSION(:), POINTER :: lri_coef
376 TYPE(pw_c1d_gs_type), DIMENSION(:), POINTER :: rho_g
377 TYPE(pw_r3d_rs_type), DIMENSION(:), POINTER :: rho_r
378
379 CALL timeset(routinen, handle)
380
381 nspin = SIZE(pmatrix, 1)
382 CALL contract12_o3c(lri_env%o3c, pmatrix)
383 CALL calculate_tvec_ri(lri_env, lri_env%o3c, para_env)
384 CALL calculate_avec_ri(lri_env, pmatrix)
385 !
386 CALL get_qs_env(qs_env, lri_density=lri_density, natom=natom)
387 IF (ASSOCIATED(lri_density)) THEN
388 CALL lri_density_release(lri_density)
389 ELSE
390 ALLOCATE (lri_density)
391 END IF
392 CALL lri_density_create(lri_density)
393 lri_density%nspin = nspin
394 ! allocate the arrays to hold RI expansion coefficients lri_coefs
395 CALL allocate_lri_coefs(lri_env, lri_density, atomic_kind_set)
396 ! reassign avec
397 avec => lri_env%ri_fit%avec
398 bas_ptr => lri_env%ri_fit%bas_ptr
399 CALL get_atomic_kind_set(atomic_kind_set, atom_of_kind=atom_of_kind, kind_of=kind_of)
400 DO ispin = 1, nspin
401 lri_coef => lri_density%lri_coefs(ispin)%lri_kinds
402 DO iatom = 1, natom
403 ikind = kind_of(iatom)
404 atom_a = atom_of_kind(iatom)
405 i1 = bas_ptr(1, iatom)
406 i2 = bas_ptr(2, iatom)
407 n = i2 - i1 + 1
408 lri_coef(ikind)%acoef(atom_a, 1:n) = avec(i1:i2, ispin)
409 END DO
410 END DO
411 CALL set_qs_env(qs_env, lri_density=lri_density)
412 DEALLOCATE (atom_of_kind, kind_of)
413 !
414 CALL qs_rho_get(lri_rho_struct, rho_r=rho_r, rho_g=rho_g, tot_rho_r=tot_rho_r)
415 DO ispin = 1, nspin
416 lri_coef => lri_density%lri_coefs(ispin)%lri_kinds
417 CALL calculate_lri_rho_elec(rho_g(ispin), rho_r(ispin), qs_env, &
418 lri_coef, tot_rho_r(ispin), "RI_HXC", .false.)
419 END DO
420
421 CALL timestop(handle)
422
423 END SUBROUTINE calculate_ri_densities
424
425! **************************************************************************************************
426!> \brief assembles the global vector t=<P.T>
427!> \param lri_env the lri environment
428!> \param o3c ...
429!> \param para_env ...
430! **************************************************************************************************
431 SUBROUTINE calculate_tvec_ri(lri_env, o3c, para_env)
432
433 TYPE(lri_environment_type), POINTER :: lri_env
434 TYPE(o3c_container_type), POINTER :: o3c
435 TYPE(mp_para_env_type), POINTER :: para_env
436
437 CHARACTER(LEN=*), PARAMETER :: routinen = 'calculate_tvec_ri'
438 INTEGER, PARAMETER :: msweep = 32
439
440 INTEGER :: handle, i1, i2, ibl, ibu, il, ispin, &
441 isweep, it, iu, katom, m, ma, mba, &
442 mepos, natom, nspin, nsweep, nthread
443 INTEGER, DIMENSION(2, msweep) :: nlimit
444 INTEGER, DIMENSION(:, :), POINTER :: bas_ptr
445 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: ta
446 REAL(kind=dp), DIMENSION(:, :), POINTER :: rm1t, tvec, tvl
447 TYPE(o3c_iterator_type) :: o3c_iterator
448
449 CALL timeset(routinen, handle)
450
451 nspin = SIZE(lri_env%ri_fit%tvec, 2)
452 bas_ptr => lri_env%ri_fit%bas_ptr
453 tvec => lri_env%ri_fit%tvec
454 tvec = 0.0_dp
455 natom = SIZE(bas_ptr, 2)
456 nthread = 1
457!$ nthread = omp_get_max_threads()
458
459 IF (natom < 1000) THEN
460 nsweep = 1
461 ELSE
462 nsweep = min(nthread, msweep)
463 END IF
464
465 nlimit = 0
466 DO isweep = 1, nsweep
467 nlimit(1:2, isweep) = get_limit(natom, nsweep, isweep - 1)
468 END DO
469
470 DO ispin = 1, nspin
471 DO isweep = 1, nsweep
472 il = nlimit(1, isweep)
473 iu = nlimit(2, isweep)
474 ma = iu - il + 1
475 IF (ma < 1) cycle
476 ibl = bas_ptr(1, il)
477 ibu = bas_ptr(2, iu)
478 mba = ibu - ibl + 1
479 ALLOCATE (ta(mba, nthread))
480 ta = 0.0_dp
481
482 CALL o3c_iterator_create(o3c, o3c_iterator, nthread=nthread)
483
484!$OMP PARALLEL DEFAULT(NONE)&
485!$OMP SHARED (nthread,o3c_iterator,ispin,il,iu,ibl,ta,bas_ptr)&
486!$OMP PRIVATE (mepos,katom,tvl,i1,i2,m)
487
488 mepos = 0
489!$ mepos = omp_get_thread_num()
490
491 DO WHILE (o3c_iterate(o3c_iterator, mepos=mepos) == 0)
492 CALL get_o3c_iterator_info(o3c_iterator, mepos=mepos, katom=katom, tvec=tvl)
493 IF (katom < il .OR. katom > iu) cycle
494 i1 = bas_ptr(1, katom) - ibl + 1
495 i2 = bas_ptr(2, katom) - ibl + 1
496 m = i2 - i1 + 1
497 ta(i1:i2, mepos + 1) = ta(i1:i2, mepos + 1) + tvl(1:m, ispin)
498 END DO
499!$OMP END PARALLEL
500 CALL o3c_iterator_release(o3c_iterator)
501
502 ! sum over threads
503 DO it = 1, nthread
504 tvec(ibl:ibu, ispin) = tvec(ibl:ibu, ispin) + ta(1:mba, it)
505 END DO
506 DEALLOCATE (ta)
507 END DO
508 END DO
509
510 ! global summation
511 CALL para_env%sum(tvec)
512
513 rm1t => lri_env%ri_fit%rm1t
514 DO ispin = 1, nspin
515 ! solve Rx=t
516 CALL ri_metric_solver(mat=lri_env%ri_smat(1)%matrix, &
517 vecr=tvec(:, ispin), &
518 vecx=rm1t(:, ispin), &
519 matp=lri_env%ri_sinv(1)%matrix, &
520 solver=lri_env%ri_sinv_app, &
521 ptr=bas_ptr)
522 END DO
523
524 CALL timestop(handle)
525
526 END SUBROUTINE calculate_tvec_ri
527! **************************************************************************************************
528!> \brief performs the fitting of the density in the RI method
529!> \param lri_env the lri environment
530!> lri_density the environment for the fitting
531!> pmatrix density matrix
532!> \param pmatrix ...
533! **************************************************************************************************
534 SUBROUTINE calculate_avec_ri(lri_env, pmatrix)
535
536 TYPE(lri_environment_type), POINTER :: lri_env
537 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: pmatrix
538
539 CHARACTER(LEN=*), PARAMETER :: routinen = 'calculate_avec_ri'
540
541 INTEGER :: handle, ispin, nspin
542 REAL(kind=dp) :: etr, nelec, nrm1t
543 REAL(kind=dp), DIMENSION(:), POINTER :: nvec, rm1n
544 REAL(kind=dp), DIMENSION(:, :), POINTER :: avec, rm1t, tvec
545 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: smatrix
546
547 CALL timeset(routinen, handle)
548
549 nspin = SIZE(pmatrix)
550 ! number of electrons
551 smatrix => lri_env%ob_smat
552 DO ispin = 1, nspin
553 etr = 0.0_dp
554 CALL dbcsr_dot(smatrix(1)%matrix, pmatrix(ispin)%matrix, etr)
555 lri_env%ri_fit%echarge(ispin) = etr
556 END DO
557
558 tvec => lri_env%ri_fit%tvec
559 rm1t => lri_env%ri_fit%rm1t
560 nvec => lri_env%ri_fit%nvec
561 rm1n => lri_env%ri_fit%rm1n
562
563 ! calculate lambda
564 DO ispin = 1, nspin
565 nelec = lri_env%ri_fit%echarge(ispin)
566 nrm1t = sum(nvec(:)*rm1t(:, ispin))
567 lri_env%ri_fit%lambda(ispin) = 2.0_dp*(nrm1t - nelec)/lri_env%ri_fit%ntrm1n
568 END DO
569
570 ! calculate avec = rm1t - lambda/2 * rm1n
571 avec => lri_env%ri_fit%avec
572 DO ispin = 1, nspin
573 avec(:, ispin) = rm1t(:, ispin) - 0.5_dp*lri_env%ri_fit%lambda(ispin)*rm1n(:)
574 END DO
575
576 CALL timestop(handle)
577
578 END SUBROUTINE calculate_avec_ri
579
580! **************************************************************************************************
581!> \brief Mutiplies a replicated vector with a DBCSR matrix: vo = mat*vi
582!> \param mat ...
583!> \param vi ...
584!> \param vo ...
585!> \param ptr ...
586! **************************************************************************************************
587 SUBROUTINE ri_matvec(mat, vi, vo, ptr)
588
589 TYPE(dbcsr_type) :: mat
590 REAL(kind=dp), DIMENSION(:), INTENT(IN) :: vi
591 REAL(kind=dp), DIMENSION(:), INTENT(OUT) :: vo
592 INTEGER, DIMENSION(:, :), INTENT(IN) :: ptr
593
594 CHARACTER(LEN=*), PARAMETER :: routinen = 'ri_matvec'
595
596 CHARACTER :: matrix_type
597 INTEGER :: handle, iatom, jatom, m1, m2, mb, n1, &
598 n2, nb
599 LOGICAL :: symm
600 REAL(kind=dp), DIMENSION(:, :), POINTER :: block
601 TYPE(dbcsr_iterator_type) :: iter
602 TYPE(mp_comm_type) :: group
603
604 CALL timeset(routinen, handle)
605
606 CALL dbcsr_get_info(mat, matrix_type=matrix_type, group=group)
607
608 SELECT CASE (matrix_type)
609 CASE (dbcsr_type_no_symmetry)
610 symm = .false.
611 CASE (dbcsr_type_symmetric)
612 symm = .true.
613 CASE (dbcsr_type_antisymmetric)
614 cpabort("NYI, antisymmetric matrix not permitted")
615 CASE DEFAULT
616 cpabort("Unknown matrix type, ...")
617 END SELECT
618
619 vo(:) = 0.0_dp
620 CALL dbcsr_iterator_start(iter, mat)
621 DO WHILE (dbcsr_iterator_blocks_left(iter))
622 CALL dbcsr_iterator_next_block(iter, iatom, jatom, block)
623 n1 = ptr(1, iatom)
624 n2 = ptr(2, iatom)
625 nb = n2 - n1 + 1
626 m1 = ptr(1, jatom)
627 m2 = ptr(2, jatom)
628 mb = m2 - m1 + 1
629 cpassert(nb == SIZE(block, 1))
630 cpassert(mb == SIZE(block, 2))
631 vo(n1:n2) = vo(n1:n2) + matmul(block, vi(m1:m2))
632 IF (symm .AND. (iatom /= jatom)) THEN
633 vo(m1:m2) = vo(m1:m2) + matmul(transpose(block), vi(n1:n2))
634 END IF
635 END DO
636 CALL dbcsr_iterator_stop(iter)
637
638 CALL group%sum(vo)
639
640 CALL timestop(handle)
641
642 END SUBROUTINE ri_matvec
643
644END MODULE ri_environment_methods
arnoldi iteration using dbcsr
Definition arnoldi_api.F:16
subroutine, public arnoldi_conjugate_gradient(matrix_a, vec_x, matrix_p, converged, threshold, max_iter)
Wrapper for conjugated gradient algorithm for Ax=b.
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind_set(atomic_kind_set, atom_of_kind, kind_of, natom_of_kind, maxatom, natom, nshell, fist_potential_present, shell_present, shell_adiabatic, shell_check_distance, damping_present)
Get attributes of an atomic kind set.
methods related to the blacs parallel environment
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_iterator_next_block(iterator, row, column, block, block_number_argument_has_been_removed, row_size, col_size, row_offset, col_offset)
...
logical function, public dbcsr_iterator_blocks_left(iterator)
...
subroutine, public dbcsr_iterator_stop(iterator)
...
subroutine, public dbcsr_multiply(transa, transb, alpha, matrix_a, matrix_b, beta, matrix_c, first_row, last_row, first_column, last_column, first_k, last_k, retain_sparsity, filter_eps, flop)
...
subroutine, public dbcsr_get_info(matrix, nblkrows_total, nblkcols_total, nfullrows_total, nfullcols_total, nblkrows_local, nblkcols_local, nfullrows_local, nfullcols_local, my_prow, my_pcol, local_rows, local_cols, proc_row_dist, proc_col_dist, row_blk_size, col_blk_size, row_blk_offset, col_blk_offset, distribution, name, matrix_type, group)
...
subroutine, public dbcsr_iterator_start(iterator, matrix, shared, dynamic, dynamic_byrows)
...
subroutine, public dbcsr_set(matrix, alpha)
...
subroutine, public dbcsr_release(matrix)
...
subroutine, public dbcsr_add(matrix_a, matrix_b, alpha_scalar, beta_scalar)
...
subroutine, public dbcsr_dot(matrix_a, matrix_b, trace)
Computes the dot product of two matrices, also known as the trace of their matrix product.
subroutine, public dbcsr_scale_by_vector(matrix, alpha, side)
Scales the rows/columns of given matrix.
Interface to (sca)lapack for the Cholesky based procedures.
subroutine, public cp_dbcsr_syevd(matrix, eigenvectors, eigenvalues, para_env, blacs_env)
...
DBCSR operations in CP2K.
Routines useful for iterative matrix calculations.
subroutine, public invert_hotelling(matrix_inverse, matrix, threshold, use_inv_as_guess, norm_convergence, filter_eps, accelerator_order, max_iter_lanczos, eps_lanczos, silent)
invert a symmetric positive definite matrix by Hotelling's method explicit symmetrization makes this ...
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
contains the types and subroutines for dealing with the lri_env lri : local resolution of the identit...
subroutine, public allocate_lri_coefs(lri_env, lri_density, atomic_kind_set)
creates and initializes lri_coefs
subroutine, public lri_density_release(lri_density)
releases the given lri_density
subroutine, public lri_density_create(lri_density)
creates and initializes an lri_density environment
Interface to the message passing library MPI.
Calculate the plane wave density by collocating the primitive Gaussian functions (pgf).
subroutine, public calculate_lri_rho_elec(lri_rho_g, lri_rho_r, qs_env, lri_coef, total_rho, basis_type, exact_1c_terms, pmat, atomlist)
Collocates the fitted lri density on a grid.
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, mimic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, sab_cneo, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, xcint_weights, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs, do_rixs, tb_tblite)
Get the QUICKSTEP environment.
subroutine, public set_qs_env(qs_env, super_cell, mos, qmmm, qmmm_periodic, mimic, ewald_env, ewald_pw, mpools, rho_external, external_vxc, mask, scf_control, rel_control, qs_charges, ks_env, ks_qmmm_env, wf_history, scf_env, active_space, input, oce, rho_atom_set, rho0_atom_set, rho0_mpole, run_rtp, rtp, rhoz_set, rhoz_tot, ecoul_1c, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, efield, rhoz_cneo_set, linres_control, xas_env, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, ls_scf_env, do_transport, transport_env, lri_env, lri_density, exstate_env, ec_env, dispersion_env, harris_env, gcp_env, mp2_env, bs_env, kg_env, force, kpoints, wanniercentres, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs, do_rixs, tb_tblite)
Set the QUICKSTEP environment.
Methods used with 3-center overlap type integrals containers.
subroutine, public calculate_o3c_integrals(o3c, calculate_forces, matrix_p)
...
subroutine, public contract12_o3c(o3c, matrix_p)
Contraction of 3-tensor over indices 1 and 2 (assuming symmetry) t(k) = sum_ij (ijk)*p(ij)
3-center overlap type integrals containers
subroutine, public get_o3c_iterator_info(o3c_iterator, mepos, iatom, jatom, katom, ikind, jkind, kkind, rij, rik, cellj, cellk, integral, tvec, force_i, force_j, force_k)
...
subroutine, public o3c_iterator_create(o3c, o3c_iterator, nthread)
...
subroutine, public release_o3c_container(o3c_container)
...
subroutine, public o3c_iterator_release(o3c_iterator)
...
integer function, public o3c_iterate(o3c_iterator, mepos)
...
subroutine, public init_o3c_container(o3c, nspin, basis_set_list_a, basis_set_list_b, basis_set_list_c, sab_nl, sac_nl, only_bc_same_center)
...
Calculation of overlap matrix, its derivatives and forces.
Definition qs_overlap.F:19
subroutine, public build_overlap_matrix(ks_env, matrix_s, matrixkp_s, matrix_name, nderivative, basis_type_a, basis_type_b, sab_nl, calculate_forces, matrix_p, matrixkp_p)
Calculation of the overlap matrix over Cartesian Gaussian functions.
Definition qs_overlap.F:120
superstucture that hold various representations of the density and keeps track of which ones are vali...
subroutine, public qs_rho_get(rho_struct, rho_ao, rho_ao_im, rho_ao_kp, rho_ao_im_kp, rho_r, drho_r, rho_g, drho_g, tau_r, tau_g, rho_r_valid, drho_r_valid, rho_g_valid, drho_g_valid, tau_r_valid, tau_g_valid, tot_rho_r, tot_rho_g, rho_r_sccs, soft_valid, complex_rho_ao)
returns info about the density described by this object. If some representation is not available an e...
Calculates integral matrices for RIGPW method.
subroutine, public build_ri_matrices(lri_env, qs_env, calculate_forces)
creates and initializes an lri_env
subroutine, public ri_metric_solver(mat, vecr, vecx, matp, solver, ptr)
solver for RI systems (R*x=n)
subroutine, public calculate_ri_densities(lri_env, qs_env, pmatrix, lri_rho_struct, atomic_kind_set, para_env)
performs the fitting of the density and distributes the fitted density on the grid
All kind of helpful little routines.
Definition util.F:14
pure integer function, dimension(2), public get_limit(m, n, me)
divide m entries into n parts, return size of part me
Definition util.F:333
Provides all information about an atomic kind.
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
stores all the informations relevant to an mpi environment
calculation environment to calculate the ks matrix, holds all the needed vars. assumes that the core ...
keeps the density in various representations, keeping track of which ones are valid.