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