(git:374b731)
Loading...
Searching...
No Matches
se_fock_matrix_dbg.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
9 USE dbcsr_api, ONLY: dbcsr_dot,&
10 dbcsr_p_type,&
11 dbcsr_set
12 USE kinds, ONLY: dp
18#include "./base/base_uses.f90"
19
20 IMPLICIT NONE
21 PRIVATE
22
23 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'se_fock_matrix_dbg'
24 LOGICAL, PARAMETER, PRIVATE :: debug_this_module = .false.
25
26 PUBLIC :: dbg_energy_coulomb_lr
27
28CONTAINS
29
30! **************************************************************************************************
31!> \brief Debug routine for long-range energy (debug value of EWALD vs VALUE KS)
32!> \param energy ...
33!> \param ks_matrix ...
34!> \param nspins ...
35!> \param qs_env ...
36!> \param matrix_p ...
37!> \param calculate_forces ...
38!> \param store_int_env ...
39!> \author Teodoro Laino [tlaino] - 04.2009
40! **************************************************************************************************
41 SUBROUTINE dbg_energy_coulomb_lr(energy, ks_matrix, nspins, qs_env, matrix_p, &
42 calculate_forces, store_int_env)
43 TYPE(qs_energy_type), POINTER :: energy
44 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: ks_matrix
45 INTEGER, INTENT(IN) :: nspins
46 TYPE(qs_environment_type), POINTER :: qs_env
47 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_p
48 LOGICAL, INTENT(IN) :: calculate_forces
49 TYPE(semi_empirical_si_type), POINTER :: store_int_env
50
51 INTEGER :: ispin
52 REAL(kind=dp) :: ecoul
53
54! Zero structures only for debugging purpose
55
56 CALL init_qs_energy(energy)
57 DO ispin = 1, nspins
58 CALL dbcsr_set(ks_matrix(ispin)%matrix, 0.0_dp)
59 END DO
60
61 ! Evaluate Coulomb Long-Range
62 CALL build_fock_matrix_coulomb_lr(qs_env, ks_matrix, matrix_p, energy, calculate_forces, &
63 store_int_env)
64
65 ! Compute the Hartree energy
66 DO ispin = 1, nspins
67 CALL dbcsr_dot(ks_matrix(ispin)%matrix, matrix_p(ispin)%matrix, ecoul)
68 energy%hartree = energy%hartree + ecoul
69
70 WRITE (*, *) ispin, "ECOUL ", ecoul
71 END DO
72 WRITE (*, *) "ENUC in DBG:", energy%core_overlap
73
74 ! Debug statements
75 WRITE (*, *) "TOTAL ENE", 0.5_dp*energy%hartree + energy%core_overlap
76 cpabort("Debug energy for Coulomb Long-Range")
77
78 END SUBROUTINE dbg_energy_coulomb_lr
79
80END MODULE se_fock_matrix_dbg
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
subroutine, public init_qs_energy(qs_energy)
Initialise a Quickstep energy data structure.
Module that collects all Coulomb parts of the fock matrix construction.
subroutine, public build_fock_matrix_coulomb_lr(qs_env, ks_matrix, matrix_p, energy, calculate_forces, store_int_env)
Long-Range part for SE Coulomb interactions.
subroutine, public dbg_energy_coulomb_lr(energy, ks_matrix, nspins, qs_env, matrix_p, calculate_forces, store_int_env)
Debug routine for long-range energy (debug value of EWALD vs VALUE KS)
Type to store integrals for semi-empirical calculations.