(git:ccc2433)
lri_environment_types.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 contains the types and subroutines for dealing with the lri_env
10 !> lri : local resolution of the identity
11 !> \par History
12 !> created JGH [08.2012]
13 !> Dorothea Golze [02.2014] (1) extended, re-structured, cleaned
14 !> (2) debugged
15 !> \authors JGH
16 !> Dorothea Golze
17 ! **************************************************************************************************
19  USE atomic_kind_types, ONLY: atomic_kind_type,&
22  gto_basis_set_p_type,&
23  gto_basis_set_type
25  USE dbcsr_api, ONLY: dbcsr_p_type
26  USE kinds, ONLY: int_8,&
27  dp,&
28  sp
29  USE mathlib, ONLY: pswitch
33  neighbor_list_iterator_p_type,&
35  neighbor_list_set_p_type,&
37  USE qs_o3c_types, ONLY: o3c_container_type,&
39 #include "./base/base_uses.f90"
40 
41  IMPLICIT NONE
42 
43  PRIVATE
44 
45 ! **************************************************************************************************
46 ! Integral container
47  TYPE carray
48  INTEGER :: compression
49  REAL(KIND=dp), DIMENSION(:), POINTER :: cdp
50  REAL(KIND=sp), DIMENSION(:), POINTER :: csp
51  INTEGER(INT_8), DIMENSION(:), POINTER :: cip
52  END TYPE carray
53  TYPE int_container
54  INTEGER :: na, nb, nc
55  TYPE(carray), DIMENSION(:), POINTER :: ca
56  END TYPE int_container
57 ! **************************************************************************************************
58  TYPE lri_rhoab_type
59  ! number of spherical basis functions (a)
60  INTEGER :: nba
61  ! number of spherical basis functions (b)
62  INTEGER :: nbb
63  ! number of spherical fit basis functions (ai)
64  INTEGER :: nfa
65  ! number of spherical fit basis functions (bi)
66  INTEGER :: nfb
67  ! expansion coeffs for RI density
68  REAL(KIND=dp), DIMENSION(:), POINTER :: avec
69  REAL(KIND=dp), DIMENSION(:), POINTER :: aveca
70  REAL(KIND=dp), DIMENSION(:), POINTER :: avecb
71  ! projection coeffs for RI density: SUM_ab (ab,i)*Pab
72  REAL(KIND=dp), DIMENSION(:), POINTER :: tvec
73  REAL(KIND=dp), DIMENSION(:), POINTER :: tveca
74  REAL(KIND=dp), DIMENSION(:), POINTER :: tvecb
75  ! integral (ai) * sinv * tvec
76  REAL(KIND=dp) :: nst
77  REAL(KIND=dp) :: nsta
78  REAL(KIND=dp) :: nstb
79  ! Lagrange parameter
80  REAL(KIND=dp) :: lambda
81  REAL(KIND=dp) :: lambdaa
82  REAL(KIND=dp) :: lambdab
83  ! Charge of pair density
84  REAL(KIND=dp) :: charge
85  REAL(KIND=dp) :: chargea
86  REAL(KIND=dp) :: chargeb
87  END TYPE lri_rhoab_type
88 
89 ! **************************************************************************************************
90 
91  TYPE lri_int_type
92  ! whether to calculate force for pair
93  LOGICAL :: calc_force_pair
94  ! number of spherical basis functions (a)
95  INTEGER :: nba
96  ! number of spherical basis functions (b)
97  INTEGER :: nbb
98  ! number of spherical fit basis functions (ai)
99  INTEGER :: nfa
100  ! number of spherical fit basis functions (bi)
101  INTEGER :: nfb
102  ! condition number of overlap matrix
103  REAL(KIND=dp) :: cond_num
104  ! integrals (a,b,ai)
105  REAL(KIND=dp), DIMENSION(:, :, :), ALLOCATABLE :: abaint
106  REAL(KIND=dp), DIMENSION(:), ALLOCATABLE :: abascr
107  ! integrals (a,b,b)
108  REAL(KIND=dp), DIMENSION(:, :, :), ALLOCATABLE :: abbint
109  REAL(KIND=dp), DIMENSION(:), ALLOCATABLE :: abbscr
110  ! compressed aba integrals
111  TYPE(int_container) :: cabai
112  ! compressed abb integrals
113  TYPE(int_container) :: cabbi
114  ! integrals (da/dA,b,dai/dA)
115  REAL(KIND=dp), DIMENSION(:, :, :, :), ALLOCATABLE :: dabdaint
116  ! integrals (da/dA,b,bi)
117  REAL(KIND=dp), DIMENSION(:, :, :, :), ALLOCATABLE :: dabbint
118  ! integrals (a,b)
119  REAL(KIND=dp), DIMENSION(:, :), ALLOCATABLE :: soo
120  ! derivative d(a,b)/dA
121  REAL(KIND=dp), DIMENSION(:, :, :), ALLOCATABLE :: dsoo
122  ! integrals (ai,bi)
123  REAL(KIND=dp), DIMENSION(:, :), ALLOCATABLE :: sab
124  ! derivative d(ai,bi)/dA
125  REAL(KIND=dp), DIMENSION(:, :, :), ALLOCATABLE :: dsab
126  ! inverse of integrals (ai,bi)
127  REAL(KIND=dp), DIMENSION(:, :), POINTER :: sinv
128  ! integral (ai) / (bi), dim(1..nfa,nfa+1..nfa+nfb)
129  REAL(KIND=dp), DIMENSION(:), POINTER :: n
130  ! sinv * (ai)
131  REAL(KIND=dp), DIMENSION(:), POINTER :: sn
132  ! (ai) * sinv * (ai)
133  REAL(KIND=dp) :: nsn
134  ! distant pair approximation
135  LOGICAL :: lrisr
136  LOGICAL :: lriff
137  REAL(KIND=dp) :: wsr, wff, dwsr, dwff
138  REAL(KIND=dp), DIMENSION(:, :), POINTER :: asinv, bsinv
139  REAL(KIND=dp), DIMENSION(:), POINTER :: na, nb
140  REAL(KIND=dp), DIMENSION(:), POINTER :: sna, snb
141  REAL(KIND=dp) :: nsna, nsnb
142  !
143  ! dmax: max deviation for integrals of primitive gtos; for debugging
144  ! dmax for overlap integrals (ai,bi); fit bas
145  REAL(KIND=dp) :: dmax_ab
146  ! dmax for overlap integrals (a,b); orb bas
147  REAL(KIND=dp) :: dmax_oo
148  ! dmax for integrals (a,b,ai)
149  REAL(KIND=dp) :: dmax_aba
150  ! dmax for integrals (a,b,bi)
151  REAL(KIND=dp) :: dmax_abb
152  END TYPE lri_int_type
153 
154  TYPE lri_int_rho_type
155  ! integrals (aa,bb), orb basis
156  REAL(KIND=dp), DIMENSION(:, :, :, :), POINTER :: soaabb
157  ! dmax for (aa,bb) integrals; for debugging
158  REAL(KIND=dp) :: dmax_aabb
159  END TYPE lri_int_rho_type
160 
161  TYPE lri_node_type
162  INTEGER :: nnode = 0
163  TYPE(lri_int_type), DIMENSION(:), POINTER :: lri_int => null()
164  TYPE(lri_int_rho_type), DIMENSION(:), POINTER :: lri_int_rho => null()
165  TYPE(lri_rhoab_type), DIMENSION(:), POINTER :: lri_rhoab => null()
166  END TYPE lri_node_type
167 
168  TYPE lri_atom_type
169  INTEGER :: natom = 0
170  TYPE(lri_node_type), DIMENSION(:), POINTER :: lri_node => null()
171  END TYPE lri_atom_type
172 
173  TYPE lri_list_type
174  INTEGER :: nkind = 0
175  TYPE(lri_atom_type), DIMENSION(:), POINTER :: lri_atom => null()
176  END TYPE lri_list_type
177 
178  TYPE lri_list_p_type
179  TYPE(lri_list_type), POINTER :: lri_list => null()
180  END TYPE lri_list_p_type
181 
182 ! **************************************************************************************************
183 
184  TYPE lri_bas_type
185  INTEGER, DIMENSION(:, :, :), POINTER :: orb_index
186  INTEGER, DIMENSION(:, :, :), POINTER :: ri_index
187  ! integral of ri basis fbas
188  REAL(KIND=dp), DIMENSION(:), POINTER :: int_fbas
189  ! self overlap ri basis
190  REAL(KIND=dp), DIMENSION(:, :), POINTER :: ri_ovlp
191  ! inverse of self overlap ri basis
192  REAL(KIND=dp), DIMENSION(:, :), POINTER :: ri_ovlp_inv
193  ! self overlap orb basis
194  REAL(KIND=dp), DIMENSION(:, :), POINTER :: orb_ovlp
195  ! self overlap (a,a,fa)
196  REAL(KIND=dp), DIMENSION(:, :, :), ALLOCATABLE :: ovlp3
197  ! contraction matrix for SHG integrals ri basis
198  REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: scon_ri
199  ! contraction matrix for SHG integrals orb basis
200  REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: scon_orb
201  ! contraction matrix for SHG integrals aba/abb
202  REAL(KIND=dp), DIMENSION(:, :, :, :), POINTER :: scon_mix
203  END TYPE lri_bas_type
204 
205 ! **************************************************************************************************
206 
207  TYPE lri_clebsch_gordon_type
208  ! Clebsch-Gordon (CG) coefficients
209  REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: cg_coeff
210  ! list of non-zero CG coefficients
211  INTEGER, DIMENSION(:, :, :), POINTER :: cg_none0_list
212  ! number of non-zero CG coefficients
213  INTEGER, DIMENSION(:, :), POINTER :: ncg_none0
214  END TYPE lri_clebsch_gordon_type
215 
216 ! **************************************************************************************************
217 
218  TYPE lri_ppl_type
219  ! integrals Vppl*fbas (potential*fit basis) dim(natom,nsgf)
220  REAL(KIND=dp), DIMENSION(:, :), POINTER :: v_int
221  END TYPE lri_ppl_type
222 
223  TYPE lri_ppl_int_type
224  TYPE(lri_ppl_type), DIMENSION(:), POINTER :: lri_ppl
225  REAL(KIND=dp) :: ecore_pp_ri
226  END TYPE lri_ppl_int_type
227 
228 ! **************************************************************************************************
229 
230  TYPE ri_fit_type
231  INTEGER, DIMENSION(:, :), POINTER :: bas_ptr
232  REAL(KIND=dp), DIMENSION(:), POINTER :: nvec
233  REAL(KIND=dp), DIMENSION(:), POINTER :: rm1n
234  REAL(KIND=dp), DIMENSION(:, :), POINTER :: tvec
235  REAL(KIND=dp), DIMENSION(:, :), POINTER :: rm1t
236  REAL(KIND=dp), DIMENSION(:, :), POINTER :: avec
237  REAL(KIND=dp), DIMENSION(:, :), POINTER :: fout
238  REAL(KIND=dp) :: ntrm1n
239  REAL(KIND=dp), DIMENSION(2) :: ftrm1n
240  REAL(KIND=dp), DIMENSION(2) :: echarge
241  REAL(KIND=dp), DIMENSION(2) :: lambda
242  END TYPE ri_fit_type
243 
244 ! **************************************************************************************************
245  TYPE wmat_type
246  REAL(KIND=dp), DIMENSION(:, :), POINTER :: mat
247  END TYPE wmat_type
248  TYPE wbas_type
249  REAL(KIND=dp), DIMENSION(:), POINTER :: vec
250  END TYPE wbas_type
251 ! **************************************************************************************************
252  TYPE stat_type
253  REAL(KIND=dp) :: pairs_tt
254  REAL(KIND=dp) :: pairs_sr
255  REAL(KIND=dp) :: pairs_ff
256  REAL(KIND=dp) :: overlap_error
257  REAL(KIND=dp) :: rho_tt
258  REAL(KIND=dp) :: rho_sr
259  REAL(KIND=dp) :: rho_ff
260  REAL(KIND=dp) :: rho_1c
261  REAL(KIND=dp) :: coef_mem
262  REAL(KIND=dp) :: oint_mem
263  REAL(KIND=dp) :: rhos_mem
264  REAL(KIND=dp) :: abai_mem
265  REAL(KIND=dp) :: ppli_mem
266  END TYPE stat_type
267 ! **************************************************************************************************
268 
269  TYPE lri_environment_type
270  ! parameter for (pseudo)inverse of overlap
271  INTEGER :: lri_overlap_inv
272  ! flag for debugging lri integrals
273  LOGICAL :: debug
274  ! flag for shg (solid haromonic Gaussian) integrals
275  LOGICAL :: use_shg_integrals
276  ! parameter for inversion (autoselect); maximal condition
277  ! number up to where inversion is legal
278  REAL(KIND=dp) :: cond_max
279  ! parameter for checking distance between atom pairs
280  REAL(KIND=dp) :: delta
281  ! threshold for aba and abb integrals
282  REAL(KIND=dp) :: eps_o3_int
283  ! orbital basis set
284  TYPE(gto_basis_set_p_type), DIMENSION(:), POINTER :: orb_basis
285  ! lri (fit) basis set
286  TYPE(gto_basis_set_p_type), DIMENSION(:), POINTER :: ri_basis
287  ! orb_basis neighborlist
288  TYPE(neighbor_list_set_p_type), DIMENSION(:), POINTER :: soo_list
289  TYPE(neighbor_list_set_p_type), DIMENSION(:), POINTER :: saa_list
290  TYPE(neighbor_list_set_p_type), DIMENSION(:), POINTER :: soa_list
291  ! local RI integrals
292  TYPE(lri_list_type), POINTER :: lri_ints
293  ! local Vppl integrals
294  TYPE(lri_ppl_int_type), POINTER :: lri_ppl_ints
295  ! local integral of rho**2; for optimization
296  TYPE(lri_list_type), POINTER :: lri_ints_rho
297  ! properties of orb and aux basis
298  TYPE(lri_bas_type), DIMENSION(:), POINTER :: bas_prop
299  ! Clebsch-Gordon for solid harmonics
300  TYPE(lri_clebsch_gordon_type), POINTER :: cg_shg
301  ! orbital basis overlap
302  TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: ob_smat
303  ! statistics
304  LOGICAL :: statistics
305  TYPE(stat_type) :: stat
306  ! exact one-center terms
307  LOGICAL :: exact_1c_terms
308  ! use RI for local pp
309  LOGICAL :: ppl_ri
310  ! store integrals (needed for basis optimization)
311  LOGICAL :: store_integrals
312  ! distant pair approximation
313  LOGICAL :: distant_pair_approximation
314  CHARACTER(len=10) :: distant_pair_method
315  REAL(KIND=dp) :: r_in
316  REAL(KIND=dp) :: r_out
317  REAL(KIND=dp), DIMENSION(:), POINTER :: aradius
318  TYPE(wbas_type), DIMENSION(:), POINTER :: wbas
319  TYPE(wmat_type), DIMENSION(:, :), POINTER :: wmat
320  ! RI overlap and inverse
321  TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: ri_smat, &
322  ri_sinv
323  TYPE(ri_fit_type), POINTER :: ri_fit
324  CHARACTER(len=10) :: ri_sinv_app
325  TYPE(o3c_container_type), POINTER :: o3c
326  END TYPE lri_environment_type
327 
328 ! **************************************************************************************************
329 
330  TYPE lri_kind_type
331  ! expansion coeff for lri density dim(natom,nsgf)
332  REAL(KIND=dp), DIMENSION(:, :), POINTER :: acoef
333  ! integrals V*fbas (potential*fit basis) dim(natom,nsgf)
334  REAL(KIND=dp), DIMENSION(:, :), POINTER :: v_int
335  ! SUM_i integral(V*fbas_i)*davec/dR dim(natom,3)
336  REAL(KIND=dp), DIMENSION(:, :), POINTER :: v_dadr
337  ! integrals V*dfbas/dR
338  REAL(KIND=dp), DIMENSION(:, :), POINTER :: v_dfdr
339  END TYPE lri_kind_type
340 
341  TYPE lri_spin_type
342  TYPE(lri_kind_type), DIMENSION(:), POINTER :: lri_kinds
343  END TYPE lri_spin_type
344 
345 ! **************************************************************************************************
346 
347  TYPE lri_force_type
348  REAL(KIND=dp), DIMENSION(:), POINTER :: st
349  REAL(KIND=dp), DIMENSION(:, :), POINTER :: dssn, &
350  dsst
351  ! derivative dtvec/dR
352  REAL(KIND=dp), DIMENSION(:, :), POINTER :: dtvec
353  END TYPE lri_force_type
354 
355 ! **************************************************************************************************
356 
357  TYPE lri_density_type
358  INTEGER :: nspin = 0
359  ! pair density expansion (nspin)
360  TYPE(lri_list_p_type), DIMENSION(:), POINTER :: lri_rhos => null()
361  ! coefficients of RI expansion and gradients (nspin)
362  TYPE(lri_spin_type), DIMENSION(:), POINTER :: lri_coefs => null()
363  END TYPE lri_density_type
364 
365 ! **************************************************************************************************
366 
367  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'lri_environment_types'
368 
369  PUBLIC :: lri_environment_type, &
370  lri_force_type, lri_list_type, &
371  lri_int_type, lri_int_rho_type, lri_density_type, &
372  lri_kind_type, lri_rhoab_type
373  PUBLIC :: int_container, carray
380 
381 ! **************************************************************************************************
382 
383 CONTAINS
384 
385 ! **************************************************************************************************
386 !> \brief creates and initializes an lri_env
387 !> \param lri_env the lri_environment you want to create
388 ! **************************************************************************************************
389  SUBROUTINE lri_env_create(lri_env)
390 
391  TYPE(lri_environment_type), INTENT(OUT) :: lri_env
392 
393  lri_env%debug = .false.
394  lri_env%delta = 1.e-6_dp
395 
396  lri_env%store_integrals = .false.
397 
398  NULLIFY (lri_env%orb_basis)
399  NULLIFY (lri_env%ri_basis)
400 
401  NULLIFY (lri_env%soo_list)
402  NULLIFY (lri_env%saa_list)
403  NULLIFY (lri_env%soa_list)
404  NULLIFY (lri_env%lri_ints)
405  NULLIFY (lri_env%lri_ppl_ints)
406  NULLIFY (lri_env%lri_ints_rho)
407  NULLIFY (lri_env%bas_prop)
408 
409  NULLIFY (lri_env%ob_smat)
410  NULLIFY (lri_env%ri_smat)
411  NULLIFY (lri_env%ri_sinv)
412  NULLIFY (lri_env%ri_fit)
413  NULLIFY (lri_env%o3c)
414  NULLIFY (lri_env%aradius)
415  NULLIFY (lri_env%wmat)
416  NULLIFY (lri_env%wbas)
417 
418  NULLIFY (lri_env%cg_shg)
419  ALLOCATE (lri_env%cg_shg)
420  NULLIFY (lri_env%cg_shg%cg_coeff)
421  NULLIFY (lri_env%cg_shg%cg_none0_list)
422  NULLIFY (lri_env%cg_shg%ncg_none0)
423 
424  END SUBROUTINE lri_env_create
425 
426 ! **************************************************************************************************
427 !> \brief releases the given lri_env
428 !> \param lri_env the lri environment to release
429 ! **************************************************************************************************
430  SUBROUTINE lri_env_release(lri_env)
431 
432  TYPE(lri_environment_type), INTENT(INOUT) :: lri_env
433 
434  INTEGER :: i, ikind, j, nkind
435 
436  ! deallocate basis sets
437  IF (ASSOCIATED(lri_env%orb_basis)) THEN
438  nkind = SIZE(lri_env%orb_basis)
439  DO ikind = 1, nkind
440  CALL deallocate_gto_basis_set(lri_env%orb_basis(ikind)%gto_basis_set)
441  END DO
442  DEALLOCATE (lri_env%orb_basis)
443  END IF
444  IF (ASSOCIATED(lri_env%ri_basis)) THEN
445  nkind = SIZE(lri_env%ri_basis)
446  DO ikind = 1, nkind
447  CALL deallocate_gto_basis_set(lri_env%ri_basis(ikind)%gto_basis_set)
448  END DO
449  DEALLOCATE (lri_env%ri_basis)
450  END IF
451  CALL release_neighbor_list_sets(lri_env%soo_list)
452  CALL release_neighbor_list_sets(lri_env%saa_list)
453  CALL release_neighbor_list_sets(lri_env%soa_list)
454  IF (ASSOCIATED(lri_env%lri_ints)) THEN
455  CALL deallocate_lri_ints(lri_env%lri_ints)
456  END IF
457  IF (ASSOCIATED(lri_env%lri_ppl_ints)) THEN
458  CALL deallocate_lri_ppl_ints(lri_env%lri_ppl_ints)
459  END IF
460  IF (ASSOCIATED(lri_env%lri_ints_rho)) THEN
461  CALL deallocate_lri_ints_rho(lri_env%lri_ints_rho)
462  END IF
463  CALL deallocate_bas_properties(lri_env)
464  IF (ASSOCIATED(lri_env%aradius)) THEN
465  DEALLOCATE (lri_env%aradius)
466  END IF
467  IF (ASSOCIATED(lri_env%wmat)) THEN
468  DO i = 1, SIZE(lri_env%wmat, 1)
469  DO j = 1, SIZE(lri_env%wmat, 2)
470  IF (ASSOCIATED(lri_env%wmat(i, j)%mat)) THEN
471  DEALLOCATE (lri_env%wmat(i, j)%mat)
472  END IF
473  END DO
474  END DO
475  DEALLOCATE (lri_env%wmat)
476  END IF
477  IF (ASSOCIATED(lri_env%wbas)) THEN
478  DO i = 1, SIZE(lri_env%wbas, 1)
479  IF (ASSOCIATED(lri_env%wbas(i)%vec)) THEN
480  DEALLOCATE (lri_env%wbas(i)%vec)
481  END IF
482  END DO
483  DEALLOCATE (lri_env%wbas)
484  END IF
485  IF (ASSOCIATED(lri_env%cg_shg)) THEN
486  IF (ASSOCIATED(lri_env%cg_shg%cg_coeff)) THEN
487  DEALLOCATE (lri_env%cg_shg%cg_coeff)
488  END IF
489  IF (ASSOCIATED(lri_env%cg_shg%cg_none0_list)) THEN
490  DEALLOCATE (lri_env%cg_shg%cg_none0_list)
491  END IF
492  IF (ASSOCIATED(lri_env%cg_shg%ncg_none0)) THEN
493  DEALLOCATE (lri_env%cg_shg%ncg_none0)
494  END IF
495  DEALLOCATE (lri_env%cg_shg)
496  END IF
497  ! RI
498  IF (ASSOCIATED(lri_env%ob_smat)) CALL dbcsr_deallocate_matrix_set(lri_env%ob_smat)
499  IF (ASSOCIATED(lri_env%ri_smat)) CALL dbcsr_deallocate_matrix_set(lri_env%ri_smat)
500  IF (ASSOCIATED(lri_env%ri_sinv)) CALL dbcsr_deallocate_matrix_set(lri_env%ri_sinv)
501  IF (ASSOCIATED(lri_env%ri_fit)) THEN
502  IF (ASSOCIATED(lri_env%ri_fit%nvec)) THEN
503  DEALLOCATE (lri_env%ri_fit%nvec)
504  END IF
505  IF (ASSOCIATED(lri_env%ri_fit%rm1n)) THEN
506  DEALLOCATE (lri_env%ri_fit%rm1n)
507  END IF
508  IF (ASSOCIATED(lri_env%ri_fit%tvec)) THEN
509  DEALLOCATE (lri_env%ri_fit%tvec)
510  END IF
511  IF (ASSOCIATED(lri_env%ri_fit%rm1t)) THEN
512  DEALLOCATE (lri_env%ri_fit%rm1t)
513  END IF
514  IF (ASSOCIATED(lri_env%ri_fit%avec)) THEN
515  DEALLOCATE (lri_env%ri_fit%avec)
516  END IF
517  IF (ASSOCIATED(lri_env%ri_fit%fout)) THEN
518  DEALLOCATE (lri_env%ri_fit%fout)
519  END IF
520  IF (ASSOCIATED(lri_env%ri_fit%bas_ptr)) THEN
521  DEALLOCATE (lri_env%ri_fit%bas_ptr)
522  END IF
523  DEALLOCATE (lri_env%ri_fit)
524  END IF
525  IF (ASSOCIATED(lri_env%o3c)) THEN
526  CALL release_o3c_container(lri_env%o3c)
527  DEALLOCATE (lri_env%o3c)
528  END IF
529 
530  END SUBROUTINE lri_env_release
531 
532 ! **************************************************************************************************
533 !> \brief creates and initializes an lri_density environment
534 !> \param lri_density the lri_density environment you want to create
535 ! **************************************************************************************************
536  SUBROUTINE lri_density_create(lri_density)
537 
538  TYPE(lri_density_type), INTENT(OUT) :: lri_density
539 
540  lri_density%nspin = 0
541 
542  NULLIFY (lri_density%lri_rhos)
543  NULLIFY (lri_density%lri_coefs)
544 
545  END SUBROUTINE lri_density_create
546 
547 ! **************************************************************************************************
548 !> \brief releases the given lri_density
549 !> \param lri_density the lri_density to release
550 ! **************************************************************************************************
551  SUBROUTINE lri_density_release(lri_density)
552 
553  TYPE(lri_density_type), INTENT(INOUT) :: lri_density
554 
555  CALL deallocate_lri_rhos(lri_density%lri_rhos)
556  CALL deallocate_lri_coefs(lri_density%lri_coefs)
557 
558  END SUBROUTINE lri_density_release
559 
560 ! **************************************************************************************************
561 !> \brief allocate lri_ints, matrices that store LRI integrals
562 !> \param lri_env ...
563 !> \param lri_ints structure storing the LRI integrals
564 !> \param nkind number of atom kinds
565 ! **************************************************************************************************
566  SUBROUTINE allocate_lri_ints(lri_env, lri_ints, nkind)
567 
568  TYPE(lri_environment_type), POINTER :: lri_env
569  TYPE(lri_list_type), POINTER :: lri_ints
570  INTEGER, INTENT(IN) :: nkind
571 
572  INTEGER :: i, iac, iatom, ikind, ilist, jatom, &
573  jkind, jneighbor, nba, nbb, nfa, nfb, &
574  nlist, nn, nneighbor
575  LOGICAL :: dpa, e1c
576  REAL(kind=dp) :: dab, ra, rab(3), rb
577  TYPE(gto_basis_set_type), POINTER :: fbasa, fbasb, obasa, obasb
578  TYPE(lri_int_type), POINTER :: lrii
579  TYPE(neighbor_list_iterator_p_type), &
580  DIMENSION(:), POINTER :: nl_iterator
581 
582  cpassert(ASSOCIATED(lri_env))
583 
584  NULLIFY (fbasa, fbasb, lrii, nl_iterator, obasa, obasb)
585 
586  ALLOCATE (lri_ints)
587 
588  dpa = lri_env%distant_pair_approximation
589  ra = lri_env%r_in
590  rb = lri_env%r_out
591  lri_env%stat%oint_mem = 0.0_dp
592 
593  lri_ints%nkind = nkind
594  ALLOCATE (lri_ints%lri_atom(nkind*nkind))
595 
596  DO i = 1, nkind*nkind
597  NULLIFY (lri_ints%lri_atom(i)%lri_node)
598  lri_ints%lri_atom(i)%natom = 0
599  END DO
600 
601  CALL neighbor_list_iterator_create(nl_iterator, lri_env%soo_list)
602  DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
603 
604  CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, &
605  nlist=nlist, ilist=ilist, nnode=nneighbor, inode=jneighbor, &
606  iatom=iatom, jatom=jatom, r=rab)
607 
608  iac = ikind + nkind*(jkind - 1)
609  dab = sqrt(sum(rab*rab))
610 
611  obasa => lri_env%orb_basis(ikind)%gto_basis_set
612  obasb => lri_env%orb_basis(jkind)%gto_basis_set
613  fbasa => lri_env%ri_basis(ikind)%gto_basis_set
614  fbasb => lri_env%ri_basis(jkind)%gto_basis_set
615 
616  IF (.NOT. ASSOCIATED(obasa)) cycle
617  IF (.NOT. ASSOCIATED(obasb)) cycle
618 
619  IF (.NOT. ASSOCIATED(lri_ints%lri_atom(iac)%lri_node)) THEN
620  lri_ints%lri_atom(iac)%natom = nlist
621  ALLOCATE (lri_ints%lri_atom(iac)%lri_node(nlist))
622  DO i = 1, nlist
623  NULLIFY (lri_ints%lri_atom(iac)%lri_node(i)%lri_int)
624  lri_ints%lri_atom(iac)%lri_node(i)%nnode = 0
625  END DO
626  END IF
627  IF (.NOT. ASSOCIATED(lri_ints%lri_atom(iac)%lri_node(ilist)%lri_int)) THEN
628  lri_ints%lri_atom(iac)%lri_node(ilist)%nnode = nneighbor
629  ALLOCATE (lri_ints%lri_atom(iac)%lri_node(ilist)%lri_int(nneighbor))
630  END IF
631 
632  lrii => lri_ints%lri_atom(iac)%lri_node(ilist)%lri_int(jneighbor)
633 
634  nba = obasa%nsgf
635  nbb = obasb%nsgf
636  nfa = fbasa%nsgf
637  nfb = fbasb%nsgf
638  nn = nfa + nfb
639 
640  IF (iatom == jatom .AND. dab < lri_env%delta) THEN
641  e1c = lri_env%exact_1c_terms
642  ELSE
643  e1c = .false.
644  END IF
645 
646  IF (.NOT. e1c) THEN
647  ALLOCATE (lrii%abascr(nfa))
648  lrii%abascr = 0._dp
649  ALLOCATE (lrii%abbscr(nfb))
650  lrii%abbscr = 0._dp
651  lri_env%stat%oint_mem = lri_env%stat%oint_mem + nfa + nfb
652  END IF
653 
654  IF (dpa) THEN
655  lrii%wsr = pswitch(dab, ra, rb, 0)
656  lrii%wff = 1.0_dp - lrii%wsr
657  lrii%dwsr = pswitch(dab, ra, rb, 1)
658  lrii%dwff = -lrii%dwsr
659  lrii%lrisr = (lrii%wsr > 0.0_dp)
660  lrii%lriff = (lrii%wff > 0.0_dp)
661  NULLIFY (lrii%asinv, lrii%bsinv)
662  ELSE
663  lrii%lrisr = .true.
664  lrii%lriff = .false.
665  lrii%wsr = 1.0_dp
666  lrii%wff = 0.0_dp
667  lrii%dwsr = 0.0_dp
668  lrii%dwff = 0.0_dp
669  NULLIFY (lrii%asinv, lrii%bsinv)
670  END IF
671 
672  ! compressed storage
673  NULLIFY (lrii%cabai%ca, lrii%cabbi%ca)
674 
675  ! full LRI method term
676  IF (lrii%lrisr) THEN
677  IF (e1c) THEN
678  NULLIFY (lrii%n, lrii%sn)
679  NULLIFY (lrii%sinv)
680  ELSE
681  ALLOCATE (lrii%soo(nba, nbb))
682  lri_env%stat%oint_mem = lri_env%stat%oint_mem + nba*nbb
683  lrii%soo = 0._dp
684 
685  IF (iatom == jatom .AND. dab < lri_env%delta) THEN
686  ALLOCATE (lrii%sinv(nfa, nfa))
687  lri_env%stat%oint_mem = lri_env%stat%oint_mem + nfa*nfa
688  ELSE
689  ALLOCATE (lrii%sinv(nn, nn))
690  lri_env%stat%oint_mem = lri_env%stat%oint_mem + nn*nn
691  END IF
692  lrii%sinv = 0._dp
693 
694  IF (iatom == jatom .AND. dab < lri_env%delta) THEN
695  ALLOCATE (lrii%n(nfa), lrii%sn(nfa))
696  lri_env%stat%oint_mem = lri_env%stat%oint_mem + 2.*nfa
697  ELSE
698  ALLOCATE (lrii%n(nn), lrii%sn(nn))
699  lri_env%stat%oint_mem = lri_env%stat%oint_mem + 2.*nn
700  END IF
701  lrii%n = 0._dp
702  lrii%sn = 0._dp
703  END IF
704  ELSE
705  NULLIFY (lrii%n, lrii%sn)
706  NULLIFY (lrii%sinv)
707  END IF
708 
709  ! far field approximation
710  IF (lrii%lriff) THEN
711  lrii%asinv => lri_env%bas_prop(ikind)%ri_ovlp_inv
712  lrii%bsinv => lri_env%bas_prop(jkind)%ri_ovlp_inv
713  ALLOCATE (lrii%na(nfa), lrii%sna(nfa))
714  lri_env%stat%oint_mem = lri_env%stat%oint_mem + 2.*nfa
715  lrii%na = 0._dp
716  lrii%sna = 0._dp
717  ALLOCATE (lrii%nb(nfb), lrii%snb(nfb))
718  lri_env%stat%oint_mem = lri_env%stat%oint_mem + 2.*nfb
719  lrii%nb = 0._dp
720  lrii%snb = 0._dp
721  IF (.NOT. ALLOCATED(lrii%soo)) THEN
722  ALLOCATE (lrii%soo(nba, nbb))
723  lri_env%stat%oint_mem = lri_env%stat%oint_mem + nba*nbb
724  lrii%soo = 0._dp
725  ELSE
726  cpassert(SIZE(lrii%soo, 1) == nba .AND. SIZE(lrii%soo, 2) == nbb)
727  END IF
728  ELSE
729  NULLIFY (lrii%na, lrii%sna)
730  NULLIFY (lrii%nb, lrii%snb)
731  END IF
732 
733  lrii%dmax_ab = 0._dp
734  lrii%dmax_oo = 0._dp
735  lrii%dmax_aba = 0._dp
736  lrii%dmax_abb = 0._dp
737 
738  lrii%calc_force_pair = .false.
739 
740  END DO
741 
742  CALL neighbor_list_iterator_release(nl_iterator)
743 
744  END SUBROUTINE allocate_lri_ints
745 
746 ! **************************************************************************************************
747 !> \brief allocate lri_ppl_ints, matrices that store LRI integrals
748 !> \param lri_env ...
749 !> \param lri_ppl_ints structure storing the LRI ppl integrals
750 !> \param atomic_kind_set ...
751 ! **************************************************************************************************
752  SUBROUTINE allocate_lri_ppl_ints(lri_env, lri_ppl_ints, atomic_kind_set)
753 
754  TYPE(lri_environment_type), POINTER :: lri_env
755  TYPE(lri_ppl_int_type), POINTER :: lri_ppl_ints
756  TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
757 
758  INTEGER :: ikind, natom, nfa, nkind
759  TYPE(atomic_kind_type), POINTER :: atomic_kind
760  TYPE(gto_basis_set_type), POINTER :: fbasa
761 
762  cpassert(ASSOCIATED(lri_env))
763 
764  lri_env%stat%ppli_mem = 0.0_dp
765  nkind = SIZE(atomic_kind_set)
766  ALLOCATE (lri_ppl_ints)
767  ALLOCATE (lri_ppl_ints%lri_ppl(nkind))
768  DO ikind = 1, nkind
769  fbasa => lri_env%ri_basis(ikind)%gto_basis_set
770  nfa = fbasa%nsgf
771  atomic_kind => atomic_kind_set(ikind)
772  CALL get_atomic_kind(atomic_kind=atomic_kind, natom=natom)
773  ALLOCATE (lri_ppl_ints%lri_ppl(ikind)%v_int(natom, nfa))
774  lri_env%stat%ppli_mem = lri_env%stat%ppli_mem + natom*nfa
775  END DO
776 
777  END SUBROUTINE allocate_lri_ppl_ints
778 
779 ! **************************************************************************************************
780 !> \brief allocate lri_ints_rho, storing integral for the exact density
781 !> \param lri_env ...
782 !> \param lri_ints_rho structure storing the integrals (aa,bb)
783 !> \param nkind number of atom kinds
784 ! **************************************************************************************************
785  SUBROUTINE allocate_lri_ints_rho(lri_env, lri_ints_rho, nkind)
786 
787  TYPE(lri_environment_type), POINTER :: lri_env
788  TYPE(lri_list_type), POINTER :: lri_ints_rho
789  INTEGER, INTENT(IN) :: nkind
790 
791  INTEGER :: i, iac, iatom, ikind, ilist, jatom, &
792  jkind, jneighbor, nba, nbb, nlist, &
793  nneighbor
794  TYPE(gto_basis_set_type), POINTER :: obasa, obasb
795  TYPE(lri_int_rho_type), POINTER :: lriir
796  TYPE(neighbor_list_iterator_p_type), &
797  DIMENSION(:), POINTER :: nl_iterator
798 
799  cpassert(ASSOCIATED(lri_env))
800 
801  ALLOCATE (lri_ints_rho)
802 
803  lri_ints_rho%nkind = nkind
804  ALLOCATE (lri_ints_rho%lri_atom(nkind*nkind))
805 
806  DO i = 1, nkind*nkind
807  NULLIFY (lri_ints_rho%lri_atom(i)%lri_node)
808  lri_ints_rho%lri_atom(i)%natom = 0
809  END DO
810 
811  CALL neighbor_list_iterator_create(nl_iterator, lri_env%soo_list)
812  DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
813 
814  CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, &
815  nlist=nlist, ilist=ilist, nnode=nneighbor, inode=jneighbor, &
816  iatom=iatom, jatom=jatom)
817 
818  iac = ikind + nkind*(jkind - 1)
819 
820  obasa => lri_env%orb_basis(ikind)%gto_basis_set
821  obasb => lri_env%orb_basis(jkind)%gto_basis_set
822 
823  IF (.NOT. ASSOCIATED(obasa)) cycle
824  IF (.NOT. ASSOCIATED(obasb)) cycle
825 
826  IF (.NOT. ASSOCIATED(lri_ints_rho%lri_atom(iac)%lri_node)) THEN
827  lri_ints_rho%lri_atom(iac)%natom = nlist
828  ALLOCATE (lri_ints_rho%lri_atom(iac)%lri_node(nlist))
829  DO i = 1, nlist
830  NULLIFY (lri_ints_rho%lri_atom(iac)%lri_node(i)%lri_int_rho)
831  lri_ints_rho%lri_atom(iac)%lri_node(i)%nnode = 0
832  END DO
833  END IF
834  IF (.NOT. ASSOCIATED(lri_ints_rho%lri_atom(iac)%lri_node(ilist)%lri_int_rho)) THEN
835  lri_ints_rho%lri_atom(iac)%lri_node(ilist)%nnode = nneighbor
836  ALLOCATE (lri_ints_rho%lri_atom(iac)%lri_node(ilist)%lri_int_rho(nneighbor))
837  END IF
838 
839  lriir => lri_ints_rho%lri_atom(iac)%lri_node(ilist)%lri_int_rho(jneighbor)
840 
841  nba = obasa%nsgf
842  nbb = obasb%nsgf
843 
844  ALLOCATE (lriir%soaabb(nba, nba, nbb, nbb))
845  lriir%soaabb = 0._dp
846  lriir%dmax_aabb = 0._dp
847 
848  END DO
849 
850  CALL neighbor_list_iterator_release(nl_iterator)
851 
852  END SUBROUTINE allocate_lri_ints_rho
853 
854 ! **************************************************************************************************
855 !> \brief creates and initializes lri_rhos
856 !> \param lri_env ...
857 !> \param lri_rhos structure storing tvec and avec
858 !> \param nspin ...
859 !> \param nkind number of atom kinds
860 ! **************************************************************************************************
861  SUBROUTINE allocate_lri_rhos(lri_env, lri_rhos, nspin, nkind)
862 
863  TYPE(lri_environment_type), POINTER :: lri_env
864  TYPE(lri_list_p_type), DIMENSION(:), POINTER :: lri_rhos
865  INTEGER, INTENT(IN) :: nspin, nkind
866 
867  INTEGER :: i, iac, iatom, ikind, ilist, ispin, &
868  jatom, jkind, jneighbor, nfa, nfb, &
869  nlist, nn, nneighbor
870  REAL(kind=dp) :: dab, rab(3)
871  TYPE(lri_int_type), POINTER :: lrii
872  TYPE(lri_list_type), POINTER :: lri_rho
873  TYPE(lri_rhoab_type), POINTER :: lrho
874  TYPE(neighbor_list_iterator_p_type), &
875  DIMENSION(:), POINTER :: nl_iterator
876 
877  cpassert(ASSOCIATED(lri_env))
878 
879  NULLIFY (lri_rho, lrho, lrii, nl_iterator)
880 
881  ALLOCATE (lri_rhos(nspin))
882 
883  lri_env%stat%rhos_mem = 0.0_dp
884 
885  DO ispin = 1, nspin
886 
887  ALLOCATE (lri_rhos(ispin)%lri_list)
888 
889  lri_rhos(ispin)%lri_list%nkind = nkind
890  ALLOCATE (lri_rhos(ispin)%lri_list%lri_atom(nkind*nkind))
891 
892  DO i = 1, nkind*nkind
893  NULLIFY (lri_rhos(ispin)%lri_list%lri_atom(i)%lri_node)
894  lri_rhos(ispin)%lri_list%lri_atom(i)%natom = 0
895  END DO
896 
897  lri_rho => lri_rhos(ispin)%lri_list
898 
899  CALL neighbor_list_iterator_create(nl_iterator, lri_env%soo_list)
900  DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
901  CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, &
902  iatom=iatom, jatom=jatom, nlist=nlist, ilist=ilist, &
903  nnode=nneighbor, inode=jneighbor, r=rab)
904 
905  iac = ikind + nkind*(jkind - 1)
906  dab = sqrt(sum(rab*rab))
907 
908  IF (.NOT. ASSOCIATED(lri_env%lri_ints%lri_atom(iac)%lri_node)) cycle
909 
910  IF (.NOT. ASSOCIATED(lri_rho%lri_atom(iac)%lri_node)) THEN
911  lri_rho%lri_atom(iac)%natom = nlist
912  ALLOCATE (lri_rho%lri_atom(iac)%lri_node(nlist))
913  DO i = 1, nlist
914  NULLIFY (lri_rho%lri_atom(iac)%lri_node(i)%lri_rhoab)
915  lri_rho%lri_atom(iac)%lri_node(i)%nnode = 0
916  END DO
917  END IF
918  IF (.NOT. ASSOCIATED(lri_rho%lri_atom(iac)%lri_node(ilist)%lri_rhoab)) THEN
919  lri_rho%lri_atom(iac)%lri_node(ilist)%nnode = nneighbor
920  ALLOCATE (lri_rho%lri_atom(iac)%lri_node(ilist)%lri_rhoab(nneighbor))
921  END IF
922 
923  lrho => lri_rho%lri_atom(iac)%lri_node(ilist)%lri_rhoab(jneighbor)
924  lrii => lri_env%lri_ints%lri_atom(iac)%lri_node(ilist)%lri_int(jneighbor)
925 
926  lrho%nba = lrii%nba
927  lrho%nbb = lrii%nbb
928  lrho%nfa = lrii%nfa
929  lrho%nfb = lrii%nfb
930 
931  nfa = lrho%nfa
932  nfb = lrho%nfb
933  nn = nfa + nfb
934 
935  NULLIFY (lrho%avec, lrho%tvec)
936  IF (lrii%lrisr) THEN
937  IF (iatom == jatom .AND. dab < lri_env%delta) THEN
938  IF (.NOT. lri_env%exact_1c_terms) THEN
939  ALLOCATE (lrho%avec(nfa))
940  ALLOCATE (lrho%tvec(nfa))
941  lri_env%stat%rhos_mem = lri_env%stat%rhos_mem + 2*nfa
942  END IF
943  ELSE
944  ALLOCATE (lrho%avec(nn))
945  ALLOCATE (lrho%tvec(nn))
946  lri_env%stat%rhos_mem = lri_env%stat%rhos_mem + 2*nn
947  END IF
948  END IF
949  NULLIFY (lrho%aveca, lrho%tveca)
950  NULLIFY (lrho%avecb, lrho%tvecb)
951  IF (lrii%lriff) THEN
952  ALLOCATE (lrho%aveca(nfa))
953  ALLOCATE (lrho%avecb(nfb))
954  ALLOCATE (lrho%tveca(nfa))
955  ALLOCATE (lrho%tvecb(nfb))
956  lri_env%stat%rhos_mem = lri_env%stat%rhos_mem + 2*(nfa + nfb)
957  END IF
958 
959  END DO
960 
961  CALL neighbor_list_iterator_release(nl_iterator)
962 
963  END DO
964 
965  END SUBROUTINE allocate_lri_rhos
966 
967 ! **************************************************************************************************
968 !> \brief creates and initializes lri_coefs
969 !> \param lri_env ...
970 !> \param lri_density ...
971 !> \param atomic_kind_set ...
972 ! **************************************************************************************************
973  SUBROUTINE allocate_lri_coefs(lri_env, lri_density, atomic_kind_set)
974 
975  TYPE(lri_environment_type), POINTER :: lri_env
976  TYPE(lri_density_type), POINTER :: lri_density
977  TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
978 
979  INTEGER :: ikind, ispin, natom, nkind, nsgf, nspin
980  TYPE(atomic_kind_type), POINTER :: atomic_kind
981  TYPE(gto_basis_set_type), POINTER :: fbas
982  TYPE(lri_spin_type), DIMENSION(:), POINTER :: lri_coefs
983 
984  cpassert(ASSOCIATED(lri_density))
985 
986  NULLIFY (atomic_kind, fbas, lri_coefs)
987 
988  nkind = SIZE(atomic_kind_set)
989  nspin = lri_density%nspin
990 
991  lri_env%stat%coef_mem = 0.0_dp
992 
993  ALLOCATE (lri_density%lri_coefs(nspin))
994  lri_coefs => lri_density%lri_coefs
995 
996  DO ispin = 1, nspin
997  ALLOCATE (lri_coefs(ispin)%lri_kinds(nkind))
998  DO ikind = 1, nkind
999  NULLIFY (lri_coefs(ispin)%lri_kinds(ikind)%acoef)
1000  NULLIFY (lri_coefs(ispin)%lri_kinds(ikind)%v_int)
1001  NULLIFY (lri_coefs(ispin)%lri_kinds(ikind)%v_dadr)
1002  NULLIFY (lri_coefs(ispin)%lri_kinds(ikind)%v_dfdr)
1003  atomic_kind => atomic_kind_set(ikind)
1004  CALL get_atomic_kind(atomic_kind=atomic_kind, natom=natom)
1005  fbas => lri_env%ri_basis(ikind)%gto_basis_set
1006  nsgf = fbas%nsgf
1007  ALLOCATE (lri_coefs(ispin)%lri_kinds(ikind)%acoef(natom, nsgf))
1008  lri_coefs(ispin)%lri_kinds(ikind)%acoef = 0._dp
1009  ALLOCATE (lri_coefs(ispin)%lri_kinds(ikind)%v_int(natom, nsgf))
1010  lri_coefs(ispin)%lri_kinds(ikind)%v_int = 0._dp
1011  ALLOCATE (lri_coefs(ispin)%lri_kinds(ikind)%v_dadr(natom, 3))
1012  lri_coefs(ispin)%lri_kinds(ikind)%v_dadr = 0._dp
1013  ALLOCATE (lri_coefs(ispin)%lri_kinds(ikind)%v_dfdr(natom, 3))
1014  lri_coefs(ispin)%lri_kinds(ikind)%v_dfdr = 0._dp
1015  !
1016  lri_env%stat%coef_mem = lri_env%stat%coef_mem + 2._dp*natom*(nsgf + 3)
1017  END DO
1018  END DO
1019 
1020  END SUBROUTINE allocate_lri_coefs
1021 
1022 ! **************************************************************************************************
1023 !> \brief creates and initializes lri_force
1024 !> \param lri_force ...
1025 !> \param nfa and nfb number of fit functions on a/b
1026 !> \param nfb ...
1027 ! **************************************************************************************************
1028  SUBROUTINE allocate_lri_force_components(lri_force, nfa, nfb)
1029 
1030  TYPE(lri_force_type), POINTER :: lri_force
1031  INTEGER, INTENT(IN) :: nfa, nfb
1032 
1033  INTEGER :: nn
1034 
1035  nn = nfa + nfb
1036 
1037  cpassert(.NOT. ASSOCIATED(lri_force))
1038 
1039  ALLOCATE (lri_force)
1040 
1041  ALLOCATE (lri_force%st(nn))
1042  lri_force%st = 0._dp
1043  ALLOCATE (lri_force%dsst(nn, 3))
1044  lri_force%dsst = 0._dp
1045  ALLOCATE (lri_force%dssn(nn, 3))
1046  lri_force%dssn = 0._dp
1047  ALLOCATE (lri_force%dtvec(nn, 3))
1048  lri_force%dtvec = 0._dp
1049 
1050  END SUBROUTINE allocate_lri_force_components
1051 
1052 ! **************************************************************************************************
1053 !> \brief deallocates one-center overlap integrals, integral of ri basis
1054 !> and scon matrices
1055 !> \param lri_env ...
1056 ! **************************************************************************************************
1057  SUBROUTINE deallocate_bas_properties(lri_env)
1058 
1059  TYPE(lri_environment_type), INTENT(INOUT) :: lri_env
1060 
1061  INTEGER :: i
1062 
1063  IF (ASSOCIATED(lri_env%bas_prop)) THEN
1064  DO i = 1, SIZE(lri_env%bas_prop)
1065  IF (ASSOCIATED(lri_env%bas_prop(i)%int_fbas)) THEN
1066  DEALLOCATE (lri_env%bas_prop(i)%int_fbas)
1067  END IF
1068  IF (ASSOCIATED(lri_env%bas_prop(i)%ri_ovlp)) THEN
1069  DEALLOCATE (lri_env%bas_prop(i)%ri_ovlp)
1070  END IF
1071  IF (ASSOCIATED(lri_env%bas_prop(i)%ri_ovlp_inv)) THEN
1072  DEALLOCATE (lri_env%bas_prop(i)%ri_ovlp_inv)
1073  END IF
1074  IF (ASSOCIATED(lri_env%bas_prop(i)%orb_ovlp)) THEN
1075  DEALLOCATE (lri_env%bas_prop(i)%orb_ovlp)
1076  END IF
1077  IF (ALLOCATED(lri_env%bas_prop(i)%ovlp3)) THEN
1078  DEALLOCATE (lri_env%bas_prop(i)%ovlp3)
1079  END IF
1080  IF (ASSOCIATED(lri_env%bas_prop(i)%scon_ri)) THEN
1081  DEALLOCATE (lri_env%bas_prop(i)%scon_ri)
1082  END IF
1083  IF (ASSOCIATED(lri_env%bas_prop(i)%scon_orb)) THEN
1084  DEALLOCATE (lri_env%bas_prop(i)%scon_orb)
1085  END IF
1086  IF (ASSOCIATED(lri_env%bas_prop(i)%orb_index)) THEN
1087  DEALLOCATE (lri_env%bas_prop(i)%orb_index)
1088  END IF
1089  IF (ASSOCIATED(lri_env%bas_prop(i)%ri_index)) THEN
1090  DEALLOCATE (lri_env%bas_prop(i)%ri_index)
1091  END IF
1092  IF (ASSOCIATED(lri_env%bas_prop(i)%scon_mix)) THEN
1093  DEALLOCATE (lri_env%bas_prop(i)%scon_mix)
1094  END IF
1095  END DO
1096  DEALLOCATE (lri_env%bas_prop)
1097  END IF
1098 
1099  END SUBROUTINE deallocate_bas_properties
1100 
1101 ! **************************************************************************************************
1102 !> \brief deallocates the given lri_ints
1103 !> \param lri_ints ...
1104 ! **************************************************************************************************
1105  SUBROUTINE deallocate_lri_ints(lri_ints)
1106 
1107  TYPE(lri_list_type), POINTER :: lri_ints
1108 
1109  INTEGER :: i, iatom, ijkind, inode, natom, nkind, &
1110  nnode
1111  TYPE(lri_int_type), POINTER :: lri_int
1112 
1113  cpassert(ASSOCIATED(lri_ints))
1114  nkind = lri_ints%nkind
1115 
1116  IF (nkind > 0) THEN
1117  DO ijkind = 1, SIZE(lri_ints%lri_atom)
1118  natom = lri_ints%lri_atom(ijkind)%natom
1119  IF (natom > 0) THEN
1120  DO iatom = 1, natom
1121  nnode = lri_ints%lri_atom(ijkind)%lri_node(iatom)%nnode
1122  IF (nnode > 0) THEN
1123  IF (ASSOCIATED(lri_ints%lri_atom(ijkind)%lri_node(iatom)%lri_int)) THEN
1124  DO inode = 1, nnode
1125  lri_int => lri_ints%lri_atom(ijkind)%lri_node(iatom)%lri_int(inode)
1126  IF (ALLOCATED(lri_int%sab)) THEN
1127  DEALLOCATE (lri_int%sab)
1128  END IF
1129  IF (ALLOCATED(lri_int%soo)) THEN
1130  DEALLOCATE (lri_int%soo)
1131  END IF
1132  IF (ALLOCATED(lri_int%abaint)) THEN
1133  DEALLOCATE (lri_int%abaint)
1134  END IF
1135  IF (ALLOCATED(lri_int%abascr)) THEN
1136  DEALLOCATE (lri_int%abascr)
1137  END IF
1138  IF (ALLOCATED(lri_int%abbint)) THEN
1139  DEALLOCATE (lri_int%abbint)
1140  END IF
1141  IF (ALLOCATED(lri_int%abbscr)) THEN
1142  DEALLOCATE (lri_int%abbscr)
1143  END IF
1144  IF (ALLOCATED(lri_int%dsab)) THEN
1145  DEALLOCATE (lri_int%dsab)
1146  END IF
1147  IF (ALLOCATED(lri_int%dsoo)) THEN
1148  DEALLOCATE (lri_int%dsoo)
1149  END IF
1150  IF (ALLOCATED(lri_int%dabbint)) THEN
1151  DEALLOCATE (lri_int%dabbint)
1152  END IF
1153  IF (ALLOCATED(lri_int%dabdaint)) THEN
1154  DEALLOCATE (lri_int%dabdaint)
1155  END IF
1156  IF (ASSOCIATED(lri_int%sinv)) THEN
1157  DEALLOCATE (lri_int%sinv)
1158  END IF
1159  IF (ASSOCIATED(lri_int%n)) THEN
1160  DEALLOCATE (lri_int%n)
1161  END IF
1162  IF (ASSOCIATED(lri_int%sn)) THEN
1163  DEALLOCATE (lri_int%sn)
1164  END IF
1165  IF (ASSOCIATED(lri_int%na)) THEN
1166  DEALLOCATE (lri_int%na)
1167  END IF
1168  IF (ASSOCIATED(lri_int%nb)) THEN
1169  DEALLOCATE (lri_int%nb)
1170  END IF
1171  IF (ASSOCIATED(lri_int%sna)) THEN
1172  DEALLOCATE (lri_int%sna)
1173  END IF
1174  IF (ASSOCIATED(lri_int%snb)) THEN
1175  DEALLOCATE (lri_int%snb)
1176  END IF
1177  ! integral container
1178  IF (ASSOCIATED(lri_int%cabai%ca)) THEN
1179  DO i = 1, SIZE(lri_int%cabai%ca)
1180  IF (ASSOCIATED(lri_int%cabai%ca(i)%cdp)) DEALLOCATE (lri_int%cabai%ca(i)%cdp)
1181  IF (ASSOCIATED(lri_int%cabai%ca(i)%csp)) DEALLOCATE (lri_int%cabai%ca(i)%csp)
1182  IF (ASSOCIATED(lri_int%cabai%ca(i)%cip)) DEALLOCATE (lri_int%cabai%ca(i)%cip)
1183  END DO
1184  DEALLOCATE (lri_int%cabai%ca)
1185  END IF
1186  IF (ASSOCIATED(lri_int%cabbi%ca)) THEN
1187  DO i = 1, SIZE(lri_int%cabbi%ca)
1188  IF (ASSOCIATED(lri_int%cabbi%ca(i)%cdp)) DEALLOCATE (lri_int%cabbi%ca(i)%cdp)
1189  IF (ASSOCIATED(lri_int%cabbi%ca(i)%csp)) DEALLOCATE (lri_int%cabbi%ca(i)%csp)
1190  IF (ASSOCIATED(lri_int%cabbi%ca(i)%cip)) DEALLOCATE (lri_int%cabbi%ca(i)%cip)
1191  END DO
1192  DEALLOCATE (lri_int%cabbi%ca)
1193  END IF
1194  END DO
1195  DEALLOCATE (lri_ints%lri_atom(ijkind)%lri_node(iatom)%lri_int)
1196  END IF
1197  END IF
1198  END DO
1199  DEALLOCATE (lri_ints%lri_atom(ijkind)%lri_node)
1200  END IF
1201  END DO
1202  DEALLOCATE (lri_ints%lri_atom)
1203  END IF
1204  DEALLOCATE (lri_ints)
1205 
1206  END SUBROUTINE deallocate_lri_ints
1207 
1208 ! **************************************************************************************************
1209 !> \brief deallocates the given lri_ppl_ints
1210 !> \param lri_ppl_ints ...
1211 ! **************************************************************************************************
1212  SUBROUTINE deallocate_lri_ppl_ints(lri_ppl_ints)
1213 
1214  TYPE(lri_ppl_int_type), POINTER :: lri_ppl_ints
1215 
1216  INTEGER :: ikind, nkind
1217 
1218  cpassert(ASSOCIATED(lri_ppl_ints))
1219  IF (ASSOCIATED(lri_ppl_ints%lri_ppl)) THEN
1220  nkind = SIZE(lri_ppl_ints%lri_ppl)
1221  DO ikind = 1, nkind
1222  IF (ASSOCIATED(lri_ppl_ints%lri_ppl(ikind)%v_int)) THEN
1223  DEALLOCATE (lri_ppl_ints%lri_ppl(ikind)%v_int)
1224  END IF
1225  END DO
1226  DEALLOCATE (lri_ppl_ints%lri_ppl)
1227  END IF
1228  DEALLOCATE (lri_ppl_ints)
1229 
1230  END SUBROUTINE deallocate_lri_ppl_ints
1231 
1232 ! **************************************************************************************************
1233 !> \brief deallocates the given lri_ints_rho
1234 !> \param lri_ints_rho ...
1235 ! **************************************************************************************************
1236  SUBROUTINE deallocate_lri_ints_rho(lri_ints_rho)
1237 
1238  TYPE(lri_list_type), POINTER :: lri_ints_rho
1239 
1240  INTEGER :: iatom, ijkind, inode, natom, nkind, nnode
1241 
1242  cpassert(ASSOCIATED(lri_ints_rho))
1243  nkind = lri_ints_rho%nkind
1244 
1245  IF (nkind > 0) THEN
1246  DO ijkind = 1, SIZE(lri_ints_rho%lri_atom)
1247  natom = lri_ints_rho%lri_atom(ijkind)%natom
1248  IF (natom > 0) THEN
1249  DO iatom = 1, natom
1250  nnode = lri_ints_rho%lri_atom(ijkind)%lri_node(iatom)%nnode
1251  IF (nnode > 0) THEN
1252  IF (ASSOCIATED(lri_ints_rho%lri_atom(ijkind)%lri_node(iatom)%lri_int_rho)) THEN
1253  DO inode = 1, nnode
1254  IF (ASSOCIATED(lri_ints_rho%lri_atom(ijkind)%lri_node(iatom)%lri_int_rho(inode)%soaabb)) THEN
1255  DEALLOCATE (lri_ints_rho%lri_atom(ijkind)%lri_node(iatom)%lri_int_rho(inode)%soaabb)
1256  END IF
1257  END DO
1258  DEALLOCATE (lri_ints_rho%lri_atom(ijkind)%lri_node(iatom)%lri_int_rho)
1259  END IF
1260  END IF
1261  END DO
1262  DEALLOCATE (lri_ints_rho%lri_atom(ijkind)%lri_node)
1263  END IF
1264  END DO
1265  DEALLOCATE (lri_ints_rho%lri_atom)
1266  END IF
1267  DEALLOCATE (lri_ints_rho)
1268 
1269  END SUBROUTINE deallocate_lri_ints_rho
1270 
1271 ! **************************************************************************************************
1272 !> \brief deallocates the given lri_rhos
1273 !> \param lri_rhos ...
1274 ! **************************************************************************************************
1275  SUBROUTINE deallocate_lri_rhos(lri_rhos)
1276 
1277  TYPE(lri_list_p_type), DIMENSION(:), POINTER :: lri_rhos
1278 
1279  INTEGER :: i, iatom, ijkind, inode, natom, nkind, &
1280  nnode
1281  TYPE(lri_list_type), POINTER :: lri_rho
1282  TYPE(lri_rhoab_type), POINTER :: lri_rhoab
1283 
1284  IF (ASSOCIATED(lri_rhos)) THEN
1285  DO i = 1, SIZE(lri_rhos)
1286  lri_rho => lri_rhos(i)%lri_list
1287  cpassert(ASSOCIATED(lri_rho))
1288  nkind = lri_rho%nkind
1289  IF (nkind > 0) THEN
1290  cpassert(ASSOCIATED(lri_rho%lri_atom))
1291  DO ijkind = 1, SIZE(lri_rho%lri_atom)
1292  natom = lri_rho%lri_atom(ijkind)%natom
1293  IF (natom > 0) THEN
1294  cpassert(ASSOCIATED(lri_rho%lri_atom(ijkind)%lri_node))
1295  DO iatom = 1, natom
1296  nnode = lri_rho%lri_atom(ijkind)%lri_node(iatom)%nnode
1297  IF (nnode > 0) THEN
1298  IF (ASSOCIATED(lri_rho%lri_atom(ijkind)%lri_node(iatom)%lri_rhoab)) THEN
1299  DO inode = 1, nnode
1300  lri_rhoab => lri_rho%lri_atom(ijkind)%lri_node(iatom)%lri_rhoab(inode)
1301  IF (ASSOCIATED(lri_rhoab%avec)) DEALLOCATE (lri_rhoab%avec)
1302  IF (ASSOCIATED(lri_rhoab%tvec)) DEALLOCATE (lri_rhoab%tvec)
1303  IF (ASSOCIATED(lri_rhoab%aveca)) DEALLOCATE (lri_rhoab%aveca)
1304  IF (ASSOCIATED(lri_rhoab%tveca)) DEALLOCATE (lri_rhoab%tveca)
1305  IF (ASSOCIATED(lri_rhoab%avecb)) DEALLOCATE (lri_rhoab%avecb)
1306  IF (ASSOCIATED(lri_rhoab%tvecb)) DEALLOCATE (lri_rhoab%tvecb)
1307  END DO
1308  DEALLOCATE (lri_rho%lri_atom(ijkind)%lri_node(iatom)%lri_rhoab)
1309  END IF
1310  END IF
1311  END DO
1312  DEALLOCATE (lri_rho%lri_atom(ijkind)%lri_node)
1313  END IF
1314  END DO
1315  DEALLOCATE (lri_rho%lri_atom)
1316  END IF
1317  DEALLOCATE (lri_rho)
1318  END DO
1319 
1320  DEALLOCATE (lri_rhos)
1321 
1322  END IF
1323  NULLIFY (lri_rhos)
1324 
1325  END SUBROUTINE deallocate_lri_rhos
1326 
1327 ! **************************************************************************************************
1328 !> \brief releases the given lri_coefs
1329 !> \param lri_coefs the integral storage environment that is released
1330 ! **************************************************************************************************
1331  SUBROUTINE deallocate_lri_coefs(lri_coefs)
1332  TYPE(lri_spin_type), DIMENSION(:), POINTER :: lri_coefs
1333 
1334  INTEGER :: i, j
1335 
1336  IF (ASSOCIATED(lri_coefs)) THEN
1337  DO i = 1, SIZE(lri_coefs)
1338  DO j = 1, SIZE(lri_coefs(i)%lri_kinds)
1339  IF (ASSOCIATED(lri_coefs(i)%lri_kinds(j)%acoef)) THEN
1340  DEALLOCATE (lri_coefs(i)%lri_kinds(j)%acoef)
1341  END IF
1342  IF (ASSOCIATED(lri_coefs(i)%lri_kinds(j)%v_int)) THEN
1343  DEALLOCATE (lri_coefs(i)%lri_kinds(j)%v_int)
1344  END IF
1345  IF (ASSOCIATED(lri_coefs(i)%lri_kinds(j)%v_dadr)) THEN
1346  DEALLOCATE (lri_coefs(i)%lri_kinds(j)%v_dadr)
1347  END IF
1348  IF (ASSOCIATED(lri_coefs(i)%lri_kinds(j)%v_dfdr)) THEN
1349  DEALLOCATE (lri_coefs(i)%lri_kinds(j)%v_dfdr)
1350  END IF
1351  END DO
1352  DEALLOCATE (lri_coefs(i)%lri_kinds)
1353  END DO
1354  DEALLOCATE (lri_coefs)
1355  END IF
1356 
1357  NULLIFY (lri_coefs)
1358 
1359  END SUBROUTINE deallocate_lri_coefs
1360 
1361 ! **************************************************************************************************
1362 !> \brief releases the given lri_force_type
1363 !> \param lri_force the integral storage environment that is released
1364 ! **************************************************************************************************
1365  SUBROUTINE deallocate_lri_force_components(lri_force)
1366 
1367  TYPE(lri_force_type), POINTER :: lri_force
1368 
1369  IF (ASSOCIATED(lri_force)) THEN
1370  IF (ASSOCIATED(lri_force%st)) THEN
1371  DEALLOCATE (lri_force%st)
1372  END IF
1373  IF (ASSOCIATED(lri_force%dssn)) THEN
1374  DEALLOCATE (lri_force%dssn)
1375  END IF
1376  IF (ASSOCIATED(lri_force%dsst)) THEN
1377  DEALLOCATE (lri_force%dsst)
1378  END IF
1379  IF (ASSOCIATED(lri_force%dtvec)) THEN
1380  DEALLOCATE (lri_force%dtvec)
1381  END IF
1382  DEALLOCATE (lri_force)
1383  END IF
1384 
1385  NULLIFY (lri_force)
1386 
1387  END SUBROUTINE deallocate_lri_force_components
1388 
1389 END MODULE lri_environment_types
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
subroutine, public deallocate_gto_basis_set(gto_basis_set)
...
DBCSR operations in CP2K.
Defines the basic variable types.
Definition: kinds.F:23
integer, parameter, public int_8
Definition: kinds.F:54
integer, parameter, public dp
Definition: kinds.F:34
integer, parameter, public sp
Definition: kinds.F:33
contains the types and subroutines for dealing with the lri_env lri : local resolution of the identit...
subroutine, public deallocate_bas_properties(lri_env)
deallocates one-center overlap integrals, integral of ri basis and scon matrices
subroutine, public allocate_lri_ints_rho(lri_env, lri_ints_rho, nkind)
allocate lri_ints_rho, storing integral for the exact density
subroutine, public lri_env_create(lri_env)
creates and initializes an lri_env
subroutine, public deallocate_lri_ppl_ints(lri_ppl_ints)
deallocates the given lri_ppl_ints
subroutine, public allocate_lri_ppl_ints(lri_env, lri_ppl_ints, atomic_kind_set)
allocate lri_ppl_ints, matrices that store LRI integrals
subroutine, public deallocate_lri_force_components(lri_force)
releases the given lri_force_type
subroutine, public allocate_lri_coefs(lri_env, lri_density, atomic_kind_set)
creates and initializes lri_coefs
subroutine, public allocate_lri_force_components(lri_force, nfa, nfb)
creates and initializes lri_force
subroutine, public lri_density_release(lri_density)
releases the given lri_density
subroutine, public allocate_lri_ints(lri_env, lri_ints, nkind)
allocate lri_ints, matrices that store LRI integrals
subroutine, public deallocate_lri_ints(lri_ints)
deallocates the given lri_ints
subroutine, public deallocate_lri_ints_rho(lri_ints_rho)
deallocates the given lri_ints_rho
subroutine, public lri_env_release(lri_env)
releases the given lri_env
subroutine, public lri_density_create(lri_density)
creates and initializes an lri_density environment
subroutine, public allocate_lri_rhos(lri_env, lri_rhos, nspin, nkind)
creates and initializes lri_rhos
Collection of simple mathematical functions and subroutines.
Definition: mathlib.F:15
real(kind=dp) function, public pswitch(x, a, b, order)
Polynomial (5th degree) switching function f(a) = 1 .... f(b) = 0 with f'(a) = f"(a) = f'(b) = f"(b) ...
Definition: mathlib.F:99
Define the neighbor list data types and the corresponding functionality.
subroutine, public release_neighbor_list_sets(nlists)
releases an array of neighbor_list_sets
subroutine, public neighbor_list_iterator_create(iterator_set, nl, search, nthread)
Neighbor list iterator functions.
subroutine, public neighbor_list_iterator_release(iterator_set)
...
integer function, public neighbor_list_iterate(iterator_set, mepos)
...
subroutine, public get_iterator_info(iterator_set, mepos, ikind, jkind, nkind, ilist, nlist, inode, nnode, iatom, jatom, r, cell)
...
3-center overlap type integrals containers
Definition: qs_o3c_types.F:12
subroutine, public release_o3c_container(o3c_container)
...
Definition: qs_o3c_types.F:225