Octopus
oct_exchange.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
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 debug_oct_m
25 use global_oct_m
26 use grid_oct_m
27 use mesh_oct_m
34 use xc_oct_m
36
37 implicit none
38
39 private
40 public :: &
47
49 private
50 logical :: oct_exchange = .false.
51 type(states_elec_t), pointer :: oct_st => null()
52 real(real64), allocatable :: oct_fxc(:, :, :)
53 real(real64), allocatable :: oct_pot(:, :)
54 real(real64), allocatable :: oct_rho(:, :)
55 end type oct_exchange_t
56
57contains
58
59 ! ---------------------------------------------------------
60 logical function oct_exchange_enabled(this) result(oct_exchange)
61 type(oct_exchange_t), intent(in) :: this
62
63 push_sub(oct_exchange_enabled)
64
65 oct_exchange = this%oct_exchange
66
68 end function oct_exchange_enabled
69
70
71 ! ---------------------------------------------------------
72 subroutine oct_exchange_set(this, st, mesh)
73 type(oct_exchange_t), intent(inout) :: this
74 type(states_elec_t), target, intent(in) :: st
75 class(mesh_t), intent(in) :: mesh
76
77 integer :: np, nspin
78
79 push_sub(oct_exchange_set)
80
81 ! In this release, no non-local part for the QOCT Hamiltonian.
82 nullify(this%oct_st)
83
84 this%oct_st => st
85 this%oct_exchange = .true.
86 np = mesh%np
87 nspin = this%oct_st%d%nspin
88
89 safe_allocate(this%oct_fxc(1:np, 1:nspin, 1:nspin))
90 safe_allocate(this%oct_pot(1:np, 1:nspin))
91 safe_allocate(this%oct_rho(1:np, 1:nspin))
92
93 this%oct_fxc = m_zero
94 this%oct_pot = m_zero
95 this%oct_rho = m_zero
96
97 pop_sub(oct_exchange_set)
98 end subroutine oct_exchange_set
99
100
101 ! ---------------------------------------------------------
102 subroutine oct_exchange_prepare(this, gr, psi, xc, psolver, namespace)
103 type(oct_exchange_t), intent(inout) :: this
104 type(grid_t), intent(in) :: gr
105 complex(real64), intent(in) :: psi(:, :, :, :)
106 type(xc_t), intent(in) :: xc
107 type(poisson_t), intent(in) :: psolver
108 type(namespace_t), intent(in) :: namespace
109
110 integer :: jst, ip, ik
111 complex(real64), allocatable :: psi2(:, :)
112
113 push_sub(oct_exchange_prepare)
114
115 safe_allocate(psi2(1:gr%np, 1:this%oct_st%d%dim))
117 select case (this%oct_st%d%ispin)
118 case (unpolarized)
119 assert(this%oct_st%nik == 1)
120
121 this%oct_pot = m_zero
122 this%oct_rho = m_zero
123 do jst = 1, this%oct_st%nst
124 call states_elec_get_state(this%oct_st, gr, jst, 1, psi2)
125 do ip = 1, gr%np
126 this%oct_rho(ip, 1) = this%oct_rho(ip, 1) + this%oct_st%occ(jst, 1)*aimag(conjg(psi2(ip, 1))*psi(ip, 1, jst, 1))
127 end do
128 end do
129 call dpoisson_solve(psolver, namespace, this%oct_pot(:, 1), this%oct_rho(:, 1), all_nodes = .false.)
130
131 case (spin_polarized)
132 assert(this%oct_st%nik == 2)
133
134 this%oct_pot = m_zero
135 this%oct_rho = m_zero
136 do ik = 1, 2
137 do jst = 1, this%oct_st%nst
138 call states_elec_get_state(this%oct_st, gr, jst, ik, psi2)
139 do ip = 1, gr%np
140 this%oct_rho(ip, ik) = this%oct_rho(ip, ik) + this%oct_st%occ(jst, ik) * aimag(conjg(psi2(ip, 1))*psi(ip, 1, jst, ik))
141 end do
142 end do
143 end do
144
145 do ik = 1, 2
146 call dpoisson_solve(psolver, namespace, this%oct_pot(:, ik), this%oct_rho(:, ik), all_nodes = .false.)
147 end do
149 end select
150
151 call xc_get_fxc(xc, gr, namespace, this%oct_st%rho, this%oct_st%d%ispin, this%oct_fxc)
152
153 safe_deallocate_a(psi2)
154 pop_sub(oct_exchange_prepare)
155 end subroutine oct_exchange_prepare
156
157
158 ! ---------------------------------------------------------
159 subroutine oct_exchange_remove(this)
160 type(oct_exchange_t), intent(inout) :: this
161
162 push_sub(oct_exchange_remove)
163
164 nullify(this%oct_st)
165 this%oct_exchange = .false.
166 safe_deallocate_a(this%oct_fxc)
167 safe_deallocate_a(this%oct_pot)
168 safe_deallocate_a(this%oct_rho)
169
170 pop_sub(oct_exchange_remove)
171 end subroutine oct_exchange_remove
172
173 ! ---------------------------------------------------------
174 subroutine oct_exchange_operator(this, namespace, mesh, hpsi, ist, ik)
175 type(oct_exchange_t), intent(in) :: this
176 type(namespace_t), intent(in) :: namespace
177 class(mesh_t), intent(in) :: mesh
178 complex(real64), intent(inout) :: hpsi(:, :)
179 integer, intent(in) :: ist
180 integer, intent(in) :: ik
181
182 integer :: ik2
183 complex(real64), allocatable :: psi(:, :), psi2(:, :)
184 integer :: ip
185
186 push_sub(oct_exchange_operator)
187
188 safe_allocate(psi(1:mesh%np, 1:this%oct_st%d%dim))
189 safe_allocate(psi2(1:mesh%np, 1:this%oct_st%d%dim))
190
191 select case (this%oct_st%d%ispin)
192 case (unpolarized)
193 assert(this%oct_st%nik == 1)
194 call states_elec_get_state(this%oct_st, mesh, ist, 1, psi2)
195 do ip = 1, mesh%np
196 hpsi(ip, 1) = hpsi(ip, 1) + m_two*m_zi*psi2(ip, 1)*(this%oct_pot(ip, 1) + this%oct_fxc(ip, 1, 1)*this%oct_rho(ip, 1))
197 end do
198
199 case (spin_polarized)
200 assert(this%oct_st%nik == 2)
201
202 call states_elec_get_state(this%oct_st, mesh, ist, ik, psi2)
203
204 do ik2 = 1, 2
205 do ip = 1, mesh%np
206 hpsi(ip, 1) = hpsi(ip, 1) + m_two * m_zi * this%oct_st%occ(ist, ik) * &
207 psi2(ip, 1) * (this%oct_pot(ip, ik2) + this%oct_fxc(ip, ik, ik2)*this%oct_rho(ip, ik2))
208 end do
209 end do
210
211 case (spinors)
212 call messages_not_implemented("Function oct_exchange_operator for spinors", namespace=namespace)
213 end select
214
215 safe_deallocate_a(psi)
216 safe_deallocate_a(psi2)
217 pop_sub(oct_exchange_operator)
218 end subroutine oct_exchange_operator
219
220end module oct_exchange_oct_m
221
222!! Local Variables:
223!! mode: f90
224!! coding: utf-8
225!! End:
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
integer, parameter, public spin_polarized
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
complex(real64), parameter, public m_zi
Definition: global.F90:205
This module implements the underlying real-space grid.
Definition: grid.F90:119
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1067
subroutine, public oct_exchange_prepare(this, gr, psi, xc, psolver, namespace)
logical function, public oct_exchange_enabled(this)
subroutine, public oct_exchange_remove(this)
subroutine, public oct_exchange_set(this, st, mesh)
subroutine, public oct_exchange_operator(this, namespace, mesh, hpsi, ist, ik)
subroutine, public dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel, reset)
Calculates the Poisson equation. Given the density returns the corresponding potential.
Definition: poisson.F90:872
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public xc_get_fxc(xcs, gr, namespace, rho, ispin, fxc, fxc_grad, fxc_grad_spin)
Returns the exchange-correlation kernel.
Definition: xc_kernel.F90:172
Definition: xc.F90:116
Describes mesh distribution to nodes.
Definition: mesh.F90:187
int true(void)