Octopus
xc_interaction.F90
Go to the documentation of this file.
1!! Copyright (C) 2021 Nicolas Tancogne-Dejean
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22 use comm_oct_m
23 use debug_oct_m
26 use global_oct_m
29 use, intrinsic :: iso_fortran_env
30 use mesh_oct_m
37 use space_oct_m
39 use xc_cam_oct_m
40 use xc_f03_lib_m
42
43 implicit none
44
45 private
46 public :: &
51
53 private
54
55 contains
56 procedure :: init => xc_interaction_init
57 procedure :: calculate => xc_interaction_calculate
58 procedure :: calculate_energy => xc_interaction_calculate_energy
59 procedure :: end => xc_interaction_end
61 end type xc_interaction_t
62
63 interface xc_interaction_t
64 module procedure xc_interaction_constructor
65 end interface xc_interaction_t
66
67 integer, public, parameter :: &
68 FUNC_X = 1, &
69 func_c = 2
70
71
72contains
73
74 ! ---------------------------------------------------------
75 function xc_interaction_constructor(partner) result(this)
76 class(interaction_partner_t), target, intent(inout) :: partner
77 class(xc_interaction_t), pointer :: this
78
80
81 allocate(this)
82
83 this%label = "exchange-correlation interaction"
84
85 this%partner => partner
86
89
90 ! ---------------------------------------------------------
91 subroutine xc_interaction_init(this)
92 class(xc_interaction_t), intent(inout) :: this
93
94 push_sub(xc_interaction_init)
95
96 pop_sub(xc_interaction_init)
97 end subroutine xc_interaction_init
98
99 ! ---------------------------------------------------------
100 subroutine xc_interaction_calculate(this)
101 class(xc_interaction_t), intent(inout) :: this
102
104
105 this%density = m_zero
106
108 end subroutine xc_interaction_calculate
109
110 ! ---------------------------------------------------------
111 subroutine xc_interaction_calculate_energy(this)
112 class(xc_interaction_t), intent(inout) :: this
113
115
116 this%energy = m_zero
117
120
121
122 ! ---------------------------------------------------------
123 subroutine xc_interaction_compute(this)
124 class(xc_interaction_t), intent(inout) :: this
125 push_sub(xc_interaction_compute)
126
127 this%density = m_zero
128
130 end subroutine xc_interaction_compute
131
132 ! ---------------------------------------------------------
133 subroutine xc_interaction_end(this)
134 class(xc_interaction_t), intent(inout) :: this
135
136 push_sub(xc_interaction_end)
137
138 safe_deallocate_a(this%density)
139
140 call interaction_end(this)
141
142 pop_sub(xc_interaction_end)
143 end subroutine xc_interaction_end
144
146 ! ---------------------------------------------------------
147 subroutine xc_interaction_finalize(this)
148 type(xc_interaction_t), intent(inout) :: this
155 end subroutine xc_interaction_finalize
156
157 ! -----------------------------------------------------
158 subroutine calc_tb09_c(mesh, space, functl, dens, gdens, ispin, rcell_volume)
159 class(mesh_t), intent(in) :: mesh
160 class(space_t), intent(in) :: space
161 type(xc_functional_t), intent(inout) :: functl(:)
162 real(real64), intent(in) :: dens(:,:)
163 real(real64), intent(in) :: gdens(:,:,:)
164 integer, intent(in) :: ispin
165 real(real64), intent(in) :: rcell_volume
166
167 real(real64), allocatable :: gnon(:)
168 real(real64) :: gn(space%dim), n, parameters(1)
169 integer :: ii
170
171 push_sub(calc_tb09_c)
172
173 safe_allocate(gnon(1:mesh%np))
174
175 do ii = 1, mesh%np
176 if (ispin == unpolarized) then
177 n = dens(ii, 1)
178 gn(1:space%dim) = gdens(ii, 1:space%dim, 1)
179 else
180 n = dens(ii, 1) + dens(ii, 2)
181 gn(1:space%dim) = gdens(ii, 1:space%dim, 1) + gdens(ii, 1:space%dim, 2)
182 end if
183
184 if (n <= 1e-7_real64) then
185 gnon(ii) = m_zero
186 else
187 gnon(ii) = sqrt(sum((gn(1:space%dim)/n)**2))
188 end if
189 end do
190
191 parameters(1) = -0.012_real64 + 1.023_real64*sqrt(dmf_integrate(mesh, gnon)/rcell_volume)
192
193 call xc_f03_func_set_ext_params(functl(1)%conf, parameters)
194
195 safe_deallocate_a(gnon)
196
197 pop_sub(calc_tb09_c)
198 end subroutine calc_tb09_c
199
200 ! ---------------------------------------------------------
201 subroutine calc_mvorb_alpha(mesh, namespace, space, functl, dens, gdens, ispin, rcell_volume, cam)
202 class(mesh_t), intent(in) :: mesh
203 type(namespace_t), intent(in) :: namespace
204 class(space_t), intent(in) :: space
205 type(xc_functional_t), intent(inout) :: functl(:)
206 real(real64), intent(in) :: dens(:,:)
207 real(real64), intent(in) :: gdens(:,:,:)
208 integer, intent(in) :: ispin
209 real(real64), intent(in) :: rcell_volume
210 type(xc_cam_t), intent(inout) :: cam
211
212 real(real64), allocatable :: gnon(:)
213 real(real64) :: tb09_c, alpha
214 real(real64) :: gn(space%dim), n
215 integer :: ii
216 real(real64) :: parameters(3)
217 real(real64), parameter :: tol_small_dens = 1e-11_real64
218
219 push_sub(calc_mvorb_alpha)
220
221 safe_allocate(gnon(1:mesh%np))
222
223 do ii = 1, mesh%np
224 if (ispin == unpolarized) then
225 n = dens(ii, 1)
226 gn(1:space%dim) = gdens(ii, 1:space%dim, 1)
227 else
228 n = dens(ii, 1) + dens(ii, 2)
229 gn(1:space%dim) = gdens(ii, 1:space%dim, 1) + gdens(ii, 1:space%dim, 2)
230 end if
231
232 if (n <= tol_small_dens) then
233 gnon(ii) = m_zero
234 else
235 gnon(ii) = sqrt(sum((gn(1:space%dim)/n)**2))
236 gnon(ii) = sqrt(gnon(ii))
237 end if
238 end do
239
240 tb09_c = dmf_integrate(mesh, gnon)/rcell_volume
241
242 safe_deallocate_a(gnon)
243
244 select case (functl(func_c)%id)
246 alpha = 0.121983_real64+0.130711_real64*tb09_c**4
247
248 if (alpha > 1) then
249 write(message(1), '(a,f6.3,a)') 'MVORB mixing parameter bigger than one (' , alpha ,').'
250 call messages_warning(1, namespace=namespace)
251 alpha = 0.25_real64
252 else
253 write(message(1), '(a,f6.3,a)') 'MVORB mixing parameter is (' , alpha ,').'
254 call messages_info(1, namespace=namespace, debug_only=.true.)
255 end if
256
257
258 parameters(1) = alpha
259 parameters(2) = cam%omega
260 parameters(3) = cam%omega
261 call xc_f03_func_set_ext_params(functl(func_c)%conf, parameters)
262 !The name is confusing. Here alpha is the beta of hybrids in functionals,
263 !but is called alpha in the original paper.
264 cam%beta = alpha
265
267 alpha = -1.00778_real64+1.10507_real64*tb09_c
268
269 if (alpha > 1) then
270 write(message(1), '(a,f6.3,a)') 'MVORB mixing parameter bigger than one (' , alpha ,').'
271 call messages_warning(1, namespace=namespace)
272 alpha = 0.25_real64
273 else if (alpha < 0) then
274 write(message(1), '(a,f6.3,a)') 'MVORB mixing parameter smaller than zero (' , alpha ,').'
275 call messages_warning(1, namespace=namespace)
276 alpha = 0.25_real64
277 else
278 write(message(1), '(a,f6.3,a)') 'MVORB mixing parameter is (' , alpha ,').'
279 call messages_info(1, namespace=namespace, debug_only=.true.)
280 end if
281
282 parameters(1) = alpha
283 call xc_f03_func_set_ext_params(functl(func_c)%conf, parameters)
284 cam%alpha = alpha
285 case default
286 call messages_not_implemented("MVORB density-based mixing for functionals other than PBE0 and HSE06", namespace=namespace)
287 end select
288
289 pop_sub(calc_mvorb_alpha)
290 end subroutine calc_mvorb_alpha
291
292
293end module xc_interaction_oct_m
295!! Local Variables:
296!! mode: f90
297!! coding: utf-8
298!! End:
double sqrt(double __x) __attribute__((__nothrow__
integer, parameter, public unpolarized
Parameters...
real(real64), parameter, public m_zero
Definition: global.F90:188
This module defines the abstract interaction_t class, and some auxiliary classes for interactions.
subroutine, public interaction_end(this)
This module defines classes and functions for interaction partners.
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1113
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:537
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
Definition: quantity.F90:138
This module handles spin dimensions of the states and the k-point distribution.
integer, parameter, public xc_hyb_gga_xc_mvorb_pbeh
Density-based mixing parameter of PBE0.
integer, parameter, public xc_hyb_gga_xc_mvorb_hse06
Density-based mixing parameter of HSE06.
subroutine xc_interaction_finalize(this)
class(xc_interaction_t) function, pointer xc_interaction_constructor(partner)
subroutine, public xc_interaction_compute(this)
subroutine xc_interaction_end(this)
subroutine, public calc_tb09_c(mesh, space, functl, dens, gdens, ispin, rcell_volume)
integer, parameter, public func_c
subroutine xc_interaction_init(this)
subroutine xc_interaction_calculate(this)
subroutine xc_interaction_calculate_energy(this)
subroutine, public calc_mvorb_alpha(mesh, namespace, space, functl, dens, gdens, ispin, rcell_volume, cam)
Describes mesh distribution to nodes.
Definition: mesh.F90:186
Coulomb-attenuating method parameters, used in the partitioning of the Coulomb potential into a short...
Definition: xc_cam.F90:139
int true(void)