(git:ed6f26b)
Loading...
Searching...
No Matches
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-2025 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,&
21 weber2009,&
22 cite_reference
24 USE cp_dbcsr_api, ONLY: dbcsr_p_type
31 use_qmmm,&
33 USE input_constants, ONLY: lr_current,&
34 lr_none,&
45 USE kinds, ONLY: dp
46 USE qs_dcdr, ONLY: apt_dr,&
65 epr_g_so,&
66 epr_g_soo,&
67 epr_g_zke,&
73 issc_issc,&
90 USE qs_linres_types, ONLY: &
93 USE qs_mfp, ONLY: mfp_aat,&
97 USE qs_mo_types, ONLY: mo_set_type
100 USE qs_p_env_types, ONLY: p_env_release,&
103 USE qs_rho_types, ONLY: qs_rho_get,&
105 USE qs_vcd, ONLY: aat_dv,&
106 apt_dv,&
110 USE qs_vcd_utils, ONLY: vcd_env_cleanup,&
113#include "./base/base_uses.f90"
114
115 IMPLICIT NONE
116
117 PRIVATE
119
120 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_linres_module'
121
122CONTAINS
123! *****************************************************************************
124!> \brief Calculates the derivatives of the MO coefficients dC/dV^lambda_beta
125!> wrt to nuclear velocities. The derivative is indexed by `beta`, the
126!> electric dipole operator by `alpha`.
127!> Calculates the APT and AAT in velocity form
128!> P^lambda_alpha,beta = d< mu_alpha >/dV^lambda_beta
129!> M^lambda_alpha,beta = d< m_alpha >/dV^lambda_beta
130!> \param qs_env ...
131!> \param p_env ...
132!> \author Edward Ditler
133! **************************************************************************************************
134 SUBROUTINE vcd_linres(qs_env, p_env)
135 TYPE(qs_environment_type), POINTER :: qs_env
136 TYPE(qs_p_env_type) :: p_env
137
138 INTEGER :: beta, i, latom
139 LOGICAL :: mfp_is_done, mfp_repeat
140 TYPE(vcd_env_type) :: vcd_env
141
142 CALL cite_reference(ditler2022)
143
144 ! We need the position perturbation for the velocity perturbation operator
145 CALL vcd_env_init(vcd_env, qs_env)
146
147 mfp_repeat = vcd_env%distributed_origin
148 mfp_is_done = .false.
149
150 qs_env%linres_control%linres_restart = .true.
151
152 ! Iterate over the list of atoms for which we want to calculate the APTs/AATs
153 ! default is all atoms.
154 DO latom = 1, SIZE(vcd_env%dcdr_env%list_of_atoms)
155 vcd_env%dcdr_env%lambda = vcd_env%dcdr_env%list_of_atoms(latom)
156
157 CALL prepare_per_atom(vcd_env%dcdr_env, qs_env)
158 CALL prepare_per_atom_vcd(vcd_env, qs_env)
159
160 DO beta = 1, 3 ! in every direction
161
162 vcd_env%dcdr_env%beta = beta
163 vcd_env%dcdr_env%deltaR(vcd_env%dcdr_env%beta, vcd_env%dcdr_env%lambda) = 1._dp
164
165 ! Since we do the heavy lifting anyways, we might also calculate the length form APTs here
166 CALL dcdr_build_op_dr(vcd_env%dcdr_env, qs_env)
167 CALL dcdr_response_dr(vcd_env%dcdr_env, p_env, qs_env)
168 CALL apt_dr(qs_env, vcd_env%dcdr_env)
169
170 ! And with the position perturbation ready, we can calculate the NVP
171 CALL vcd_build_op_dv(vcd_env, qs_env)
172 CALL vcd_response_dv(vcd_env, p_env, qs_env)
173
174 CALL apt_dv(vcd_env, qs_env)
175 CALL aat_dv(vcd_env, qs_env)
176
177 IF (vcd_env%do_mfp) THEN
178 ! Since we came so far, we might as well calculate the MFP AATs
179 ! If we use a distributed origin we need to compute the MFP response again for each
180 ! atom, because the reference point changes.
181 IF (.NOT. mfp_is_done .OR. mfp_repeat) THEN
182 DO i = 1, 3
183 IF (vcd_env%origin_dependent_op_mfp) THEN
184 cpwarn("Using the origin dependent MFP operator")
185 CALL mfp_build_operator_gauge_dependent(vcd_env, qs_env, i)
186 ELSE
187 CALL mfp_build_operator_gauge_independent(vcd_env, qs_env, i)
188 END IF
189 CALL mfp_response(vcd_env, p_env, qs_env, i)
190 END DO
191 mfp_is_done = .true.
192 END IF
193
194 CALL mfp_aat(vcd_env, qs_env)
195 END IF
196 END DO ! beta
197
198 vcd_env%dcdr_env%apt_total_dcdr(:, :, vcd_env%dcdr_env%lambda) = &
199 vcd_env%dcdr_env%apt_el_dcdr(:, :, vcd_env%dcdr_env%lambda) &
200 + vcd_env%dcdr_env%apt_nuc_dcdr(:, :, vcd_env%dcdr_env%lambda)
201
202 vcd_env%apt_total_nvpt(:, :, vcd_env%dcdr_env%lambda) = &
203 vcd_env%apt_el_nvpt(:, :, vcd_env%dcdr_env%lambda) + vcd_env%apt_nuc_nvpt(:, :, vcd_env%dcdr_env%lambda)
204
205 IF (vcd_env%do_mfp) &
206 vcd_env%aat_atom_mfp(:, :, vcd_env%dcdr_env%lambda) = vcd_env%aat_atom_mfp(:, :, vcd_env%dcdr_env%lambda)*4._dp
207
208 END DO !lambda
209
210 CALL vcd_print(vcd_env, qs_env)
211 CALL vcd_env_cleanup(qs_env, vcd_env)
212
213 END SUBROUTINE vcd_linres
214
215! **************************************************************************************************
216!> \brief Calculates the derivatives of the MO coefficients dC/dR^lambda_beta
217!> wrt to nuclear coordinates. The derivative is index by `beta`, the
218!> electric dipole operator by `alpha`.
219!> Also calculates the APT
220!> P^lambda_alpha,beta = d< mu_alpha >/dR^lambda_beta
221!> and calculates the sum rules for the APT elements.
222!> \param qs_env ...
223!> \param p_env ...
224! **************************************************************************************************
225 SUBROUTINE dcdr_linres(qs_env, p_env)
226 TYPE(qs_environment_type), POINTER :: qs_env
227 TYPE(qs_p_env_type) :: p_env
228
229 INTEGER :: beta, latom
230 TYPE(dcdr_env_type) :: dcdr_env
231 TYPE(polar_env_type), POINTER :: polar_env
232
233 CALL cite_reference(ditler2021)
234 CALL dcdr_env_init(dcdr_env, qs_env)
235
236 IF (.NOT. dcdr_env%z_matrix_method) THEN
237
238 DO latom = 1, SIZE(dcdr_env%list_of_atoms)
239 dcdr_env%lambda = dcdr_env%list_of_atoms(latom)
240 CALL prepare_per_atom(dcdr_env, qs_env)
241
242 DO beta = 1, 3 ! in every direction
243 dcdr_env%beta = beta
244 dcdr_env%deltaR(dcdr_env%beta, dcdr_env%lambda) = 1._dp
245
246 CALL dcdr_build_op_dr(dcdr_env, qs_env)
247 CALL dcdr_response_dr(dcdr_env, p_env, qs_env)
248
249 IF (.NOT. dcdr_env%localized_psi0) THEN
250 CALL apt_dr(qs_env, dcdr_env)
251 ELSE IF (dcdr_env%localized_psi0) THEN
252 CALL apt_dr_localization(qs_env, dcdr_env)
253 END IF
254
255 END DO !beta
256
257 dcdr_env%apt_total_dcdr(:, :, dcdr_env%lambda) = &
258 dcdr_env%apt_el_dcdr(:, :, dcdr_env%lambda) + dcdr_env%apt_nuc_dcdr(:, :, dcdr_env%lambda)
259 END DO !lambda
260
261 ELSE
262
263 CALL polar_env_init(qs_env)
264 CALL get_qs_env(qs_env=qs_env, polar_env=polar_env)
265 CALL get_polar_env(polar_env=polar_env)
266
267 IF (.NOT. dcdr_env%localized_psi0) THEN
268 CALL polar_operators_local(qs_env)
269 ELSE
270 CALL polar_operators_local_wannier(qs_env, dcdr_env)
271 END IF
272
273 polar_env%do_periodic = .false.
274 CALL polar_response(p_env, qs_env)
275
276 DO latom = 1, SIZE(dcdr_env%list_of_atoms)
277 dcdr_env%lambda = dcdr_env%list_of_atoms(latom)
278 CALL prepare_per_atom(dcdr_env, qs_env)
279
280 DO beta = 1, 3 ! in every direction
281 dcdr_env%beta = beta
282 dcdr_env%deltaR(dcdr_env%beta, dcdr_env%lambda) = 1._dp
283
284 CALL dcdr_build_op_dr(dcdr_env, qs_env)
285 IF (.NOT. dcdr_env%localized_psi0) THEN
286 CALL apt_dr(qs_env, dcdr_env)
287 ELSE
288 CALL apt_dr_localization(qs_env, dcdr_env)
289 END IF
290 END DO !beta
291
292 dcdr_env%apt_total_dcdr(:, :, dcdr_env%lambda) = &
293 dcdr_env%apt_el_dcdr(:, :, dcdr_env%lambda) + dcdr_env%apt_nuc_dcdr(:, :, dcdr_env%lambda)
294 END DO !lambda
295
296 END IF
297
298 CALL dcdr_print(dcdr_env, qs_env)
299 CALL dcdr_env_cleanup(qs_env, dcdr_env)
300 END SUBROUTINE dcdr_linres
301
302! **************************************************************************************************
303!> \brief Driver for the linear response calculatios
304!> \param force_env ...
305!> \par History
306!> 06.2005 created [MI]
307!> \author MI
308! **************************************************************************************************
309 SUBROUTINE linres_calculation(force_env)
310
311 TYPE(force_env_type), POINTER :: force_env
312
313 CHARACTER(LEN=*), PARAMETER :: routinen = 'linres_calculation'
314
315 INTEGER :: handle
316 TYPE(qs_environment_type), POINTER :: qs_env
317
318 CALL timeset(routinen, handle)
319
320 NULLIFY (qs_env)
321
322 cpassert(ASSOCIATED(force_env))
323 cpassert(force_env%ref_count > 0)
324
325 SELECT CASE (force_env%in_use)
326 CASE (use_qs_force)
327 CALL force_env_get(force_env, qs_env=qs_env)
328 CASE (use_qmmm)
329 qs_env => force_env%qmmm_env%qs_env
330 CASE DEFAULT
331 cpabort("Does not recognize this force_env")
332 END SELECT
333
334 qs_env%linres_run = .true.
335
336 CALL linres_calculation_low(qs_env)
337
338 CALL timestop(handle)
339
340 END SUBROUTINE linres_calculation
341
342! **************************************************************************************************
343!> \brief Linear response can be called as run type or as post scf calculation
344!> Initialize the perturbation environment
345!> Define which properties is to be calculated
346!> Start up the optimization of the response density and wfn
347!> \param qs_env ...
348!> \par History
349!> 06.2005 created [MI]
350!> 02.2013 added polarizability section [SL]
351!> \author MI
352! **************************************************************************************************
353 SUBROUTINE linres_calculation_low(qs_env)
354
355 TYPE(qs_environment_type), POINTER :: qs_env
356
357 CHARACTER(LEN=*), PARAMETER :: routinen = 'linres_calculation_low'
358
359 INTEGER :: every_n_step, handle, iounit
360 LOGICAL :: dcdr_present, epr_present, issc_present, &
361 lr_calculation, nmr_present, &
362 polar_present, vcd_present
363 TYPE(cp_logger_type), POINTER :: logger
364 TYPE(dft_control_type), POINTER :: dft_control
365 TYPE(linres_control_type), POINTER :: linres_control
366 TYPE(qs_p_env_type) :: p_env
367 TYPE(section_vals_type), POINTER :: lr_section, prop_section
368
369 CALL timeset(routinen, handle)
370
371 lr_calculation = .false.
372 nmr_present = .false.
373 epr_present = .false.
374 issc_present = .false.
375 polar_present = .false.
376 dcdr_present = .false.
377
378 NULLIFY (dft_control, linres_control, logger, prop_section, lr_section)
379 logger => cp_get_default_logger()
380
381 lr_section => section_vals_get_subs_vals(qs_env%input, "PROPERTIES%LINRES")
382 CALL section_vals_get(lr_section, explicit=lr_calculation)
383
384 CALL section_vals_val_get(lr_section, "EVERY_N_STEP", i_val=every_n_step)
385
386 IF (lr_calculation .AND. modulo(qs_env%sim_step, every_n_step) == 0) THEN
387 CALL linres_init(lr_section, p_env, qs_env)
388 iounit = cp_print_key_unit_nr(logger, lr_section, "PRINT%PROGRAM_RUN_INFO", &
389 extension=".linresLog")
390 CALL get_qs_env(qs_env=qs_env, dft_control=dft_control, &
391 linres_control=linres_control)
392
393 ! The type of perturbation has not been defined yet
394 linres_control%property = lr_none
395
396 ! We do NMR or EPR, then compute the current response
397 prop_section => section_vals_get_subs_vals(lr_section, "NMR")
398 CALL section_vals_get(prop_section, explicit=nmr_present)
399 prop_section => section_vals_get_subs_vals(lr_section, "EPR")
400 CALL section_vals_get(prop_section, explicit=epr_present)
401
402 IF (nmr_present .OR. epr_present) THEN
403 CALL nmr_epr_linres(linres_control, qs_env, p_env, dft_control, &
404 nmr_present, epr_present, iounit)
405 END IF
406
407 ! We do the indirect spin-spin coupling calculation
408 prop_section => section_vals_get_subs_vals(lr_section, "SPINSPIN")
409 CALL section_vals_get(prop_section, explicit=issc_present)
410
411 IF (issc_present) THEN
412 CALL issc_linres(linres_control, qs_env, p_env, dft_control)
413 END IF
414
415 ! We do the polarizability calculation
416 prop_section => section_vals_get_subs_vals(lr_section, "POLAR")
417 CALL section_vals_get(prop_section, explicit=polar_present)
418 IF (polar_present) THEN
419 CALL polar_linres(qs_env, p_env)
420 END IF
421
422 ! Nuclear Position Perturbation
423 prop_section => section_vals_get_subs_vals(lr_section, "dcdr")
424 CALL section_vals_get(prop_section, explicit=dcdr_present)
425
426 IF (dcdr_present) THEN
427 CALL dcdr_linres(qs_env, p_env)
428 END IF
429
430 ! VCD
431 prop_section => section_vals_get_subs_vals(lr_section, "VCD")
432 CALL section_vals_get(prop_section, explicit=vcd_present)
433
434 IF (vcd_present) THEN
435 CALL vcd_linres(qs_env, p_env)
436 END IF
437
438 ! Other possible LR calculations can be introduced here
439
440 CALL p_env_release(p_env)
441
442 IF (iounit > 0) THEN
443 WRITE (unit=iounit, fmt="(/,T2,A,/,T25,A,/,T2,A,/)") &
444 repeat("=", 79), &
445 "ENDED LINRES CALCULATION", &
446 repeat("=", 79)
447 END IF
448 CALL cp_print_key_finished_output(iounit, logger, lr_section, &
449 "PRINT%PROGRAM_RUN_INFO")
450 END IF
451
452 CALL timestop(handle)
453
454 END SUBROUTINE linres_calculation_low
455
456! **************************************************************************************************
457!> \brief Initialize some general settings like the p_env
458!> Localize the psi0 if required
459!> \param lr_section ...
460!> \param p_env ...
461!> \param qs_env ...
462!> \par History
463!> 06.2005 created [MI]
464!> \author MI
465!> \note
466!> - The localization should probably be always for all the occupied states
467! **************************************************************************************************
468 SUBROUTINE linres_init(lr_section, p_env, qs_env)
469
470 TYPE(section_vals_type), POINTER :: lr_section
471 TYPE(qs_p_env_type), INTENT(OUT) :: p_env
472 TYPE(qs_environment_type), POINTER :: qs_env
473
474 INTEGER :: iounit, ispin
475 LOGICAL :: do_it
476 TYPE(cp_logger_type), POINTER :: logger
477 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_ks, rho_ao
478 TYPE(dft_control_type), POINTER :: dft_control
479 TYPE(linres_control_type), POINTER :: linres_control
480 TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
481 TYPE(qs_rho_type), POINTER :: rho
482 TYPE(section_vals_type), POINTER :: loc_section
483
484 NULLIFY (logger)
485 logger => cp_get_default_logger()
486 iounit = cp_print_key_unit_nr(logger, lr_section, "PRINT%PROGRAM_RUN_INFO", &
487 extension=".linresLog")
488 NULLIFY (dft_control, linres_control, loc_section, rho, mos, matrix_ks, rho_ao)
489
490 ALLOCATE (linres_control)
491 CALL set_qs_env(qs_env=qs_env, linres_control=linres_control)
492 CALL get_qs_env(qs_env=qs_env, &
493 dft_control=dft_control, matrix_ks=matrix_ks, mos=mos, rho=rho)
494 CALL qs_rho_get(rho, rho_ao=rho_ao)
495
496 ! Localized Psi0 are required when the position operator has to be defined (nmr)
497 loc_section => section_vals_get_subs_vals(lr_section, "LOCALIZE")
498 CALL section_vals_val_get(loc_section, "_SECTION_PARAMETERS_", &
499 l_val=linres_control%localized_psi0)
500 IF (linres_control%localized_psi0) THEN
501 IF (iounit > 0) THEN
502 WRITE (unit=iounit, fmt="(/,T3,A,A)") &
503 "Localization of ground state orbitals", &
504 " before starting linear response calculation"
505 END IF
506
507 CALL linres_localize(qs_env, linres_control, dft_control%nspins)
508
509 DO ispin = 1, dft_control%nspins
510 CALL calculate_density_matrix(mos(ispin), rho_ao(ispin)%matrix)
511 END DO
512 ! ** update qs_env%rho
513 CALL qs_rho_update_rho(rho, qs_env=qs_env)
514 END IF
515
516 CALL section_vals_val_get(lr_section, "RESTART", l_val=linres_control%linres_restart)
517 CALL section_vals_val_get(lr_section, "MAX_ITER", i_val=linres_control%max_iter)
518 CALL section_vals_val_get(lr_section, "EPS", r_val=linres_control%eps)
519 CALL section_vals_val_get(lr_section, "EPS_FILTER", r_val=linres_control%eps_filter)
520 CALL section_vals_val_get(lr_section, "RESTART_EVERY", i_val=linres_control%restart_every)
521 CALL section_vals_val_get(lr_section, "PRECONDITIONER", i_val=linres_control%preconditioner_type)
522 CALL section_vals_val_get(lr_section, "ENERGY_GAP", r_val=linres_control%energy_gap)
523
524 IF (iounit > 0) THEN
525 WRITE (unit=iounit, fmt="(/,T2,A,/,T25,A,/,T2,A,/)") &
526 repeat("=", 79), &
527 "START LINRES CALCULATION", &
528 repeat("=", 79)
529
530 WRITE (unit=iounit, fmt="(T2,A)") &
531 "LINRES| Properties to be calculated:"
532 CALL section_vals_val_get(lr_section, "NMR%_SECTION_PARAMETERS_", l_val=do_it)
533 IF (do_it) WRITE (unit=iounit, fmt="(T62,A)") "NMR Chemical Shift"
534 CALL section_vals_val_get(lr_section, "EPR%_SECTION_PARAMETERS_", l_val=do_it)
535 IF (do_it) WRITE (unit=iounit, fmt="(T68,A)") "EPR g Tensor"
536 CALL section_vals_val_get(lr_section, "SPINSPIN%_SECTION_PARAMETERS_", l_val=do_it)
537 IF (do_it) WRITE (unit=iounit, fmt="(T43,A)") "Indirect spin-spin coupling constants"
538 CALL section_vals_val_get(lr_section, "POLAR%_SECTION_PARAMETERS_", l_val=do_it)
539 IF (do_it) WRITE (unit=iounit, fmt="(T57,A)") "Electric Polarizability"
540
541 IF (linres_control%localized_psi0) WRITE (unit=iounit, fmt="(T2,A,T65,A)") &
542 "LINRES|", " LOCALIZED PSI0"
543
544 WRITE (unit=iounit, fmt="(T2,A,T60,A)") &
545 "LINRES| Optimization algorithm", " Conjugate Gradients"
546
547 SELECT CASE (linres_control%preconditioner_type)
548 CASE (ot_precond_none)
549 WRITE (unit=iounit, fmt="(T2,A,T60,A)") &
550 "LINRES| Preconditioner", " NONE"
552 WRITE (unit=iounit, fmt="(T2,A,T60,A)") &
553 "LINRES| Preconditioner", " FULL_SINGLE"
555 WRITE (unit=iounit, fmt="(T2,A,T60,A)") &
556 "LINRES| Preconditioner", " FULL_KINETIC"
558 WRITE (unit=iounit, fmt="(T2,A,T60,A)") &
559 "LINRES| Preconditioner", " FULL_S_INVERSE"
561 WRITE (unit=iounit, fmt="(T2,A,T60,A)") &
562 "LINRES| Preconditioner", " FULL_SINGLE_INVERSE"
564 WRITE (unit=iounit, fmt="(T2,A,T60,A)") &
565 "LINRES| Preconditioner", " FULL_ALL"
566 CASE DEFAULT
567 cpabort("Preconditioner NYI")
568 END SELECT
569
570 WRITE (unit=iounit, fmt="(T2,A,T72,ES8.1)") &
571 "LINRES| EPS", linres_control%eps
572 WRITE (unit=iounit, fmt="(T2,A,T72,I8)") &
573 "LINRES| MAX_ITER", linres_control%max_iter
574 END IF
575
576 !------------------!
577 ! create the p_env !
578 !------------------!
579 CALL p_env_create(p_env, qs_env, orthogonal_orbitals=.true., linres_control=linres_control)
580
581 ! update the m_epsilon matrix
582 CALL p_env_psi0_changed(p_env, qs_env)
583
584 p_env%new_preconditioner = .true.
585 CALL cp_print_key_finished_output(iounit, logger, lr_section, &
586 "PRINT%PROGRAM_RUN_INFO")
587
588 END SUBROUTINE linres_init
589
590! **************************************************************************************************
591!> \brief ...
592!> \param linres_control ...
593!> \param qs_env ...
594!> \param p_env ...
595!> \param dft_control ...
596!> \param nmr_present ...
597!> \param epr_present ...
598!> \param iounit ...
599! **************************************************************************************************
600 SUBROUTINE nmr_epr_linres(linres_control, qs_env, p_env, dft_control, nmr_present, epr_present, iounit)
601
602 TYPE(linres_control_type), POINTER :: linres_control
603 TYPE(qs_environment_type), POINTER :: qs_env
604 TYPE(qs_p_env_type) :: p_env
605 TYPE(dft_control_type), POINTER :: dft_control
606 LOGICAL :: nmr_present, epr_present
607 INTEGER :: iounit
608
609 INTEGER :: ib
610 LOGICAL :: do_qmmm
611 TYPE(current_env_type) :: current_env
612 TYPE(epr_env_type) :: epr_env
613 TYPE(nmr_env_type) :: nmr_env
614
615 linres_control%property = lr_current
616
617 CALL cite_reference(weber2009)
618
619 IF (.NOT. linres_control%localized_psi0) THEN
620 CALL cp_abort(__location__, &
621 "Are you sure that you want to calculate the chemical "// &
622 "shift without localized psi0?")
623 CALL linres_localize(qs_env, linres_control, &
624 dft_control%nspins, centers_only=.true.)
625 END IF
626 IF (dft_control%nspins /= 2 .AND. epr_present) THEN
627 cpabort("LSD is needed to perform a g tensor calculation!")
628 END IF
629 !
630 !Initialize the current environment
631 do_qmmm = .false.
632 IF (qs_env%qmmm) do_qmmm = .true.
633 current_env%do_qmmm = do_qmmm
634 !current_env%prop='nmr'
635 CALL current_env_init(current_env, qs_env)
636 CALL current_operators(current_env, qs_env)
637 CALL current_response(current_env, p_env, qs_env)
638 !
639 IF (current_env%all_pert_op_done) THEN
640 !Initialize the nmr environment
641 IF (nmr_present) THEN
642 CALL nmr_env_init(nmr_env, qs_env)
643 END IF
644 !
645 !Initialize the epr environment
646 IF (epr_present) THEN
647 CALL epr_env_init(epr_env, qs_env)
648 CALL epr_g_zke(epr_env, qs_env)
649 CALL epr_nablavks(epr_env, qs_env)
650 END IF
651 !
652 ! Build the rs_gauge if needed
653 !CALL current_set_gauge(current_env,qs_env)
654 !
655 ! Loop over field direction
656 DO ib = 1, 3
657 !
658 ! Build current response and succeptibility
659 CALL current_build_current(current_env, qs_env, ib)
660 CALL current_build_chi(current_env, qs_env, ib)
661 !
662 ! Compute NMR shift
663 IF (nmr_present) THEN
664 CALL nmr_shift(nmr_env, current_env, qs_env, ib)
665 END IF
666 !
667 ! Compute EPR
668 IF (epr_present) THEN
669 CALL epr_ind_magnetic_field(epr_env, current_env, qs_env, ib)
670 CALL epr_g_so(epr_env, current_env, qs_env, ib)
671 CALL epr_g_soo(epr_env, current_env, qs_env, ib)
672 END IF
673 END DO
674 !
675 ! Finalized the nmr environment
676 IF (nmr_present) THEN
677 CALL nmr_shift_print(nmr_env, current_env, qs_env)
678 CALL nmr_env_cleanup(nmr_env)
679 END IF
680 !
681 ! Finalized the epr environment
682 IF (epr_present) THEN
683 CALL epr_g_print(epr_env, qs_env)
684 CALL epr_env_cleanup(epr_env)
685 END IF
686 !
687 ELSE
688 IF (iounit > 0) THEN
689 WRITE (iounit, "(T10,A,/T20,A,/)") &
690 "CURRENT: Not all responses to perturbation operators could be calculated.", &
691 " Hence: NO nmr and NO epr possible."
692 END IF
693 END IF
694 ! Finalized the current environment
695 CALL current_env_cleanup(current_env)
696
697 END SUBROUTINE nmr_epr_linres
698
699! **************************************************************************************************
700!> \brief ...
701!> \param linres_control ...
702!> \param qs_env ...
703!> \param p_env ...
704!> \param dft_control ...
705! **************************************************************************************************
706 SUBROUTINE issc_linres(linres_control, qs_env, p_env, dft_control)
707
708 TYPE(linres_control_type), POINTER :: linres_control
709 TYPE(qs_environment_type), POINTER :: qs_env
710 TYPE(qs_p_env_type) :: p_env
711 TYPE(dft_control_type), POINTER :: dft_control
712
713 INTEGER :: iatom
714 LOGICAL :: do_qmmm
715 TYPE(current_env_type) :: current_env
716 TYPE(issc_env_type) :: issc_env
717
718 linres_control%property = lr_current
719 IF (.NOT. linres_control%localized_psi0) THEN
720 CALL cp_abort(__location__, &
721 "Are you sure that you want to calculate the chemical "// &
722 "shift without localized psi0?")
723 CALL linres_localize(qs_env, linres_control, &
724 dft_control%nspins, centers_only=.true.)
725 END IF
726 !
727 !Initialize the current environment
728 do_qmmm = .false.
729 IF (qs_env%qmmm) do_qmmm = .true.
730 current_env%do_qmmm = do_qmmm
731 !current_env%prop='issc'
732 !CALL current_env_init(current_env,qs_env)
733 !CALL current_response(current_env,p_env,qs_env)
734 !
735 !Initialize the issc environment
736 CALL issc_env_init(issc_env, qs_env)
737 !
738 ! Loop over atoms
739 DO iatom = 1, issc_env%issc_natms
740 CALL issc_operators(issc_env, qs_env, iatom)
741 CALL issc_response(issc_env, p_env, qs_env)
742 CALL issc_issc(issc_env, qs_env, iatom)
743 END DO
744 !
745 ! Finalized the issc environment
746 CALL issc_print(issc_env, qs_env)
747 CALL issc_env_cleanup(issc_env)
748
749 END SUBROUTINE issc_linres
750
751! **************************************************************************************************
752!> \brief ...
753!> \param qs_env ...
754!> \param p_env ...
755!> \par History
756!> 06.2018 polar_env integrated into qs_env (MK)
757! **************************************************************************************************
758 SUBROUTINE polar_linres(qs_env, p_env)
759
760 TYPE(qs_environment_type), POINTER :: qs_env
761 TYPE(qs_p_env_type) :: p_env
762
763 CALL polar_env_init(qs_env)
764 CALL polar_operators(qs_env)
765 CALL polar_response(p_env, qs_env)
766 CALL polar_polar(qs_env)
767 CALL polar_print(qs_env)
768
769 END SUBROUTINE polar_linres
770
771END MODULE qs_linres_module
static GRID_HOST_DEVICE int modulo(int a, int m)
Equivalent of Fortran's MODULO, which always return a positive number. https://gcc....
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public weber2009
integer, save, public ditler2021
integer, save, public ditler2022
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, ipi_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 do_qmmm
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.
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:90
subroutine, public apt_dr_localization(qs_env, dcdr_env)
Calculate atomic polar tensor using the localized dipole operator.
Definition qs_dcdr.F:519
subroutine, public apt_dr(qs_env, dcdr_env)
Calculate atomic polar tensor.
Definition qs_dcdr.F:376
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:269
subroutine, public dcdr_build_op_dr(dcdr_env, qs_env)
Build the operator for the position perturbation.
Definition qs_dcdr.F:194
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_pp, 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, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs)
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, harris_env, gcp_env, mp2_env, bs_env, kg_env, force, kpoints, wanniercentres, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs)
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_soo(epr_env, current_env, qs_env, ib)
Calculates g_soo (soft part only for now)
subroutine, public epr_g_so(epr_env, current_env, qs_env, ib)
Calculates g_so.
subroutine, public epr_g_zke(epr_env, qs_env)
Calculate zke part of the g tensor.
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_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 ...
subroutine, public nmr_shift(nmr_env, current_env, qs_env, ib)
...
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...
subroutine, public polar_operators_local_wannier(qs_env, dcdr_env)
Calculate the dipole operator referenced at the Wannier centers in the MO basis.
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...
subroutine, public polar_operators(qs_env)
Calculate the dipole operator in the AO basis and its derivative wrt to MOs.
subroutine, public polar_operators_local(qs_env)
Calculate the Berry phase operator in the AO basis and then the derivative of the Berry phase operato...
subroutine, public issc_operators(issc_env, qs_env, iatom)
...
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.
subroutine, public get_polar_env(polar_env, do_raman, do_periodic, dberry_psi0, polar, psi1_dberry, run_stopped)
...
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...
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...
subroutine, public vcd_print(vcd_env, qs_env)
Print the APTs, AATs, and sum rules.
subroutine, public vcd_env_cleanup(qs_env, vcd_env)
Deallocate the vcd environment.
subroutine, public vcd_env_init(vcd_env, qs_env)
Initialize the vcd environment.
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
type of a logger, at the moment it contains just a print level starting at which level it should be l...
wrapper to abstract the force evaluation of the various methods
General settings for linear response calculations.
Represent a qs system that is perturbed. Can calculate the linear operator and the rhs of the system ...
keeps the density in various representations, keeping track of which ones are valid.