(git:374b731)
Loading...
Searching...
No Matches
lri_optimize_ri_basis.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 Optimizes exponents and contraction coefficients of the lri auxiliary
10!> basis sets using the UOBYQA minimizer
11!> lri : local resolution of the identity
12!> \par History
13!> created Dorothea Golze [05.2014]
14!> \authors Dorothea Golze
15! **************************************************************************************************
17
22 USE cell_types, ONLY: cell_type
25 USE cp_output_handling, ONLY: cp_p_file,&
30 USE dbcsr_api, ONLY: dbcsr_get_block_p,&
31 dbcsr_p_type,&
32 dbcsr_type
41 USE kinds, ONLY: default_path_length,&
42 dp
62 USE powell, ONLY: opt_state_type,&
73 USE qs_rho_types, ONLY: qs_rho_get,&
75
76!$ USE OMP_LIB, ONLY: omp_get_max_threads, omp_get_thread_num
77#include "./base/base_uses.f90"
78
79 IMPLICIT NONE
80
81 PRIVATE
82
83! **************************************************************************************************
84
85 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'lri_optimize_ri_basis'
86
87 PUBLIC :: optimize_lri_basis, &
89
90! **************************************************************************************************
91
92CONTAINS
93
94! **************************************************************************************************
95!> \brief optimizes the lri basis set
96!> \param qs_env qs environment
97! **************************************************************************************************
98 SUBROUTINE optimize_lri_basis(qs_env)
99
100 TYPE(qs_environment_type), POINTER :: qs_env
101
102 INTEGER :: iunit, nkind
103 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
104 TYPE(cp_logger_type), POINTER :: logger
105 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: pmatrix
106 TYPE(lri_density_type), POINTER :: lri_density
107 TYPE(lri_environment_type), POINTER :: lri_env
108 TYPE(lri_opt_type), POINTER :: lri_opt
109 TYPE(mp_para_env_type), POINTER :: para_env
110 TYPE(opt_state_type) :: opt_state
111 TYPE(qs_rho_type), POINTER :: rho_struct
112 TYPE(section_vals_type), POINTER :: dft_section, input, lri_optbas_section
113
114 NULLIFY (atomic_kind_set, dft_section, lri_density, lri_env, &
115 lri_opt, lri_optbas_section, rho_struct)
116 NULLIFY (input, logger, para_env)
117
118 CALL get_qs_env(qs_env, atomic_kind_set=atomic_kind_set, input=input, &
119 lri_env=lri_env, lri_density=lri_density, nkind=nkind, &
120 para_env=para_env, rho=rho_struct)
121
122 ! density matrix
123 CALL qs_rho_get(rho_struct, rho_ao_kp=pmatrix)
124
125 logger => cp_get_default_logger()
126 dft_section => section_vals_get_subs_vals(input, "DFT")
127 lri_optbas_section => section_vals_get_subs_vals(input, &
128 "DFT%QS%OPTIMIZE_LRI_BASIS")
129 iunit = cp_print_key_unit_nr(logger, input, "PRINT%PROGRAM_RUN_INFO", &
130 extension=".opt")
131
132 IF (iunit > 0) THEN
133 WRITE (iunit, '(/," POWELL| Start optimization procedure")')
134 END IF
135
136 ! *** initialization
137 CALL create_lri_opt(lri_opt)
138 CALL init_optimization(lri_env, lri_opt, lri_optbas_section, &
139 opt_state, lri_opt%x, lri_opt%zet_init, nkind, iunit)
140
141 CALL calculate_lri_overlap_aabb(lri_env, qs_env)
142
143 ! *** ======================= START optimization =====================
144 opt_state%state = 0
145 DO
146 IF (opt_state%state == 2) THEN
147 CALL calc_lri_integrals_get_objective(lri_env, lri_density, qs_env, &
148 lri_opt, opt_state, pmatrix, para_env, &
149 nkind)
150 ! lri_density has been re-initialized!
151 CALL set_qs_env(qs_env, lri_density=lri_density)
152 END IF
153
154 IF (opt_state%state == -1) EXIT
155
156 CALL powell_optimize(opt_state%nvar, lri_opt%x, opt_state)
157 CALL update_exponents(lri_env, lri_opt, lri_opt%x, lri_opt%zet_init, nkind)
158 CALL print_optimization_update(opt_state, lri_opt, iunit)
159 END DO
160 ! *** ======================= END optimization =======================
161
162 ! *** get final optimized parameters
163 opt_state%state = 8
164 CALL powell_optimize(opt_state%nvar, lri_opt%x, opt_state)
165 CALL update_exponents(lri_env, lri_opt, lri_opt%x, lri_opt%zet_init, nkind)
166
167 CALL write_optimized_lri_basis(lri_env, dft_section, nkind, lri_opt, &
168 atomic_kind_set)
169
170 IF (iunit > 0) THEN
171 WRITE (iunit, '(" POWELL| Number of function evaluations",T71,I10)') opt_state%nf
172 WRITE (iunit, '(" POWELL| Final value of function",T61,F20.10)') opt_state%fopt
173 WRITE (iunit, '(/," Printed optimized lri basis set to file")')
174 END IF
175
176 CALL cp_print_key_finished_output(iunit, logger, input, &
177 "PRINT%PROGRAM_RUN_INFO")
178
179 CALL deallocate_lri_opt(lri_opt)
180
181 END SUBROUTINE optimize_lri_basis
182
183! **************************************************************************************************
184!> \brief calculates overlap integrals (aabb) of the orbital basis set,
185!> required for LRI basis set optimization
186!> \param lri_env ...
187!> \param qs_env ...
188! **************************************************************************************************
189 SUBROUTINE calculate_lri_overlap_aabb(lri_env, qs_env)
190
191 TYPE(lri_environment_type), POINTER :: lri_env
192 TYPE(qs_environment_type), POINTER :: qs_env
193
194 CHARACTER(LEN=*), PARAMETER :: routinen = 'calculate_lri_overlap_aabb'
195
196 INTEGER :: handle, iac, iatom, ikind, ilist, jatom, &
197 jkind, jneighbor, nba, nbb, nkind, &
198 nlist, nneighbor
199 REAL(kind=dp) :: dab
200 REAL(kind=dp), DIMENSION(3) :: rab
201 TYPE(cell_type), POINTER :: cell
202 TYPE(gto_basis_set_type), POINTER :: obasa, obasb
203 TYPE(lri_int_rho_type), POINTER :: lriir
204 TYPE(lri_list_type), POINTER :: lri_ints_rho
206 DIMENSION(:), POINTER :: nl_iterator
207 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
208 POINTER :: soo_list
209 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
210
211 CALL timeset(routinen, handle)
212 NULLIFY (cell, lriir, lri_ints_rho, nl_iterator, obasa, obasb, &
213 particle_set, soo_list)
214
215 IF (ASSOCIATED(lri_env%soo_list)) THEN
216 soo_list => lri_env%soo_list
217
218 CALL get_qs_env(qs_env=qs_env, nkind=nkind, particle_set=particle_set, &
219 cell=cell)
220
221 IF (ASSOCIATED(lri_env%lri_ints_rho)) THEN
222 CALL deallocate_lri_ints_rho(lri_env%lri_ints_rho)
223 END IF
224
225 CALL allocate_lri_ints_rho(lri_env, lri_env%lri_ints_rho, nkind)
226 lri_ints_rho => lri_env%lri_ints_rho
227
228 CALL neighbor_list_iterator_create(nl_iterator, soo_list)
229 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
230
231 CALL get_iterator_info(nl_iterator, ikind=ikind, jkind=jkind, &
232 nlist=nlist, ilist=ilist, nnode=nneighbor, inode=jneighbor, &
233 iatom=iatom, jatom=jatom, r=rab)
234
235 iac = ikind + nkind*(jkind - 1)
236 dab = sqrt(sum(rab*rab))
237
238 obasa => lri_env%orb_basis(ikind)%gto_basis_set
239 obasb => lri_env%orb_basis(jkind)%gto_basis_set
240 IF (.NOT. ASSOCIATED(obasa)) cycle
241 IF (.NOT. ASSOCIATED(obasb)) cycle
242
243 lriir => lri_ints_rho%lri_atom(iac)%lri_node(ilist)%lri_int_rho(jneighbor)
244
245 nba = obasa%nsgf
246 nbb = obasb%nsgf
247
248 ! calculate integrals (aa,bb)
249 CALL int_overlap_aabb_os(lriir%soaabb, obasa, obasb, rab, lri_env%debug, &
250 lriir%dmax_aabb)
251
252 END DO
253
254 CALL neighbor_list_iterator_release(nl_iterator)
255
256 END IF
257
258 CALL timestop(handle)
259
260 END SUBROUTINE calculate_lri_overlap_aabb
261
262! **************************************************************************************************
263!> \brief initialize optimization parameter
264!> \param lri_env lri environment
265!> \param lri_opt optimization environment
266!> \param lri_optbas_section ...
267!> \param opt_state state of the optimizer
268!> \param x parameters to be optimized, i.e. exponents and contraction coeffs
269!> of the lri basis set
270!> \param zet_init initial values of the exponents
271!> \param nkind number of atom kinds
272!> \param iunit output unit
273! **************************************************************************************************
274 SUBROUTINE init_optimization(lri_env, lri_opt, lri_optbas_section, opt_state, &
275 x, zet_init, nkind, iunit)
276
277 TYPE(lri_environment_type), POINTER :: lri_env
278 TYPE(lri_opt_type), POINTER :: lri_opt
279 TYPE(section_vals_type), POINTER :: lri_optbas_section
280 TYPE(opt_state_type) :: opt_state
281 REAL(kind=dp), DIMENSION(:), POINTER :: x, zet_init
282 INTEGER, INTENT(IN) :: nkind, iunit
283
284 INTEGER :: ikind, iset, ishell, n, nset
285 INTEGER, DIMENSION(:), POINTER :: npgf, nshell
286 REAL(kind=dp), DIMENSION(:, :), POINTER :: zet
287 REAL(kind=dp), DIMENSION(:, :, :), POINTER :: gcc_orig
288 TYPE(gto_basis_set_type), POINTER :: fbas
289
290 NULLIFY (fbas, gcc_orig, npgf, nshell, zet)
291
292 ALLOCATE (lri_opt%ri_gcc_orig(nkind))
293
294 ! *** get parameters
295 CALL get_optimization_parameter(lri_opt, lri_optbas_section, &
296 opt_state)
297
298 opt_state%nvar = 0
299 opt_state%nf = 0
300 opt_state%iprint = 1
301 opt_state%unit = iunit
302
303 ! *** init exponents
304 IF (lri_opt%opt_exps) THEN
305 n = 0
306 DO ikind = 1, nkind
307 fbas => lri_env%ri_basis(ikind)%gto_basis_set
308 CALL get_gto_basis_set(gto_basis_set=fbas, &
309 npgf=npgf, nset=nset, zet=zet)
310 DO iset = 1, nset
311 IF (lri_opt%use_geometric_seq .AND. npgf(iset) > 2) THEN
312 opt_state%nvar = opt_state%nvar + 2
313 CALL reallocate(x, 1, opt_state%nvar)
314 x(n + 1) = maxval(zet(1:npgf(iset), iset))
315 x(n + 2) = minval(zet(1:npgf(iset), iset))
316 n = n + 2
317 ELSE
318 opt_state%nvar = opt_state%nvar + npgf(iset)
319 CALL reallocate(x, 1, opt_state%nvar)
320 x(n + 1:n + npgf(iset)) = zet(1:npgf(iset), iset)
321 n = n + npgf(iset)
322 END IF
323 lri_opt%nexp = lri_opt%nexp + npgf(iset)
324 END DO
325 END DO
326
327 ! *** constraints on exponents
328 IF (lri_opt%use_constraints) THEN
329 ALLOCATE (zet_init(SIZE(x)))
330 zet_init(:) = x
331 ELSE
332 x(:) = sqrt(x)
333 END IF
334 END IF
335
336 ! *** get the original gcc without normalization factor
337 DO ikind = 1, nkind
338 fbas => lri_env%ri_basis(ikind)%gto_basis_set
339 CALL get_original_gcc(lri_opt%ri_gcc_orig(ikind)%gcc_orig, fbas, &
340 lri_opt)
341 END DO
342
343 ! *** init coefficients
344 IF (lri_opt%opt_coeffs) THEN
345 DO ikind = 1, nkind
346 fbas => lri_env%ri_basis(ikind)%gto_basis_set
347 gcc_orig => lri_opt%ri_gcc_orig(ikind)%gcc_orig
348 CALL get_gto_basis_set(gto_basis_set=fbas, &
349 npgf=npgf, nset=nset, nshell=nshell, zet=zet)
350 ! *** Gram Schmidt orthonormalization
351 CALL orthonormalize_gcc(gcc_orig, fbas, lri_opt)
352 n = opt_state%nvar
353 DO iset = 1, nset
354 DO ishell = 1, nshell(iset)
355 opt_state%nvar = opt_state%nvar + npgf(iset)
356 CALL reallocate(x, 1, opt_state%nvar)
357 x(n + 1:n + npgf(iset)) = gcc_orig(1:npgf(iset), ishell, iset)
358 lri_opt%ncoeff = lri_opt%ncoeff + npgf(iset)
359 n = n + npgf(iset)
360 END DO
361 END DO
362 END DO
363 END IF
364
365 IF (iunit > 0) THEN
366 WRITE (iunit, '(/," POWELL| Accuracy",T69,ES12.5)') opt_state%rhoend
367 WRITE (iunit, '(" POWELL| Initial step size",T69,ES12.5)') opt_state%rhobeg
368 WRITE (iunit, '(" POWELL| Maximum number of evaluations",T71,I10)') &
369 opt_state%maxfun
370 WRITE (iunit, '(" POWELL| Total number of parameters",T71,I10)') &
371 opt_state%nvar
372 END IF
373
374 END SUBROUTINE init_optimization
375
376! **************************************************************************************************
377!> \brief read input for optimization
378!> \param lri_opt optimization environment
379!> \param lri_optbas_section ...
380!> \param opt_state state of the optimizer
381! **************************************************************************************************
382 SUBROUTINE get_optimization_parameter(lri_opt, lri_optbas_section, &
383 opt_state)
384
385 TYPE(lri_opt_type), POINTER :: lri_opt
386 TYPE(section_vals_type), POINTER :: lri_optbas_section
387 TYPE(opt_state_type) :: opt_state
388
389 INTEGER :: degree_freedom
390 TYPE(section_vals_type), POINTER :: constrain_exp_section
391
392 NULLIFY (constrain_exp_section)
393
394 ! *** parameter for POWELL optimizer
395 CALL section_vals_val_get(lri_optbas_section, "ACCURACY", &
396 r_val=opt_state%rhoend)
397 CALL section_vals_val_get(lri_optbas_section, "STEP_SIZE", &
398 r_val=opt_state%rhobeg)
399 CALL section_vals_val_get(lri_optbas_section, "MAX_FUN", &
400 i_val=opt_state%maxfun)
401
402 ! *** parameters which are optimized, i.e. exps or coeff or both
403 CALL section_vals_val_get(lri_optbas_section, "DEGREES_OF_FREEDOM", &
404 i_val=degree_freedom)
405
406 SELECT CASE (degree_freedom)
407 CASE (do_lri_opt_all)
408 lri_opt%opt_coeffs = .true.
409 lri_opt%opt_exps = .true.
410 CASE (do_lri_opt_coeff)
411 lri_opt%opt_coeffs = .true.
412 CASE (do_lri_opt_exps)
413 lri_opt%opt_exps = .true.
414 CASE DEFAULT
415 cpabort("No initialization available?????")
416 END SELECT
417
418 ! *** restraint
419 CALL section_vals_val_get(lri_optbas_section, "USE_CONDITION_NUMBER", &
420 l_val=lri_opt%use_condition_number)
421 CALL section_vals_val_get(lri_optbas_section, "CONDITION_WEIGHT", &
422 r_val=lri_opt%cond_weight)
423 CALL section_vals_val_get(lri_optbas_section, "GEOMETRIC_SEQUENCE", &
424 l_val=lri_opt%use_geometric_seq)
425
426 ! *** get constraint info
427 constrain_exp_section => section_vals_get_subs_vals(lri_optbas_section, &
428 "CONSTRAIN_EXPONENTS")
429 CALL section_vals_get(constrain_exp_section, explicit=lri_opt%use_constraints)
430
431 IF (lri_opt%use_constraints) THEN
432 CALL section_vals_val_get(constrain_exp_section, "SCALE", &
433 r_val=lri_opt%scale_exp)
434 CALL section_vals_val_get(constrain_exp_section, "FERMI_EXP", &
435 r_val=lri_opt%fermi_exp)
436 END IF
437
438 END SUBROUTINE get_optimization_parameter
439
440! **************************************************************************************************
441!> \brief update exponents after optimization step
442!> \param lri_env lri environment
443!> \param lri_opt optimization environment
444!> \param x optimization parameters
445!> \param zet_init initial values of the exponents
446!> \param nkind number of atomic kinds
447! **************************************************************************************************
448 SUBROUTINE update_exponents(lri_env, lri_opt, x, zet_init, nkind)
449
450 TYPE(lri_environment_type), POINTER :: lri_env
451 TYPE(lri_opt_type), POINTER :: lri_opt
452 REAL(kind=dp), DIMENSION(:), POINTER :: x, zet_init
453 INTEGER, INTENT(IN) :: nkind
454
455 INTEGER :: ikind, iset, ishell, n, nset, nvar_exp
456 INTEGER, DIMENSION(:), POINTER :: npgf, nshell
457 REAL(kind=dp) :: zet_max, zet_min
458 REAL(kind=dp), DIMENSION(:), POINTER :: zet, zet_trans
459 REAL(kind=dp), DIMENSION(:, :, :), POINTER :: gcc_orig
460 TYPE(gto_basis_set_type), POINTER :: fbas
461
462 NULLIFY (fbas, gcc_orig, npgf, nshell, zet_trans, zet)
463
464 ! nvar_exp: number of exponents that are variables
465 nvar_exp = SIZE(x) - lri_opt%ncoeff
466 ALLOCATE (zet_trans(nvar_exp))
467
468 ! *** update exponents
469 IF (lri_opt%opt_exps) THEN
470 IF (lri_opt%use_constraints) THEN
471 zet => x(1:nvar_exp)
472 CALL transfer_exp(lri_opt, zet, zet_init, zet_trans, nvar_exp)
473 ELSE
474 zet_trans(:) = x(1:nvar_exp)**2.0_dp
475 END IF
476 n = 0
477 DO ikind = 1, nkind
478 fbas => lri_env%ri_basis(ikind)%gto_basis_set
479 CALL get_gto_basis_set(gto_basis_set=fbas, npgf=npgf, nset=nset)
480 DO iset = 1, nset
481 IF (lri_opt%use_geometric_seq .AND. npgf(iset) > 2) THEN
482 zet_max = maxval(zet_trans(n + 1:n + 2))
483 zet_min = minval(zet_trans(n + 1:n + 2))
484 zet => fbas%zet(1:npgf(iset), iset)
485 CALL geometric_progression(zet, zet_max, zet_min, npgf(iset))
486 n = n + 2
487 ELSE
488 fbas%zet(1:npgf(iset), iset) = zet_trans(n + 1:n + npgf(iset))
489 n = n + npgf(iset)
490 END IF
491 END DO
492 END DO
493 END IF
494
495 ! *** update coefficients
496 IF (lri_opt%opt_coeffs) THEN
497 n = nvar_exp
498 DO ikind = 1, nkind
499 fbas => lri_env%ri_basis(ikind)%gto_basis_set
500 gcc_orig => lri_opt%ri_gcc_orig(ikind)%gcc_orig
501 CALL get_gto_basis_set(gto_basis_set=fbas, &
502 nshell=nshell, npgf=npgf, nset=nset)
503 DO iset = 1, nset
504 DO ishell = 1, nshell(iset)
505 gcc_orig(1:npgf(iset), ishell, iset) = x(n + 1:n + npgf(iset))
506 n = n + npgf(iset)
507 END DO
508 END DO
509 ! *** Gram Schmidt orthonormalization
510 CALL orthonormalize_gcc(gcc_orig, fbas, lri_opt)
511 END DO
512 END IF
513
514 DEALLOCATE (zet_trans)
515 END SUBROUTINE update_exponents
516
517! **************************************************************************************************
518!> \brief employ Fermi constraint, transfer exponents
519!> \param lri_opt optimization environment
520!> \param zet untransferred exponents
521!> \param zet_init initial value of the exponents
522!> \param zet_trans transferred exponents
523!> \param nvar number of optimized exponents
524! **************************************************************************************************
525 SUBROUTINE transfer_exp(lri_opt, zet, zet_init, zet_trans, nvar)
526
527 TYPE(lri_opt_type), POINTER :: lri_opt
528 REAL(kind=dp), DIMENSION(:), POINTER :: zet, zet_init, zet_trans
529 INTEGER, INTENT(IN) :: nvar
530
531 REAL(kind=dp) :: a
532 REAL(kind=dp), DIMENSION(:), POINTER :: zet_max, zet_min
533
534 ALLOCATE (zet_max(nvar), zet_min(nvar))
535
536 zet_min(:) = zet_init(:)*(1.0_dp - lri_opt%scale_exp)
537 zet_max(:) = zet_init(:)*(1.0_dp + lri_opt%scale_exp)
538
539 a = lri_opt%fermi_exp
540
541 zet_trans = zet_min + (zet_max - zet_min)/(1 + exp(-a*(zet - zet_init)))
542
543 DEALLOCATE (zet_max, zet_min)
544
545 END SUBROUTINE transfer_exp
546
547! **************************************************************************************************
548!> \brief complete geometric sequence
549!> \param zet all exponents of the set
550!> \param zet_max maximal exponent of the set
551!> \param zet_min minimal exponent of the set
552!> \param nexp number of exponents of the set
553! **************************************************************************************************
554 SUBROUTINE geometric_progression(zet, zet_max, zet_min, nexp)
555
556 REAL(kind=dp), DIMENSION(:), POINTER :: zet
557 REAL(kind=dp), INTENT(IN) :: zet_max, zet_min
558 INTEGER, INTENT(IN) :: nexp
559
560 INTEGER :: i, n
561 REAL(kind=dp) :: q
562
563 n = nexp - 1
564
565 q = (zet_min/zet_max)**(1._dp/real(n, dp))
566
567 DO i = 1, nexp
568 zet(i) = zet_max*q**(i - 1)
569 END DO
570
571 END SUBROUTINE geometric_progression
572
573! **************************************************************************************************
574!> \brief calculates the lri integrals and coefficients with the new exponents
575!> of the lri basis sets and calculates the objective function
576!> \param lri_env lri environment
577!> \param lri_density ...
578!> \param qs_env ...
579!> \param lri_opt optimization environment
580!> \param opt_state state of the optimizer
581!> \param pmatrix density matrix
582!> \param para_env ...
583!> \param nkind number of atomic kinds
584! **************************************************************************************************
585 SUBROUTINE calc_lri_integrals_get_objective(lri_env, lri_density, qs_env, &
586 lri_opt, opt_state, pmatrix, para_env, &
587 nkind)
588
589 TYPE(lri_environment_type), POINTER :: lri_env
590 TYPE(lri_density_type), POINTER :: lri_density
591 TYPE(qs_environment_type), POINTER :: qs_env
592 TYPE(lri_opt_type), POINTER :: lri_opt
593 TYPE(opt_state_type) :: opt_state
594 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: pmatrix
595 TYPE(mp_para_env_type), POINTER :: para_env
596 INTEGER, INTENT(IN) :: nkind
597
598 INTEGER :: ikind, nset
599 INTEGER, DIMENSION(:), POINTER :: npgf
600 INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index
601 TYPE(gto_basis_set_type), POINTER :: fbas
602
603 NULLIFY (fbas, npgf)
604
605 !*** build new transformation matrices sphi with new exponents
606 lri_env%store_integrals = .true.
607 DO ikind = 1, nkind
608 fbas => lri_env%ri_basis(ikind)%gto_basis_set
609 CALL get_gto_basis_set(gto_basis_set=fbas, npgf=npgf, nset=nset)
610 !build new sphi
611 fbas%gcc = lri_opt%ri_gcc_orig(ikind)%gcc_orig
612 CALL init_orb_basis_set(fbas)
613 END DO
614 CALL lri_basis_init(lri_env)
615 CALL calculate_lri_integrals(lri_env, qs_env)
616 CALL calculate_avec_lri(lri_env, lri_density, pmatrix, cell_to_index)
617 IF (lri_opt%use_condition_number) THEN
619 END IF
620 CALL calculate_objective(lri_env, lri_density, lri_opt, pmatrix, para_env, &
621 opt_state%f)
622
623 END SUBROUTINE calc_lri_integrals_get_objective
624
625! **************************************************************************************************
626!> \brief calculates the objective function defined as integral of the square
627!> of rhoexact - rhofit, i.e. integral[(rhoexact-rhofit)**2]
628!> rhoexact is the exact pair density and rhofit the lri pair density
629!> \param lri_env lri environment
630!> \param lri_density ...
631!> \param lri_opt optimization environment
632!> \param pmatrix density matrix
633!> \param para_env ...
634!> \param fobj objective function
635! **************************************************************************************************
636 SUBROUTINE calculate_objective(lri_env, lri_density, lri_opt, pmatrix, para_env, &
637 fobj)
638
639 TYPE(lri_environment_type), POINTER :: lri_env
640 TYPE(lri_density_type), POINTER :: lri_density
641 TYPE(lri_opt_type), POINTER :: lri_opt
642 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: pmatrix
643 TYPE(mp_para_env_type), POINTER :: para_env
644 REAL(kind=dp), INTENT(OUT) :: fobj
645
646 CHARACTER(LEN=*), PARAMETER :: routinen = 'calculate_objective'
647
648 INTEGER :: handle, iac, iatom, ikind, ilist, isgfa, ispin, jatom, jkind, jneighbor, jsgfa, &
649 ksgfb, lsgfb, mepos, nba, nbb, nfa, nfb, nkind, nlist, nn, nneighbor, nspin, nthread
650 LOGICAL :: found, trans
651 REAL(kind=dp) :: obj_ab, rhoexact_sq, rhofit_sq, rhomix
652 REAL(kind=dp), DIMENSION(:, :), POINTER :: pbij
653 TYPE(dbcsr_type), POINTER :: pmat
654 TYPE(lri_int_rho_type), POINTER :: lriir
655 TYPE(lri_int_type), POINTER :: lrii
656 TYPE(lri_list_type), POINTER :: lri_rho
657 TYPE(lri_rhoab_type), POINTER :: lrho
659 DIMENSION(:), POINTER :: nl_iterator
660 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
661 POINTER :: soo_list
662
663 CALL timeset(routinen, handle)
664 NULLIFY (lrii, lriir, lri_rho, lrho, nl_iterator, pmat, soo_list)
665
666 IF (ASSOCIATED(lri_env%soo_list)) THEN
667 soo_list => lri_env%soo_list
668
669 nkind = lri_env%lri_ints%nkind
670 nspin = SIZE(pmatrix, 1)
671 cpassert(SIZE(pmatrix, 2) == 1)
672 nthread = 1
673!$ nthread = omp_get_max_threads()
674
675 fobj = 0._dp
676 lri_opt%rho_diff = 0._dp
677
678 DO ispin = 1, nspin
679
680 pmat => pmatrix(ispin, 1)%matrix
681 lri_rho => lri_density%lri_rhos(ispin)%lri_list
682
683 CALL neighbor_list_iterator_create(nl_iterator, soo_list, nthread=nthread)
684!$OMP PARALLEL DEFAULT(NONE)&
685!$OMP SHARED (nthread,nl_iterator,pmat,nkind,fobj,lri_env,lri_opt,lri_rho)&
686!$OMP PRIVATE (mepos,ikind,jkind,iatom,jatom,nlist,ilist,nneighbor,jneighbor,&
687!$OMP iac,lrii,lriir,lrho,nfa,nfb,nba,nbb,nn,rhoexact_sq,rhomix,rhofit_sq,&
688!$OMP obj_ab,pbij,trans,found,isgfa,jsgfa,ksgfb,lsgfb)
689
690 mepos = 0
691!$ mepos = omp_get_thread_num()
692
693 DO WHILE (neighbor_list_iterate(nl_iterator, mepos) == 0)
694 CALL get_iterator_info(nl_iterator, mepos=mepos, ikind=ikind, jkind=jkind, iatom=iatom, &
695 jatom=jatom, nlist=nlist, ilist=ilist, nnode=nneighbor, inode=jneighbor)
696
697 iac = ikind + nkind*(jkind - 1)
698
699 IF (.NOT. ASSOCIATED(lri_env%lri_ints%lri_atom(iac)%lri_node)) cycle
700
701 lrii => lri_env%lri_ints%lri_atom(iac)%lri_node(ilist)%lri_int(jneighbor)
702 lriir => lri_env%lri_ints_rho%lri_atom(iac)%lri_node(ilist)%lri_int_rho(jneighbor)
703 lrho => lri_rho%lri_atom(iac)%lri_node(ilist)%lri_rhoab(jneighbor)
704 nfa = lrii%nfa
705 nfb = lrii%nfb
706 nba = lrii%nba
707 nbb = lrii%nbb
708 nn = nfa + nfb
709
710 rhoexact_sq = 0._dp
711 rhomix = 0._dp
712 rhofit_sq = 0._dp
713 obj_ab = 0._dp
714
715 NULLIFY (pbij)
716 IF (iatom <= jatom) THEN
717 CALL dbcsr_get_block_p(matrix=pmat, row=iatom, col=jatom, block=pbij, found=found)
718 trans = .false.
719 ELSE
720 CALL dbcsr_get_block_p(matrix=pmat, row=jatom, col=iatom, block=pbij, found=found)
721 trans = .true.
722 END IF
723 cpassert(found)
724
725 ! *** calculate integral of the square of exact density rhoexact_sq
726 IF (trans) THEN
727 DO isgfa = 1, nba
728 DO jsgfa = 1, nba
729 DO ksgfb = 1, nbb
730 DO lsgfb = 1, nbb
731 rhoexact_sq = rhoexact_sq + pbij(ksgfb, isgfa)*pbij(lsgfb, jsgfa) &
732 *lriir%soaabb(isgfa, jsgfa, ksgfb, lsgfb)
733 END DO
734 END DO
735 END DO
736 END DO
737 ELSE
738 DO isgfa = 1, nba
739 DO jsgfa = 1, nba
740 DO ksgfb = 1, nbb
741 DO lsgfb = 1, nbb
742 rhoexact_sq = rhoexact_sq + pbij(isgfa, ksgfb)*pbij(jsgfa, lsgfb) &
743 *lriir%soaabb(isgfa, jsgfa, ksgfb, lsgfb)
744 END DO
745 END DO
746 END DO
747 END DO
748 END IF
749
750 ! *** calculate integral of the square of the fitted density rhofit_sq
751 DO isgfa = 1, nfa
752 DO jsgfa = 1, nfa
753 rhofit_sq = rhofit_sq + lrho%avec(isgfa)*lrho%avec(jsgfa) &
754 *lri_env%bas_prop(ikind)%ri_ovlp(isgfa, jsgfa)
755 END DO
756 END DO
757 IF (iatom /= jatom) THEN
758 DO ksgfb = 1, nfb
759 DO lsgfb = 1, nfb
760 rhofit_sq = rhofit_sq + lrho%avec(nfa + ksgfb)*lrho%avec(nfa + lsgfb) &
761 *lri_env%bas_prop(jkind)%ri_ovlp(ksgfb, lsgfb)
762 END DO
763 END DO
764 DO isgfa = 1, nfa
765 DO ksgfb = 1, nfb
766 rhofit_sq = rhofit_sq + 2._dp*lrho%avec(isgfa)*lrho%avec(nfa + ksgfb) &
767 *lrii%sab(isgfa, ksgfb)
768 END DO
769 END DO
770 END IF
771
772 ! *** and integral of the product of exact and fitted density rhomix
773 IF (iatom == jatom) THEN
774 rhomix = sum(lrho%avec(1:nfa)*lrho%tvec(1:nfa))
775 ELSE
776 rhomix = sum(lrho%avec(1:nn)*lrho%tvec(1:nn))
777 END IF
778
779 ! *** calculate contribution to the objective function for pair ab
780 ! *** taking density matrix symmetry in account, double-count for off-diagonal blocks
781 IF (iatom == jatom) THEN
782 obj_ab = rhoexact_sq - 2._dp*rhomix + rhofit_sq
783 ELSE
784 obj_ab = 2.0_dp*(rhoexact_sq - 2._dp*rhomix + rhofit_sq)
785 END IF
786
787!$OMP CRITICAL(addfun)
788 IF (lri_opt%use_condition_number) THEN
789 fobj = fobj + obj_ab + lri_opt%cond_weight*log(lrii%cond_num)
790 lri_opt%rho_diff = lri_opt%rho_diff + obj_ab
791 ELSE
792 fobj = fobj + obj_ab
793 END IF
794!$OMP END CRITICAL(addfun)
795
796 END DO
797!$OMP END PARALLEL
798
799 CALL neighbor_list_iterator_release(nl_iterator)
800
801 END DO
802 CALL para_env%sum(fobj)
803
804 END IF
805
806 CALL timestop(handle)
807
808 END SUBROUTINE calculate_objective
809
810! **************************************************************************************************
811!> \brief get condition number of overlap matrix
812!> \param lri_env lri environment
813! **************************************************************************************************
815
816 TYPE(lri_environment_type), POINTER :: lri_env
817
818 CHARACTER(LEN=*), PARAMETER :: routinen = 'get_condition_number_of_overlap'
819
820 INTEGER :: handle, iac, iatom, ikind, ilist, info, &
821 jatom, jkind, jneighbor, lwork, mepos, &
822 nfa, nfb, nkind, nlist, nn, nneighbor, &
823 nthread
824 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: diag, off_diag, tau
825 REAL(kind=dp), DIMENSION(:), POINTER :: work
826 REAL(kind=dp), DIMENSION(:, :), POINTER :: smat
827 TYPE(lri_int_type), POINTER :: lrii
829 DIMENSION(:), POINTER :: nl_iterator
830 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
831 POINTER :: soo_list
832
833 CALL timeset(routinen, handle)
834 NULLIFY (lrii, nl_iterator, smat, soo_list)
835
836 soo_list => lri_env%soo_list
837
838 nkind = lri_env%lri_ints%nkind
839 nthread = 1
840!$ nthread = omp_get_max_threads()
841
842 CALL neighbor_list_iterator_create(nl_iterator, soo_list, nthread=nthread)
843!$OMP PARALLEL DEFAULT(NONE)&
844!$OMP SHARED (nthread,nl_iterator,nkind,lri_env)&
845!$OMP PRIVATE (mepos,ikind,jkind,iatom,jatom,nlist,ilist,nneighbor,jneighbor,&
846!$OMP diag,off_diag,smat,tau,work,iac,lrii,nfa,nfb,nn,info,lwork)
847
848 mepos = 0
849!$ mepos = omp_get_thread_num()
850
851 DO WHILE (neighbor_list_iterate(nl_iterator, mepos) == 0)
852 CALL get_iterator_info(nl_iterator, mepos=mepos, ikind=ikind, jkind=jkind, iatom=iatom, &
853 jatom=jatom, nlist=nlist, ilist=ilist, nnode=nneighbor, inode=jneighbor)
854
855 iac = ikind + nkind*(jkind - 1)
856 IF (.NOT. ASSOCIATED(lri_env%lri_ints%lri_atom(iac)%lri_node)) cycle
857 lrii => lri_env%lri_ints%lri_atom(iac)%lri_node(ilist)%lri_int(jneighbor)
858
859 nfa = lrii%nfa
860 nfb = lrii%nfb
861 nn = nfa + nfb
862
863 ! build the overlap matrix
864 IF (iatom == jatom) THEN
865 ALLOCATE (smat(nfa, nfa))
866 ELSE
867 ALLOCATE (smat(nn, nn))
868 END IF
869 smat(1:nfa, 1:nfa) = lri_env%bas_prop(ikind)%ri_ovlp(1:nfa, 1:nfa)
870 IF (iatom /= jatom) THEN
871 nn = nfa + nfb
872 smat(1:nfa, nfa + 1:nn) = lrii%sab(1:nfa, 1:nfb)
873 smat(nfa + 1:nn, 1:nfa) = transpose(lrii%sab(1:nfa, 1:nfb))
874 smat(nfa + 1:nn, nfa + 1:nn) = lri_env%bas_prop(jkind)%ri_ovlp(1:nfb, 1:nfb)
875 END IF
876
877 IF (iatom == jatom) nn = nfa
878 ALLOCATE (diag(nn), off_diag(nn - 1), tau(nn - 1), work(1))
879 diag = 0.0_dp
880 off_diag = 0.0_dp
881 tau = 0.0_dp
882 work = 0.0_dp
883 lwork = -1
884 ! get lwork
885 CALL dsytrd('U', nn, smat, nn, diag, off_diag, tau, work, lwork, info)
886 lwork = int(work(1))
887 CALL reallocate(work, 1, lwork)
888 ! get the eigenvalues
889 CALL dsytrd('U', nn, smat, nn, diag, off_diag, tau, work, lwork, info)
890 CALL dsterf(nn, diag, off_diag, info)
891
892 lrii%cond_num = maxval(abs(diag))/minval(abs(diag))
893
894 DEALLOCATE (diag, off_diag, smat, tau, work)
895 END DO
896!$OMP END PARALLEL
897
898 CALL neighbor_list_iterator_release(nl_iterator)
899
900 CALL timestop(handle)
901
903
904! **************************************************************************************************
905!> \brief print recent information on optimization
906!> \param opt_state state of the optimizer
907!> \param lri_opt optimization environment
908!> \param iunit ...
909! **************************************************************************************************
910 SUBROUTINE print_optimization_update(opt_state, lri_opt, iunit)
911
912 TYPE(opt_state_type) :: opt_state
913 TYPE(lri_opt_type), POINTER :: lri_opt
914 INTEGER, INTENT(IN) :: iunit
915
916 INTEGER :: n10
917
918 n10 = max(opt_state%maxfun/100, 1)
919
920 IF (opt_state%nf == 2 .AND. opt_state%state == 2 .AND. iunit > 0) THEN
921 WRITE (iunit, '(/," POWELL| Initial value of function",T61,F20.10)') opt_state%f
922 END IF
923 IF (mod(opt_state%nf, n10) == 0 .AND. opt_state%nf > 1 .AND. iunit > 0) THEN
924 WRITE (iunit, '(" POWELL| Reached",i4,"% of maximal function calls",T61,F20.10)') &
925 int(real(opt_state%nf, dp)/real(opt_state%maxfun, dp)*100._dp), opt_state%fopt
926 END IF
927 IF (lri_opt%use_condition_number) THEN
928 IF (mod(opt_state%nf, n10) == 0 .AND. opt_state%nf > 1 .AND. iunit > 0) THEN
929 WRITE (iunit, '(" POWELL| Recent value of function without condition nr.",T61,F20.10)') &
930 lri_opt%rho_diff
931 END IF
932 END IF
933
934 END SUBROUTINE print_optimization_update
935
936! **************************************************************************************************
937!> \brief write optimized LRI basis set to file
938!> \param lri_env ...
939!> \param dft_section ...
940!> \param nkind ...
941!> \param lri_opt ...
942!> \param atomic_kind_set ...
943! **************************************************************************************************
944 SUBROUTINE write_optimized_lri_basis(lri_env, dft_section, nkind, lri_opt, &
945 atomic_kind_set)
946
947 TYPE(lri_environment_type), POINTER :: lri_env
948 TYPE(section_vals_type), POINTER :: dft_section
949 INTEGER, INTENT(IN) :: nkind
950 TYPE(lri_opt_type), POINTER :: lri_opt
951 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
952
953 CHARACTER(LEN=default_path_length) :: filename
954 INTEGER :: cc_l, ikind, ipgf, iset, ishell, nset, &
955 output_file
956 INTEGER, DIMENSION(:), POINTER :: lmax, lmin, npgf, nshell
957 INTEGER, DIMENSION(:, :), POINTER :: l
958 REAL(kind=dp), DIMENSION(:, :), POINTER :: zet
959 REAL(kind=dp), DIMENSION(:, :, :), POINTER :: gcc_orig
960 TYPE(cp_logger_type), POINTER :: logger
961 TYPE(gto_basis_set_type), POINTER :: fbas
962 TYPE(section_vals_type), POINTER :: print_key
963
964 NULLIFY (fbas, gcc_orig, l, lmax, lmin, logger, npgf, nshell, print_key, zet)
965
966 !*** do the printing
967 print_key => section_vals_get_subs_vals(dft_section, &
968 "PRINT%OPTIMIZE_LRI_BASIS")
969 logger => cp_get_default_logger()
970 IF (btest(cp_print_key_should_output(logger%iter_info, &
971 dft_section, "PRINT%OPTIMIZE_LRI_BASIS"), &
972 cp_p_file)) THEN
973 output_file = cp_print_key_unit_nr(logger, dft_section, &
974 "PRINT%OPTIMIZE_LRI_BASIS", &
975 extension=".opt", &
976 file_status="REPLACE", &
977 file_action="WRITE", &
978 file_form="FORMATTED")
979
980 IF (output_file > 0) THEN
981
982 filename = cp_print_key_generate_filename(logger, &
983 print_key, extension=".opt", &
984 my_local=.true.)
985
986 DO ikind = 1, nkind
987 fbas => lri_env%ri_basis(ikind)%gto_basis_set
988 gcc_orig => lri_opt%ri_gcc_orig(ikind)%gcc_orig
989 CALL get_gto_basis_set(gto_basis_set=fbas, &
990 l=l, lmax=lmax, lmin=lmin, &
991 npgf=npgf, nshell=nshell, &
992 nset=nset, zet=zet)
993 WRITE (output_file, '(T1,A2,T5,A)') trim(atomic_kind_set(ikind)%name), &
994 trim(fbas%name)
995 WRITE (output_file, '(T1,I4)') nset
996 DO iset = 1, nset
997 WRITE (output_file, '(4(1X,I0))', advance='no') 2, lmin(iset), &
998 lmax(iset), npgf(iset)
999 cc_l = 1
1000 DO ishell = 1, nshell(iset)
1001 IF (ishell /= nshell(iset)) THEN
1002 IF (l(ishell, iset) == l(ishell + 1, iset)) THEN
1003 cc_l = cc_l + 1
1004 ELSE
1005 WRITE (output_file, '(1X,I0)', advance='no') cc_l
1006 cc_l = 1
1007 END IF
1008 ELSE
1009 WRITE (output_file, '(1X,I0)') cc_l
1010 END IF
1011 END DO
1012 DO ipgf = 1, npgf(iset)
1013 WRITE (output_file, '(F18.12)', advance='no') zet(ipgf, iset)
1014 DO ishell = 1, nshell(iset)
1015 IF (ishell == nshell(iset)) THEN
1016 WRITE (output_file, '(T5,F18.12)') gcc_orig(ipgf, ishell, iset)
1017 ELSE
1018 WRITE (output_file, '(T5,F18.12)', advance='no') gcc_orig(ipgf, ishell, iset)
1019 END IF
1020 END DO
1021 END DO
1022 END DO
1023 END DO
1024
1025 END IF
1026
1027 CALL cp_print_key_finished_output(output_file, logger, dft_section, &
1028 "PRINT%OPTIMIZE_LRI_BASIS")
1029 END IF
1030
1031 END SUBROUTINE write_optimized_lri_basis
1032
1033END MODULE lri_optimize_ri_basis
Define the atomic kind types and their sub types.
subroutine, public init_orb_basis_set(gto_basis_set)
Initialise a Gaussian-type orbital (GTO) basis set data set.
subroutine, public get_gto_basis_set(gto_basis_set, name, aliases, norm_type, kind_radius, ncgf, nset, nsgf, cgf_symbol, sgf_symbol, norm_cgf, set_radius, lmax, lmin, lx, ly, lz, m, ncgf_set, npgf, nsgf_set, nshell, cphi, pgf_radius, sphi, scon, zet, first_cgf, first_sgf, l, last_cgf, last_sgf, n, gcc, maxco, maxl, maxpgf, maxsgf_set, maxshell, maxso, nco_sum, npgf_sum, nshell_sum, maxder, short_kind_radius)
...
Handles all functions related to the CELL.
Definition cell_types.F:15
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
character(len=default_path_length) function, public cp_print_key_generate_filename(logger, print_key, middle_name, extension, my_local)
Utility function that returns a unit number to write the print key. Might open a file with a unique f...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
integer, parameter, public cp_p_file
integer function, public cp_print_key_should_output(iteration_info, basis_section, print_key_path, used_print_key, first_time)
returns what should be done with the given property if btest(res,cp_p_store) then the property should...
Calculation of contracted, spherical Gaussian integrals using the (OS) integral scheme....
subroutine, public int_overlap_aabb_os(saabb, oba, obb, rab, debug, dmax)
calculate overlap integrals (aa,bb)
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public do_lri_opt_coeff
integer, parameter, public do_lri_opt_all
integer, parameter, public do_lri_opt_exps
objects that represent the structure of input sections and the data contained in an input section
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
subroutine, public section_vals_get(section_vals, ref_count, n_repetition, n_subs_vals_rep, section, explicit)
returns various attributes about the section_vals
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_path_length
Definition kinds.F:58
initializes the environment for lri lri : local resolution of the identity
subroutine, public lri_basis_init(lri_env)
initializes the lri basis: calculates the norm, self-overlap and integral of the ri basis
Calculates integral matrices for LRIGPW method lri : local resolution of the identity.
subroutine, public calculate_lri_integrals(lri_env, qs_env)
calculates integrals needed for the LRI density fitting, integrals are calculated once,...
subroutine, public calculate_avec_lri(lri_env, lri_density, pmatrix, cell_to_index, response_density)
performs the fitting of the density; solves the linear system of equations; yield the expansion coeff...
contains the types and subroutines for dealing with the lri_env lri : local resolution of the identit...
subroutine, public allocate_lri_ints_rho(lri_env, lri_ints_rho, nkind)
allocate lri_ints_rho, storing integral for the exact density
subroutine, public deallocate_lri_ints_rho(lri_ints_rho)
deallocates the given lri_ints_rho
sets the environment for optimization of exponents and contraction coefficients of the lri auxiliary ...
subroutine, public orthonormalize_gcc(gcc, gto_basis_set, lri_opt)
orthonormalize contraction coefficients using Gram-Schmidt
subroutine, public create_lri_opt(lri_opt)
creates lri_opt
subroutine, public deallocate_lri_opt(lri_opt)
deallocates lri_opt
subroutine, public get_original_gcc(gcc_orig, gto_basis_set, lri_opt)
primitive Cartesian Gaussian functions are normalized. The normalization factor is included in the Ga...
Optimizes exponents and contraction coefficients of the lri auxiliary basis sets using the UOBYQA min...
subroutine, public get_condition_number_of_overlap(lri_env)
get condition number of overlap matrix
subroutine, public optimize_lri_basis(qs_env)
optimizes the lri basis set
Utility routines for the memory handling.
Interface to the message passing library MPI.
Define the data structure for the particle information.
Definition powell.F:9
subroutine, public powell_optimize(n, x, optstate)
...
Definition powell.F:55
subroutine, public set_qs_env(qs_env, super_cell, mos, qmmm, qmmm_periodic, ewald_env, ewald_pw, mpools, rho_external, external_vxc, mask, scf_control, rel_control, qs_charges, ks_env, ks_qmmm_env, wf_history, scf_env, active_space, input, oce, rho_atom_set, rho0_atom_set, rho0_mpole, run_rtp, rtp, rhoz_set, rhoz_tot, ecoul_1c, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, efield, linres_control, xas_env, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, ls_scf_env, do_transport, transport_env, lri_env, lri_density, exstate_env, ec_env, dispersion_env, gcp_env, mp2_env, bs_env, kg_env, force, kpoints, wanniercentres, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, rhs)
Set the QUICKSTEP environment.
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, ecoul_1c, rho0_s_rs, rho0_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, rhs)
Get the QUICKSTEP environment.
Define the neighbor list data types and the corresponding functionality.
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)
...
superstucture that hold various representations of the density and keeps track of which ones are vali...
subroutine, public qs_rho_get(rho_struct, rho_ao, rho_ao_im, rho_ao_kp, rho_ao_im_kp, rho_r, drho_r, rho_g, drho_g, tau_r, tau_g, rho_r_valid, drho_r_valid, rho_g_valid, drho_g_valid, tau_r_valid, tau_g_valid, tot_rho_r, tot_rho_g, rho_r_sccs, soft_valid, complex_rho_ao)
returns info about the density described by this object. If some representation is not available an e...
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
type of a logger, at the moment it contains just a print level starting at which level it should be l...
stores all the informations relevant to an mpi environment
keeps the density in various representations, keeping track of which ones are valid.