(git:d18deda)
Loading...
Searching...
No Matches
kg_environment_types.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Types needed for a Kim-Gordon-like partitioning into molecular
10!> subunits
11!> \par History
12!> 2012.07 created [Martin Haeufel]
13!> \author Martin Haeufel
14! **************************************************************************************************
16 USE cp_dbcsr_api, ONLY: dbcsr_p_type
21 USE kinds, ONLY: dp
34#include "./base/base_uses.f90"
35
36 IMPLICIT NONE
37
38 PRIVATE
39
40 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'kg_environment_types'
41
43
44 TYPE subset_type
45 TYPE(neighbor_list_set_p_type), DIMENSION(:), POINTER :: sab_orb => null()
46 TYPE(task_list_type), POINTER :: task_list => null()
47 END TYPE subset_type
48
49! *****************************************************************************
50!> \brief Contains information on the energy correction functional for KG
51!> \par History
52!> 03.2014 created
53!> \author JGH
54! *****************************************************************************
56 CHARACTER(len=20) :: ec_name = ""
57 INTEGER :: energy_functional = -1
58 INTEGER :: ks_solver = -1
59 INTEGER :: factorization = -1
60 REAL(kind=dp) :: eps_default = 0.0_dp
61 ! basis set
62 CHARACTER(len=20) :: basis = ""
63 LOGICAL :: mao = .false.
64 INTEGER :: mao_max_iter = -1
65 REAL(kind=dp) :: mao_eps_grad = 0.0_dp
66 ! energy components
67 REAL(kind=dp) :: etotal = 0.0_dp
68 REAL(kind=dp) :: eband = 0.0_dp, exc = 0.0_dp, ehartree = 0.0_dp, vhxc = 0.0_dp
69 REAL(kind=dp) :: edispersion = 0.0_dp
70 ! full neighbor lists and corresponding task list
72 DIMENSION(:), POINTER :: sab_orb => null(), sac_ppl => null(), sap_ppnl => null()
73 TYPE(task_list_type), POINTER :: task_list => null()
74 ! the XC function to be used for the correction, dispersion info
75 TYPE(section_vals_type), POINTER :: xc_section => null()
76 TYPE(qs_dispersion_type), POINTER :: dispersion_env => null()
77 ! matrices in complete basis
78 ! KS: Kohn-Sham; H: Core; S: overlap; T: kinetic energy;
79 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_ks => null()
80 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_h => null()
81 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_s => null()
82 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_t => null()
83 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_p => null()
84 ! reduce basis
85 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: mao_coef => null()
87
88! **************************************************************************************************
89!> \brief Contains all the info needed for KG runs...
90!> \param xc_section_kg: XC section with only the KE functional
91!> \param molecule_set: set of molecular entities as in qs_env
92!> \param sab_orb_full: full neighborlist (build with molecular=.FALSE.)
93!> needed for the coloring
94!> \param subset_of_mol: ith entry contains the index of the subset, the ith
95!> molecule belongs to
96!> \param subset: task list and neighbor list of each subset of molecules
97!> \param nsubsets: number of subsets
98!> \par History
99!> 2012.07 created [Martin Haeufel]
100!> \author Martin Haeufel
101! **************************************************************************************************
103 INTEGER :: nspins = -1
104 INTEGER :: natom = -1
105 TYPE(section_vals_type), POINTER :: xc_section_kg => null()
106 INTEGER, ALLOCATABLE, DIMENSION(:) :: atom_to_molecule
107 TYPE(molecule_type), DIMENSION(:), POINTER :: molecule_set => null()
108 INTEGER :: tnadd_method = -1
110 DIMENSION(:), POINTER :: sab_orb_full => null(), sac_kin => null()
111 !
112 INTEGER, DIMENSION(:), POINTER :: subset_of_mol => null()
113 TYPE(subset_type), DIMENSION(:), POINTER :: subset => null()
114 INTEGER :: nsubsets = -1
115 INTEGER :: maxdegree = -1
116 INTEGER :: coloring_method = -1
117 !
118 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: tnadd_mat => null()
119 ! LRI
120 TYPE(lri_environment_type), POINTER :: lri_env => null(), lri_env1 => null()
121 TYPE(lri_density_type), POINTER :: lri_density => null(), lri_rho1 => null()
122 ! atomic grid
123 TYPE(atom_integration_grid_type), POINTER :: int_grid_atom => null()
124 TYPE(integration_grid_type), POINTER :: int_grid_molecules => null()
125 TYPE(integration_grid_type), POINTER :: int_grid_full => null()
126 END TYPE kg_environment_type
127
128CONTAINS
129
130! **************************************************************************************************
131!> \brief ...
132!> \param kg_env ...
133! **************************************************************************************************
134 SUBROUTINE kg_env_release(kg_env)
135 TYPE(kg_environment_type), POINTER :: kg_env
136
137 CHARACTER(LEN=*), PARAMETER :: routinen = 'kg_env_release'
138
139 INTEGER :: handle, isub
140
141 CALL timeset(routinen, handle)
142
143 cpassert(ASSOCIATED(kg_env))
144
145 CALL release_neighbor_list_sets(kg_env%sab_orb_full)
146 CALL release_neighbor_list_sets(kg_env%sac_kin)
147
148 IF (ASSOCIATED(kg_env%tnadd_mat)) THEN
149 CALL dbcsr_deallocate_matrix_set(kg_env%tnadd_mat)
150 END IF
151
152 DO isub = 1, kg_env%nsubsets
153 CALL release_neighbor_list_sets(kg_env%subset(isub)%sab_orb)
154 CALL deallocate_task_list(kg_env%subset(isub)%task_list)
155 END DO
156
157 IF (ASSOCIATED(kg_env%subset_of_mol)) DEALLOCATE (kg_env%subset_of_mol)
158 IF (ASSOCIATED(kg_env%subset)) DEALLOCATE (kg_env%subset)
159
160 IF (ALLOCATED(kg_env%atom_to_molecule)) DEALLOCATE (kg_env%atom_to_molecule)
161
162 ! LRI
163 IF (ASSOCIATED(kg_env%lri_env)) THEN
164 CALL lri_env_release(kg_env%lri_env)
165 DEALLOCATE (kg_env%lri_env)
166 END IF
167 IF (ASSOCIATED(kg_env%lri_density)) THEN
168 CALL lri_density_release(kg_env%lri_density)
169 DEALLOCATE (kg_env%lri_density)
170 END IF
171 IF (ASSOCIATED(kg_env%lri_env1)) THEN
172 CALL lri_env_release(kg_env%lri_env1)
173 DEALLOCATE (kg_env%lri_env1)
174 END IF
175 IF (ASSOCIATED(kg_env%lri_rho1)) THEN
176 CALL lri_density_release(kg_env%lri_rho1)
177 DEALLOCATE (kg_env%lri_rho1)
178 END IF
179 ! atom grids
180 IF (ASSOCIATED(kg_env%int_grid_atom)) THEN
181 CALL deallocate_atom_int_grid(kg_env%int_grid_atom)
182 END IF
183 IF (ASSOCIATED(kg_env%int_grid_molecules)) THEN
184 CALL deallocate_intgrid(kg_env%int_grid_molecules)
185 END IF
186 IF (ASSOCIATED(kg_env%int_grid_full)) THEN
187 CALL deallocate_intgrid(kg_env%int_grid_full)
188 END IF
189
190 DEALLOCATE (kg_env)
191
192 CALL timestop(handle)
193
194 END SUBROUTINE kg_env_release
195
196END MODULE kg_environment_types
DBCSR operations in CP2K.
objects that represent the structure of input sections and the data contained in an input section
subroutine, public deallocate_intgrid(int_grid)
Deallocate integration_grid_type.
Types needed for a Kim-Gordon-like partitioning into molecular subunits.
subroutine, public kg_env_release(kg_env)
...
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
contains the types and subroutines for dealing with the lri_env lri : local resolution of the identit...
subroutine, public lri_density_release(lri_density)
releases the given lri_density
subroutine, public lri_env_release(lri_env)
releases the given lri_env
Define the data structure for the molecule information.
Definition of disperson types for DFT calculations.
subroutine, public deallocate_atom_int_grid(int_grid)
...
Define the neighbor list data types and the corresponding functionality.
subroutine, public release_neighbor_list_sets(nlists)
releases an array of neighbor_list_sets
types for task lists
subroutine, public deallocate_task_list(task_list)
deallocates the components and the object itself
Contains information on the energy correction functional for KG.
Contains all the info needed for KG runs...