(git:618e178)
Loading...
Searching...
No Matches
cp2k_debug.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2026 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! **************************************************************************************************
27 USE cell_types, ONLY: cell_type,&
51 USE kinds, ONLY: default_string_length,&
52 dp
60 USE virial_types, ONLY: virial_set,&
62#include "./base/base_uses.f90"
63
64 IMPLICIT NONE
65 PRIVATE
66 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'cp2k_debug'
67
69
70CONTAINS
71
72! **************************************************************************************************
73!> \brief ...
74!> \param force_env ...
75! **************************************************************************************************
76 SUBROUTINE cp2k_debug_energy_and_forces(force_env)
77
78 TYPE(force_env_type), POINTER :: force_env
79
80 CHARACTER(LEN=3) :: cval1
81 CHARACTER(LEN=3*default_string_length) :: message
82 CHARACTER(LEN=60) :: line
83 CHARACTER(LEN=80), DIMENSION(:), POINTER :: cval2
84 CHARACTER(LEN=default_string_length) :: description
85 INTEGER :: i, ip, irep, iw, j, k, n_periodic, np, &
86 nrep, stress_tensor
87 INTEGER, DIMENSION(3) :: periodic
88 LOGICAL :: check_failed, debug_dipole, &
89 debug_forces, debug_polar, &
90 debug_stress_tensor, skip, &
91 stop_on_mismatch
92 LOGICAL, ALLOCATABLE, DIMENSION(:, :) :: do_dof_atom_coor
93 LOGICAL, DIMENSION(3) :: do_dof_dipole
94 REAL(kind=dp) :: amplitude, dd, de, derr, difference, dx, eps_no_error_check, errmax, &
95 maxerr, periodic_stress_sum, std_value, sum_of_differences
96 REAL(kind=dp), DIMENSION(2) :: numer_energy
97 REAL(kind=dp), DIMENSION(3) :: dipole_moment, dipole_numer, err, &
98 my_maxerr, poldir
99 REAL(kind=dp), DIMENSION(3, 2) :: dipn
100 REAL(kind=dp), DIMENSION(3, 3) :: polar_analytic, polar_numeric
101 REAL(kind=dp), DIMENSION(9) :: pvals
102 REAL(kind=dp), DIMENSION(:, :), POINTER :: analyt_forces, numer_forces
103 TYPE(cell_type), POINTER :: cell
104 TYPE(cp_logger_type), POINTER :: logger
105 TYPE(cp_result_type), POINTER :: results
106 TYPE(cp_subsys_type), POINTER :: subsys
107 TYPE(dft_control_type), POINTER :: dft_control
108 TYPE(mp_para_env_type), POINTER :: para_env
109 TYPE(particle_type), DIMENSION(:), POINTER :: particles
110 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
111 TYPE(section_vals_type), POINTER :: root_section, subsys_section
112
113 NULLIFY (analyt_forces, numer_forces, subsys, particles)
114
115 root_section => force_env%root_section
116
117 CALL force_env_get(force_env, para_env=para_env, subsys=subsys, cell=cell)
118 subsys_section => section_vals_get_subs_vals(force_env%force_env_section, &
119 "SUBSYS")
120
121 CALL section_vals_val_get(root_section, "DEBUG%DEBUG_STRESS_TENSOR", &
122 l_val=debug_stress_tensor)
123 CALL section_vals_val_get(root_section, "DEBUG%DEBUG_FORCES", &
124 l_val=debug_forces)
125 CALL section_vals_val_get(root_section, "DEBUG%DEBUG_DIPOLE", &
126 l_val=debug_dipole)
127 CALL section_vals_val_get(root_section, "DEBUG%DEBUG_POLARIZABILITY", &
128 l_val=debug_polar)
129 CALL section_vals_val_get(root_section, "DEBUG%DX", &
130 r_val=dx)
131 CALL section_vals_val_get(root_section, "DEBUG%DE", &
132 r_val=de)
133 CALL section_vals_val_get(root_section, "DEBUG%CHECK_DIPOLE_DIRS", &
134 c_val=cval1)
135 dx = abs(dx)
136 CALL section_vals_val_get(root_section, "DEBUG%MAX_RELATIVE_ERROR", &
137 r_val=maxerr)
138 CALL section_vals_val_get(root_section, "DEBUG%EPS_NO_ERROR_CHECK", &
139 r_val=eps_no_error_check)
140 eps_no_error_check = max(eps_no_error_check, epsilon(0.0_dp))
141 CALL section_vals_val_get(root_section, "DEBUG%STOP_ON_MISMATCH", &
142 l_val=stop_on_mismatch)
143
144 ! set active DOF
145 CALL uppercase(cval1)
146 do_dof_dipole(1) = (index(cval1, "X") /= 0)
147 do_dof_dipole(2) = (index(cval1, "Y") /= 0)
148 do_dof_dipole(3) = (index(cval1, "Z") /= 0)
149 NULLIFY (cval2)
150 IF (debug_forces) THEN
151 np = subsys%particles%n_els
152 ALLOCATE (do_dof_atom_coor(3, np))
153 CALL section_vals_val_get(root_section, "DEBUG%CHECK_ATOM_FORCE", n_rep_val=nrep)
154 IF (nrep > 0) THEN
155 do_dof_atom_coor = .false.
156 DO irep = 1, nrep
157 CALL section_vals_val_get(root_section, "DEBUG%CHECK_ATOM_FORCE", i_rep_val=irep, &
158 c_vals=cval2)
159 READ (cval2(1), fmt="(I10)") k
160 CALL uppercase(cval2(2))
161 do_dof_atom_coor(1, k) = (index(cval2(2), "X") /= 0)
162 do_dof_atom_coor(2, k) = (index(cval2(2), "Y") /= 0)
163 do_dof_atom_coor(3, k) = (index(cval2(2), "Z") /= 0)
164 END DO
165 ELSE
166 do_dof_atom_coor = .true.
167 END IF
168 END IF
169
170 logger => cp_get_default_logger()
171 iw = cp_print_key_unit_nr(logger, root_section, "DEBUG%PROGRAM_RUN_INFO", &
172 extension=".log")
173 IF (debug_stress_tensor) THEN
174 IF (sum(cell%perd) == 0) THEN
175 CALL cp_warn(__location__, &
176 "DEBUG_STRESS_TENSOR requested for PERIODIC NONE. "// &
177 "The isolated-system virial is useful for finite-difference diagnostics, "// &
178 "but it is not a physically meaningful bulk stress.")
179 END IF
180 ! To debug stress tensor the stress tensor calculation must be
181 ! first enabled..
182 CALL section_vals_val_get(force_env%force_env_section, "STRESS_TENSOR", &
183 i_val=stress_tensor)
184 skip = .false.
185 SELECT CASE (stress_tensor)
187 ! OK
189 ! Nothing to check
190 CALL cp_warn(__location__, "Numerical stress tensor was requested in "// &
191 "the FORCE_EVAL section. Nothing to debug!")
192 skip = .true.
193 CASE DEFAULT
194 CALL cp_warn(__location__, "Stress tensor calculation was not enabled in "// &
195 "the FORCE_EVAL section. Nothing to debug!")
196 skip = .true.
197 END SELECT
198
199 IF (.NOT. skip) THEN
200
201 block
202 TYPE(virial_type) :: virial_analytical, virial_numerical
203 TYPE(virial_type), POINTER :: virial
204
205 ! Compute the analytical stress tensor
206 CALL cp_subsys_get(subsys, virial=virial)
207 CALL virial_set(virial, pv_numer=.false.)
208 CALL force_env_calc_energy_force(force_env, &
209 calc_force=.true., &
210 calc_stress_tensor=.true.)
211
212 ! Retrieve the analytical virial
213 virial_analytical = virial
214
215 ! Debug stress tensor (numerical vs analytical)
216 CALL virial_set(virial, pv_numer=.true.)
217 CALL force_env_calc_num_pressure(force_env, dx=dx)
218
219 ! Retrieve the numerical virial
220 CALL cp_subsys_get(subsys, virial=virial)
221 virial_numerical = virial
222
223 ! Numerical diagonal stress checks only perturb diagonal cell elements.
224 IF (virial_analytical%pv_diagonal .OR. virial_numerical%pv_diagonal) THEN
225 DO i = 1, 3
226 DO k = 1, 3
227 IF (i /= k) THEN
228 virial_analytical%pv_virial(i, k) = 0.0_dp
229 virial_numerical%pv_virial(i, k) = 0.0_dp
230 END IF
231 END DO
232 END DO
233 END IF
234
235 ! Print results
236 IF (iw > 0) THEN
237 WRITE (unit=iw, fmt="((T2,A))") &
238 "DEBUG| Numerical pv_virial [a.u.]"
239 WRITE (unit=iw, fmt="((T2,A,T21,3(1X,F19.12)))") &
240 ("DEBUG|", virial_numerical%pv_virial(i, 1:3), i=1, 3)
241 WRITE (unit=iw, fmt="(/,(T2,A))") &
242 "DEBUG| Analytical pv_virial [a.u.]"
243 WRITE (unit=iw, fmt="((T2,A,T21,3(1X,F19.12)))") &
244 ("DEBUG|", virial_analytical%pv_virial(i, 1:3), i=1, 3)
245 WRITE (unit=iw, fmt="(/,(T2,A))") &
246 "DEBUG| Difference pv_virial [a.u.]"
247 WRITE (unit=iw, fmt="((T2,A,T21,3(1X,F19.12)))") &
248 ("DEBUG|", virial_numerical%pv_virial(i, 1:3) - virial_analytical%pv_virial(i, 1:3), i=1, 3)
249 WRITE (unit=iw, fmt="(T2,A,T61,F20.12)") &
250 "DEBUG| Sum of differences", &
251 sum(abs(virial_numerical%pv_virial(:, :) - virial_analytical%pv_virial(:, :)))
252 CALL get_cell(cell=cell, periodic=periodic)
253 n_periodic = count(periodic /= 0)
254 IF (n_periodic > 0 .AND. n_periodic < 3) THEN
255 periodic_stress_sum = 0.0_dp
256 DO i = 1, 3
257 DO k = 1, 3
258 IF (periodic(i) /= 0 .AND. periodic(k) /= 0) THEN
259 periodic_stress_sum = periodic_stress_sum + &
260 abs(virial_numerical%pv_virial(i, k) - &
261 virial_analytical%pv_virial(i, k))
262 END IF
263 END DO
264 END DO
265 WRITE (unit=iw, fmt="(T2,A,T61,F20.12)") &
266 "DEBUG| Periodic-subspace sum of differences", periodic_stress_sum
267 END IF
268 END IF
269
270 ! Check and abort on failure
271 check_failed = .false.
272 IF (iw > 0) THEN
273 WRITE (unit=iw, fmt="(/T2,A)") &
274 "DEBUG| Relative error pv_virial"
275 WRITE (unit=iw, fmt="(T2,A,T61,ES20.8)") &
276 "DEBUG| Threshold value for error check [a.u.]", eps_no_error_check
277 END IF
278 DO i = 1, 3
279 err(:) = 0.0_dp
280 DO k = 1, 3
281 IF (abs(virial_analytical%pv_virial(i, k)) >= eps_no_error_check) THEN
282 err(k) = 100.0_dp*(virial_numerical%pv_virial(i, k) - virial_analytical%pv_virial(i, k))/ &
283 virial_analytical%pv_virial(i, k)
284 WRITE (unit=line(20*(k - 1) + 1:20*k), fmt="(1X,F17.2,A2)") err(k), " %"
285 ELSE
286 WRITE (unit=line(20*(k - 1) + 1:20*k), fmt="(17X,A3)") "- %"
287 END IF
288 END DO
289 IF (iw > 0) THEN
290 WRITE (unit=iw, fmt="(T2,A,T21,A60)") &
291 "DEBUG|", line
292 END IF
293 IF (any(abs(err(1:3)) > maxerr)) check_failed = .true.
294 END DO
295 IF (iw > 0) THEN
296 WRITE (unit=iw, fmt="(T2,A,T61,F18.2,A2)") &
297 "DEBUG| Maximum accepted error", maxerr, " %"
298 END IF
299 IF (check_failed) THEN
300 message = "A mismatch between the analytical and the numerical "// &
301 "stress tensor has been detected. Check the implementation "// &
302 "of the stress tensor"
303 IF (stop_on_mismatch) THEN
304 cpabort(trim(message))
305 ELSE
306 cpwarn(trim(message))
307 END IF
308 END IF
309 END block
310 END IF
311 END IF
312
313 IF (debug_forces) THEN
314 ! Debug forces (numerical vs analytical)
315 particles => subsys%particles%els
316 SELECT CASE (force_env%in_use)
317 CASE (use_qs_force)
318 CALL get_qs_env(force_env%qs_env, qs_kind_set=qs_kind_set)
319 CALL write_qs_particle_coordinates(particles, qs_kind_set, subsys_section, "DEBUG")
320 CASE DEFAULT
321 CALL write_fist_particle_coordinates(particles, subsys_section)
322 END SELECT
323 ! First evaluate energy and forces
324 CALL force_env_calc_energy_force(force_env, &
325 calc_force=.true., &
326 calc_stress_tensor=.false.)
327 ! Copy forces in array and start the numerical calculation
328 IF (ASSOCIATED(analyt_forces)) DEALLOCATE (analyt_forces)
329 np = subsys%particles%n_els
330 ALLOCATE (analyt_forces(np, 3))
331 DO ip = 1, np
332 analyt_forces(ip, 1:3) = particles(ip)%f(1:3)
333 END DO
334 ! Loop on atoms and coordinates
335 IF (ASSOCIATED(numer_forces)) DEALLOCATE (numer_forces)
336 ALLOCATE (numer_forces(subsys%particles%n_els, 3))
337 atom: DO ip = 1, np
338 coord: DO k = 1, 3
339 IF (do_dof_atom_coor(k, ip)) THEN
340 numer_energy = 0.0_dp
341 std_value = particles(ip)%r(k)
342 DO j = 1, 2
343 particles(ip)%r(k) = std_value - (-1.0_dp)**j*dx
344 SELECT CASE (force_env%in_use)
345 CASE (use_qs_force)
346 CALL get_qs_env(force_env%qs_env, qs_kind_set=qs_kind_set)
347 CALL write_qs_particle_coordinates(particles, qs_kind_set, subsys_section, "DEBUG")
348 CASE DEFAULT
349 CALL write_fist_particle_coordinates(particles, subsys_section)
350 END SELECT
351 ! Compute energy
352 CALL force_env_calc_energy_force(force_env, &
353 calc_force=.false., &
354 calc_stress_tensor=.false., &
355 consistent_energies=.true.)
356 CALL force_env_get(force_env, potential_energy=numer_energy(j))
357 END DO
358 particles(ip)%r(k) = std_value
359 numer_forces(ip, k) = -0.5_dp*(numer_energy(1) - numer_energy(2))/dx
360 IF (iw > 0) THEN
361 WRITE (unit=iw, fmt="(/,T2,A,T17,A,F7.4,A,T34,A,F7.4,A,T52,A,T68,A)") &
362 "DEBUG| Atom", "E("//achar(119 + k)//" +", dx, ")", &
363 "E("//achar(119 + k)//" -", dx, ")", &
364 "f(numerical)", "f(analytical)"
365 WRITE (unit=iw, fmt="(T2,A,I5,4(1X,F16.8))") &
366 "DEBUG|", ip, numer_energy(1:2), numer_forces(ip, k), analyt_forces(ip, k)
367 END IF
368 ELSE
369 numer_forces(ip, k) = 0.0_dp
370 END IF
371 END DO coord
372 ! Check analytical forces using the numerical forces as reference
373 my_maxerr = maxerr
374 err(1:3) = 0.0_dp
375 DO k = 1, 3
376 IF (do_dof_atom_coor(k, ip)) THEN
377 ! Calculate percentage but ignore very small force values
378 IF (abs(analyt_forces(ip, k)) >= eps_no_error_check) THEN
379 err(k) = 100.0_dp*(numer_forces(ip, k) - analyt_forces(ip, k))/analyt_forces(ip, k)
380 END IF
381 ! Increase error tolerance for small force values
382 IF (abs(analyt_forces(ip, k)) <= 0.0001_dp) my_maxerr(k) = 5.0_dp*my_maxerr(k)
383 ELSE
384 err(k) = 0.0_dp
385 END IF
386 END DO
387 IF (iw > 0) THEN
388 IF (any(do_dof_atom_coor(1:3, ip))) THEN
389 WRITE (unit=iw, fmt="(/,T2,A)") &
390 "DEBUG| Atom Coordinate f(numerical) f(analytical) Difference Error [%]"
391 END IF
392 DO k = 1, 3
393 IF (do_dof_atom_coor(k, ip)) THEN
394 difference = analyt_forces(ip, k) - numer_forces(ip, k)
395 IF (abs(analyt_forces(ip, k)) >= eps_no_error_check) THEN
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 END IF
403 END DO
404 END IF
405 IF (any(abs(err(1:3)) > my_maxerr(1:3))) THEN
406 message = "A mismatch between analytical and numerical forces "// &
407 "has been detected. Check the implementation of the "// &
408 "analytical force calculation"
409 IF (stop_on_mismatch) THEN
410 cpabort(message)
411 ELSE
412 cpwarn(message)
413 END IF
414 END IF
415 END DO atom
416 ! Print summary
417 IF (iw > 0) THEN
418 WRITE (unit=iw, fmt="(/,(T2,A))") &
419 "DEBUG|======================== BEGIN OF SUMMARY ===============================", &
420 "DEBUG| Atom Coordinate f(numerical) f(analytical) Difference Error [%]"
421 sum_of_differences = 0.0_dp
422 errmax = 0.0_dp
423 DO ip = 1, np
424 err(1:3) = 0.0_dp
425 DO k = 1, 3
426 IF (do_dof_atom_coor(k, ip)) THEN
427 difference = analyt_forces(ip, k) - numer_forces(ip, k)
428 IF (abs(analyt_forces(ip, k)) >= eps_no_error_check) THEN
429 err(k) = 100_dp*(numer_forces(ip, k) - analyt_forces(ip, k))/analyt_forces(ip, k)
430 errmax = max(errmax, abs(err(k)))
431 WRITE (unit=iw, fmt="(T2,A,I5,T19,A1,T26,F14.8,T42,F14.8,T57,F12.8,T71,F10.2)") &
432 "DEBUG|", ip, achar(119 + k), numer_forces(ip, k), analyt_forces(ip, k), difference, err(k)
433 ELSE
434 WRITE (unit=iw, fmt="(T2,A,I5,T19,A1,T26,F14.8,T42,F14.8,T57,F12.8,T80,A1)") &
435 "DEBUG|", ip, achar(119 + k), numer_forces(ip, k), analyt_forces(ip, k), difference, "-"
436 END IF
437 sum_of_differences = sum_of_differences + abs(difference)
438 END IF
439 END DO
440 END DO
441 WRITE (unit=iw, fmt="(T2,A,T57,F12.8,T71,F10.2)") &
442 "DEBUG| Sum of differences:", sum_of_differences, errmax
443 WRITE (unit=iw, fmt="(T2,A)") &
444 "DEBUG|======================== END OF SUMMARY ================================="
445 END IF
446 ! Release work storage
447 IF (ASSOCIATED(analyt_forces)) DEALLOCATE (analyt_forces)
448 IF (ASSOCIATED(numer_forces)) DEALLOCATE (numer_forces)
449 DEALLOCATE (do_dof_atom_coor)
450 END IF
451
452 IF (debug_dipole) THEN
453 IF (force_env%in_use == use_qs_force) THEN
454 CALL get_qs_env(force_env%qs_env, dft_control=dft_control)
455 poldir = [0.0_dp, 0.0_dp, 1.0_dp]
456 amplitude = 0.0_dp
457 CALL set_efield(dft_control, amplitude, poldir)
458 CALL force_env_calc_energy_force(force_env, calc_force=.true.)
459 CALL get_qs_env(force_env%qs_env, results=results)
460 description = "[DIPOLE]"
461 IF (test_for_result(results, description=description)) THEN
462 CALL get_results(results, description=description, values=dipole_moment)
463 ELSE
464 CALL cp_warn(__location__, "Debug of dipole moments needs DFT/PRINT/MOMENTS section enabled")
465 cpabort("DEBUG")
466 END IF
467 amplitude = de
468 DO k = 1, 3
469 IF (do_dof_dipole(k)) THEN
470 poldir = 0.0_dp
471 poldir(k) = 1.0_dp
472 DO j = 1, 2
473 poldir = -1.0_dp*poldir
474 CALL set_efield(dft_control, amplitude, poldir)
475 CALL force_env_calc_energy_force(force_env, calc_force=.false.)
476 CALL force_env_get(force_env, potential_energy=numer_energy(j))
477 END DO
478 dipole_numer(k) = 0.5_dp*(numer_energy(1) - numer_energy(2))/de
479 ELSE
480 dipole_numer(k) = 0.0_dp
481 END IF
482 END DO
483 IF (iw > 0) THEN
484 WRITE (unit=iw, fmt="(/,(T2,A))") &
485 "DEBUG|========================= DIPOLE MOMENTS ================================", &
486 "DEBUG| Coordinate D(numerical) D(analytical) Difference Error [%]"
487 err(1:3) = 0.0_dp
488 DO k = 1, 3
489 IF (do_dof_dipole(k)) THEN
490 dd = dipole_moment(k) - dipole_numer(k)
491 IF (abs(dipole_moment(k)) > eps_no_error_check) THEN
492 derr = 100._dp*dd/dipole_moment(k)
493 WRITE (unit=iw, fmt="(T2,A,T13,A1,T21,F16.8,T38,F16.8,T56,G12.3,T72,F9.3)") &
494 "DEBUG|", achar(119 + k), dipole_numer(k), dipole_moment(k), dd, derr
495 ELSE
496 derr = 0.0_dp
497 WRITE (unit=iw, fmt="(T2,A,T13,A1,T21,F16.8,T38,F16.8,T56,G12.3)") &
498 "DEBUG|", achar(119 + k), dipole_numer(k), dipole_moment(k), dd
499 END IF
500 err(k) = derr
501 ELSE
502 WRITE (unit=iw, fmt="(T2,A,T13,A1,T21,A16,T38,F16.8)") &
503 "DEBUG|", achar(119 + k), " skipped", dipole_moment(k)
504 END IF
505 END DO
506 WRITE (unit=iw, fmt="((T2,A))") &
507 "DEBUG|========================================================================="
508 WRITE (unit=iw, fmt="(T2,A,T61,E20.12)") 'DIPOLE : CheckSum =', sum(dipole_moment)
509 IF (any(abs(err(1:3)) > maxerr)) THEN
510 message = "A mismatch between analytical and numerical dipoles "// &
511 "has been detected. Check the implementation of the "// &
512 "analytical dipole calculation"
513 IF (stop_on_mismatch) THEN
514 cpabort(message)
515 ELSE
516 cpwarn(message)
517 END IF
518 END IF
519 END IF
520
521 ELSE
522 CALL cp_warn(__location__, "Debug of dipole moments only for Quickstep code available")
523 END IF
524 END IF
525
526 IF (debug_polar) THEN
527 IF (force_env%in_use == use_qs_force) THEN
528 CALL get_qs_env(force_env%qs_env, dft_control=dft_control)
529 poldir = [0.0_dp, 0.0_dp, 1.0_dp]
530 amplitude = 0.0_dp
531 CALL set_efield(dft_control, amplitude, poldir)
532 CALL force_env_calc_energy_force(force_env, calc_force=.false.)
533 CALL get_qs_env(force_env%qs_env, results=results)
534 description = "[POLAR]"
535 IF (test_for_result(results, description=description)) THEN
536 CALL get_results(results, description=description, values=pvals)
537 polar_analytic(1:3, 1:3) = reshape(pvals(1:9), [3, 3])
538 ELSE
539 CALL cp_warn(__location__, "Debug of polarizabilities needs PROPERTIES/LINRES/POLAR section enabled")
540 cpabort("DEBUG")
541 END IF
542 description = "[DIPOLE]"
543 IF (.NOT. test_for_result(results, description=description)) THEN
544 CALL cp_warn(__location__, "Debug of polarizabilities need DFT/PRINT/MOMENTS section enabled")
545 cpabort("DEBUG")
546 END IF
547 amplitude = de
548 DO k = 1, 3
549 poldir = 0.0_dp
550 poldir(k) = 1.0_dp
551 DO j = 1, 2
552 poldir = -1.0_dp*poldir
553 CALL set_efield(dft_control, amplitude, poldir)
554 CALL force_env_calc_energy_force(force_env, calc_force=.false., linres=.true.)
555 CALL get_results(results, description=description, values=dipn(1:3, j))
556 END DO
557 polar_numeric(1:3, k) = 0.5_dp*(dipn(1:3, 2) - dipn(1:3, 1))/de
558 END DO
559 IF (iw > 0) THEN
560 WRITE (unit=iw, fmt="(/,(T2,A))") &
561 "DEBUG|========================= POLARIZABILITY ================================", &
562 "DEBUG| Coordinates P(numerical) P(analytical) Difference Error [%]"
563 DO k = 1, 3
564 DO j = 1, 3
565 dd = polar_analytic(k, j) - polar_numeric(k, j)
566 IF (abs(polar_analytic(k, j)) > eps_no_error_check) THEN
567 derr = 100._dp*dd/polar_analytic(k, j)
568 WRITE (unit=iw, fmt="(T2,A,T12,A1,A1,T21,F16.8,T38,F16.8,T56,G12.3,T72,F9.3)") &
569 "DEBUG|", achar(119 + k), achar(119 + j), polar_numeric(k, j), polar_analytic(k, j), dd, derr
570 ELSE
571 WRITE (unit=iw, fmt="(T2,A,T12,A1,A1,T21,F16.8,T38,F16.8,T56,G12.3)") &
572 "DEBUG|", achar(119 + k), achar(119 + j), polar_numeric(k, j), polar_analytic(k, j), dd
573 END IF
574 END DO
575 END DO
576 WRITE (unit=iw, fmt="((T2,A))") &
577 "DEBUG|========================================================================="
578 WRITE (unit=iw, fmt="(T2,A,T61,E20.12)") ' POLAR : CheckSum =', sum(polar_analytic)
579 END IF
580 ELSE
581 CALL cp_warn(__location__, "Debug of polarizabilities only for Quickstep code available")
582 END IF
583 END IF
584
585 CALL cp_print_key_finished_output(iw, logger, root_section, "DEBUG%PROGRAM_RUN_INFO")
586
587 END SUBROUTINE cp2k_debug_energy_and_forces
588
589! **************************************************************************************************
590!> \brief ...
591!> \param dft_control ...
592!> \param amplitude ...
593!> \param poldir ...
594! **************************************************************************************************
595 SUBROUTINE set_efield(dft_control, amplitude, poldir)
596 TYPE(dft_control_type), POINTER :: dft_control
597 REAL(kind=dp), INTENT(IN) :: amplitude
598 REAL(kind=dp), DIMENSION(3), INTENT(IN) :: poldir
599
600 IF (dft_control%apply_efield) THEN
601 dft_control%efield_fields(1)%efield%strength = amplitude
602 dft_control%efield_fields(1)%efield%polarisation(1:3) = poldir(1:3)
603 ELSEIF (dft_control%apply_period_efield) THEN
604 dft_control%period_efield%strength = amplitude
605 dft_control%period_efield%polarisation(1:3) = poldir(1:3)
606 ELSE
607 cpabort("No EFIELD definition available")
608 END IF
609
610 END SUBROUTINE set_efield
611
612END MODULE cp2k_debug
Definition atom.F:9
Handles all functions related to the CELL.
Definition cell_types.F:15
subroutine, public get_cell(cell, alpha, beta, gamma, deth, orthorhombic, abc, periodic, h, h_inv, symmetry_id, tag)
Get informations about a simulation cell.
Definition cell_types.F:200
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:77
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, cell_ref, use_ref_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, 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 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, mimic, 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, sab_cneo, 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, xcint_weights, 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, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_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, do_rixs, tb_tblite)
Get the QUICKSTEP environment.
Define the quickstep kind type and their sub types.
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)
...
Type defining parameters related to the simulation cell.
Definition cell_types.F:60
type of a logger, at the moment it contains just a print level starting at which level it should be l...
contains arbitrary information which need to be stored
represents a system: atoms, molecules, their pos,vel,...
wrapper to abstract the force evaluation of the various methods
stores all the informations relevant to an mpi environment
Provides all information about a quickstep kind.