(git:34ef472)
cp2k_debug.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 Debug energy and derivatives w.r.t. finite differences
10 !> \note
11 !> Use INTERPOLATION USE_GUESS, in order to perform force and energy
12 !> calculations with the same density. This is not compulsory when iterating
13 !> to selfconsistency, but essential in the non-selfconsistent case [08.2005,TdK].
14 !> \par History
15 !> 12.2004 created [tlaino]
16 !> 08.2005 consistent_energies option added, to allow FD calculations
17 !> with the correct energies in the non-selfconsistent case, but
18 !> keep in mind, that the QS energies and forces are then NOT
19 !> consistent to each other [TdK].
20 !> 08.2005 In case the Harris functional is used, consistent_energies is
21 !> et to .FALSE., otherwise the QS energies are spuriously used [TdK].
22 !> 01.2015 Remove Harris functional option
23 !> - Revised (20.11.2013,MK)
24 !> \author Teodoro Laino
25 ! **************************************************************************************************
26 MODULE cp2k_debug
27  USE cell_types, ONLY: cell_type
28  USE cp_control_types, ONLY: dft_control_type
30  cp_logger_type
33  USE cp_result_methods, ONLY: get_results,&
35  USE cp_result_types, ONLY: cp_result_type
36  USE cp_subsys_types, ONLY: cp_subsys_get,&
37  cp_subsys_type
40  USE force_env_types, ONLY: force_env_get,&
41  force_env_type,&
48  section_vals_type,&
50  USE kinds, ONLY: default_string_length,&
51  dp
52  USE message_passing, ONLY: mp_para_env_type
55  USE particle_types, ONLY: particle_type
57  USE qs_kind_types, ONLY: qs_kind_type
58  USE string_utilities, ONLY: uppercase
59  USE virial_types, ONLY: virial_set,&
60  virial_type
61 #include "./base/base_uses.f90"
62 
63  IMPLICIT NONE
64  PRIVATE
65  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'cp2k_debug'
66 
68 
69 CONTAINS
70 
71 ! **************************************************************************************************
72 !> \brief ...
73 !> \param force_env ...
74 ! **************************************************************************************************
75  SUBROUTINE cp2k_debug_energy_and_forces(force_env)
76 
77  TYPE(force_env_type), POINTER :: force_env
78 
79  CHARACTER(LEN=3) :: cval1
80  CHARACTER(LEN=3*default_string_length) :: message
81  CHARACTER(LEN=60) :: line
82  CHARACTER(LEN=80), DIMENSION(:), POINTER :: cval2
83  CHARACTER(LEN=default_string_length) :: description
84  INTEGER :: i, ip, irep, iw, j, k, np, nrep, &
85  stress_tensor
86  LOGICAL :: check_failed, debug_dipole, &
87  debug_forces, debug_polar, &
88  debug_stress_tensor, skip, &
89  stop_on_mismatch
90  LOGICAL, ALLOCATABLE, DIMENSION(:, :) :: do_dof_atom_coor
91  LOGICAL, DIMENSION(3) :: do_dof_dipole
92  REAL(kind=dp) :: amplitude, dd, de, derr, difference, dx, &
93  eps_no_error_check, errmax, maxerr, &
94  std_value, sum_of_differences
95  REAL(kind=dp), DIMENSION(2) :: numer_energy
96  REAL(kind=dp), DIMENSION(3) :: dipole_moment, dipole_numer, err, &
97  my_maxerr, poldir
98  REAL(kind=dp), DIMENSION(3, 2) :: dipn
99  REAL(kind=dp), DIMENSION(3, 3) :: polar_analytic, polar_numeric
100  REAL(kind=dp), DIMENSION(9) :: pvals
101  REAL(kind=dp), DIMENSION(:, :), POINTER :: analyt_forces, numer_forces
102  TYPE(cell_type), POINTER :: cell
103  TYPE(cp_logger_type), POINTER :: logger
104  TYPE(cp_result_type), POINTER :: results
105  TYPE(cp_subsys_type), POINTER :: subsys
106  TYPE(dft_control_type), POINTER :: dft_control
107  TYPE(mp_para_env_type), POINTER :: para_env
108  TYPE(particle_type), DIMENSION(:), POINTER :: particles
109  TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
110  TYPE(section_vals_type), POINTER :: root_section, subsys_section
111 
112  NULLIFY (analyt_forces, numer_forces, subsys, particles)
113 
114  root_section => force_env%root_section
115 
116  CALL force_env_get(force_env, para_env=para_env, subsys=subsys, cell=cell)
117  subsys_section => section_vals_get_subs_vals(force_env%force_env_section, &
118  "SUBSYS")
119 
120  CALL section_vals_val_get(root_section, "DEBUG%DEBUG_STRESS_TENSOR", &
121  l_val=debug_stress_tensor)
122  CALL section_vals_val_get(root_section, "DEBUG%DEBUG_FORCES", &
123  l_val=debug_forces)
124  CALL section_vals_val_get(root_section, "DEBUG%DEBUG_DIPOLE", &
125  l_val=debug_dipole)
126  CALL section_vals_val_get(root_section, "DEBUG%DEBUG_POLARIZABILITY", &
127  l_val=debug_polar)
128  CALL section_vals_val_get(root_section, "DEBUG%DX", &
129  r_val=dx)
130  CALL section_vals_val_get(root_section, "DEBUG%DE", &
131  r_val=de)
132  CALL section_vals_val_get(root_section, "DEBUG%CHECK_DIPOLE_DIRS", &
133  c_val=cval1)
134  dx = abs(dx)
135  CALL section_vals_val_get(root_section, "DEBUG%MAX_RELATIVE_ERROR", &
136  r_val=maxerr)
137  CALL section_vals_val_get(root_section, "DEBUG%EPS_NO_ERROR_CHECK", &
138  r_val=eps_no_error_check)
139  eps_no_error_check = max(eps_no_error_check, epsilon(0.0_dp))
140  CALL section_vals_val_get(root_section, "DEBUG%STOP_ON_MISMATCH", &
141  l_val=stop_on_mismatch)
142 
143  ! set active DOF
144  CALL uppercase(cval1)
145  do_dof_dipole(1) = (index(cval1, "X") /= 0)
146  do_dof_dipole(2) = (index(cval1, "Y") /= 0)
147  do_dof_dipole(3) = (index(cval1, "Z") /= 0)
148  NULLIFY (cval2)
149  IF (debug_forces) THEN
150  np = subsys%particles%n_els
151  ALLOCATE (do_dof_atom_coor(3, np))
152  CALL section_vals_val_get(root_section, "DEBUG%CHECK_ATOM_FORCE", n_rep_val=nrep)
153  IF (nrep > 0) THEN
154  do_dof_atom_coor = .false.
155  DO irep = 1, nrep
156  CALL section_vals_val_get(root_section, "DEBUG%CHECK_ATOM_FORCE", i_rep_val=irep, &
157  c_vals=cval2)
158  READ (cval2(1), fmt="(I10)") k
159  CALL uppercase(cval2(2))
160  do_dof_atom_coor(1, k) = (index(cval2(2), "X") /= 0)
161  do_dof_atom_coor(2, k) = (index(cval2(2), "Y") /= 0)
162  do_dof_atom_coor(3, k) = (index(cval2(2), "Z") /= 0)
163  END DO
164  ELSE
165  do_dof_atom_coor = .true.
166  END IF
167  END IF
168 
169  logger => cp_get_default_logger()
170  iw = cp_print_key_unit_nr(logger, root_section, "DEBUG%PROGRAM_RUN_INFO", &
171  extension=".log")
172  IF (debug_stress_tensor) THEN
173  ! To debug stress tensor the stress tensor calculation must be
174  ! first enabled..
175  CALL section_vals_val_get(force_env%force_env_section, "STRESS_TENSOR", &
176  i_val=stress_tensor)
177  skip = .false.
178  SELECT CASE (stress_tensor)
180  ! OK
182  ! Nothing to check
183  CALL cp_warn(__location__, "Numerical stress tensor was requested in "// &
184  "the FORCE_EVAL section. Nothing to debug!")
185  skip = .true.
186  CASE DEFAULT
187  CALL cp_warn(__location__, "Stress tensor calculation was not enabled in "// &
188  "the FORCE_EVAL section. Nothing to debug!")
189  skip = .true.
190  END SELECT
191 
192  IF (.NOT. skip) THEN
193 
194  block
195  TYPE(virial_type) :: virial_analytical, virial_numerical
196  TYPE(virial_type), POINTER :: virial
197 
198  ! Compute the analytical stress tensor
199  CALL cp_subsys_get(subsys, virial=virial)
200  CALL virial_set(virial, pv_numer=.false.)
201  CALL force_env_calc_energy_force(force_env, &
202  calc_force=.true., &
203  calc_stress_tensor=.true.)
204 
205  ! Retrieve the analytical virial
206  virial_analytical = virial
207 
208  ! Debug stress tensor (numerical vs analytical)
209  CALL virial_set(virial, pv_numer=.true.)
210  CALL force_env_calc_num_pressure(force_env, dx=dx)
211 
212  ! Retrieve the numerical virial
213  CALL cp_subsys_get(subsys, virial=virial)
214  virial_numerical = virial
215 
216  ! Print results
217  IF (iw > 0) THEN
218  WRITE (unit=iw, fmt="((T2,A))") &
219  "DEBUG| Numerical pv_virial [a.u.]"
220  WRITE (unit=iw, fmt="((T2,A,T21,3(1X,F19.12)))") &
221  ("DEBUG|", virial_numerical%pv_virial(i, 1:3), i=1, 3)
222  WRITE (unit=iw, fmt="(/,(T2,A))") &
223  "DEBUG| Analytical pv_virial [a.u.]"
224  WRITE (unit=iw, fmt="((T2,A,T21,3(1X,F19.12)))") &
225  ("DEBUG|", virial_analytical%pv_virial(i, 1:3), i=1, 3)
226  WRITE (unit=iw, fmt="(/,(T2,A))") &
227  "DEBUG| Difference pv_virial [a.u.]"
228  WRITE (unit=iw, fmt="((T2,A,T21,3(1X,F19.12)))") &
229  ("DEBUG|", virial_numerical%pv_virial(i, 1:3) - virial_analytical%pv_virial(i, 1:3), i=1, 3)
230  WRITE (unit=iw, fmt="(T2,A,T61,F20.12)") &
231  "DEBUG| Sum of differences", &
232  sum(abs(virial_numerical%pv_virial(:, :) - virial_analytical%pv_virial(:, :)))
233  END IF
234 
235  ! Check and abort on failure
236  check_failed = .false.
237  IF (iw > 0) THEN
238  WRITE (unit=iw, fmt="(/T2,A)") &
239  "DEBUG| Relative error pv_virial"
240  WRITE (unit=iw, fmt="(T2,A,T61,ES20.8)") &
241  "DEBUG| Threshold value for error check [a.u.]", eps_no_error_check
242  END IF
243  DO i = 1, 3
244  err(:) = 0.0_dp
245  DO k = 1, 3
246  IF (abs(virial_analytical%pv_virial(i, k)) >= eps_no_error_check) THEN
247  err(k) = 100.0_dp*(virial_numerical%pv_virial(i, k) - virial_analytical%pv_virial(i, k))/ &
248  virial_analytical%pv_virial(i, k)
249  WRITE (unit=line(20*(k - 1) + 1:20*k), fmt="(1X,F17.2,A2)") err(k), " %"
250  ELSE
251  WRITE (unit=line(20*(k - 1) + 1:20*k), fmt="(17X,A3)") "- %"
252  END IF
253  END DO
254  IF (iw > 0) THEN
255  WRITE (unit=iw, fmt="(T2,A,T21,A60)") &
256  "DEBUG|", line
257  END IF
258  IF (any(abs(err(1:3)) > maxerr)) check_failed = .true.
259  END DO
260  IF (iw > 0) THEN
261  WRITE (unit=iw, fmt="(T2,A,T61,F18.2,A2)") &
262  "DEBUG| Maximum accepted error", maxerr, " %"
263  END IF
264  IF (check_failed) THEN
265  message = "A mismatch between the analytical and the numerical "// &
266  "stress tensor has been detected. Check the implementation "// &
267  "of the stress tensor"
268  IF (stop_on_mismatch) THEN
269  cpabort(trim(message))
270  ELSE
271  cpwarn(trim(message))
272  END IF
273  END IF
274  END block
275  END IF
276  END IF
277 
278  IF (debug_forces) THEN
279  ! Debug forces (numerical vs analytical)
280  particles => subsys%particles%els
281  SELECT CASE (force_env%in_use)
282  CASE (use_qs_force)
283  CALL get_qs_env(force_env%qs_env, qs_kind_set=qs_kind_set)
284  CALL write_qs_particle_coordinates(particles, qs_kind_set, subsys_section, "DEBUG")
285  CASE DEFAULT
286  CALL write_fist_particle_coordinates(particles, subsys_section)
287  END SELECT
288  ! First evaluate energy and forces
289  CALL force_env_calc_energy_force(force_env, &
290  calc_force=.true., &
291  calc_stress_tensor=.false.)
292  ! Copy forces in array and start the numerical calculation
293  IF (ASSOCIATED(analyt_forces)) DEALLOCATE (analyt_forces)
294  np = subsys%particles%n_els
295  ALLOCATE (analyt_forces(np, 3))
296  DO ip = 1, np
297  analyt_forces(ip, 1:3) = particles(ip)%f(1:3)
298  END DO
299  ! Loop on atoms and coordinates
300  IF (ASSOCIATED(numer_forces)) DEALLOCATE (numer_forces)
301  ALLOCATE (numer_forces(subsys%particles%n_els, 3))
302  atom: DO ip = 1, np
303  coord: DO k = 1, 3
304  IF (do_dof_atom_coor(k, ip)) THEN
305  numer_energy = 0.0_dp
306  std_value = particles(ip)%r(k)
307  DO j = 1, 2
308  particles(ip)%r(k) = std_value - (-1.0_dp)**j*dx
309  SELECT CASE (force_env%in_use)
310  CASE (use_qs_force)
311  CALL get_qs_env(force_env%qs_env, qs_kind_set=qs_kind_set)
312  CALL write_qs_particle_coordinates(particles, qs_kind_set, subsys_section, "DEBUG")
313  CASE DEFAULT
314  CALL write_fist_particle_coordinates(particles, subsys_section)
315  END SELECT
316  ! Compute energy
317  CALL force_env_calc_energy_force(force_env, &
318  calc_force=.false., &
319  calc_stress_tensor=.false., &
320  consistent_energies=.true.)
321  CALL force_env_get(force_env, potential_energy=numer_energy(j))
322  END DO
323  particles(ip)%r(k) = std_value
324  numer_forces(ip, k) = -0.5_dp*(numer_energy(1) - numer_energy(2))/dx
325  IF (iw > 0) THEN
326  WRITE (unit=iw, fmt="(/,T2,A,T17,A,F7.4,A,T34,A,F7.4,A,T52,A,T68,A)") &
327  "DEBUG| Atom", "E("//achar(119 + k)//" +", dx, ")", &
328  "E("//achar(119 + k)//" -", dx, ")", &
329  "f(numerical)", "f(analytical)"
330  WRITE (unit=iw, fmt="(T2,A,I5,4(1X,F16.8))") &
331  "DEBUG|", ip, numer_energy(1:2), numer_forces(ip, k), analyt_forces(ip, k)
332  END IF
333  ELSE
334  numer_forces(ip, k) = 0.0_dp
335  END IF
336  END DO coord
337  ! Check analytical forces using the numerical forces as reference
338  my_maxerr = maxerr
339  err(1:3) = 0.0_dp
340  DO k = 1, 3
341  IF (do_dof_atom_coor(k, ip)) THEN
342  ! Calculate percentage but ignore very small force values
343  IF (abs(analyt_forces(ip, k)) >= eps_no_error_check) THEN
344  err(k) = 100.0_dp*(numer_forces(ip, k) - analyt_forces(ip, k))/analyt_forces(ip, k)
345  END IF
346  ! Increase error tolerance for small force values
347  IF (abs(analyt_forces(ip, k)) <= 0.0001_dp) my_maxerr(k) = 5.0_dp*my_maxerr(k)
348  ELSE
349  err(k) = 0.0_dp
350  END IF
351  END DO
352  IF (iw > 0) THEN
353  IF (any(do_dof_atom_coor(1:3, ip))) THEN
354  WRITE (unit=iw, fmt="(/,T2,A)") &
355  "DEBUG| Atom Coordinate f(numerical) f(analytical) Difference Error [%]"
356  END IF
357  DO k = 1, 3
358  IF (do_dof_atom_coor(k, ip)) THEN
359  difference = analyt_forces(ip, k) - numer_forces(ip, k)
360  IF (abs(analyt_forces(ip, k)) >= eps_no_error_check) THEN
361  WRITE (unit=iw, fmt="(T2,A,I5,T19,A1,T26,F14.8,T42,F14.8,T57,F12.8,T71,F10.2)") &
362  "DEBUG|", ip, achar(119 + k), numer_forces(ip, k), analyt_forces(ip, k), difference, err(k)
363  ELSE
364  WRITE (unit=iw, fmt="(T2,A,I5,T19,A1,T26,F14.8,T42,F14.8,T57,F12.8,T80,A1)") &
365  "DEBUG|", ip, achar(119 + k), numer_forces(ip, k), analyt_forces(ip, k), difference, "-"
366  END IF
367  END IF
368  END DO
369  END IF
370  IF (any(abs(err(1:3)) > my_maxerr(1:3))) THEN
371  message = "A mismatch between analytical and numerical forces "// &
372  "has been detected. Check the implementation of the "// &
373  "analytical force calculation"
374  IF (stop_on_mismatch) THEN
375  cpabort(message)
376  ELSE
377  cpwarn(message)
378  END IF
379  END IF
380  END DO atom
381  ! Print summary
382  IF (iw > 0) THEN
383  WRITE (unit=iw, fmt="(/,(T2,A))") &
384  "DEBUG|======================== BEGIN OF SUMMARY ===============================", &
385  "DEBUG| Atom Coordinate f(numerical) f(analytical) Difference Error [%]"
386  sum_of_differences = 0.0_dp
387  errmax = 0.0_dp
388  DO ip = 1, np
389  err(1:3) = 0.0_dp
390  DO k = 1, 3
391  IF (do_dof_atom_coor(k, ip)) THEN
392  difference = analyt_forces(ip, k) - numer_forces(ip, k)
393  IF (abs(analyt_forces(ip, k)) >= eps_no_error_check) THEN
394  err(k) = 100_dp*(numer_forces(ip, k) - analyt_forces(ip, k))/analyt_forces(ip, k)
395  errmax = max(errmax, abs(err(k)))
396  WRITE (unit=iw, fmt="(T2,A,I5,T19,A1,T26,F14.8,T42,F14.8,T57,F12.8,T71,F10.2)") &
397  "DEBUG|", ip, achar(119 + k), numer_forces(ip, k), analyt_forces(ip, k), difference, err(k)
398  ELSE
399  WRITE (unit=iw, fmt="(T2,A,I5,T19,A1,T26,F14.8,T42,F14.8,T57,F12.8,T80,A1)") &
400  "DEBUG|", ip, achar(119 + k), numer_forces(ip, k), analyt_forces(ip, k), difference, "-"
401  END IF
402  sum_of_differences = sum_of_differences + abs(difference)
403  END IF
404  END DO
405  END DO
406  WRITE (unit=iw, fmt="(T2,A,T57,F12.8,T71,F10.2)") &
407  "DEBUG| Sum of differences:", sum_of_differences, errmax
408  WRITE (unit=iw, fmt="(T2,A)") &
409  "DEBUG|======================== END OF SUMMARY ================================="
410  END IF
411  ! Release work storage
412  IF (ASSOCIATED(analyt_forces)) DEALLOCATE (analyt_forces)
413  IF (ASSOCIATED(numer_forces)) DEALLOCATE (numer_forces)
414  DEALLOCATE (do_dof_atom_coor)
415  END IF
416 
417  IF (debug_dipole) THEN
418  IF (force_env%in_use == use_qs_force) THEN
419  CALL get_qs_env(force_env%qs_env, dft_control=dft_control)
420  poldir = (/0.0_dp, 0.0_dp, 1.0_dp/)
421  amplitude = 0.0_dp
422  CALL set_efield(dft_control, amplitude, poldir)
423  CALL force_env_calc_energy_force(force_env, calc_force=.true.)
424  CALL get_qs_env(force_env%qs_env, results=results)
425  description = "[DIPOLE]"
426  IF (test_for_result(results, description=description)) THEN
427  CALL get_results(results, description=description, values=dipole_moment)
428  ELSE
429  CALL cp_warn(__location__, "Debug of dipole moments needs DFT/PRINT/MOMENTS section enabled")
430  cpabort("DEBUG")
431  END IF
432  amplitude = de
433  DO k = 1, 3
434  IF (do_dof_dipole(k)) THEN
435  poldir = 0.0_dp
436  poldir(k) = 1.0_dp
437  DO j = 1, 2
438  poldir = -1.0_dp*poldir
439  CALL set_efield(dft_control, amplitude, poldir)
440  CALL force_env_calc_energy_force(force_env, calc_force=.false.)
441  CALL force_env_get(force_env, potential_energy=numer_energy(j))
442  END DO
443  dipole_numer(k) = 0.5_dp*(numer_energy(1) - numer_energy(2))/de
444  ELSE
445  dipole_numer(k) = 0.0_dp
446  END IF
447  END DO
448  IF (iw > 0) THEN
449  WRITE (unit=iw, fmt="(/,(T2,A))") &
450  "DEBUG|========================= DIPOLE MOMENTS ================================", &
451  "DEBUG| Coordinate D(numerical) D(analytical) Difference Error [%]"
452  err(1:3) = 0.0_dp
453  DO k = 1, 3
454  IF (do_dof_dipole(k)) THEN
455  dd = dipole_moment(k) - dipole_numer(k)
456  IF (abs(dipole_moment(k)) > eps_no_error_check) THEN
457  derr = 100._dp*dd/dipole_moment(k)
458  WRITE (unit=iw, fmt="(T2,A,T13,A1,T21,F16.8,T38,F16.8,T56,G12.3,T72,F9.3)") &
459  "DEBUG|", achar(119 + k), dipole_numer(k), dipole_moment(k), dd, derr
460  ELSE
461  derr = 0.0_dp
462  WRITE (unit=iw, fmt="(T2,A,T13,A1,T21,F16.8,T38,F16.8,T56,G12.3)") &
463  "DEBUG|", achar(119 + k), dipole_numer(k), dipole_moment(k), dd
464  END IF
465  err(k) = derr
466  ELSE
467  WRITE (unit=iw, fmt="(T2,A,T13,A1,T21,A16,T38,F16.8)") &
468  "DEBUG|", achar(119 + k), " skipped", dipole_moment(k)
469  END IF
470  END DO
471  WRITE (unit=iw, fmt="((T2,A))") &
472  "DEBUG|========================================================================="
473  WRITE (unit=iw, fmt="(T2,A,T61,E20.12)") 'DIPOLE : CheckSum =', sum(dipole_moment)
474  IF (any(abs(err(1:3)) > maxerr)) THEN
475  message = "A mismatch between analytical and numerical dipoles "// &
476  "has been detected. Check the implementation of the "// &
477  "analytical dipole calculation"
478  IF (stop_on_mismatch) THEN
479  cpabort(message)
480  ELSE
481  cpwarn(message)
482  END IF
483  END IF
484  END IF
485 
486  ELSE
487  CALL cp_warn(__location__, "Debug of dipole moments only for Quickstep code available")
488  END IF
489  END IF
490 
491  IF (debug_polar) THEN
492  IF (force_env%in_use == use_qs_force) THEN
493  CALL get_qs_env(force_env%qs_env, dft_control=dft_control)
494  poldir = (/0.0_dp, 0.0_dp, 1.0_dp/)
495  amplitude = 0.0_dp
496  CALL set_efield(dft_control, amplitude, poldir)
497  CALL force_env_calc_energy_force(force_env, calc_force=.false.)
498  CALL get_qs_env(force_env%qs_env, results=results)
499  description = "[POLAR]"
500  IF (test_for_result(results, description=description)) THEN
501  CALL get_results(results, description=description, values=pvals)
502  polar_analytic(1:3, 1:3) = reshape(pvals(1:9), (/3, 3/))
503  ELSE
504  CALL cp_warn(__location__, "Debug of polarizabilities needs PROPERTIES/LINRES/POLAR section enabled")
505  cpabort("DEBUG")
506  END IF
507  description = "[DIPOLE]"
508  IF (.NOT. test_for_result(results, description=description)) THEN
509  CALL cp_warn(__location__, "Debug of polarizabilities need DFT/PRINT/MOMENTS section enabled")
510  cpabort("DEBUG")
511  END IF
512  amplitude = de
513  DO k = 1, 3
514  poldir = 0.0_dp
515  poldir(k) = 1.0_dp
516  DO j = 1, 2
517  poldir = -1.0_dp*poldir
518  CALL set_efield(dft_control, amplitude, poldir)
519  CALL force_env_calc_energy_force(force_env, calc_force=.false., linres=.true.)
520  CALL get_results(results, description=description, values=dipn(1:3, j))
521  END DO
522  polar_numeric(1:3, k) = 0.5_dp*(dipn(1:3, 2) - dipn(1:3, 1))/de
523  END DO
524  IF (iw > 0) THEN
525  WRITE (unit=iw, fmt="(/,(T2,A))") &
526  "DEBUG|========================= POLARIZABILITY ================================", &
527  "DEBUG| Coordinates P(numerical) P(analytical) Difference Error [%]"
528  DO k = 1, 3
529  DO j = 1, 3
530  dd = polar_analytic(k, j) - polar_numeric(k, j)
531  IF (abs(polar_analytic(k, j)) > eps_no_error_check) THEN
532  derr = 100._dp*dd/polar_analytic(k, j)
533  WRITE (unit=iw, fmt="(T2,A,T12,A1,A1,T21,F16.8,T38,F16.8,T56,G12.3,T72,F9.3)") &
534  "DEBUG|", achar(119 + k), achar(119 + j), polar_numeric(k, j), polar_analytic(k, j), dd, derr
535  ELSE
536  WRITE (unit=iw, fmt="(T2,A,T12,A1,A1,T21,F16.8,T38,F16.8,T56,G12.3)") &
537  "DEBUG|", achar(119 + k), achar(119 + j), polar_numeric(k, j), polar_analytic(k, j), dd
538  END IF
539  END DO
540  END DO
541  WRITE (unit=iw, fmt="((T2,A))") &
542  "DEBUG|========================================================================="
543  WRITE (unit=iw, fmt="(T2,A,T61,E20.12)") ' POLAR : CheckSum =', sum(polar_analytic)
544  END IF
545  ELSE
546  CALL cp_warn(__location__, "Debug of polarizabilities only for Quickstep code available")
547  END IF
548  END IF
549 
550  CALL cp_print_key_finished_output(iw, logger, root_section, "DEBUG%PROGRAM_RUN_INFO")
551 
552  END SUBROUTINE cp2k_debug_energy_and_forces
553 
554 ! **************************************************************************************************
555 !> \brief ...
556 !> \param dft_control ...
557 !> \param amplitude ...
558 !> \param poldir ...
559 ! **************************************************************************************************
560  SUBROUTINE set_efield(dft_control, amplitude, poldir)
561  TYPE(dft_control_type), POINTER :: dft_control
562  REAL(kind=dp), INTENT(IN) :: amplitude
563  REAL(kind=dp), DIMENSION(3), INTENT(IN) :: poldir
564 
565  IF (dft_control%apply_efield) THEN
566  dft_control%efield_fields(1)%efield%strength = amplitude
567  dft_control%efield_fields(1)%efield%polarisation(1:3) = poldir(1:3)
568  ELSEIF (dft_control%apply_period_efield) THEN
569  dft_control%period_efield%strength = amplitude
570  dft_control%period_efield%polarisation(1:3) = poldir(1:3)
571  ELSE
572  cpabort("No EFIELD definition available")
573  END IF
574 
575  END SUBROUTINE set_efield
576 
577 END MODULE cp2k_debug
Definition: atom.F:9
Handles all functions related to the CELL.
Definition: cell_types.F:15
Debug energy and derivatives w.r.t. finite differences.
Definition: cp2k_debug.F:26
subroutine, public cp2k_debug_energy_and_forces(force_env)
...
Definition: cp2k_debug.F:76
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,...
set of type/routines to handle the storage of results in force_envs
logical function, public test_for_result(results, description)
test for a certain result in the result_list
set of type/routines to handle the storage of results in force_envs
types that represent a subsys, i.e. a part of the system
subroutine, public cp_subsys_get(subsys, ref_count, atomic_kinds, atomic_kind_set, particles, particle_set, local_particles, molecules, molecule_set, molecule_kinds, molecule_kind_set, local_molecules, para_env, colvar_p, shell_particles, core_particles, gci, multipoles, natom, nparticle, ncore, nshell, nkind, atprop, virial, results, cell)
returns information about various attributes of the given subsys
Interface for the force calculations.
recursive subroutine, public force_env_calc_energy_force(force_env, calc_force, consistent_energies, skip_external_control, eval_energy_forces, require_consistent_energy_force, linres, calc_stress_tensor)
Interface routine for force and energy calculations.
subroutine, public force_env_calc_num_pressure(force_env, dx)
Evaluates the stress tensor and pressure numerically.
Interface for the force calculations.
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 do_stress_analytical
integer, parameter, public do_stress_diagonal_anal
integer, parameter, public do_stress_diagonal_numer
integer, parameter, public do_stress_numerical
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_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_string_length
Definition: kinds.F:57
Interface to the message passing library MPI.
Define methods related to particle_type.
subroutine, public write_qs_particle_coordinates(particle_set, qs_kind_set, subsys_section, label)
Write the atomic coordinates to the output unit.
subroutine, public write_fist_particle_coordinates(particle_set, subsys_section, charges)
Write the atomic coordinates to the output unit.
Define the data structure for the particle information.
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 quickstep kind type and their sub types.
Definition: qs_kind_types.F:23
Utilities for string manipulations.
elemental subroutine, public uppercase(string)
Convert all lower case characters in a string to upper case.
subroutine, public virial_set(virial, pv_total, pv_kinetic, pv_virial, pv_xc, pv_fock_4c, pv_constraint, pv_overlap, pv_ekinetic, pv_ppl, pv_ppnl, pv_ecore_overlap, pv_ehartree, pv_exc, pv_exx, pv_vdw, pv_mp2, pv_nlcc, pv_gapw, pv_lrigpw, pv_availability, pv_calculate, pv_numer, pv_diagonal)
...
Definition: virial_types.F:192