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 :: &
66
67 type orbitalset_t
68 ! Components are public by default
69 integer :: nn, ll, ii
70 real(real64) :: jj
71 integer :: norbs
72 integer :: ndim
73 integer :: iatom
74 type(submesh_t) :: sphere
75 complex(real64), allocatable :: phase(:,:)
76 ! !< if the sphere cross the border of the box
77 real(real64) :: Ueff
78 real(real64) :: Ubar, Jbar
79 real(real64) :: alpha
80 integer :: nneighbors
81 ! !< interaction is considered
82 real(real64), allocatable :: V_ij(:,:)
83 real(real64), allocatable :: coulomb_IIJJ(:,:,:,:,:)
84 complex(real64), allocatable :: zcoulomb_IIJJ(:,:,:,:,:,:,:)
85 integer, allocatable:: map_os(:)
86 complex(real64), allocatable :: phase_shift(:,:)
87
88 real(real64) :: radius
89 class(species_t), pointer :: spec => null()
90 integer :: spec_index
91
92 real(real64), allocatable :: dorb(:,:,:)
93 complex(real64), allocatable :: zorb(:,:,:)
94 complex(real64), allocatable :: eorb_submesh(:,:,:,:)
95 ! !! (for isolated system with TD phase)
96 complex(real64), allocatable :: eorb_mesh(:,:,:,:)
97 ! !! (for periodic systems GS and TD)
98 integer :: ldorbs, ldorbs_eorb
99 type(accel_mem_t) :: dbuff_orb, zbuff_orb
100 type(accel_mem_t), allocatable :: buff_eorb (:)
101
102 logical :: use_submesh
103 logical :: allocated_on_mesh
104 ! This stores how the localized orbitals (without phase) were allocated.
105 ! This is different from "use_submesh", which dictate how the orbitals (with phase) will be applied.
106
107 type(poisson_t) :: poisson
108 contains
109 procedure :: set_species_index => orbitalset_set_species_index
110 end type orbitalset_t
111
112contains
113
114 subroutine orbitalset_init(this)
115 type(orbitalset_t), intent(inout) :: this
117 push_sub(orbitalset_init)
118
119 this%iatom = -1
120 this%nneighbors = 0
121 this%nn = 0
122 this%ll = 0
123 this%jj = m_one
124 this%ii = 0
125 this%iatom = 0
126 this%ndim = 1
127 this%spec_index = 0
128 this%norbs = 0
129
130 this%Ueff = m_zero
131 this%Ubar = m_zero
132 this%Jbar = m_zero
133 this%alpha = m_zero
134 this%radius = m_zero
135
136 pop_sub(orbitalset_init)
137 end subroutine orbitalset_init
138
139
140 subroutine orbitalset_end(this)
141 type(orbitalset_t), intent(inout) :: this
142
143 integer :: ik
144
145 push_sub(orbitalset_end)
146
147 safe_deallocate_a(this%phase)
148 safe_deallocate_a(this%dorb)
149 safe_deallocate_a(this%zorb)
150 safe_deallocate_a(this%eorb_submesh)
151 safe_deallocate_a(this%eorb_mesh)
152 nullify(this%spec)
153 call submesh_end(this%sphere)
154
155 safe_deallocate_a(this%V_ij)
156 safe_deallocate_a(this%coulomb_IIJJ)
157 safe_deallocate_a(this%map_os)
158 safe_deallocate_a(this%phase_shift)
159
160 if(accel_is_enabled()) then
161 call accel_release_buffer(this%dbuff_orb)
162 call accel_release_buffer(this%zbuff_orb)
163 if(allocated(this%buff_eorb)) then
164 do ik = lbound(this%buff_eorb, dim=1), ubound(this%buff_eorb, dim=1)
165 call accel_release_buffer(this%buff_eorb(ik))
166 end do
167 safe_deallocate_a(this%buff_eorb)
168 end if
169 end if
171 pop_sub(orbitalset_end)
172 end subroutine orbitalset_end
174 subroutine orbitalset_set_jln(this, jj, ll, nn)
175 type(orbitalset_t), intent(inout) :: this
176 real(real64), intent(in) :: jj
177 integer, intent(in) :: ll, nn
181 this%jj = jj
182 this%ll = ll
183 this%nn = nn
186 end subroutine orbitalset_set_jln
190 subroutine orbitalset_update_phase(os, dim, kpt, kpoints, spin_polarized, vec_pot, vec_pot_var, kpt_max)
191 type(orbitalset_t), intent(inout) :: os
192 integer, intent(in) :: dim
193 type(distributed_t), intent(in) :: kpt
194 type(kpoints_t), intent(in) :: kpoints
195 logical, intent(in) :: spin_polarized
196 real(real64), optional, allocatable, intent(in) :: vec_pot(:)
197 real(real64), optional, allocatable, intent(in) :: vec_pot_var(:, :)
198 integer, optional, intent(in) :: kpt_max
199
200 integer :: ns, iq, is, ikpoint, im, idim, kpt_end
201 real(real64) :: kr, kpoint(1:dim), dx(1:dim)
202 integer :: iorb
203
205
206 ns = os%sphere%np
207
208 kpt_end = kpt%end
209 if (present(kpt_max)) kpt_end = min(kpt_max, kpt_end)
210
211 do iq = kpt%start, kpt_end
212 !This is durty but avoids to refer to states_get_kpoint_index
213 if (spin_polarized) then
214 ikpoint = 1 + (iq - 1)/2
215 else
216 ikpoint = iq
217 end if
218
219 ! if this fails, it probably means that sb is not compatible with std
220 assert(ikpoint <= kpoints_number(kpoints))
221
222 kpoint(1:dim) = kpoints%get_point(ikpoint)
223
224 !$omp parallel do private(dx, kr)
225 do is = 1, ns
226 ! this is only the correction to the global phase, that can
227 ! appear if the sphere crossed the boundary of the cell.
228 dx = os%sphere%rel_x(1:dim, is) - os%sphere%mesh%x(1:dim, os%sphere%map(is)) + os%sphere%center(1:dim)
229 kr = dot_product(kpoint, dx)
230 if (present(vec_pot)) then
231 if (allocated(vec_pot)) kr = kr + dot_product(vec_pot, dx)
232 end if
233
234 if (present(vec_pot_var)) then
235 if (allocated(vec_pot_var)) kr = kr + sum(vec_pot_var(1:dim, os%sphere%map(is)) &
236 *(os%sphere%rel_x(1:dim, is)+os%sphere%center))
237 end if
238
239 os%phase(is, iq) = exp(m_zi*kr)
240 end do
241
242 if (.not. os%use_submesh) then
243 ! We now compute the so-called Bloch sum of the localized orbitals
244 do idim = 1, os%ndim
245 do im = 1, os%norbs
246 os%eorb_mesh(:, im, idim, iq) = m_z0
247 do is = 1, ns
248 os%eorb_mesh(os%sphere%map(is),im,idim,iq) = os%eorb_mesh(os%sphere%map(is),im,idim,iq) &
249 + os%zorb(is, idim, im) * os%phase(is, iq)
250 end do
251 end do
252 end do
253 else ! In the case of the isolated system, we still use the submesh
254 !$omp parallel private(idim, is)
255 do im = 1, os%norbs
256 do idim = 1, os%ndim
257 !$omp do simd
258 do is = 1, ns
259 os%eorb_submesh(is,idim,im,iq) = os%zorb(is,idim,im)*os%phase(is, iq)
260 end do
261 !$omp end do simd
262 end do
263 end do
264 !$omp end parallel
265 end if
266
267 if(accel_is_enabled() .and. os%ndim == 1) then
268 if(os%use_submesh) then
269 do iorb = 1, os%norbs
270 call accel_write_buffer(os%buff_eorb(iq), ns, os%eorb_submesh(:, 1, iorb, iq), offset = (iorb - 1)*os%ldorbs_eorb)
271 end do
272 else
273 do iorb = 1, os%norbs
274 call accel_write_buffer(os%buff_eorb(iq), os%sphere%mesh%np, &
275 os%eorb_mesh(:, iorb, 1, iq), offset = (iorb - 1)*os%ldorbs_eorb)
276 end do
277 end if
278 end if
279 end do
280
282 end subroutine orbitalset_update_phase
283
284
286 subroutine orbitalset_update_phase_shift(os, dim, kpt, kpoints, spin_polarized, vec_pot, vec_pot_var, kpt_max)
287 type(orbitalset_t), intent(inout) :: os
288 integer, intent(in) :: dim
289 type(distributed_t), intent(in) :: kpt
290 type(kpoints_t), intent(in) :: kpoints
291 logical, intent(in) :: spin_polarized
292 real(real64), optional, allocatable, intent(in) :: vec_pot(:)
293 real(real64), optional, allocatable, intent(in) :: vec_pot_var(:, :)
294 integer, optional, intent(in) :: kpt_max
295
296 integer :: iq, ikpoint
297 real(real64) :: kr, kpoint(dim), dx(dim)
298 integer :: inn, kpt_end
299
301
302 kpt_end = kpt%end
303 if(present(kpt_max)) kpt_end = min(kpt_max, kpt_end)
304
305 do iq = kpt%start, kpt_end
306 !This is durty but avoids to refer to states_get_kpoint_index
307 if(spin_polarized) then
308 ikpoint = 1 + (iq - 1)/2
309 else
310 ikpoint = iq
311 end if
312
313 ! if this fails, it probably means that sb is not compatible with std
314 assert(ikpoint <= kpoints_number(kpoints))
315
316 kpoint(1:dim) = kpoints%get_point(ikpoint)
317
318 if (os%nneighbors > 0) then
319 do inn = 1, os%nneighbors
320 dx(1:dim) = os%V_ij(inn,1:dim)
321 kr = sum(kpoint(1:dim)*dx(1:dim))
322 if (present(vec_pot)) then
323 if (allocated(vec_pot)) kr = kr + sum(vec_pot(1:dim)*dx(1:dim))
324 end if
325
326 !At the moment the uniform vector potential is in vec_pot_var
327 if (present(vec_pot_var)) then
328 if (allocated(vec_pot_var)) kr = kr + sum(vec_pot_var(1:dim, 1)*dx(1:dim))
329 end if
330
331 !The sign is different as this is applied on the wavefunction and not the orbitals
332 os%phase_shift(inn, iq) = exp(-m_zi*kr)
333 end do
334 end if
335 end do
336
337
339 end subroutine orbitalset_update_phase_shift
340
341 ! ---------------------------------------------------------
347 subroutine orbitalset_set_species_index(os, ions)
348 class(orbitalset_t), intent(inout) :: os
349 type(ions_t), intent(in) :: ions
350
351 integer :: ja
352
354
355 os%spec_index = os%spec%get_index()
356 do ja = 1, ions%natoms
357 if(ions%atom(ja)%species == os%spec) then
358 os%spec_index = min(os%spec_index, ions%atom(ja)%species%get_index())
359 end if
360 end do
361
363 end subroutine orbitalset_set_species_index
364
365#include "undef.F90"
366#include "real.F90"
367#include "orbitalset_inc.F90"
368
369#include "undef.F90"
370#include "complex.F90"
371#include "orbitalset_inc.F90"
372
373end module orbitalset_oct_m
double exp(double __x) __attribute__((__nothrow__
subroutine, public accel_release_buffer(this, async)
Definition: accel.F90:879
pure logical function, public accel_is_enabled()
Definition: accel.F90:392
This module implements batches of mesh functions.
Definition: batch.F90:135
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
This module contains interfaces for BLAS routines You should not use these routines directly....
Definition: blas.F90:120
real(real64), parameter, public m_zero
Definition: global.F90:191
complex(real64), parameter, public m_z0
Definition: global.F90:201
complex(real64), parameter, public m_zi
Definition: global.F90:205
real(real64), parameter, public m_one
Definition: global.F90:192
integer pure function, public kpoints_number(this)
Definition: kpoints.F90:1126
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
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:382
subroutine orbitalset_set_species_index(os, ions)
Set the species index for a given orbital set.
Definition: orbitalset.F90:443
subroutine, public orbitalset_init(this)
Definition: orbitalset.F90:210
subroutine, public orbitalset_end(this)
Definition: orbitalset.F90:236
subroutine, public dorbitalset_transfer_to_device(os, kpt, use_mesh)
Allocate and transfer the orbitals to the device.
subroutine, public dorbitalset_add_to_batch(os, ndim, psib, weight)
Definition: orbitalset.F90:865
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:286
subroutine, public dorbitalset_get_coefficients(os, ndim, psi, ik, has_phase, dot, reduce)
Definition: orbitalset.F90:528
subroutine, public dorbitalset_get_coeff_batch(os, ndim, psib, dot, reduce)
Definition: orbitalset.F90:590
subroutine, public orbitalset_set_jln(this, jj, ll, nn)
Definition: orbitalset.F90:270
subroutine, public zorbitalset_get_coeff_batch(os, ndim, psib, dot, reduce)
subroutine, public zorbitalset_transfer_to_device(os, kpt, use_mesh)
Allocate and transfer the orbitals to the device.
subroutine, public zorbitalset_get_coefficients(os, ndim, psi, ik, has_phase, dot, reduce)
subroutine, public submesh_end(this)
Definition: submesh.F90:733
Distribution of N instances over mpi_grpsize processes, for the local rank mpi_grprank....