(git:e7e05ae)
qs_linres_module.F
Go to the documentation of this file.
1 !--------------------------------------------------------------------------------------------------!
2 ! CP2K: A general program to perform molecular dynamics simulations !
3 ! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4 ! !
5 ! SPDX-License-Identifier: GPL-2.0-or-later !
6 !--------------------------------------------------------------------------------------------------!
7 
8 ! **************************************************************************************************
9 !> \brief Contains the setup for the calculation of properties by linear response
10 !> by the application of second order density functional perturbation theory.
11 !> The knowledge of the ground state energy, density and wavefunctions is assumed.
12 !> Uses the self consistent approach.
13 !> Properties that can be calculated : none
14 !> \par History
15 !> created 06-2005 [MI]
16 !> \author MI
17 ! **************************************************************************************************
19  USE bibliography, ONLY: ditler2021,&
20  ditler2022,&
21  weber2009,&
22  cite_reference
23  USE cp_control_types, ONLY: dft_control_type
25  cp_logger_type
28  USE dbcsr_api, ONLY: dbcsr_p_type
29  USE force_env_types, ONLY: force_env_get,&
30  force_env_type,&
31  use_qmmm,&
33  USE input_constants, ONLY: lr_current,&
34  lr_none,&
43  section_vals_type,&
45  USE kinds, ONLY: dp
46  USE qs_dcdr, ONLY: apt_dr,&
51  USE qs_dcdr_utils, ONLY: dcdr_env_cleanup,&
54  USE qs_density_matrices, ONLY: calculate_density_matrix
55  USE qs_environment_types, ONLY: get_qs_env,&
56  qs_environment_type,&
65  epr_g_so,&
66  epr_g_soo,&
67  epr_g_zke,&
73  issc_issc,&
74  issc_print,&
77  USE qs_linres_nmr_shift, ONLY: nmr_shift,&
81  USE qs_linres_op, ONLY: current_operators,&
85  polar_polar,&
86  polar_print,&
88  USE qs_linres_types, ONLY: current_env_type,&
89  dcdr_env_type,&
90  epr_env_type,&
91  issc_env_type,&
92  linres_control_type,&
93  nmr_env_type,&
94  vcd_env_type
95  USE qs_mfp, ONLY: mfp_aat,&
99  USE qs_mo_types, ONLY: mo_set_type
100  USE qs_p_env_methods, ONLY: p_env_create,&
102  USE qs_p_env_types, ONLY: p_env_release,&
103  qs_p_env_type
105  USE qs_rho_types, ONLY: qs_rho_get,&
106  qs_rho_type
107  USE qs_vcd, ONLY: aat_dv,&
108  apt_dv,&
112  USE qs_vcd_utils, ONLY: vcd_env_cleanup,&
113  vcd_env_init,&
114  vcd_print
115 #include "./base/base_uses.f90"
116 
117  IMPLICIT NONE
118 
119  PRIVATE
121 
122  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_linres_module'
123 
124 CONTAINS
125 
126 ! *****************************************************************************
127 !> \brief Calculates the derivatives of the MO coefficients dC/dV^lambda_beta
128 !> wrt to nuclear velocities. The derivative is indexed by `beta`, the
129 !> electric dipole operator by `alpha`.
130 !> Calculates the APT and AAT in velocity form
131 !> P^lambda_alpha,beta = d< mu_alpha >/dV^lambda_beta
132 !> M^lambda_alpha,beta = d< m_alpha >/dV^lambda_beta
133 !> \param qs_env ...
134 !> \param p_env ...
135 !> \author Edward Ditler
136 ! **************************************************************************************************
137  SUBROUTINE vcd_linres(qs_env, p_env)
138  TYPE(qs_environment_type), POINTER :: qs_env
139  TYPE(qs_p_env_type) :: p_env
140 
141  INTEGER :: beta, i, latom
142  LOGICAL :: mfp_is_done, mfp_repeat
143  TYPE(vcd_env_type) :: vcd_env
144 
145  CALL cite_reference(ditler2022)
146 
147  ! We need the position perturbation for the velocity perturbation operator
148  CALL vcd_env_init(vcd_env, qs_env)
149 
150  mfp_repeat = vcd_env%distributed_origin
151  mfp_is_done = .false.
152 
153  qs_env%linres_control%linres_restart = .true.
154 
155  ! Iterate over the list of atoms for which we want to calculate the APTs/AATs
156  ! default is all atoms.
157  DO latom = 1, SIZE(vcd_env%dcdr_env%list_of_atoms)
158  vcd_env%dcdr_env%lambda = vcd_env%dcdr_env%list_of_atoms(latom)
159 
160  CALL prepare_per_atom(vcd_env%dcdr_env, qs_env)
161  CALL prepare_per_atom_vcd(vcd_env, qs_env)
162 
163  DO beta = 1, 3 ! in every direction
164 
165  vcd_env%dcdr_env%beta = beta
166  vcd_env%dcdr_env%deltaR(vcd_env%dcdr_env%beta, vcd_env%dcdr_env%lambda) = 1._dp
167 
168  ! Since we do the heavy lifting anyways, we might also calculate the length form APTs here
169  CALL dcdr_build_op_dr(vcd_env%dcdr_env, qs_env)
170  CALL dcdr_response_dr(vcd_env%dcdr_env, p_env, qs_env)
171  CALL apt_dr(qs_env, vcd_env%dcdr_env)
172 
173  ! And with the position perturbation ready, we can calculate the NVP
174  CALL vcd_build_op_dv(vcd_env, qs_env)
175  CALL vcd_response_dv(vcd_env, p_env, qs_env)
176 
177  CALL apt_dv(vcd_env, qs_env)
178  CALL aat_dv(vcd_env, qs_env)
179 
180  IF (vcd_env%do_mfp) THEN
181  ! Since we came so far, we might as well calculate the MFP AATs
182  ! If we use a distributed origin we need to compute the MFP response again for each
183  ! atom, because the reference point changes.
184  IF (.NOT. mfp_is_done .OR. mfp_repeat) THEN
185  DO i = 1, 3
186  IF (vcd_env%origin_dependent_op_mfp) THEN
187  cpwarn("Using the origin dependent MFP operator")
188  CALL mfp_build_operator_gauge_dependent(vcd_env, qs_env, i)
189  ELSE
190  CALL mfp_build_operator_gauge_independent(vcd_env, qs_env, i)
191  END IF
192  CALL mfp_response(vcd_env, p_env, qs_env, i)
193  END DO
194  mfp_is_done = .true.
195  END IF
196 
197  CALL mfp_aat(vcd_env, qs_env)
198  END IF
199  END DO ! beta
200 
201  vcd_env%dcdr_env%apt_total_dcdr(:, :, vcd_env%dcdr_env%lambda) = &
202  vcd_env%dcdr_env%apt_el_dcdr(:, :, vcd_env%dcdr_env%lambda) &
203  + vcd_env%dcdr_env%apt_nuc_dcdr(:, :, vcd_env%dcdr_env%lambda)
204 
205  vcd_env%apt_total_nvpt(:, :, vcd_env%dcdr_env%lambda) = &
206  vcd_env%apt_el_nvpt(:, :, vcd_env%dcdr_env%lambda) + vcd_env%apt_nuc_nvpt(:, :, vcd_env%dcdr_env%lambda)
207 
208  IF (vcd_env%do_mfp) &
209  vcd_env%aat_atom_mfp(:, :, vcd_env%dcdr_env%lambda) = vcd_env%aat_atom_mfp(:, :, vcd_env%dcdr_env%lambda)*4._dp
210 
211  END DO !lambda
212 
213  CALL vcd_print(vcd_env, qs_env)
214  CALL vcd_env_cleanup(qs_env, vcd_env)
215 
216  END SUBROUTINE vcd_linres
217 
218 ! **************************************************************************************************
219 !> \brief Calculates the derivatives of the MO coefficients dC/dR^lambda_beta
220 !> wrt to nuclear coordinates. The derivative is index by `beta`, the
221 !> electric dipole operator by `alpha`.
222 !> Also calculates the APT
223 !> P^lambda_alpha,beta = d< mu_alpha >/dR^lambda_beta
224 !> and calculates the sum rules for the APT elements.
225 !> \param qs_env ...
226 !> \param p_env ...
227 ! **************************************************************************************************
228  SUBROUTINE dcdr_linres(qs_env, p_env)
229  TYPE(qs_environment_type), POINTER :: qs_env
230  TYPE(qs_p_env_type) :: p_env
231 
232  INTEGER :: beta, latom
233  TYPE(dcdr_env_type) :: dcdr_env
234 
235  CALL cite_reference(ditler2021)
236  CALL dcdr_env_init(dcdr_env, qs_env)
237  DO latom = 1, SIZE(dcdr_env%list_of_atoms)
238  dcdr_env%lambda = dcdr_env%list_of_atoms(latom)
239  CALL prepare_per_atom(dcdr_env, qs_env)
240 
241  DO beta = 1, 3 ! in every direction
242  dcdr_env%beta = beta
243  dcdr_env%deltaR(dcdr_env%beta, dcdr_env%lambda) = 1._dp
244 
245  CALL dcdr_build_op_dr(dcdr_env, qs_env)
246  CALL dcdr_response_dr(dcdr_env, p_env, qs_env)
247 
248  IF (.NOT. dcdr_env%localized_psi0) THEN
249  CALL apt_dr(qs_env, dcdr_env)
250  ELSE IF (dcdr_env%localized_psi0) THEN
251  CALL apt_dr_localization(qs_env, dcdr_env)
252  END IF
253 
254  END DO !beta
255 
256  dcdr_env%apt_total_dcdr(:, :, dcdr_env%lambda) = &
257  dcdr_env%apt_el_dcdr(:, :, dcdr_env%lambda) + dcdr_env%apt_nuc_dcdr(:, :, dcdr_env%lambda)
258  END DO !lambda
259 
260  CALL dcdr_print(dcdr_env, qs_env)
261  CALL dcdr_env_cleanup(qs_env, dcdr_env)
262  END SUBROUTINE dcdr_linres
263 
264 ! **************************************************************************************************
265 !> \brief Driver for the linear response calculatios
266 !> \param force_env ...
267 !> \par History
268 !> 06.2005 created [MI]
269 !> \author MI
270 ! **************************************************************************************************
271  SUBROUTINE linres_calculation(force_env)
272 
273  TYPE(force_env_type), POINTER :: force_env
274 
275  CHARACTER(LEN=*), PARAMETER :: routinen = 'linres_calculation'
276 
277  INTEGER :: handle
278  TYPE(qs_environment_type), POINTER :: qs_env
279 
280  CALL timeset(routinen, handle)
281 
282  NULLIFY (qs_env)
283 
284  cpassert(ASSOCIATED(force_env))
285  cpassert(force_env%ref_count > 0)
286 
287  SELECT CASE (force_env%in_use)
288  CASE (use_qs_force)
289  CALL force_env_get(force_env, qs_env=qs_env)
290  CASE (use_qmmm)
291  qs_env => force_env%qmmm_env%qs_env
292  CASE DEFAULT
293  cpabort("Does not recognize this force_env")
294  END SELECT
295 
296  qs_env%linres_run = .true.
297 
298  CALL linres_calculation_low(qs_env)
299 
300  CALL timestop(handle)
301 
302  END SUBROUTINE linres_calculation
303 
304 ! **************************************************************************************************
305 !> \brief Linear response can be called as run type or as post scf calculation
306 !> Initialize the perturbation environment
307 !> Define which properties is to be calculated
308 !> Start up the optimization of the response density and wfn
309 !> \param qs_env ...
310 !> \par History
311 !> 06.2005 created [MI]
312 !> 02.2013 added polarizability section [SL]
313 !> \author MI
314 ! **************************************************************************************************
315  SUBROUTINE linres_calculation_low(qs_env)
316 
317  TYPE(qs_environment_type), POINTER :: qs_env
318 
319  CHARACTER(LEN=*), PARAMETER :: routinen = 'linres_calculation_low'
320 
321  INTEGER :: handle, iounit
322  LOGICAL :: dcdr_present, epr_present, issc_present, &
323  lr_calculation, nmr_present, &
324  polar_present, vcd_present
325  TYPE(cp_logger_type), POINTER :: logger
326  TYPE(dft_control_type), POINTER :: dft_control
327  TYPE(linres_control_type), POINTER :: linres_control
328  TYPE(qs_p_env_type) :: p_env
329  TYPE(section_vals_type), POINTER :: lr_section, prop_section
330 
331  CALL timeset(routinen, handle)
332 
333  lr_calculation = .false.
334  nmr_present = .false.
335  epr_present = .false.
336  issc_present = .false.
337  polar_present = .false.
338  dcdr_present = .false.
339 
340  NULLIFY (dft_control, linres_control, logger, prop_section, lr_section)
341  logger => cp_get_default_logger()
342 
343  lr_section => section_vals_get_subs_vals(qs_env%input, "PROPERTIES%LINRES")
344  CALL section_vals_get(lr_section, explicit=lr_calculation)
345 
346  IF (lr_calculation) THEN
347  CALL linres_init(lr_section, p_env, qs_env)
348  iounit = cp_print_key_unit_nr(logger, lr_section, "PRINT%PROGRAM_RUN_INFO", &
349  extension=".linresLog")
350  CALL get_qs_env(qs_env=qs_env, dft_control=dft_control, &
351  linres_control=linres_control)
352 
353  ! The type of perturbation has not been defined yet
354  linres_control%property = lr_none
355 
356  ! We do NMR or EPR, then compute the current response
357  prop_section => section_vals_get_subs_vals(lr_section, "NMR")
358  CALL section_vals_get(prop_section, explicit=nmr_present)
359  prop_section => section_vals_get_subs_vals(lr_section, "EPR")
360  CALL section_vals_get(prop_section, explicit=epr_present)
361 
362  IF (nmr_present .OR. epr_present) THEN
363  CALL nmr_epr_linres(linres_control, qs_env, p_env, dft_control, &
364  nmr_present, epr_present, iounit)
365  END IF
366 
367  ! We do the indirect spin-spin coupling calculation
368  prop_section => section_vals_get_subs_vals(lr_section, "SPINSPIN")
369  CALL section_vals_get(prop_section, explicit=issc_present)
370 
371  IF (issc_present) THEN
372  CALL issc_linres(linres_control, qs_env, p_env, dft_control)
373  END IF
374 
375  ! We do the polarizability calculation
376  prop_section => section_vals_get_subs_vals(lr_section, "POLAR")
377  CALL section_vals_get(prop_section, explicit=polar_present)
378  IF (polar_present) THEN
379  CALL polar_linres(qs_env, p_env)
380  END IF
381 
382  ! Nuclear Position Perturbation
383  prop_section => section_vals_get_subs_vals(lr_section, "dcdr")
384  CALL section_vals_get(prop_section, explicit=dcdr_present)
385 
386  IF (dcdr_present) THEN
387  CALL dcdr_linres(qs_env, p_env)
388  END IF
389 
390  ! VCD
391  prop_section => section_vals_get_subs_vals(lr_section, "VCD")
392  CALL section_vals_get(prop_section, explicit=vcd_present)
393 
394  IF (vcd_present) THEN
395  CALL vcd_linres(qs_env, p_env)
396  END IF
397 
398  ! Other possible LR calculations can be introduced here
399 
400  CALL p_env_release(p_env)
401 
402  IF (iounit > 0) THEN
403  WRITE (unit=iounit, fmt="(/,T2,A,/,T25,A,/,T2,A,/)") &
404  repeat("=", 79), &
405  "ENDED LINRES CALCULATION", &
406  repeat("=", 79)
407  END IF
408  CALL cp_print_key_finished_output(iounit, logger, lr_section, &
409  "PRINT%PROGRAM_RUN_INFO")
410  END IF
411 
412  CALL timestop(handle)
413 
414  END SUBROUTINE linres_calculation_low
415 
416 ! **************************************************************************************************
417 !> \brief Initialize some general settings like the p_env
418 !> Localize the psi0 if required
419 !> \param lr_section ...
420 !> \param p_env ...
421 !> \param qs_env ...
422 !> \par History
423 !> 06.2005 created [MI]
424 !> \author MI
425 !> \note
426 !> - The localization should probably be always for all the occupied states
427 ! **************************************************************************************************
428  SUBROUTINE linres_init(lr_section, p_env, qs_env)
429 
430  TYPE(section_vals_type), POINTER :: lr_section
431  TYPE(qs_p_env_type), INTENT(OUT) :: p_env
432  TYPE(qs_environment_type), POINTER :: qs_env
433 
434  INTEGER :: iounit, ispin
435  LOGICAL :: do_it
436  TYPE(cp_logger_type), POINTER :: logger
437  TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_ks, rho_ao
438  TYPE(dft_control_type), POINTER :: dft_control
439  TYPE(linres_control_type), POINTER :: linres_control
440  TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
441  TYPE(qs_rho_type), POINTER :: rho
442  TYPE(section_vals_type), POINTER :: loc_section
443 
444  NULLIFY (logger)
445  logger => cp_get_default_logger()
446  iounit = cp_print_key_unit_nr(logger, lr_section, "PRINT%PROGRAM_RUN_INFO", &
447  extension=".linresLog")
448  NULLIFY (dft_control, linres_control, loc_section, rho, mos, matrix_ks, rho_ao)
449 
450  ALLOCATE (linres_control)
451  CALL set_qs_env(qs_env=qs_env, linres_control=linres_control)
452  CALL get_qs_env(qs_env=qs_env, &
453  dft_control=dft_control, matrix_ks=matrix_ks, mos=mos, rho=rho)
454  CALL qs_rho_get(rho, rho_ao=rho_ao)
455 
456  ! Localized Psi0 are required when the position operator has to be defined (nmr)
457  loc_section => section_vals_get_subs_vals(lr_section, "LOCALIZE")
458  CALL section_vals_val_get(loc_section, "_SECTION_PARAMETERS_", &
459  l_val=linres_control%localized_psi0)
460  IF (linres_control%localized_psi0) THEN
461  IF (iounit > 0) THEN
462  WRITE (unit=iounit, fmt="(/,T3,A,A)") &
463  "Localization of ground state orbitals", &
464  " before starting linear response calculation"
465  END IF
466 
467  CALL linres_localize(qs_env, linres_control, dft_control%nspins)
468 
469  DO ispin = 1, dft_control%nspins
470  CALL calculate_density_matrix(mos(ispin), rho_ao(ispin)%matrix)
471  END DO
472  ! ** update qs_env%rho
473  CALL qs_rho_update_rho(rho, qs_env=qs_env)
474  END IF
475 
476  CALL section_vals_val_get(lr_section, "RESTART", l_val=linres_control%linres_restart)
477  CALL section_vals_val_get(lr_section, "MAX_ITER", i_val=linres_control%max_iter)
478  CALL section_vals_val_get(lr_section, "EPS", r_val=linres_control%eps)
479  CALL section_vals_val_get(lr_section, "EPS_FILTER", r_val=linres_control%eps_filter)
480  CALL section_vals_val_get(lr_section, "RESTART_EVERY", i_val=linres_control%restart_every)
481  CALL section_vals_val_get(lr_section, "PRECONDITIONER", i_val=linres_control%preconditioner_type)
482  CALL section_vals_val_get(lr_section, "ENERGY_GAP", r_val=linres_control%energy_gap)
483 
484  IF (iounit > 0) THEN
485  WRITE (unit=iounit, fmt="(/,T2,A,/,T25,A,/,T2,A,/)") &
486  repeat("=", 79), &
487  "START LINRES CALCULATION", &
488  repeat("=", 79)
489 
490  WRITE (unit=iounit, fmt="(T2,A)") &
491  "LINRES| Properties to be calculated:"
492  CALL section_vals_val_get(lr_section, "NMR%_SECTION_PARAMETERS_", l_val=do_it)
493  IF (do_it) WRITE (unit=iounit, fmt="(T62,A)") "NMR Chemical Shift"
494  CALL section_vals_val_get(lr_section, "EPR%_SECTION_PARAMETERS_", l_val=do_it)
495  IF (do_it) WRITE (unit=iounit, fmt="(T68,A)") "EPR g Tensor"
496  CALL section_vals_val_get(lr_section, "SPINSPIN%_SECTION_PARAMETERS_", l_val=do_it)
497  IF (do_it) WRITE (unit=iounit, fmt="(T43,A)") "Indirect spin-spin coupling constants"
498  CALL section_vals_val_get(lr_section, "POLAR%_SECTION_PARAMETERS_", l_val=do_it)
499  IF (do_it) WRITE (unit=iounit, fmt="(T57,A)") "Electric Polarizability"
500 
501  IF (linres_control%localized_psi0) WRITE (unit=iounit, fmt="(T2,A,T65,A)") &
502  "LINRES|", " LOCALIZED PSI0"
503 
504  WRITE (unit=iounit, fmt="(T2,A,T60,A)") &
505  "LINRES| Optimization algorithm", " Conjugate Gradients"
506 
507  SELECT CASE (linres_control%preconditioner_type)
508  CASE (ot_precond_none)
509  WRITE (unit=iounit, fmt="(T2,A,T60,A)") &
510  "LINRES| Preconditioner", " NONE"
512  WRITE (unit=iounit, fmt="(T2,A,T60,A)") &
513  "LINRES| Preconditioner", " FULL_SINGLE"
515  WRITE (unit=iounit, fmt="(T2,A,T60,A)") &
516  "LINRES| Preconditioner", " FULL_KINETIC"
517  CASE (ot_precond_s_inverse)
518  WRITE (unit=iounit, fmt="(T2,A,T60,A)") &
519  "LINRES| Preconditioner", " FULL_S_INVERSE"
521  WRITE (unit=iounit, fmt="(T2,A,T60,A)") &
522  "LINRES| Preconditioner", " FULL_SINGLE_INVERSE"
523  CASE (ot_precond_full_all)
524  WRITE (unit=iounit, fmt="(T2,A,T60,A)") &
525  "LINRES| Preconditioner", " FULL_ALL"
526  CASE DEFAULT
527  cpabort("Preconditioner NYI")
528  END SELECT
529 
530  WRITE (unit=iounit, fmt="(T2,A,T72,ES8.1)") &
531  "LINRES| EPS", linres_control%eps
532  WRITE (unit=iounit, fmt="(T2,A,T72,I8)") &
533  "LINRES| MAX_ITER", linres_control%max_iter
534  END IF
535 
536  !------------------!
537  ! create the p_env !
538  !------------------!
539  CALL p_env_create(p_env, qs_env, orthogonal_orbitals=.true., linres_control=linres_control)
540 
541  ! update the m_epsilon matrix
542  CALL p_env_psi0_changed(p_env, qs_env)
543 
544  p_env%new_preconditioner = .true.
545  CALL cp_print_key_finished_output(iounit, logger, lr_section, &
546  "PRINT%PROGRAM_RUN_INFO")
547 
548  END SUBROUTINE linres_init
549 
550 ! **************************************************************************************************
551 !> \brief ...
552 !> \param linres_control ...
553 !> \param qs_env ...
554 !> \param p_env ...
555 !> \param dft_control ...
556 !> \param nmr_present ...
557 !> \param epr_present ...
558 !> \param iounit ...
559 ! **************************************************************************************************
560  SUBROUTINE nmr_epr_linres(linres_control, qs_env, p_env, dft_control, nmr_present, epr_present, iounit)
561 
562  TYPE(linres_control_type), POINTER :: linres_control
563  TYPE(qs_environment_type), POINTER :: qs_env
564  TYPE(qs_p_env_type) :: p_env
565  TYPE(dft_control_type), POINTER :: dft_control
566  LOGICAL :: nmr_present, epr_present
567  INTEGER :: iounit
568 
569  INTEGER :: ib
570  LOGICAL :: do_qmmm
571  TYPE(current_env_type) :: current_env
572  TYPE(epr_env_type) :: epr_env
573  TYPE(nmr_env_type) :: nmr_env
574 
575  linres_control%property = lr_current
576 
577  CALL cite_reference(weber2009)
578 
579  IF (.NOT. linres_control%localized_psi0) THEN
580  CALL cp_abort(__location__, &
581  "Are you sure that you want to calculate the chemical "// &
582  "shift without localized psi0?")
583  CALL linres_localize(qs_env, linres_control, &
584  dft_control%nspins, centers_only=.true.)
585  END IF
586  IF (dft_control%nspins /= 2 .AND. epr_present) THEN
587  cpabort("LSD is needed to perform a g tensor calculation!")
588  END IF
589  !
590  !Initialize the current environment
591  do_qmmm = .false.
592  IF (qs_env%qmmm) do_qmmm = .true.
593  current_env%do_qmmm = do_qmmm
594  !current_env%prop='nmr'
595  CALL current_env_init(current_env, qs_env)
596  CALL current_operators(current_env, qs_env)
597  CALL current_response(current_env, p_env, qs_env)
598  !
599  IF (current_env%all_pert_op_done) THEN
600  !Initialize the nmr environment
601  IF (nmr_present) THEN
602  CALL nmr_env_init(nmr_env, qs_env)
603  END IF
604  !
605  !Initialize the epr environment
606  IF (epr_present) THEN
607  CALL epr_env_init(epr_env, qs_env)
608  CALL epr_g_zke(epr_env, qs_env)
609  CALL epr_nablavks(epr_env, qs_env)
610  END IF
611  !
612  ! Build the rs_gauge if needed
613  !CALL current_set_gauge(current_env,qs_env)
614  !
615  ! Loop over field direction
616  DO ib = 1, 3
617  !
618  ! Build current response and succeptibility
619  CALL current_build_current(current_env, qs_env, ib)
620  CALL current_build_chi(current_env, qs_env, ib)
621  !
622  ! Compute NMR shift
623  IF (nmr_present) THEN
624  CALL nmr_shift(nmr_env, current_env, qs_env, ib)
625  END IF
626  !
627  ! Compute EPR
628  IF (epr_present) THEN
629  CALL epr_ind_magnetic_field(epr_env, current_env, qs_env, ib)
630  CALL epr_g_so(epr_env, current_env, qs_env, ib)
631  CALL epr_g_soo(epr_env, current_env, qs_env, ib)
632  END IF
633  END DO
634  !
635  ! Finalized the nmr environment
636  IF (nmr_present) THEN
637  CALL nmr_shift_print(nmr_env, current_env, qs_env)
638  CALL nmr_env_cleanup(nmr_env)
639  END IF
640  !
641  ! Finalized the epr environment
642  IF (epr_present) THEN
643  CALL epr_g_print(epr_env, qs_env)
644  CALL epr_env_cleanup(epr_env)
645  END IF
646  !
647  ELSE
648  IF (iounit > 0) THEN
649  WRITE (iounit, "(T10,A,/T20,A,/)") &
650  "CURRENT: Not all responses to perturbation operators could be calculated.", &
651  " Hence: NO nmr and NO epr possible."
652  END IF
653  END IF
654  ! Finalized the current environment
655  CALL current_env_cleanup(current_env)
656 
657  END SUBROUTINE nmr_epr_linres
658 
659 ! **************************************************************************************************
660 !> \brief ...
661 !> \param linres_control ...
662 !> \param qs_env ...
663 !> \param p_env ...
664 !> \param dft_control ...
665 ! **************************************************************************************************
666  SUBROUTINE issc_linres(linres_control, qs_env, p_env, dft_control)
667 
668  TYPE(linres_control_type), POINTER :: linres_control
669  TYPE(qs_environment_type), POINTER :: qs_env
670  TYPE(qs_p_env_type) :: p_env
671  TYPE(dft_control_type), POINTER :: dft_control
672 
673  INTEGER :: iatom
674  LOGICAL :: do_qmmm
675  TYPE(current_env_type) :: current_env
676  TYPE(issc_env_type) :: issc_env
677 
678  linres_control%property = lr_current
679  IF (.NOT. linres_control%localized_psi0) THEN
680  CALL cp_abort(__location__, &
681  "Are you sure that you want to calculate the chemical "// &
682  "shift without localized psi0?")
683  CALL linres_localize(qs_env, linres_control, &
684  dft_control%nspins, centers_only=.true.)
685  END IF
686  !
687  !Initialize the current environment
688  do_qmmm = .false.
689  IF (qs_env%qmmm) do_qmmm = .true.
690  current_env%do_qmmm = do_qmmm
691  !current_env%prop='issc'
692  !CALL current_env_init(current_env,qs_env)
693  !CALL current_response(current_env,p_env,qs_env)
694  !
695  !Initialize the issc environment
696  CALL issc_env_init(issc_env, qs_env)
697  !
698  ! Loop over atoms
699  DO iatom = 1, issc_env%issc_natms
700  CALL issc_operators(issc_env, qs_env, iatom)
701  CALL issc_response(issc_env, p_env, qs_env)
702  CALL issc_issc(issc_env, qs_env, iatom)
703  END DO
704  !
705  ! Finalized the issc environment
706  CALL issc_print(issc_env, qs_env)
707  CALL issc_env_cleanup(issc_env)
708 
709  END SUBROUTINE issc_linres
710 
711 ! **************************************************************************************************
712 !> \brief ...
713 !> \param qs_env ...
714 !> \param p_env ...
715 !> \par History
716 !> 06.2018 polar_env integrated into qs_env (MK)
717 ! **************************************************************************************************
718  SUBROUTINE polar_linres(qs_env, p_env)
719 
720  TYPE(qs_environment_type), POINTER :: qs_env
721  TYPE(qs_p_env_type) :: p_env
722 
723  CALL polar_env_init(qs_env)
724  CALL polar_operators(qs_env)
725  CALL polar_response(p_env, qs_env)
726  CALL polar_polar(qs_env)
727  CALL polar_print(qs_env)
728 
729  END SUBROUTINE polar_linres
730 
731 END MODULE qs_linres_module
collects all references to literature in CP2K as new algorithms / method are included from literature...
Definition: bibliography.F:28
integer, save, public weber2009
Definition: bibliography.F:43
integer, save, public ditler2021
Definition: bibliography.F:43
integer, save, public ditler2022
Definition: bibliography.F:43
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
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)
...
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,...
Interface for the force calculations.
integer, parameter, public use_qmmm
recursive subroutine, public force_env_get(force_env, in_use, fist_env, qs_env, meta_env, fp_env, subsys, para_env, potential_energy, additional_potential, kinetic_energy, harmonic_shell, kinetic_shell, cell, sub_force_env, qmmm_env, qmmmx_env, eip_env, pwdft_env, globenv, input, force_env_section, method_name_id, root_section, mixed_env, nnp_env, embed_env)
returns various attributes about the force environment
integer, parameter, public use_qs_force
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public lr_none
integer, parameter, public ot_precond_full_kinetic
integer, parameter, public ot_precond_full_single
integer, parameter, public ot_precond_none
integer, parameter, public ot_precond_full_single_inverse
integer, parameter, public lr_current
integer, parameter, public ot_precond_s_inverse
integer, parameter, public ot_precond_full_all
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
Calculate the derivatives of the MO coefficients wrt nuclear coordinates.
Definition: qs_dcdr_utils.F:13
subroutine, public dcdr_env_cleanup(qs_env, dcdr_env)
Deallocate the dcdr environment.
subroutine, public dcdr_print(dcdr_env, qs_env)
Print the APT and sum rules.
subroutine, public dcdr_env_init(dcdr_env, qs_env)
Initialize the dcdr environment.
Calculate the derivatives of the MO coefficients wrt nuclear coordinates.
Definition: qs_dcdr.F:13
subroutine, public prepare_per_atom(dcdr_env, qs_env)
Prepare the environment for a choice of lambda.
Definition: qs_dcdr.F:88
subroutine, public apt_dr_localization(qs_env, dcdr_env)
Calculate atomic polar tensor using the localized dipole operator.
Definition: qs_dcdr.F:472
subroutine, public apt_dr(qs_env, dcdr_env)
Calculate atomic polar tensor.
Definition: qs_dcdr.F:363
subroutine, public dcdr_response_dr(dcdr_env, p_env, qs_env)
Get the dC/dR by solving the Sternheimer equation, using the op_dR matrix.
Definition: qs_dcdr.F:259
subroutine, public dcdr_build_op_dr(dcdr_env, qs_env)
Build the operator for the position perturbation.
Definition: qs_dcdr.F:189
collects routines that calculate density matrices
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.
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.
Chemical shift calculation by dfpt Initialization of the nmr_env, creation of the special neighbor li...
subroutine, public current_response(current_env, p_env, qs_env)
...
subroutine, public current_env_init(current_env, qs_env)
...
subroutine, public current_env_cleanup(current_env)
...
given the response wavefunctions obtained by the application of the (rxp), p, and ((dk-dl)xp) operato...
subroutine, public current_build_current(current_env, qs_env, iB)
First calculate the density matrixes, for each component of the current they are 3 because of the r d...
subroutine, public current_build_chi(current_env, qs_env, iB)
...
Calculates Nabla V_KS (local part if PSP) on the different grids.
subroutine, public epr_nablavks(epr_env, qs_env)
Evaluates Nabla V_KS on the grids.
subroutine, public epr_g_zke(epr_env, qs_env)
Calculate zke part of the g tensor.
subroutine, public epr_g_so(epr_env, current_env, qs_env, iB)
Calculates g_so.
subroutine, public epr_g_soo(epr_env, current_env, qs_env, iB)
Calculates g_soo (soft part only for now)
subroutine, public epr_ind_magnetic_field(epr_env, current_env, qs_env, iB)
...
subroutine, public epr_g_print(epr_env, qs_env)
Prints the g tensor.
g tensor calculation by dfpt Initialization of the epr_env, creation of the special neighbor lists Pe...
subroutine, public epr_env_cleanup(epr_env)
Deallocate the epr environment.
subroutine, public epr_env_init(epr_env, qs_env)
Initialize the epr environment.
Chemical shift calculation by dfpt Initialization of the issc_env, creation of the special neighbor l...
subroutine, public issc_issc(issc_env, qs_env, iatom)
...
subroutine, public issc_response(issc_env, p_env, qs_env)
Initialize the issc environment.
subroutine, public issc_env_cleanup(issc_env)
Deallocate the issc environment.
subroutine, public issc_env_init(issc_env, qs_env)
Initialize the issc environment.
subroutine, public issc_print(issc_env, qs_env)
...
localize wavefunctions linear response scf
subroutine, public linres_localize(qs_env, linres_control, nspins, centers_only)
Find the centers and spreads of the wfn, if required apply a localization algorithm.
Contains the setup for the calculation of properties by linear response by the application of second ...
subroutine, public linres_calculation(force_env)
Driver for the linear response calculatios.
subroutine, public linres_calculation_low(qs_env)
Linear response can be called as run type or as post scf calculation Initialize the perturbation envi...
from the response current density calculates the shift tensor and the susceptibility
subroutine, public nmr_shift(nmr_env, current_env, qs_env, iB)
...
subroutine, public nmr_shift_print(nmr_env, current_env, qs_env)
Shielding tensor and Chi are printed into a file if required from input It is possible to print only ...
Chemical shift calculation by dfpt Initialization of the nmr_env, creation of the special neighbor li...
subroutine, public nmr_env_cleanup(nmr_env)
Deallocate the nmr environment.
subroutine, public nmr_env_init(nmr_env, qs_env)
Initialize the nmr environment.
Calculate the operators p rxp and D needed in the optimization of the different contribution of the f...
Definition: qs_linres_op.F:18
subroutine, public current_operators(current_env, qs_env)
Calculate the first order hamiltonian applied to the ao and then apply them to the ground state orbit...
Definition: qs_linres_op.F:116
subroutine, public polar_operators(qs_env)
Calculate the dipole operator in the AO basis and its derivative wrt to MOs.
Definition: qs_linres_op.F:639
subroutine, public issc_operators(issc_env, qs_env, iatom)
...
Definition: qs_linres_op.F:441
Polarizability calculation by dfpt Initialization of the polar_env, Perturbation Hamiltonian by appli...
subroutine, public polar_polar(qs_env)
...
subroutine, public polar_print(qs_env)
Print information related to the polarisability tensor.
subroutine, public polar_env_init(qs_env)
Initialize the polar environment.
subroutine, public polar_response(p_env, qs_env)
Calculate the polarisability tensor using response theory.
Type definitiona for linear response calculations.
Definition: qs_mfp.F:7
subroutine, public mfp_build_operator_gauge_independent(vcd_env, qs_env, alpha)
...
Definition: qs_mfp.F:848
subroutine, public mfp_response(vcd_env, p_env, qs_env, alpha)
Get the dC/dB using the vcd_envop_dB.
Definition: qs_mfp.F:1150
subroutine, public mfp_build_operator_gauge_dependent(vcd_env, qs_env, alpha)
...
Definition: qs_mfp.F:522
subroutine, public mfp_aat(vcd_env, qs_env)
...
Definition: qs_mfp.F:90
Definition and initialisation of the mo data type.
Definition: qs_mo_types.F:22
Utility functions for the perturbation calculations.
subroutine, public p_env_psi0_changed(p_env, qs_env)
To be called after the value of psi0 has changed. Recalculates the quantities S_psi0 and m_epsilon.
subroutine, public p_env_create(p_env, qs_env, p1_option, p1_admm_option, orthogonal_orbitals, linres_control)
allocates and initializes the perturbation environment (no setup)
basis types for the calculation of the perturbation of density theory.
subroutine, public p_env_release(p_env)
relases the given p_env (see doc/ReferenceCounting.html)
methods of the rho structure (defined in qs_rho_types)
subroutine, public qs_rho_update_rho(rho_struct, qs_env, rho_xc_external, local_rho_set, task_list_external, task_list_external_soft, pw_env_external, para_env_external)
updates rho_r and rho_g to the rhorho_ao. if use_kinetic_energy_density also computes tau_r and tau_g...
superstucture that hold various representations of the density and keeps track of which ones are vali...
Definition: qs_rho_types.F:18
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...
Definition: qs_rho_types.F:229
subroutine, public vcd_print(vcd_env, qs_env)
Print the APTs, AATs, and sum rules.
Definition: qs_vcd_utils.F:757
subroutine, public vcd_env_cleanup(qs_env, vcd_env)
Deallocate the vcd environment.
Definition: qs_vcd_utils.F:439
subroutine, public vcd_env_init(vcd_env, qs_env)
Initialize the vcd environment.
Definition: qs_vcd_utils.F:96
Definition: qs_vcd.F:7
subroutine, public vcd_build_op_dv(vcd_env, qs_env)
What we are building here is the operator for the NVPT response: H0 * C1 - S0 * E0 * C1 = - op_dV lin...
Definition: qs_vcd.F:955
subroutine, public vcd_response_dv(vcd_env, p_env, qs_env)
Get the dC/dV using the vcd_envop_dV.
Definition: qs_vcd.F:1038
subroutine, public aat_dv(vcd_env, qs_env)
Compute I_{alpha beta}^lambda = d/dV^lambda_beta <m_alpha> = d/dV^lambda_beta < r x.
Definition: qs_vcd.F:82
subroutine, public apt_dv(vcd_env, qs_env)
Compute E_{alpha beta}^lambda = d/dV^lambda_beta <\mu_alpha> = d/dV^lambda_beta <.
Definition: qs_vcd.F:593
subroutine, public prepare_per_atom_vcd(vcd_env, qs_env)
Initialize the matrices for the NVPT calculation.
Definition: qs_vcd.F:830