Line data Source code
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 Quickstep force driver routine
10 : !> \author MK (12.06.2002)
11 : ! **************************************************************************************************
12 : MODULE qs_force
13 : USE atomic_kind_types, ONLY: atomic_kind_type,&
14 : get_atomic_kind_set
15 : USE cp_control_types, ONLY: dft_control_type
16 : USE cp_dbcsr_api, ONLY: dbcsr_copy,&
17 : dbcsr_p_type,&
18 : dbcsr_set
19 : USE cp_dbcsr_operations, ONLY: dbcsr_allocate_matrix_set,&
20 : dbcsr_deallocate_matrix_set
21 : USE cp_dbcsr_output, ONLY: cp_dbcsr_write_sparse_matrix
22 : USE cp_log_handling, ONLY: cp_get_default_logger,&
23 : cp_logger_get_default_io_unit,&
24 : cp_logger_type
25 : USE cp_output_handling, ONLY: cp_p_file,&
26 : cp_print_key_finished_output,&
27 : cp_print_key_should_output,&
28 : cp_print_key_unit_nr
29 : USE dft_plus_u, ONLY: plus_u
30 : USE ec_env_types, ONLY: energy_correction_type
31 : USE efield_utils, ONLY: calculate_ecore_efield,&
32 : efield_potential_lengh_gauge
33 : USE energy_corrections, ONLY: energy_correction
34 : USE excited_states, ONLY: excited_state_energy
35 : USE hfx_exx, ONLY: calculate_exx
36 : USE input_constants, ONLY: ri_mp2_laplace,&
37 : ri_mp2_method_gpw,&
38 : ri_rpa_method_gpw
39 : USE input_section_types, ONLY: section_vals_get,&
40 : section_vals_get_subs_vals,&
41 : section_vals_type,&
42 : section_vals_val_get
43 : USE kinds, ONLY: dp
44 : USE lri_environment_types, ONLY: lri_environment_type
45 : USE message_passing, ONLY: mp_para_env_type
46 : USE mp2_cphf, ONLY: update_mp2_forces
47 : USE mulliken, ONLY: mulliken_restraint
48 : USE particle_types, ONLY: particle_type
49 : USE qs_core_energies, ONLY: calculate_ecore_overlap,&
50 : calculate_ecore_self
51 : USE qs_core_hamiltonian, ONLY: build_core_hamiltonian_matrix
52 : USE qs_dftb_dispersion, ONLY: calculate_dftb_dispersion
53 : USE qs_dftb_matrices, ONLY: build_dftb_matrices
54 : USE qs_energy, ONLY: qs_energies
55 : USE qs_energy_types, ONLY: qs_energy_type
56 : USE qs_environment_methods, ONLY: qs_env_rebuild_pw_env
57 : USE qs_environment_types, ONLY: get_qs_env,&
58 : qs_environment_type
59 : USE qs_external_potential, ONLY: external_c_potential,&
60 : external_e_potential
61 : USE qs_force_types, ONLY: allocate_qs_force,&
62 : qs_force_type,&
63 : replicate_qs_force,&
64 : zero_qs_force
65 : USE qs_ks_methods, ONLY: qs_ks_update_qs_env
66 : USE qs_ks_types, ONLY: qs_ks_env_type,&
67 : set_ks_env
68 : USE qs_rho_types, ONLY: qs_rho_get,&
69 : qs_rho_type
70 : USE qs_scf_post_scf, ONLY: qs_scf_compute_properties
71 : USE qs_subsys_types, ONLY: qs_subsys_set,&
72 : qs_subsys_type
73 : USE ri_environment_methods, ONLY: build_ri_matrices
74 : USE rt_propagation_forces, ONLY: calc_c_mat_force,&
75 : rt_admm_force
76 : USE rt_propagation_velocity_gauge, ONLY: velocity_gauge_ks_matrix,&
77 : velocity_gauge_nl_force
78 : USE se_core_core, ONLY: se_core_core_interaction
79 : USE se_core_matrix, ONLY: build_se_core_matrix
80 : USE virial_types, ONLY: symmetrize_virial,&
81 : virial_type
82 : USE xtb_matrices, ONLY: build_xtb_matrices
83 : #include "./base/base_uses.f90"
84 :
85 : IMPLICIT NONE
86 :
87 : PRIVATE
88 :
89 : ! *** Global parameters ***
90 :
91 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_force'
92 :
93 : ! *** Public subroutines ***
94 :
95 : PUBLIC :: qs_calc_energy_force
96 :
97 : CONTAINS
98 :
99 : ! **************************************************************************************************
100 : !> \brief ...
101 : !> \param qs_env ...
102 : !> \param calc_force ...
103 : !> \param consistent_energies ...
104 : !> \param linres ...
105 : ! **************************************************************************************************
106 21419 : SUBROUTINE qs_calc_energy_force(qs_env, calc_force, consistent_energies, linres)
107 : TYPE(qs_environment_type), POINTER :: qs_env
108 : LOGICAL :: calc_force, consistent_energies, linres
109 :
110 21419 : qs_env%linres_run = linres
111 21419 : IF (calc_force) THEN
112 10285 : CALL qs_forces(qs_env)
113 : ELSE
114 : CALL qs_energies(qs_env, calc_forces=.FALSE., &
115 11134 : consistent_energies=consistent_energies)
116 : END IF
117 :
118 21419 : END SUBROUTINE qs_calc_energy_force
119 :
120 : ! **************************************************************************************************
121 : !> \brief Calculate the Quickstep forces.
122 : !> \param qs_env ...
123 : !> \date 29.10.2002
124 : !> \author MK
125 : !> \version 1.0
126 : ! **************************************************************************************************
127 10285 : SUBROUTINE qs_forces(qs_env)
128 :
129 : TYPE(qs_environment_type), POINTER :: qs_env
130 :
131 : CHARACTER(len=*), PARAMETER :: routineN = 'qs_forces'
132 :
133 : INTEGER :: after, handle, i, iatom, ic, ikind, &
134 : ispin, iw, natom, nkind, nspin, &
135 : output_unit
136 10285 : INTEGER, ALLOCATABLE, DIMENSION(:) :: atom_of_kind, kind_of, natom_of_kind
137 : LOGICAL :: do_admm, do_exx, do_gw, do_im_time, &
138 : has_unit_metric, omit_headers, &
139 : perform_ec, reuse_hfx
140 : REAL(dp) :: dummy_real, dummy_real2(2)
141 10285 : TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
142 : TYPE(cp_logger_type), POINTER :: logger
143 10285 : TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_s, matrix_w, rho_ao
144 10285 : TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_w_kp
145 : TYPE(dft_control_type), POINTER :: dft_control
146 : TYPE(energy_correction_type), POINTER :: ec_env
147 : TYPE(lri_environment_type), POINTER :: lri_env
148 : TYPE(mp_para_env_type), POINTER :: para_env
149 10285 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
150 : TYPE(qs_energy_type), POINTER :: energy
151 10285 : TYPE(qs_force_type), DIMENSION(:), POINTER :: force
152 : TYPE(qs_ks_env_type), POINTER :: ks_env
153 : TYPE(qs_rho_type), POINTER :: rho
154 : TYPE(qs_subsys_type), POINTER :: subsys
155 : TYPE(section_vals_type), POINTER :: hfx_sections, print_section
156 : TYPE(virial_type), POINTER :: virial
157 :
158 10285 : CALL timeset(routineN, handle)
159 10285 : NULLIFY (logger)
160 10285 : logger => cp_get_default_logger()
161 :
162 : ! rebuild plane wave environment
163 10285 : CALL qs_env_rebuild_pw_env(qs_env)
164 :
165 : ! zero out the forces in particle set
166 10285 : CALL get_qs_env(qs_env, particle_set=particle_set)
167 10285 : natom = SIZE(particle_set)
168 76200 : DO iatom = 1, natom
169 273945 : particle_set(iatom)%f = 0.0_dp
170 : END DO
171 :
172 : ! get atom mapping
173 10285 : NULLIFY (atomic_kind_set)
174 10285 : CALL get_qs_env(qs_env, atomic_kind_set=atomic_kind_set)
175 : CALL get_atomic_kind_set(atomic_kind_set=atomic_kind_set, &
176 : atom_of_kind=atom_of_kind, &
177 10285 : kind_of=kind_of)
178 :
179 10285 : NULLIFY (force, subsys, dft_control)
180 : CALL get_qs_env(qs_env, &
181 : force=force, &
182 : subsys=subsys, &
183 10285 : dft_control=dft_control)
184 10285 : IF (.NOT. ASSOCIATED(force)) THEN
185 : ! *** Allocate the force data structure ***
186 2823 : nkind = SIZE(atomic_kind_set)
187 2823 : CALL get_atomic_kind_set(atomic_kind_set=atomic_kind_set, natom_of_kind=natom_of_kind)
188 2823 : CALL allocate_qs_force(force, natom_of_kind)
189 2823 : DEALLOCATE (natom_of_kind)
190 2823 : CALL qs_subsys_set(subsys, force=force)
191 : END IF
192 10285 : CALL zero_qs_force(force)
193 :
194 : ! Check if CDFT potential is needed and save it until forces have been calculated
195 10285 : IF (dft_control%qs_control%cdft) THEN
196 118 : dft_control%qs_control%cdft_control%save_pot = .TRUE.
197 : END IF
198 :
199 : ! recalculate energy with forces
200 10285 : CALL qs_energies(qs_env, calc_forces=.TRUE.)
201 :
202 10285 : NULLIFY (para_env)
203 : CALL get_qs_env(qs_env, &
204 10285 : para_env=para_env)
205 :
206 : ! Now we handle some special cases
207 : ! Maybe some of these would be better dealt with in qs_energies?
208 10285 : IF (qs_env%run_rtp) THEN
209 1222 : NULLIFY (matrix_w, matrix_s, ks_env)
210 : CALL get_qs_env(qs_env, &
211 : ks_env=ks_env, &
212 : matrix_w=matrix_w, &
213 1222 : matrix_s=matrix_s)
214 1222 : CALL dbcsr_allocate_matrix_set(matrix_w, dft_control%nspins)
215 2686 : DO ispin = 1, dft_control%nspins
216 1464 : ALLOCATE (matrix_w(ispin)%matrix)
217 : CALL dbcsr_copy(matrix_w(ispin)%matrix, matrix_s(1)%matrix, &
218 1464 : name="W MATRIX")
219 2686 : CALL dbcsr_set(matrix_w(ispin)%matrix, 0.0_dp)
220 : END DO
221 1222 : CALL set_ks_env(ks_env, matrix_w=matrix_w)
222 :
223 1222 : CALL calc_c_mat_force(qs_env)
224 1222 : IF (dft_control%do_admm) CALL rt_admm_force(qs_env)
225 1222 : IF (dft_control%rtp_control%velocity_gauge .AND. dft_control%rtp_control%nl_gauge_transform) &
226 22 : CALL velocity_gauge_nl_force(qs_env, particle_set)
227 : END IF
228 : ! from an eventual Mulliken restraint
229 10285 : IF (dft_control%qs_control%mulliken_restraint) THEN
230 6 : NULLIFY (matrix_w, matrix_s, rho)
231 : CALL get_qs_env(qs_env, &
232 : matrix_w=matrix_w, &
233 : matrix_s=matrix_s, &
234 6 : rho=rho)
235 6 : NULLIFY (rho_ao)
236 6 : CALL qs_rho_get(rho, rho_ao=rho_ao)
237 : CALL mulliken_restraint(dft_control%qs_control%mulliken_restraint_control, &
238 6 : para_env, matrix_s(1)%matrix, rho_ao, w_matrix=matrix_w)
239 : END IF
240 : ! Add non-Pulay contribution of DFT+U to W matrix, since it has also to be
241 : ! digested with overlap matrix derivatives
242 10285 : IF (dft_control%dft_plus_u) THEN
243 72 : NULLIFY (matrix_w)
244 72 : CALL get_qs_env(qs_env, matrix_w=matrix_w)
245 72 : CALL plus_u(qs_env=qs_env, matrix_w=matrix_w)
246 : END IF
247 :
248 : ! Write W Matrix to output (if requested)
249 10285 : CALL get_qs_env(qs_env, has_unit_metric=has_unit_metric)
250 10285 : IF (.NOT. has_unit_metric) THEN
251 7283 : NULLIFY (matrix_w_kp)
252 7283 : CALL get_qs_env(qs_env, matrix_w_kp=matrix_w_kp)
253 7283 : nspin = SIZE(matrix_w_kp, 1)
254 15566 : DO ispin = 1, nspin
255 8283 : IF (BTEST(cp_print_key_should_output(logger%iter_info, &
256 7283 : qs_env%input, "DFT%PRINT%AO_MATRICES/W_MATRIX"), cp_p_file)) THEN
257 : iw = cp_print_key_unit_nr(logger, qs_env%input, "DFT%PRINT%AO_MATRICES/W_MATRIX", &
258 8 : extension=".Log")
259 8 : CALL section_vals_val_get(qs_env%input, "DFT%PRINT%AO_MATRICES%NDIGITS", i_val=after)
260 8 : CALL section_vals_val_get(qs_env%input, "DFT%PRINT%AO_MATRICES%OMIT_HEADERS", l_val=omit_headers)
261 8 : after = MIN(MAX(after, 1), 16)
262 16 : DO ic = 1, SIZE(matrix_w_kp, 2)
263 : CALL cp_dbcsr_write_sparse_matrix(matrix_w_kp(ispin, ic)%matrix, 4, after, qs_env, &
264 16 : para_env, output_unit=iw, omit_headers=omit_headers)
265 : END DO
266 : CALL cp_print_key_finished_output(iw, logger, qs_env%input, &
267 8 : "DFT%PRINT%AO_MATRICES/W_MATRIX")
268 : END IF
269 : END DO
270 : END IF
271 :
272 : ! Check if energy correction should be skipped
273 10285 : perform_ec = .FALSE.
274 10285 : IF (qs_env%energy_correction) THEN
275 436 : CALL get_qs_env(qs_env, ec_env=ec_env)
276 436 : IF (.NOT. ec_env%do_skip) perform_ec = .TRUE.
277 : END IF
278 :
279 : ! Compute core forces (also overwrites matrix_w)
280 10285 : IF (dft_control%qs_control%semi_empirical) THEN
281 : CALL build_se_core_matrix(qs_env=qs_env, para_env=para_env, &
282 3002 : calculate_forces=.TRUE.)
283 3002 : CALL se_core_core_interaction(qs_env, para_env, calculate_forces=.TRUE.)
284 7283 : ELSEIF (dft_control%qs_control%dftb) THEN
285 : CALL build_dftb_matrices(qs_env=qs_env, para_env=para_env, &
286 724 : calculate_forces=.TRUE.)
287 : CALL calculate_dftb_dispersion(qs_env=qs_env, para_env=para_env, &
288 724 : calculate_forces=.TRUE.)
289 6559 : ELSEIF (dft_control%qs_control%xtb) THEN
290 456 : CALL build_xtb_matrices(qs_env=qs_env, calculate_forces=.TRUE.)
291 6103 : ELSEIF (perform_ec) THEN
292 : ! Calculates core and grid based forces
293 436 : CALL energy_correction(qs_env, ec_init=.FALSE., calculate_forces=.TRUE.)
294 : ELSE
295 : ! Dispersion energy and forces are calculated in qs_energy?
296 5667 : CALL build_core_hamiltonian_matrix(qs_env=qs_env, calculate_forces=.TRUE.)
297 : ! The above line reset the core H, which should be re-updated in case a TD field is applied:
298 5667 : IF (qs_env%run_rtp) THEN
299 810 : IF (dft_control%apply_efield_field) &
300 160 : CALL efield_potential_lengh_gauge(qs_env)
301 810 : IF (dft_control%rtp_control%velocity_gauge) &
302 22 : CALL velocity_gauge_ks_matrix(qs_env, subtract_nl_term=.FALSE.)
303 :
304 : END IF
305 5667 : CALL calculate_ecore_self(qs_env)
306 5667 : CALL calculate_ecore_overlap(qs_env, para_env, calculate_forces=.TRUE.)
307 5667 : CALL calculate_ecore_efield(qs_env, calculate_forces=.TRUE.)
308 : !swap external_e_potential before external_c_potential, to ensure
309 : !that external potential on grid is loaded before calculating energy of cores
310 5667 : CALL external_e_potential(qs_env)
311 5667 : IF (.NOT. dft_control%qs_control%gapw) THEN
312 5247 : CALL external_c_potential(qs_env, calculate_forces=.TRUE.)
313 : END IF
314 : ! RIGPW matrices
315 5667 : IF (dft_control%qs_control%rigpw) THEN
316 0 : CALL get_qs_env(qs_env=qs_env, lri_env=lri_env)
317 0 : CALL build_ri_matrices(lri_env, qs_env, calculate_forces=.TRUE.)
318 : END IF
319 : END IF
320 :
321 : ! MP2 Code
322 10285 : IF (ASSOCIATED(qs_env%mp2_env)) THEN
323 314 : NULLIFY (energy)
324 314 : CALL get_qs_env(qs_env, energy=energy)
325 314 : CALL qs_scf_compute_properties(qs_env, wf_type='MP2 ', do_mp2=.TRUE.)
326 314 : CALL qs_ks_update_qs_env(qs_env, just_energy=.TRUE.)
327 314 : energy%total = energy%total + energy%mp2
328 :
329 : IF ((qs_env%mp2_env%method == ri_mp2_method_gpw .OR. qs_env%mp2_env%method == ri_mp2_laplace .OR. &
330 : qs_env%mp2_env%method == ri_rpa_method_gpw) &
331 314 : .AND. .NOT. qs_env%mp2_env%do_im_time) THEN
332 264 : CALL update_mp2_forces(qs_env)
333 : END IF
334 :
335 : !RPA EXX energy and forces
336 314 : IF (qs_env%mp2_env%method == ri_rpa_method_gpw) THEN
337 : do_exx = .FALSE.
338 48 : hfx_sections => section_vals_get_subs_vals(qs_env%input, "DFT%XC%WF_CORRELATION%RI_RPA%HF")
339 48 : CALL section_vals_get(hfx_sections, explicit=do_exx)
340 48 : IF (do_exx) THEN
341 26 : do_gw = qs_env%mp2_env%ri_rpa%do_ri_g0w0
342 26 : do_admm = qs_env%mp2_env%ri_rpa%do_admm
343 26 : reuse_hfx = qs_env%mp2_env%ri_rpa%reuse_hfx
344 26 : do_im_time = qs_env%mp2_env%do_im_time
345 26 : output_unit = cp_logger_get_default_io_unit()
346 26 : dummy_real = 0.0_dp
347 :
348 : CALL calculate_exx(qs_env=qs_env, &
349 : unit_nr=output_unit, &
350 : hfx_sections=hfx_sections, &
351 : x_data=qs_env%mp2_env%ri_rpa%x_data, &
352 : do_gw=do_gw, &
353 : do_admm=do_admm, &
354 : calc_forces=.TRUE., &
355 : reuse_hfx=reuse_hfx, &
356 : do_im_time=do_im_time, &
357 : E_ex_from_GW=dummy_real, &
358 : E_admm_from_GW=dummy_real2, &
359 26 : t3=dummy_real)
360 : END IF
361 : END IF
362 9971 : ELSEIF (perform_ec) THEN
363 : ! energy correction forces postponed
364 9535 : ELSEIF (qs_env%harris_method) THEN
365 : ! Harris method forces already done in harris_energy_correction
366 : ELSE
367 : ! Compute grid-based forces
368 9531 : CALL qs_ks_update_qs_env(qs_env, calculate_forces=.TRUE.)
369 : END IF
370 :
371 : ! Excited state forces
372 10285 : CALL excited_state_energy(qs_env, calculate_forces=.TRUE.)
373 :
374 : ! replicate forces (get current pointer)
375 10285 : NULLIFY (force)
376 10285 : CALL get_qs_env(qs_env=qs_env, force=force)
377 10285 : CALL replicate_qs_force(force, para_env)
378 :
379 76200 : DO iatom = 1, natom
380 65915 : ikind = kind_of(iatom)
381 65915 : i = atom_of_kind(iatom)
382 : ! XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX
383 : ! the force is - dE/dR, what is called force is actually the gradient
384 : ! Things should have the right name
385 : ! The minus sign below is a hack
386 : ! XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX
387 527320 : force(ikind)%other(1:3, i) = -particle_set(iatom)%f(1:3) + force(ikind)%ch_pulay(1:3, i)
388 263660 : force(ikind)%total(1:3, i) = force(ikind)%total(1:3, i) + force(ikind)%other(1:3, i)
389 537605 : particle_set(iatom)%f = -force(ikind)%total(1:3, i)
390 : END DO
391 :
392 10285 : NULLIFY (virial, energy)
393 10285 : CALL get_qs_env(qs_env=qs_env, virial=virial, energy=energy)
394 10285 : IF (virial%pv_availability) THEN
395 800 : CALL para_env%sum(virial%pv_overlap)
396 800 : CALL para_env%sum(virial%pv_ekinetic)
397 800 : CALL para_env%sum(virial%pv_ppl)
398 800 : CALL para_env%sum(virial%pv_ppnl)
399 800 : CALL para_env%sum(virial%pv_ecore_overlap)
400 800 : CALL para_env%sum(virial%pv_ehartree)
401 800 : CALL para_env%sum(virial%pv_exc)
402 800 : CALL para_env%sum(virial%pv_exx)
403 800 : CALL para_env%sum(virial%pv_vdw)
404 800 : CALL para_env%sum(virial%pv_mp2)
405 800 : CALL para_env%sum(virial%pv_nlcc)
406 800 : CALL para_env%sum(virial%pv_gapw)
407 800 : CALL para_env%sum(virial%pv_lrigpw)
408 800 : CALL para_env%sum(virial%pv_virial)
409 800 : CALL symmetrize_virial(virial)
410 : ! Add the volume terms of the virial
411 800 : IF ((.NOT. virial%pv_numer) .AND. &
412 : (.NOT. (dft_control%qs_control%dftb .OR. &
413 : dft_control%qs_control%xtb .OR. &
414 : dft_control%qs_control%semi_empirical))) THEN
415 :
416 : ! Harris energy correction requires volume terms from
417 : ! 1) Harris functional contribution, and
418 : ! 2) Linear Response solver
419 600 : IF (perform_ec) THEN
420 168 : CALL get_qs_env(qs_env, ec_env=ec_env)
421 168 : energy%hartree = ec_env%ehartree
422 168 : energy%exc = ec_env%exc
423 168 : IF (dft_control%do_admm) THEN
424 38 : energy%exc_aux_fit = ec_env%exc_aux_fit
425 : END IF
426 : END IF
427 2400 : DO i = 1, 3
428 : virial%pv_ehartree(i, i) = virial%pv_ehartree(i, i) &
429 1800 : - 2.0_dp*(energy%hartree + energy%sccs_pol)
430 : virial%pv_virial(i, i) = virial%pv_virial(i, i) - energy%exc &
431 1800 : - 2.0_dp*(energy%hartree + energy%sccs_pol)
432 1800 : virial%pv_exc(i, i) = virial%pv_exc(i, i) - energy%exc
433 2400 : IF (dft_control%do_admm) THEN
434 198 : virial%pv_exc(i, i) = virial%pv_exc(i, i) - energy%exc_aux_fit
435 198 : virial%pv_virial(i, i) = virial%pv_virial(i, i) - energy%exc_aux_fit
436 : END IF
437 : ! The factor 2 is a hack. It compensates the plus sign in h_stress/pw_poisson_solve.
438 : ! The sign in pw_poisson_solve is correct for FIST, but not for QS.
439 : ! There should be a more elegant solution to that ...
440 : END DO
441 : END IF
442 : END IF
443 :
444 : output_unit = cp_print_key_unit_nr(logger, qs_env%input, "DFT%PRINT%DERIVATIVES", &
445 10285 : extension=".Log")
446 10285 : print_section => section_vals_get_subs_vals(qs_env%input, "DFT%PRINT%DERIVATIVES")
447 10285 : IF (dft_control%qs_control%semi_empirical) THEN
448 : CALL write_forces(force, atomic_kind_set, 2, output_unit=output_unit, &
449 3002 : print_section=print_section)
450 7283 : ELSE IF (dft_control%qs_control%dftb) THEN
451 : CALL write_forces(force, atomic_kind_set, 4, output_unit=output_unit, &
452 724 : print_section=print_section)
453 6559 : ELSE IF (dft_control%qs_control%xtb) THEN
454 : CALL write_forces(force, atomic_kind_set, 4, output_unit=output_unit, &
455 456 : print_section=print_section)
456 6103 : ELSE IF (dft_control%qs_control%gapw .OR. dft_control%qs_control%gapw_xc) THEN
457 : CALL write_forces(force, atomic_kind_set, 1, output_unit=output_unit, &
458 492 : print_section=print_section)
459 : ELSE
460 : CALL write_forces(force, atomic_kind_set, 0, output_unit=output_unit, &
461 5611 : print_section=print_section)
462 : END IF
463 : CALL cp_print_key_finished_output(output_unit, logger, qs_env%input, &
464 10285 : "DFT%PRINT%DERIVATIVES")
465 :
466 : ! deallocate W Matrix:
467 10285 : NULLIFY (ks_env, matrix_w_kp)
468 : CALL get_qs_env(qs_env=qs_env, &
469 : matrix_w_kp=matrix_w_kp, &
470 10285 : ks_env=ks_env)
471 10285 : CALL dbcsr_deallocate_matrix_set(matrix_w_kp)
472 10285 : NULLIFY (matrix_w_kp)
473 10285 : CALL set_ks_env(ks_env, matrix_w_kp=matrix_w_kp)
474 :
475 10285 : DEALLOCATE (atom_of_kind, kind_of)
476 :
477 10285 : CALL timestop(handle)
478 :
479 20570 : END SUBROUTINE qs_forces
480 :
481 : ! **************************************************************************************************
482 : !> \brief Write a Quickstep force data structure to output unit
483 : !> \param qs_force ...
484 : !> \param atomic_kind_set ...
485 : !> \param ftype ...
486 : !> \param output_unit ...
487 : !> \param print_section ...
488 : !> \date 05.06.2002
489 : !> \author MK
490 : !> \version 1.0
491 : ! **************************************************************************************************
492 10285 : SUBROUTINE write_forces(qs_force, atomic_kind_set, ftype, output_unit, &
493 : print_section)
494 :
495 : TYPE(qs_force_type), DIMENSION(:), POINTER :: qs_force
496 : TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
497 : INTEGER, INTENT(IN) :: ftype, output_unit
498 : TYPE(section_vals_type), POINTER :: print_section
499 :
500 : CHARACTER(LEN=13) :: fmtstr5
501 : CHARACTER(LEN=15) :: fmtstr4
502 : CHARACTER(LEN=20) :: fmtstr3
503 : CHARACTER(LEN=35) :: fmtstr2
504 : CHARACTER(LEN=48) :: fmtstr1
505 : INTEGER :: i, iatom, ikind, my_ftype, natom, ndigits
506 10285 : INTEGER, ALLOCATABLE, DIMENSION(:) :: atom_of_kind, kind_of
507 : REAL(KIND=dp), DIMENSION(3) :: grand_total
508 :
509 10285 : IF (output_unit > 0) THEN
510 :
511 154 : IF (.NOT. ASSOCIATED(qs_force)) THEN
512 : CALL cp_abort(__LOCATION__, &
513 : "The qs_force pointer is not associated "// &
514 0 : "and cannot be printed")
515 : END IF
516 :
517 : CALL get_atomic_kind_set(atomic_kind_set=atomic_kind_set, atom_of_kind=atom_of_kind, &
518 154 : kind_of=kind_of, natom=natom)
519 :
520 : ! Variable precision output of the forces
521 : CALL section_vals_val_get(print_section, "NDIGITS", &
522 154 : i_val=ndigits)
523 :
524 154 : fmtstr1 = "(/,/,T2,A,/,/,T3,A,T11,A,T23,A,T40,A1,2( X,A1))"
525 154 : WRITE (UNIT=fmtstr1(41:42), FMT="(I2)") ndigits + 5
526 :
527 154 : fmtstr2 = "(/,(T2,I5,4X,I4,T18,A,T34,3F . ))"
528 154 : WRITE (UNIT=fmtstr2(32:33), FMT="(I2)") ndigits
529 154 : WRITE (UNIT=fmtstr2(29:30), FMT="(I2)") ndigits + 6
530 :
531 154 : fmtstr3 = "(/,T3,A,T34,3F . )"
532 154 : WRITE (UNIT=fmtstr3(18:19), FMT="(I2)") ndigits
533 154 : WRITE (UNIT=fmtstr3(15:16), FMT="(I2)") ndigits + 6
534 :
535 154 : fmtstr4 = "((T34,3F . ))"
536 154 : WRITE (UNIT=fmtstr4(12:13), FMT="(I2)") ndigits
537 154 : WRITE (UNIT=fmtstr4(9:10), FMT="(I2)") ndigits + 6
538 :
539 : fmtstr5 = "(/T2,A//T3,A)"
540 :
541 : WRITE (UNIT=output_unit, FMT=fmtstr1) &
542 154 : "FORCES [a.u.]", "Atom", "Kind", "Component", "X", "Y", "Z"
543 :
544 154 : grand_total(:) = 0.0_dp
545 :
546 154 : my_ftype = ftype
547 :
548 0 : SELECT CASE (my_ftype)
549 : CASE DEFAULT
550 0 : DO iatom = 1, natom
551 0 : ikind = kind_of(iatom)
552 0 : i = atom_of_kind(iatom)
553 : WRITE (UNIT=output_unit, FMT=fmtstr2) &
554 0 : iatom, ikind, " total", qs_force(ikind)%total(1:3, i)
555 0 : grand_total(1:3) = grand_total(1:3) + qs_force(ikind)%total(1:3, i)
556 : END DO
557 : CASE (0)
558 476 : DO iatom = 1, natom
559 342 : ikind = kind_of(iatom)
560 342 : i = atom_of_kind(iatom)
561 : WRITE (UNIT=output_unit, FMT=fmtstr2) &
562 1368 : iatom, ikind, " overlap", qs_force(ikind)%overlap(1:3, i), &
563 1368 : iatom, ikind, " overlap_admm", qs_force(ikind)%overlap_admm(1:3, i), &
564 1368 : iatom, ikind, " kinetic", qs_force(ikind)%kinetic(1:3, i), &
565 1368 : iatom, ikind, " gth_ppl", qs_force(ikind)%gth_ppl(1:3, i), &
566 1368 : iatom, ikind, " gth_nlcc", qs_force(ikind)%gth_nlcc(1:3, i), &
567 1368 : iatom, ikind, " gth_ppnl", qs_force(ikind)%gth_ppnl(1:3, i), &
568 1368 : iatom, ikind, " core_overlap", qs_force(ikind)%core_overlap(1:3, i), &
569 1368 : iatom, ikind, " rho_core", qs_force(ikind)%rho_core(1:3, i), &
570 1368 : iatom, ikind, " rho_elec", qs_force(ikind)%rho_elec(1:3, i), &
571 1368 : iatom, ikind, " rho_lri_elec", qs_force(ikind)%rho_lri_elec(1:3, i), &
572 1368 : iatom, ikind, " ch_pulay", qs_force(ikind)%ch_pulay(1:3, i), &
573 1368 : iatom, ikind, " dispersion", qs_force(ikind)%dispersion(1:3, i), &
574 1368 : iatom, ikind, " gCP", qs_force(ikind)%gcp(1:3, i), &
575 1368 : iatom, ikind, " other", qs_force(ikind)%other(1:3, i), &
576 1368 : iatom, ikind, " fock_4c", qs_force(ikind)%fock_4c(1:3, i), &
577 1368 : iatom, ikind, " ehrenfest", qs_force(ikind)%ehrenfest(1:3, i), &
578 1368 : iatom, ikind, " efield", qs_force(ikind)%efield(1:3, i), &
579 1368 : iatom, ikind, " eev", qs_force(ikind)%eev(1:3, i), &
580 1368 : iatom, ikind, " mp2_non_sep", qs_force(ikind)%mp2_non_sep(1:3, i), &
581 1710 : iatom, ikind, " total", qs_force(ikind)%total(1:3, i)
582 1502 : grand_total(1:3) = grand_total(1:3) + qs_force(ikind)%total(1:3, i)
583 : END DO
584 : CASE (1)
585 0 : DO iatom = 1, natom
586 0 : ikind = kind_of(iatom)
587 0 : i = atom_of_kind(iatom)
588 : WRITE (UNIT=output_unit, FMT=fmtstr2) &
589 0 : iatom, ikind, " overlap", qs_force(ikind)%overlap(1:3, i), &
590 0 : iatom, ikind, " overlap_admm", qs_force(ikind)%overlap_admm(1:3, i), &
591 0 : iatom, ikind, " kinetic", qs_force(ikind)%kinetic(1:3, i), &
592 0 : iatom, ikind, " gth_ppl", qs_force(ikind)%gth_ppl(1:3, i), &
593 0 : iatom, ikind, " gth_nlcc", qs_force(ikind)%gth_nlcc(1:3, i), &
594 0 : iatom, ikind, " gth_ppnl", qs_force(ikind)%gth_ppnl(1:3, i), &
595 0 : iatom, ikind, " all_potential", qs_force(ikind)%all_potential(1:3, i), &
596 0 : iatom, ikind, " core_overlap", qs_force(ikind)%core_overlap(1:3, i), &
597 0 : iatom, ikind, " rho_core", qs_force(ikind)%rho_core(1:3, i), &
598 0 : iatom, ikind, " rho_elec", qs_force(ikind)%rho_elec(1:3, i), &
599 0 : iatom, ikind, " rho_lri_elec", qs_force(ikind)%rho_lri_elec(1:3, i), &
600 0 : iatom, ikind, " vhxc_atom", qs_force(ikind)%vhxc_atom(1:3, i), &
601 0 : iatom, ikind, " g0s_Vh_elec", qs_force(ikind)%g0s_Vh_elec(1:3, i), &
602 0 : iatom, ikind, " ch_pulay", qs_force(ikind)%ch_pulay(1:3, i), &
603 0 : iatom, ikind, " dispersion", qs_force(ikind)%dispersion(1:3, i), &
604 0 : iatom, ikind, " gCP", qs_force(ikind)%gcp(1:3, i), &
605 0 : iatom, ikind, " fock_4c", qs_force(ikind)%fock_4c(1:3, i), &
606 0 : iatom, ikind, " ehrenfest", qs_force(ikind)%ehrenfest(1:3, i), &
607 0 : iatom, ikind, " efield", qs_force(ikind)%efield(1:3, i), &
608 0 : iatom, ikind, " eev", qs_force(ikind)%eev(1:3, i), &
609 0 : iatom, ikind, " mp2_non_sep", qs_force(ikind)%mp2_non_sep(1:3, i), &
610 0 : iatom, ikind, " total", qs_force(ikind)%total(1:3, i)
611 0 : grand_total(1:3) = grand_total(1:3) + qs_force(ikind)%total(1:3, i)
612 : END DO
613 : CASE (2)
614 75 : DO iatom = 1, natom
615 73 : ikind = kind_of(iatom)
616 73 : i = atom_of_kind(iatom)
617 : WRITE (UNIT=output_unit, FMT=fmtstr2) &
618 292 : iatom, ikind, " all_potential", qs_force(ikind)%all_potential(1:3, i), &
619 292 : iatom, ikind, " rho_elec", qs_force(ikind)%rho_elec(1:3, i), &
620 365 : iatom, ikind, " total", qs_force(ikind)%total(1:3, i)
621 294 : grand_total(1:3) = grand_total(1:3) + qs_force(ikind)%total(1:3, i)
622 : END DO
623 : CASE (3)
624 0 : DO iatom = 1, natom
625 0 : ikind = kind_of(iatom)
626 0 : i = atom_of_kind(iatom)
627 : WRITE (UNIT=output_unit, FMT=fmtstr2) &
628 0 : iatom, ikind, " overlap", qs_force(ikind)%overlap(1:3, i), &
629 0 : iatom, ikind, "overlap_admm", qs_force(ikind)%overlap_admm(1:3, i), &
630 0 : iatom, ikind, " kinetic", qs_force(ikind)%kinetic(1:3, i), &
631 0 : iatom, ikind, " gth_ppl", qs_force(ikind)%gth_ppl(1:3, i), &
632 0 : iatom, ikind, " gth_nlcc", qs_force(ikind)%gth_nlcc(1:3, i), &
633 0 : iatom, ikind, " gth_ppnl", qs_force(ikind)%gth_ppnl(1:3, i), &
634 0 : iatom, ikind, " core_overlap", qs_force(ikind)%core_overlap(1:3, i), &
635 0 : iatom, ikind, " rho_core", qs_force(ikind)%rho_core(1:3, i), &
636 0 : iatom, ikind, " rho_elec", qs_force(ikind)%rho_elec(1:3, i), &
637 0 : iatom, ikind, " ch_pulay", qs_force(ikind)%ch_pulay(1:3, i), &
638 0 : iatom, ikind, " fock_4c", qs_force(ikind)%fock_4c(1:3, i), &
639 0 : iatom, ikind, " mp2_non_sep", qs_force(ikind)%mp2_non_sep(1:3, i), &
640 0 : iatom, ikind, " total", qs_force(ikind)%total(1:3, i)
641 0 : grand_total(1:3) = grand_total(1:3) + qs_force(ikind)%total(1:3, i)
642 : END DO
643 : CASE (4)
644 152 : DO iatom = 1, natom
645 134 : ikind = kind_of(iatom)
646 134 : i = atom_of_kind(iatom)
647 : WRITE (UNIT=output_unit, FMT=fmtstr2) &
648 536 : iatom, ikind, " all_potential", qs_force(ikind)%all_potential(1:3, i), &
649 536 : iatom, ikind, " overlap", qs_force(ikind)%overlap(1:3, i), &
650 536 : iatom, ikind, " rho_elec", qs_force(ikind)%rho_elec(1:3, i), &
651 536 : iatom, ikind, " repulsive", qs_force(ikind)%repulsive(1:3, i), &
652 536 : iatom, ikind, " dispersion", qs_force(ikind)%dispersion(1:3, i), &
653 536 : iatom, ikind, " efield", qs_force(ikind)%efield(1:3, i), &
654 536 : iatom, ikind, " ehrenfest", qs_force(ikind)%ehrenfest(1:3, i), &
655 670 : iatom, ikind, " total", qs_force(ikind)%total(1:3, i)
656 554 : grand_total(1:3) = grand_total(1:3) + qs_force(ikind)%total(1:3, i)
657 : END DO
658 : CASE (5)
659 154 : DO iatom = 1, natom
660 0 : ikind = kind_of(iatom)
661 0 : i = atom_of_kind(iatom)
662 : WRITE (UNIT=output_unit, FMT=fmtstr2) &
663 0 : iatom, ikind, " overlap", qs_force(ikind)%overlap(1:3, i), &
664 0 : iatom, ikind, " kinetic", qs_force(ikind)%kinetic(1:3, i), &
665 0 : iatom, ikind, " rho_elec", qs_force(ikind)%rho_elec(1:3, i), &
666 0 : iatom, ikind, " dispersion", qs_force(ikind)%dispersion(1:3, i), &
667 0 : iatom, ikind, " all potential", qs_force(ikind)%all_potential(1:3, i), &
668 0 : iatom, ikind, " other", qs_force(ikind)%other(1:3, i), &
669 0 : iatom, ikind, " total", qs_force(ikind)%total(1:3, i)
670 0 : grand_total(1:3) = grand_total(1:3) + qs_force(ikind)%total(1:3, i)
671 : END DO
672 : END SELECT
673 :
674 154 : WRITE (UNIT=output_unit, FMT=fmtstr3) "Sum of total", grand_total(1:3)
675 :
676 154 : DEALLOCATE (atom_of_kind)
677 154 : DEALLOCATE (kind_of)
678 :
679 : END IF
680 :
681 10285 : END SUBROUTINE write_forces
682 :
683 : END MODULE qs_force
|