Octopus
orbitalset.F90
Go to the documentation of this file.
1!! Copyright (C) 2016 N. 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 accel_oct_m
23 use batch_oct_m
25 use blas_oct_m
26 use comm_oct_m
27 use debug_oct_m
29 use global_oct_m
33 use math_oct_m
34 use mesh_oct_m
41 use types_oct_m
43
44 implicit none
45
46 private
47
48 public :: &
65
66 type orbitalset_t
67 ! Components are public by default
68 integer :: nn, ll, ii
69 real(real64) :: jj
70 integer :: norbs
71 integer :: ndim
72 integer :: iatom
73 type(submesh_t) :: sphere
74 complex(real64), allocatable :: phase(:,:)
75 ! !< if the sphere cross the border of the box
76 real(real64) :: Ueff
77 real(real64) :: Ubar, Jbar
78 real(real64) :: alpha
79 integer :: nneighbors
80 ! !< interaction is considered
81 real(real64), allocatable :: V_ij(:,:)
82 real(real64), allocatable :: coulomb_IIJJ(:,:,:,:,:)
83 integer, allocatable:: map_os(:)
84 complex(real64), allocatable :: phase_shift(:,:)
85
86 real(real64) :: radius
87 class(species_t), pointer :: spec => null()
88 integer :: spec_index
89
90 real(real64), allocatable :: dorb(:,:,:)
91 complex(real64), allocatable :: zorb(:,:,:)
92 complex(real64), allocatable :: eorb_submesh(:,:,:,:)
93 complex(real64), allocatable :: eorb_mesh(:,:,:,:)
94 integer :: ldorbs
95 type(accel_mem_t) :: dbuff_orb, zbuff_orb
96 type(accel_mem_t), allocatable :: buff_eorb (:)
97
98 logical :: use_submesh
99
100 type(poisson_t) :: poisson
101 end type orbitalset_t
102
103contains
104
105 subroutine orbitalset_init(this)
106 type(orbitalset_t), intent(inout) :: this
107
108 push_sub(orbitalset_init)
109
110 this%iatom = -1
111 this%nneighbors = 0
112 this%nn = 0
113 this%ll = 0
114 this%jj = m_one
115 this%ii = 0
116 this%iatom = 0
117 this%ndim = 1
118 this%spec_index = 0
119 this%norbs = 0
120
121 this%Ueff = m_zero
122 this%Ubar = m_zero
123 this%Jbar = m_zero
124 this%alpha = m_zero
125 this%radius = m_zero
126
127 pop_sub(orbitalset_init)
128 end subroutine orbitalset_init
129
130
131 subroutine orbitalset_end(this)
132 type(orbitalset_t), intent(inout) :: this
133
134 integer :: ik
135
136 push_sub(orbitalset_end)
137
138 safe_deallocate_a(this%phase)
139 safe_deallocate_a(this%dorb)
140 safe_deallocate_a(this%zorb)
141 safe_deallocate_a(this%eorb_submesh)
142 safe_deallocate_a(this%eorb_mesh)
143 nullify(this%spec)
144 call submesh_end(this%sphere)
145
146 safe_deallocate_a(this%V_ij)
147 safe_deallocate_a(this%coulomb_IIJJ)
148 safe_deallocate_a(this%map_os)
149 safe_deallocate_a(this%phase_shift)
150
151 if(accel_is_enabled()) then
152 call accel_release_buffer(this%dbuff_orb)
153 call accel_release_buffer(this%zbuff_orb)
154 if(allocated(this%buff_eorb)) then
155 do ik = lbound(this%buff_eorb, dim=1), ubound(this%buff_eorb, dim=1)
156 call accel_release_buffer(this%buff_eorb(ik))
157 end do
158 safe_deallocate_a(this%buff_eorb)
159 end if
160 end if
163 end subroutine orbitalset_end
165 subroutine orbitalset_set_jln(this, jj, ll, nn)
166 type(orbitalset_t), intent(inout) :: this
167 real(real64), intent(in) :: jj
168 integer, intent(in) :: ll, nn
172 this%jj = jj
173 this%ll = ll
174 this%nn = nn
177 end subroutine orbitalset_set_jln
178
181 subroutine orbitalset_update_phase(os, dim, kpt, kpoints, spin_polarized, vec_pot, vec_pot_var, kpt_max)
182 type(orbitalset_t), intent(inout) :: os
183 integer, intent(in) :: dim
184 type(distributed_t), intent(in) :: kpt
185 type(kpoints_t), intent(in) :: kpoints
186 logical, intent(in) :: spin_polarized
187 real(real64), optional, allocatable, intent(in) :: vec_pot(:)
188 real(real64), optional, allocatable, intent(in) :: vec_pot_var(:, :)
189 integer, optional, intent(in) :: kpt_max
190
191 integer :: ns, iq, is, ikpoint, im, idim, kpt_end
192 real(real64) :: kr, kpoint(1:dim), dx(1:dim)
193 integer :: iorb
194
196
197 ns = os%sphere%np
199 kpt_end = kpt%end
200 if (present(kpt_max)) kpt_end = min(kpt_max, kpt_end)
201
202 do iq = kpt%start, kpt_end
203 !This is durty but avoids to refer to states_get_kpoint_index
204 if (spin_polarized) then
205 ikpoint = 1 + (iq - 1)/2
206 else
207 ikpoint = iq
208 end if
209
210 ! if this fails, it probably means that sb is not compatible with std
211 assert(ikpoint <= kpoints_number(kpoints))
212
213 kpoint(1:dim) = kpoints%get_point(ikpoint)
214
215 !$omp parallel do private(dx, kr)
216 do is = 1, ns
217 ! this is only the correction to the global phase, that can
218 ! appear if the sphere crossed the boundary of the cell.
219 dx = os%sphere%rel_x(1:dim, is) - os%sphere%mesh%x(os%sphere%map(is), 1:dim) + os%sphere%center(1:dim)
220 kr = dot_product(kpoint, dx)
221 if (present(vec_pot)) then
222 if (allocated(vec_pot)) kr = kr + dot_product(vec_pot, dx)
223 end if
225 if (present(vec_pot_var)) then
226 if (allocated(vec_pot_var)) kr = kr + sum(vec_pot_var(1:dim, os%sphere%map(is)) &
227 *(os%sphere%rel_x(1:dim, is)+os%sphere%center))
228 end if
229
230 os%phase(is, iq) = exp(m_zi*kr)
231 end do
232
233 if (.not. os%use_submesh) then
234 ! We now compute the so-called Bloch sum of the localized orbitals
235 do idim = 1, os%ndim
236 do im = 1, os%norbs
237 os%eorb_mesh(:, im, idim, iq) = m_z0
238 do is = 1, ns
239 os%eorb_mesh(os%sphere%map(is),im,idim,iq) = os%eorb_mesh(os%sphere%map(is),im,idim,iq) &
240 + os%zorb(is, idim, im) * os%phase(is, iq)
241 end do
242 end do
243 end do
244 else ! In the case of the isolated system, we still use the submesh
245 !$omp parallel private(idim, is)
246 do im = 1, os%norbs
247 do idim = 1, os%ndim
248 !$omp do simd
249 do is = 1, ns
250 os%eorb_submesh(is,idim,im,iq) = os%zorb(is,idim,im)*os%phase(is, iq)
251 end do
252 !$omp end do simd
253 end do
254 end do
255 !$omp end parallel
256 end if
257
258 if(accel_is_enabled() .and. os%ndim == 1) then
259 if(os%use_submesh) then
260 do iorb = 1, os%norbs
261 call accel_write_buffer(os%buff_eorb(iq), ns, os%eorb_submesh(:, 1, iorb, iq), offset = (iorb - 1)*os%ldorbs)
262 end do
263 else
264 do iorb = 1, os%norbs
265 call accel_write_buffer(os%buff_eorb(iq), os%sphere%mesh%np, &
266 os%eorb_mesh(:, iorb, 1, iq), offset = (iorb - 1)*os%ldorbs)
267 end do
268 end if
269 end if
270 end do
271
273 end subroutine orbitalset_update_phase
275
277 subroutine orbitalset_update_phase_shift(os, dim, kpt, kpoints, spin_polarized, vec_pot, vec_pot_var, kpt_max)
278 type(orbitalset_t), intent(inout) :: os
279 integer, intent(in) :: dim
280 type(distributed_t), intent(in) :: kpt
281 type(kpoints_t), intent(in) :: kpoints
282 logical, intent(in) :: spin_polarized
283 real(real64), optional, allocatable, intent(in) :: vec_pot(:)
284 real(real64), optional, allocatable, intent(in) :: vec_pot_var(:, :)
285 integer, optional, intent(in) :: kpt_max
286
287 integer :: iq, ikpoint
288 real(real64) :: kr, kpoint(dim), dx(dim)
289 integer :: inn, kpt_end
290
292
293 kpt_end = kpt%end
294 if(present(kpt_max)) kpt_end = min(kpt_max, kpt_end)
295
296 do iq = kpt%start, kpt_end
297 !This is durty but avoids to refer to states_get_kpoint_index
298 if(spin_polarized) then
299 ikpoint = 1 + (iq - 1)/2
300 else
301 ikpoint = iq
302 end if
303
304 ! if this fails, it probably means that sb is not compatible with std
305 assert(ikpoint <= kpoints_number(kpoints))
306
307 kpoint(1:dim) = kpoints%get_point(ikpoint)
308
309 if (os%nneighbors > 0) then
310 do inn = 1, os%nneighbors
311 dx(1:dim) = os%V_ij(inn,1:dim)
312 kr = sum(kpoint(1:dim)*dx(1:dim))
313 if (present(vec_pot)) then
314 if (allocated(vec_pot)) kr = kr + sum(vec_pot(1:dim)*dx(1:dim))
315 end if
316
317 !At the moment the uniform vector potential is in vec_pot_var
318 if (present(vec_pot_var)) then
319 if (allocated(vec_pot_var)) kr = kr + sum(vec_pot_var(1:dim, 1)*dx(1:dim))
320 end if
321
322 !The sign is different as this is applied on the wavefunction and not the orbitals
323 os%phase_shift(inn, iq) = exp(-m_zi*kr)
324 end do
325 end if
326 end do
327
328
330 end subroutine orbitalset_update_phase_shift
331
332#include "undef.F90"
333#include "real.F90"
334#include "orbitalset_inc.F90"
335
336#include "undef.F90"
337#include "complex.F90"
338#include "orbitalset_inc.F90"
339
340end module orbitalset_oct_m
double exp(double __x) __attribute__((__nothrow__
subroutine, public accel_release_buffer(this)
Definition: accel.F90:1246
pure logical function, public accel_is_enabled()
Definition: accel.F90:400
This module implements batches of mesh functions.
Definition: batch.F90:133
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
This module contains interfaces for BLAS routines You should not use these routines directly....
Definition: blas.F90:118
real(real64), parameter, public m_zero
Definition: global.F90:187
complex(real64), parameter, public m_z0
Definition: global.F90:197
complex(real64), parameter, public m_zi
Definition: global.F90:201
real(real64), parameter, public m_one
Definition: global.F90:188
integer pure function, public kpoints_number(this)
Definition: kpoints.F90:1099
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
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 zorbitalset_get_position_matrix_elem(os, ndim, psib, idir, dot)
subroutine, public orbitalset_update_phase_shift(os, dim, kpt, kpoints, spin_polarized, vec_pot, vec_pot_var, kpt_max)
Build the phase shift for the intersite interaction.
Definition: orbitalset.F90:371
subroutine, public dorbitalset_get_coeff_batch(os, ndim, psib, dot)
Definition: orbitalset.F90:567
subroutine, public orbitalset_init(this)
Definition: orbitalset.F90:199
subroutine, public orbitalset_end(this)
Definition: orbitalset.F90:225
subroutine, public dorbitalset_add_to_batch(os, ndim, psib, weight)
Definition: orbitalset.F90:711
subroutine, public dorbitalset_get_position_matrix_elem(os, ndim, psib, idir, dot)
Definition: orbitalset.F90:965
subroutine, public zorbitalset_add_to_batch(os, ndim, psib, weight)
subroutine, public orbitalset_update_phase(os, dim, kpt, kpoints, spin_polarized, vec_pot, vec_pot_var, kpt_max)
Build the phase correction to the global phase in case the orbital crosses the border of the simulato...
Definition: orbitalset.F90:275
subroutine, public dorbitalset_get_coefficients(os, ndim, psi, ik, has_phase, dot, reduce)
Definition: orbitalset.F90:493
subroutine, public dorbitalset_get_coeff_batch_accel(os, ndim, psib, dot)
Definition: orbitalset.F90:604
subroutine, public orbitalset_set_jln(this, jj, ll, nn)
Definition: orbitalset.F90:259
subroutine, public zorbitalset_get_coeff_batch_accel(os, ndim, psib, dot)
subroutine, public zorbitalset_get_coeff_batch(os, ndim, psib, dot)
subroutine, public zorbitalset_get_coefficients(os, ndim, psi, ik, has_phase, dot, reduce)
subroutine, public submesh_end(this)
Definition: submesh.F90:737
Distribution of N instances over mpi_grpsize processes, for the local rank mpi_grprank....