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
31 use ions_oct_m
34 use math_oct_m
35 use mesh_oct_m
42 use types_oct_m
44
45 implicit none
46
47 private
48
49 public :: &
64
65 type orbitalset_t
66 ! Components are public by default
67 integer :: nn, ll, ii
68 real(real64) :: jj
69 integer :: norbs
70 integer :: ndim
71 integer :: iatom
72 type(submesh_t) :: sphere
73 complex(real64), allocatable :: phase(:,:)
74 ! !< if the sphere cross the border of the box
75 real(real64) :: Ueff
76 real(real64) :: Ubar, Jbar
77 real(real64) :: alpha
78 integer :: nneighbors
79 ! !< interaction is considered
80 real(real64), allocatable :: V_ij(:,:)
81 real(real64), allocatable :: coulomb_IIJJ(:,:,:,:,:)
82 complex(real64), allocatable :: zcoulomb_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, ldorbs_eorb
95 type(accel_mem_t) :: dbuff_orb, zbuff_orb
96 type(accel_mem_t), allocatable :: buff_eorb (:)
97
98 logical :: use_submesh
99 logical :: allocated_on_mesh
100 ! This stores how the localized orbitals (without phase) were allocated.
101 ! This is different from "use_submesh", which dictate how the orbitals (with phase) will be applied.
102
103 type(poisson_t) :: poisson
104 contains
105 procedure :: set_species_index => orbitalset_set_species_index
106 end type orbitalset_t
107
108contains
109
110 subroutine orbitalset_init(this)
111 type(orbitalset_t), intent(inout) :: this
112
113 push_sub(orbitalset_init)
115 this%iatom = -1
116 this%nneighbors = 0
117 this%nn = 0
118 this%ll = 0
119 this%jj = m_one
120 this%ii = 0
121 this%iatom = 0
122 this%ndim = 1
123 this%spec_index = 0
124 this%norbs = 0
125
126 this%Ueff = m_zero
127 this%Ubar = m_zero
128 this%Jbar = m_zero
129 this%alpha = m_zero
130 this%radius = m_zero
131
132 pop_sub(orbitalset_init)
133 end subroutine orbitalset_init
134
135
136 subroutine orbitalset_end(this)
137 type(orbitalset_t), intent(inout) :: this
138
139 integer :: ik
140
141 push_sub(orbitalset_end)
142
143 safe_deallocate_a(this%phase)
144 safe_deallocate_a(this%dorb)
145 safe_deallocate_a(this%zorb)
146 safe_deallocate_a(this%eorb_submesh)
147 safe_deallocate_a(this%eorb_mesh)
148 nullify(this%spec)
149 call submesh_end(this%sphere)
150
151 safe_deallocate_a(this%V_ij)
152 safe_deallocate_a(this%coulomb_IIJJ)
153 safe_deallocate_a(this%map_os)
154 safe_deallocate_a(this%phase_shift)
155
156 if(accel_is_enabled()) then
157 call accel_release_buffer(this%dbuff_orb)
158 call accel_release_buffer(this%zbuff_orb)
159 if(allocated(this%buff_eorb)) then
160 do ik = lbound(this%buff_eorb, dim=1), ubound(this%buff_eorb, dim=1)
161 call accel_release_buffer(this%buff_eorb(ik))
162 end do
163 safe_deallocate_a(this%buff_eorb)
164 end if
165 end if
167 pop_sub(orbitalset_end)
168 end subroutine orbitalset_end
170 subroutine orbitalset_set_jln(this, jj, ll, nn)
171 type(orbitalset_t), intent(inout) :: this
172 real(real64), intent(in) :: jj
173 integer, intent(in) :: ll, nn
177 this%jj = jj
178 this%ll = ll
179 this%nn = nn
182 end subroutine orbitalset_set_jln
186 subroutine orbitalset_update_phase(os, dim, kpt, kpoints, spin_polarized, vec_pot, vec_pot_var, kpt_max)
187 type(orbitalset_t), intent(inout) :: os
188 integer, intent(in) :: dim
189 type(distributed_t), intent(in) :: kpt
190 type(kpoints_t), intent(in) :: kpoints
191 logical, intent(in) :: spin_polarized
192 real(real64), optional, allocatable, intent(in) :: vec_pot(:)
193 real(real64), optional, allocatable, intent(in) :: vec_pot_var(:, :)
194 integer, optional, intent(in) :: kpt_max
195
196 integer :: ns, iq, is, ikpoint, im, idim, kpt_end
197 real(real64) :: kr, kpoint(1:dim), dx(1:dim)
198 integer :: iorb
199
201
202 ns = os%sphere%np
204 kpt_end = kpt%end
205 if (present(kpt_max)) kpt_end = min(kpt_max, kpt_end)
206
207 do iq = kpt%start, kpt_end
208 !This is durty but avoids to refer to states_get_kpoint_index
209 if (spin_polarized) then
210 ikpoint = 1 + (iq - 1)/2
211 else
212 ikpoint = iq
213 end if
214
215 ! if this fails, it probably means that sb is not compatible with std
216 assert(ikpoint <= kpoints_number(kpoints))
217
218 kpoint(1:dim) = kpoints%get_point(ikpoint)
219
220 !$omp parallel do private(dx, kr)
221 do is = 1, ns
222 ! this is only the correction to the global phase, that can
223 ! appear if the sphere crossed the boundary of the cell.
224 dx = os%sphere%rel_x(1:dim, is) - os%sphere%mesh%x(os%sphere%map(is), 1:dim) + os%sphere%center(1:dim)
225 kr = dot_product(kpoint, dx)
226 if (present(vec_pot)) then
227 if (allocated(vec_pot)) kr = kr + dot_product(vec_pot, dx)
228 end if
230 if (present(vec_pot_var)) then
231 if (allocated(vec_pot_var)) kr = kr + sum(vec_pot_var(1:dim, os%sphere%map(is)) &
232 *(os%sphere%rel_x(1:dim, is)+os%sphere%center))
233 end if
234
235 os%phase(is, iq) = exp(m_zi*kr)
236 end do
237
238 if (.not. os%use_submesh) then
239 ! We now compute the so-called Bloch sum of the localized orbitals
240 do idim = 1, os%ndim
241 do im = 1, os%norbs
242 os%eorb_mesh(:, im, idim, iq) = m_z0
243 do is = 1, ns
244 os%eorb_mesh(os%sphere%map(is),im,idim,iq) = os%eorb_mesh(os%sphere%map(is),im,idim,iq) &
245 + os%zorb(is, idim, im) * os%phase(is, iq)
246 end do
247 end do
248 end do
249 else ! In the case of the isolated system, we still use the submesh
250 !$omp parallel private(idim, is)
251 do im = 1, os%norbs
252 do idim = 1, os%ndim
253 !$omp do simd
254 do is = 1, ns
255 os%eorb_submesh(is,idim,im,iq) = os%zorb(is,idim,im)*os%phase(is, iq)
256 end do
257 !$omp end do simd
258 end do
259 end do
260 !$omp end parallel
261 end if
262
263 if(accel_is_enabled() .and. os%ndim == 1) then
264 if(os%use_submesh) then
265 do iorb = 1, os%norbs
266 call accel_write_buffer(os%buff_eorb(iq), ns, os%eorb_submesh(:, 1, iorb, iq), offset = (iorb - 1)*os%ldorbs_eorb)
267 end do
268 else
269 do iorb = 1, os%norbs
270 call accel_write_buffer(os%buff_eorb(iq), os%sphere%mesh%np, &
271 os%eorb_mesh(:, iorb, 1, iq), offset = (iorb - 1)*os%ldorbs_eorb)
272 end do
273 end if
274 end if
275 end do
276
278 end subroutine orbitalset_update_phase
280
282 subroutine orbitalset_update_phase_shift(os, dim, kpt, kpoints, spin_polarized, vec_pot, vec_pot_var, kpt_max)
283 type(orbitalset_t), intent(inout) :: os
284 integer, intent(in) :: dim
285 type(distributed_t), intent(in) :: kpt
286 type(kpoints_t), intent(in) :: kpoints
287 logical, intent(in) :: spin_polarized
288 real(real64), optional, allocatable, intent(in) :: vec_pot(:)
289 real(real64), optional, allocatable, intent(in) :: vec_pot_var(:, :)
290 integer, optional, intent(in) :: kpt_max
291
292 integer :: iq, ikpoint
293 real(real64) :: kr, kpoint(dim), dx(dim)
294 integer :: inn, kpt_end
295
297
298 kpt_end = kpt%end
299 if(present(kpt_max)) kpt_end = min(kpt_max, kpt_end)
300
301 do iq = kpt%start, kpt_end
302 !This is durty but avoids to refer to states_get_kpoint_index
303 if(spin_polarized) then
304 ikpoint = 1 + (iq - 1)/2
305 else
306 ikpoint = iq
307 end if
308
309 ! if this fails, it probably means that sb is not compatible with std
310 assert(ikpoint <= kpoints_number(kpoints))
311
312 kpoint(1:dim) = kpoints%get_point(ikpoint)
313
314 if (os%nneighbors > 0) then
315 do inn = 1, os%nneighbors
316 dx(1:dim) = os%V_ij(inn,1:dim)
317 kr = sum(kpoint(1:dim)*dx(1:dim))
318 if (present(vec_pot)) then
319 if (allocated(vec_pot)) kr = kr + sum(vec_pot(1:dim)*dx(1:dim))
320 end if
321
322 !At the moment the uniform vector potential is in vec_pot_var
323 if (present(vec_pot_var)) then
324 if (allocated(vec_pot_var)) kr = kr + sum(vec_pot_var(1:dim, 1)*dx(1:dim))
325 end if
326
327 !The sign is different as this is applied on the wavefunction and not the orbitals
328 os%phase_shift(inn, iq) = exp(-m_zi*kr)
329 end do
330 end if
331 end do
332
333
335 end subroutine orbitalset_update_phase_shift
336
337 ! ---------------------------------------------------------
343 subroutine orbitalset_set_species_index(os, ions)
344 class(orbitalset_t), intent(inout) :: os
345 type(ions_t), intent(in) :: ions
346
347 integer :: ja
348
350
351 os%spec_index = os%spec%get_index()
352 do ja = 1, ions%natoms
353 if(ions%atom(ja)%species == os%spec) then
354 os%spec_index = min(os%spec_index, ions%atom(ja)%species%get_index())
355 end if
356 end do
357
359 end subroutine orbitalset_set_species_index
360
361#include "undef.F90"
362#include "real.F90"
363#include "orbitalset_inc.F90"
364
365#include "undef.F90"
366#include "complex.F90"
367#include "orbitalset_inc.F90"
368
369end module orbitalset_oct_m
double exp(double __x) __attribute__((__nothrow__
subroutine, public accel_release_buffer(this)
Definition: accel.F90:1248
pure logical function, public accel_is_enabled()
Definition: accel.F90:401
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:188
complex(real64), parameter, public m_z0
Definition: global.F90:198
complex(real64), parameter, public m_zi
Definition: global.F90:202
real(real64), parameter, public m_one
Definition: global.F90:189
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:376
subroutine orbitalset_set_species_index(os, ions)
Set the species index for a given orbital set.
Definition: orbitalset.F90:437
subroutine, public orbitalset_init(this)
Definition: orbitalset.F90:204
subroutine, public orbitalset_end(this)
Definition: orbitalset.F90:230
subroutine, public dorbitalset_add_to_batch(os, ndim, psib, weight)
Definition: orbitalset.F90:831
subroutine, public dorbitalset_get_position_matrix_elem(os, ndim, psib, idir, dot)
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:280
subroutine, public dorbitalset_get_coefficients(os, ndim, psi, ik, has_phase, dot, reduce)
Definition: orbitalset.F90:522
subroutine, public dorbitalset_get_coeff_batch(os, ndim, psib, dot, reduce)
Definition: orbitalset.F90:584
subroutine, public orbitalset_set_jln(this, jj, ll, nn)
Definition: orbitalset.F90:264
subroutine, public zorbitalset_get_coeff_batch(os, ndim, psib, dot, reduce)
subroutine, public zorbitalset_get_coefficients(os, ndim, psi, ik, has_phase, dot, reduce)
subroutine, public submesh_end(this)
Definition: submesh.F90:736
Distribution of N instances over mpi_grpsize processes, for the local rank mpi_grprank....