Octopus
lda_u.F90
Go to the documentation of this file.
1!! Copyright (C) 2016-2020 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
21module lda_u_oct_m
22 use accel_oct_m
25 use batch_oct_m
27 use comm_oct_m
28 use debug_oct_m
32 use energy_oct_m
33 use global_oct_m
34 use grid_oct_m
36 use ions_oct_m
40 use loct_oct_m
43 use math_oct_m
44 use mesh_oct_m
47 use mpi_oct_m
53 use parser_oct_m
55 use phase_oct_m
58 use space_oct_m
65 use types_oct_m
68
69 implicit none
70
71 private
72
73 public :: &
74 lda_u_t, &
75 lda_u_init, &
82 lda_u_end, &
100 dlda_u_force, &
101 zlda_u_force, &
102 dlda_u_rvu, &
103 zlda_u_rvu, &
108
109
110 integer, public, parameter :: &
111 DFT_U_NONE = 0, &
112 dft_u_empirical = 1, &
113 dft_u_acbn0 = 2
114
115 integer, public, parameter :: &
116 DFT_U_FLL = 0, &
117 dft_u_amf = 1, &
118 dft_u_mix = 2
119
120
123 type lda_u_t
124 private
125 integer, public :: level = dft_u_none
126
127 ! DFT+U basic variables
128 real(real64), allocatable, public :: dn(:,:,:,:)
129 real(real64), allocatable :: dV(:,:,:,:)
130
131 ! ACBN0 variables
132 complex(real64), allocatable, public :: zn(:,:,:,:)
133 complex(real64), allocatable :: zV(:,:,:,:)
134 real(real64), allocatable, public :: dn_alt(:,:,:,:)
135 complex(real64), allocatable, public :: zn_alt(:,:,:,:)
136
137 real(real64), allocatable :: renorm_occ(:,:,:,:,:)
138
139 ! Coulomb integrales
140 real(real64), allocatable, public :: coulomb(:,:,:,:,:)
141 ! !<(for the ACBN0 functional)
142 complex(real64), allocatable, public :: zcoulomb(:,:,:,:,:,:,:)
143 ! !< (for the ACBN0 functional with spinors)
144
145 type(orbitalbasis_t), public :: basis
146 type(orbitalset_t), pointer, public :: orbsets(:) => null()
147 integer, public :: norbsets = 0
148
149 integer, public :: nspins = 0
150 integer, public :: spin_channels = 0
151 integer :: nspecies = 0
152 integer, public :: maxnorbs = 0
153 integer :: max_np = 0
154
155 logical :: useAllOrbitals = .false.
156 logical :: skipSOrbitals = .true.
157 logical :: freeze_occ = .false.
158 logical :: freeze_u = .false.
159 logical, public :: intersite = .false.
160 real(real64) :: intersite_radius = m_zero
161 logical, public :: basisfromstates = .false.
162 real(real64) :: acbn0_screening = m_one
163 integer, allocatable :: basisstates(:)
164 integer, allocatable :: basisstates_os(:)
165 ! !! the state specified in basisstate(:)
166 logical :: rot_inv = .false.
167 integer :: double_couting = dft_u_fll
168 integer :: sm_poisson = sm_poisson_direct
169 real(real64), allocatable :: dc_alpha(:)
170
171 type(lattice_vectors_t), pointer :: latt
172
173 type(distributed_t) :: orbs_dist
174
175 ! Intersite interaction variables
176 integer, public :: maxneighbors = 0
177 real(real64), allocatable, public :: dn_ij(:,:,:,:,:), dn_alt_ij(:,:,:,:,:), dn_alt_ii(:,:,:,:,:)
178 complex(real64), allocatable, public :: zn_ij(:,:,:,:,:), zn_alt_ij(:,:,:,:,:), zn_alt_ii(:,:,:,:,:)
179
180 ! Symmetrization-related variables
181 logical :: symmetrize_occ_matrices
182 integer, allocatable :: inv_map_symm(:,:)
183 integer :: nsym
184 real(real64), allocatable :: symm_weight(:,:,:,:)
185 end type lda_u_t
186
187contains
188
189 ! ---------------------------------------------------------
190 subroutine lda_u_init(this, namespace, space, level, gr, ions, st, mc, kpoints)
191 type(lda_u_t), target, intent(inout) :: this
192 type(namespace_t), intent(in) :: namespace
193 class(space_t), intent(in) :: space
194 integer, intent(in) :: level
195 type(grid_t), intent(in) :: gr
196 type(ions_t), target, intent(in) :: ions
197 type(states_elec_t), intent(in) :: st
198 type(multicomm_t), intent(in) :: mc
199 type(kpoints_t), intent(in) :: kpoints
200
201 integer :: is, ierr
202 type(block_t) :: blk
203
204 push_sub(lda_u_init)
206 assert(.not. (level == dft_u_none))
207
208 call messages_print_with_emphasis(msg="DFT+U", namespace=namespace)
209 if (gr%parallel_in_domains) call messages_experimental("dft+u parallel in domains", namespace=namespace)
210 this%level = level
211
212 this%latt => ions%latt
213
214 !%Variable DFTUBasisFromStates
215 !%Type logical
216 !%Default no
217 !%Section Hamiltonian::DFT+U
218 !%Description
219 !% If set to yes, Octopus will construct the localized basis from
220 !% user-defined states. The states are taken at the Gamma point (or the first k-point of the
221 !% states in the restart_proj folder.
222 !% The states are defined via the block DFTUBasisStates
223 !%End
224 call parse_variable(namespace, 'DFTUBasisFromStates', .false., this%basisfromstates)
225 if (this%basisfromstates) call messages_experimental("DFTUBasisFromStates", namespace=namespace)
226
227 !%Variable DFTUDoubleCounting
228 !%Type integer
229 !%Default dft_u_fll
230 !%Section Hamiltonian::DFT+U
231 !%Description
232 !% This variable selects which DFT+U
233 !% double counting term is used.
234 !%Option dft_u_fll 0
235 !% (Default) The Fully Localized Limit (FLL)
236 !%Option dft_u_amf 1
237 !% (Experimental) Around mean field double counting, as defined in PRB 44, 943 (1991) and PRB 49, 14211 (1994).
238 !%Option dft_u_mix 2
239 !% (Experimental) Mixed double countind term as introduced by Petukhov et al., PRB 67, 153106 (2003).
240 !% This recovers the FLL and AMF as limiting cases.
241 !%End
242 call parse_variable(namespace, 'DFTUDoubleCounting', dft_u_fll, this%double_couting)
243 call messages_print_var_option('DFTUDoubleCounting', this%double_couting, namespace=namespace)
244 if (this%double_couting /= dft_u_fll) call messages_experimental("DFTUDoubleCounting /= dft_u_ffl", namespace=namespace)
245 if (st%d%ispin == spinors .and. this%double_couting /= dft_u_fll) then
246 call messages_not_implemented("AMF and MIX double counting with spinors", namespace=namespace)
247 end if
249 !%Variable DFTUPoissonSolver
250 !%Type integer
251 !%Section Hamiltonian::DFT+U
252 !%Description
253 !% This variable selects which Poisson solver
254 !% is used to compute the Coulomb integrals over a submesh.
255 !% These are non-periodic Poisson solvers.
256 !% The FFT Poisson solver with spherical cutoff is used by default.
257 !%
258 !%Option dft_u_poisson_direct 0
259 !% Direct Poisson solver. Slow but working in all cases.
260 !%Option dft_u_poisson_isf 1
261 !% (Experimental) ISF Poisson solver on a submesh.
262 !% This does not work for non-orthogonal cells nor domain parallelization.
263 !%Option dft_u_poisson_psolver 2
264 !% (Experimental) PSolver Poisson solver on a submesh.
265 !% This does not work for non-orthogonal cells nor domain parallelization.
266 !% Requires the PSolver external library.
267 !%Option dft_u_poisson_fft 3
268 !% (Default) FFT Poisson solver on a submesh.
269 !% This uses the 0D periodic version of the FFT kernels.
270 !%End
271 call parse_variable(namespace, 'DFTUPoissonSolver', sm_poisson_fft, this%sm_poisson)
272 call messages_print_var_option('DFTUPoissonSolver', this%sm_poisson, namespace=namespace)
273 if (this%sm_poisson /= sm_poisson_direct .and. this%sm_poisson /= sm_poisson_fft) then
274 call messages_experimental("DFTUPoissonSolver different from dft_u_poisson_direct", namespace=namespace)
275 call messages_experimental("and dft_u_poisson_fft", namespace=namespace)
276 end if
277 if (this%sm_poisson == sm_poisson_isf) then
278 if (gr%parallel_in_domains) then
279 call messages_not_implemented("DFTUPoissonSolver=dft_u_poisson_isf with domain parallelization", namespace=namespace)
280 end if
281 if (ions%latt%nonorthogonal) then
282 call messages_not_implemented("DFTUPoissonSolver=dft_u_poisson_isf with non-orthogonal cells", namespace=namespace)
283 end if
284 end if
285 if (this%sm_poisson == sm_poisson_psolver) then
286#if !(defined HAVE_PSOLVER)
287 message(1) = "The PSolver Poisson solver cannot be used since the code was not compiled with the PSolver library."
288 call messages_fatal(1, namespace=namespace)
289#endif
290 if (gr%parallel_in_domains) then
291 call messages_not_implemented("DFTUPoissonSolver=dft_u_poisson_psolver with domain parallelization", namespace=namespace)
292 end if
293 if (ions%latt%nonorthogonal) then
294 call messages_not_implemented("DFTUPoissonSolver=dft_u_poisson_psolver with non-orthogonal cells", namespace=namespace)
295 end if
296 end if
297
298 if (this%level == dft_u_acbn0) then
299 !%Variable UseAllAtomicOrbitals
300 !%Type logical
301 !%Default no
302 !%Section Hamiltonian::DFT+U
303 !%Description
304 !% If set to yes, Octopus will determine the effective U for all atomic orbitals
305 !% from the peusopotential. Only available with ACBN0 functional.
306 !% It is strongly recommended to set AOLoewdin=yes when using the option.
307 !%End
308 call parse_variable(namespace, 'UseAllAtomicOrbitals', .false., this%useAllOrbitals)
309 if (this%useAllOrbitals) call messages_experimental("UseAllAtomicOrbitals", namespace=namespace)
310
311 !%Variable SkipSOrbitals
312 !%Type logical
313 !%Default no
314 !%Section Hamiltonian::DFT+U
315 !%Description
316 !% If set to yes, Octopus will determine the effective U for all atomic orbitals
317 !% from the peusopotential but s orbitals. Only available with ACBN0 functional.
318 !%End
319 call parse_variable(namespace, 'SkipSOrbitals', .true., this%skipSOrbitals)
320 if (.not. this%SkipSOrbitals) call messages_experimental("SkipSOrbitals", namespace=namespace)
321
322 !%Variable ACBN0Screening
323 !%Type float
324 !%Default 1.0
325 !%Section Hamiltonian::DFT+U
326 !%Description
327 !% If set to 0, no screening will be included in the ACBN0 functional, and the U
328 !% will be estimated from bare Hartree-Fock. If set to 1 (default), the full screening
329 !% of the U, as defined in the ACBN0 functional, is used.
330 !%End
331 call parse_variable(namespace, 'ACBN0Screening', m_one, this%acbn0_screening)
332 call messages_print_var_value('ACBN0Screening', this%acbn0_screening, namespace=namespace)
333
334 !%Variable ACBN0RotationallyInvariant
335 !%Type logical
336 !%Section Hamiltonian::DFT+U
337 !%Description
338 !% If set to yes, Octopus will use for U and J a formula which is rotationally invariant.
339 !% This is different from the original formula for U and J.
340 !% This is activated by default, except in the case of spinors, as this is not yet implemented in this case.
341 !%End
342 call parse_variable(namespace, 'ACBN0RotationallyInvariant', st%d%ispin /= spinors, this%rot_inv)
343 call messages_print_var_value('ACBN0RotationallyInvariant', this%rot_inv, namespace=namespace)
344 if (this%rot_inv .and. st%d%ispin == spinors) then
345 call messages_not_implemented("Rotationally invariant ACBN0 with spinors", namespace=namespace)
346 end if
347
348 !%Variable ACBN0IntersiteInteraction
349 !%Type logical
350 !%Default no
351 !%Section Hamiltonian::DFT+U
352 !%Description
353 !% If set to yes, Octopus will determine the effective intersite interaction V
354 !% Only available with ACBN0 functional.
355 !% It is strongly recommended to set AOLoewdin=yes when using the option.
356 !%End
357 call parse_variable(namespace, 'ACBN0IntersiteInteraction', .false., this%intersite)
358 call messages_print_var_value('ACBN0IntersiteInteraction', this%intersite, namespace=namespace)
359 if (this%intersite) call messages_experimental("ACBN0IntersiteInteraction", namespace=namespace)
360
361 if (this%intersite) then
362
363 !This is a non local operator. To make this working, one probably needs to apply the
364 ! symmetries to the generalized occupation matrices
365 if (kpoints%use_symmetries) then
366 call messages_not_implemented("Intersite interaction with kpoint symmetries", namespace=namespace)
367 end if
368
369 !%Variable ACBN0IntersiteCutoff
370 !%Type float
371 !%Section Hamiltonian::DFT+U
372 !%Description
373 !% The cutoff radius defining the maximal intersite distance considered.
374 !% Only available with ACBN0 functional with intersite interaction.
375 !%End
376 call parse_variable(namespace, 'ACBN0IntersiteCutoff', m_zero, this%intersite_radius, unit = units_inp%length)
377 if (abs(this%intersite_radius) < m_epsilon) then
378 call messages_write("ACBN0IntersiteCutoff must be greater than 0")
379 call messages_fatal(1, namespace=namespace)
380 end if
381
382 end if
383
384 end if
385
386 call lda_u_write_info(this, namespace=namespace)
387
388 if (.not. this%basisfromstates) then
389
390 call orbitalbasis_init(this%basis, namespace, space%periodic_dim)
391
392 if (states_are_real(st)) then
393 call dorbitalbasis_build(this%basis, namespace, ions, gr, st%d%kpt, st%d%dim, &
394 this%skipSOrbitals, this%useAllOrbitals)
395 else
396 call zorbitalbasis_build(this%basis, namespace, ions, gr, st%d%kpt, st%d%dim, &
397 this%skipSOrbitals, this%useAllOrbitals)
398 end if
399 this%orbsets => this%basis%orbsets
400 this%norbsets = this%basis%norbsets
401 this%maxnorbs = this%basis%maxnorbs
402 this%max_np = this%basis%max_np
403 this%nspins = st%d%nspin
404 this%spin_channels = st%d%spin_channels
405 this%nspecies = ions%nspecies
406
407 !We allocate the necessary ressources
408 if (states_are_real(st)) then
409 call dlda_u_allocate(this, st)
410 else
411 call zlda_u_allocate(this, st)
412 end if
413
414 call distributed_nullify(this%orbs_dist, this%norbsets)
415#ifdef HAVE_MPI
416 if (.not. gr%parallel_in_domains) then
417 call distributed_init(this%orbs_dist, this%norbsets, mpi_comm_world, "orbsets")
418 end if
419#endif
420
421 else
422
423 !%Variable DFTUBasisStates
424 !%Type block
425 !%Default none
426 !%Section Hamiltonian::DFT+U
427 !%Description
428 !% This block starts by a line containing a single integer describing the number of
429 !% orbital sets. One orbital set is a group of orbitals on which one adds a Hubbard U.
430 !% Each following line of this block contains the index of a state to be used to construct the
431 !% localized basis, followed by the index of the corresponding orbital set.
432 !% See DFTUBasisFromStates for details.
433 !%End
434 if (parse_block(namespace, 'DFTUBasisStates', blk) == 0) then
435 call parse_block_integer(blk, 0, 0, this%norbsets)
436 this%maxnorbs = parse_block_n(blk)-1
437 if (this%maxnorbs <1) then
438 write(message(1),'(a,i3,a,i3)') 'DFTUBasisStates must contains at least one state.'
439 call messages_fatal(1, namespace=namespace)
440 end if
441 safe_allocate(this%basisstates(1:this%maxnorbs))
442 safe_allocate(this%basisstates_os(1:this%maxnorbs))
443 do is = 1, this%maxnorbs
444 call parse_block_integer(blk, is, 0, this%basisstates(is))
445 call parse_block_integer(blk, is, 1, this%basisstates_os(is))
446 end do
447 call parse_block_end(blk)
448 else
449 write(message(1),'(a,i3,a,i3)') 'DFTUBasisStates must be specified if DFTUBasisFromStates=yes'
450 call messages_fatal(1, namespace=namespace)
451 end if
452
453 if (states_are_real(st)) then
454 call dorbitalbasis_build_empty(this%basis, gr, st%d%dim, this%norbsets, this%basisstates_os)
455 else
456 call zorbitalbasis_build_empty(this%basis, gr, st%d%dim, this%norbsets, this%basisstates_os)
457 end if
458
459 this%max_np = gr%np
460 this%nspins = st%d%nspin
461 this%spin_channels = st%d%spin_channels
462 this%nspecies = 1
463
464 this%orbsets => this%basis%orbsets
465
466 call distributed_nullify(this%orbs_dist, this%norbsets)
467
468 !We allocate the necessary ressources
469 if (states_are_real(st)) then
470 call dlda_u_allocate(this, st)
471 else
472 call zlda_u_allocate(this, st)
473 end if
474
475 call lda_u_loadbasis(this, namespace, space, st, gr, mc, ierr)
476 if (ierr /= 0) then
477 message(1) = "Unable to load DFT+U basis from selected states."
478 call messages_fatal(1)
479 end if
480
481 end if
482
483 ! Symmetrization of the occupation matrices
484 this%symmetrize_occ_matrices = st%symmetrize_density .or. kpoints%use_symmetries
485 if (this%basisfromstates) this%symmetrize_occ_matrices = .false.
486 if (this%symmetrize_occ_matrices) then
487 call build_symmetrization_map(this, ions, gr)
488 end if
489
490 safe_allocate(this%dc_alpha(this%norbsets))
491 this%dc_alpha = m_one
492
493 call messages_print_with_emphasis(namespace=namespace)
494
495 pop_sub(lda_u_init)
496 end subroutine lda_u_init
497
498
499 ! ---------------------------------------------------------
500 subroutine lda_u_init_coulomb_integrals(this, namespace, space, gr, st, psolver)
501 type(lda_u_t), target, intent(inout) :: this
502 type(namespace_t), intent(in) :: namespace
503 class(space_t), intent(in) :: space
504 type(grid_t), intent(in) :: gr
505 type(states_elec_t), intent(in) :: st
506 type(poisson_t), intent(in) :: psolver
507
508 logical :: complex_coulomb_integrals
509 integer :: ios
510 integer :: norbs
511
513
514 norbs = this%maxnorbs
515
516 if (.not. this%basisfromstates) then
517
518 if (this%level == dft_u_acbn0) then
519
520 complex_coulomb_integrals = .false.
521 do ios = 1, this%norbsets
522 if (this%orbsets(ios)%ndim > 1) complex_coulomb_integrals = .true.
523 end do
524
525 if (.not. complex_coulomb_integrals) then
526 write(message(1),'(a)') 'Computing the Coulomb integrals of the localized basis.'
527 call messages_info(1, namespace=namespace)
528 safe_allocate(this%coulomb(1:norbs,1:norbs,1:norbs,1:norbs, 1:this%norbsets))
529 if (states_are_real(st)) then
530 call dcompute_coulomb_integrals(this, namespace, space, gr, psolver)
531 else
532 call zcompute_coulomb_integrals(this, namespace, space, gr, psolver)
533 end if
534 else
535 assert(.not. states_are_real(st))
536 write(message(1),'(a)') 'Computing complex Coulomb integrals of the localized basis.'
537 call messages_info(1, namespace=namespace)
538 safe_allocate(this%zcoulomb(1:norbs, 1:norbs, 1:norbs, 1:norbs, 1:st%d%dim, 1:st%d%dim, 1:this%norbsets))
539 call compute_complex_coulomb_integrals(this, gr, st, psolver, namespace, space)
540 end if
541 end if
542
543 else
544 write(message(1),'(a)') 'Computing the Coulomb integrals of the localized basis.'
545 call messages_info(1, namespace=namespace)
546 safe_allocate(this%coulomb(1:norbs, 1:norbs, 1:norbs, 1:norbs, 1:this%norbsets))
547 if (states_are_real(st)) then
548 call dcompute_coulomb_integrals(this, namespace, space, gr, psolver)
549 else
550 call zcompute_coulomb_integrals(this, namespace, space, gr, psolver)
551 end if
552
553 end if
554
556
557 end subroutine lda_u_init_coulomb_integrals
558
559
560 ! ---------------------------------------------------------
561 subroutine lda_u_end(this)
562 implicit none
563 type(lda_u_t), intent(inout) :: this
564
565 push_sub(lda_u_end)
566
567 this%level = dft_u_none
568
569 safe_deallocate_a(this%dn)
570 safe_deallocate_a(this%zn)
571 safe_deallocate_a(this%dn_alt)
572 safe_deallocate_a(this%zn_alt)
573 safe_deallocate_a(this%dV)
574 safe_deallocate_a(this%zV)
575 safe_deallocate_a(this%coulomb)
576 safe_deallocate_a(this%zcoulomb)
577 safe_deallocate_a(this%renorm_occ)
578 safe_deallocate_a(this%dn_ij)
579 safe_deallocate_a(this%zn_ij)
580 safe_deallocate_a(this%dn_alt_ij)
581 safe_deallocate_a(this%zn_alt_ij)
582 safe_deallocate_a(this%dn_alt_ii)
583 safe_deallocate_a(this%zn_alt_ii)
584 safe_deallocate_a(this%basisstates)
585 safe_deallocate_a(this%basisstates_os)
586 safe_deallocate_a(this%dc_alpha)
587 safe_deallocate_a(this%inv_map_symm)
588 safe_deallocate_a(this%symm_weight)
589
590 nullify(this%orbsets)
591 call orbitalbasis_end(this%basis)
592
593 this%max_np = 0
594
595 if (.not. this%basisfromstates) then
596 call distributed_end(this%orbs_dist)
597 end if
598
599 pop_sub(lda_u_end)
600 end subroutine lda_u_end
601
602 ! ---------------------------------------------------------
607 subroutine lda_u_rebind_after_copy(this, ions)
608 type(lda_u_t), target, intent(inout) :: this
609 type(ions_t), target, intent(in) :: ions
610
611 integer :: ios, ispec
612
614
615 if (this%level == dft_u_none) then
617 return
618 end if
619
620 this%latt => ions%latt
621
622 if (allocated(this%basis%orbsets)) then
623 this%orbsets => this%basis%orbsets
624 else
625 nullify(this%orbsets)
627 return
628 end if
629
630 do ios = 1, this%norbsets
631 if (.not. associated(this%orbsets(ios)%spec)) cycle
632
633 ispec = this%orbsets(ios)%spec_index
634 if (ispec < 1 .or. ispec > ions%nspecies) then
635 call messages_not_implemented("lda_u_rebind_after_copy with invalid species index")
636 end if
637 if (.not. associated(ions%species(ispec)%s)) then
638 call messages_not_implemented("lda_u_rebind_after_copy with unassociated ion species")
639 end if
640 this%orbsets(ios)%spec => ions%species(ispec)%s
641 end do
642
644 end subroutine lda_u_rebind_after_copy
645
646 ! ---------------------------------------------------------
652 subroutine lda_u_accel_rebuild(this, kpt)
653 type(lda_u_t), intent(inout) :: this
654 type(distributed_t), intent(in) :: kpt
655
656 integer :: ios, ik
657
658 push_sub(lda_u_accel_rebuild)
659
660 if (.not. accel_is_enabled()) then
661 pop_sub(lda_u_accel_rebuild)
662 return
663 end if
664
665 if (this%level == dft_u_none .or. .not. associated(this%orbsets)) then
666 pop_sub(lda_u_accel_rebuild)
667 return
668 end if
669
670 do ios = 1, this%norbsets
671 call accel_detach_buffer(this%orbsets(ios)%dbuff_orb)
672 call accel_detach_buffer(this%orbsets(ios)%zbuff_orb)
673 call accel_detach_buffer(this%orbsets(ios)%sphere%buff_map)
674
675 if (allocated(this%orbsets(ios)%buff_eorb)) then
676 do ik = lbound(this%orbsets(ios)%buff_eorb, dim=1), ubound(this%orbsets(ios)%buff_eorb, dim=1)
677 call accel_detach_buffer(this%orbsets(ios)%buff_eorb(ik))
678 end do
679 safe_deallocate_a(this%orbsets(ios)%buff_eorb)
680 end if
681
682 if (allocated(this%orbsets(ios)%dorb)) then
683 call dorbitalset_transfer_to_device(this%orbsets(ios), kpt, this%orbsets(ios)%allocated_on_mesh)
684 else if (allocated(this%orbsets(ios)%zorb)) then
685 call zorbitalset_transfer_to_device(this%orbsets(ios), kpt, this%orbsets(ios)%allocated_on_mesh)
686 else
687 call messages_not_implemented("lda_u_accel_rebuild with missing orbital arrays")
688 end if
689 end do
690
691 pop_sub(lda_u_accel_rebuild)
692 end subroutine lda_u_accel_rebuild
693
694 ! When moving the ions, the basis must be reconstructed
695 subroutine lda_u_update_basis(this, space, gr, ions, st, psolver, namespace, kpoints, has_phase)
696 type(lda_u_t), target, intent(inout) :: this
697 class(space_t), intent(in) :: space
698 type(grid_t), intent(in) :: gr
699 type(ions_t), target, intent(in) :: ions
700 type(states_elec_t), intent(in) :: st
701 type(poisson_t), intent(in) :: psolver
702 type(namespace_t), intent(in) :: namespace
703 type(kpoints_t), intent(in) :: kpoints
704 logical, intent(in) :: has_phase
705
706 integer :: ios, maxorbs, nspin
707
708 if(this%level == dft_u_none) return
709
710 push_sub(lda_u_update_basis)
711
712 if(.not. this%basisfromstates) then
713 !We clean the orbital basis, to be able to reconstruct it
714 call orbitalbasis_end(this%basis)
715 nullify(this%orbsets)
716
717 !We now reconstruct the basis
718 if (states_are_real(st)) then
719 call dorbitalbasis_build(this%basis, namespace, ions, gr, st%d%kpt, st%d%dim, &
720 this%skipSOrbitals, this%useAllOrbitals, verbose = .false.)
721 else
722 call zorbitalbasis_build(this%basis, namespace, ions, gr, st%d%kpt, st%d%dim, &
723 this%skipSOrbitals, this%useAllOrbitals, verbose = .false.)
724 end if
725 this%orbsets => this%basis%orbsets
726 this%max_np = this%basis%max_np
727 end if
728
729 !In case of intersite interaction we need to reconstruct the basis
730 if (this%intersite) then
731 this%maxneighbors = 0
732 do ios = 1, this%norbsets
733 call orbitalset_init_intersite(this%orbsets(ios), namespace, space, ios, ions, gr%der, psolver, &
734 this%orbsets, this%norbsets, this%maxnorbs, this%intersite_radius, st%d%kpt, has_phase, &
735 this%sm_poisson, this%basisfromstates, this%basis%combine_j_orbitals)
736 this%maxneighbors = max(this%maxneighbors, this%orbsets(ios)%nneighbors)
737 end do
738
739 maxorbs = this%maxnorbs
740 nspin = this%nspins
741
742 if (states_are_real(st)) then
743 safe_deallocate_a(this%dn_ij)
744 safe_allocate(this%dn_ij(1:maxorbs,1:maxorbs,1:nspin,1:this%norbsets,1:this%maxneighbors))
745 this%dn_ij(1:maxorbs,1:maxorbs,1:nspin,1:this%norbsets,1:this%maxneighbors) = m_zero
746 safe_deallocate_a(this%dn_alt_ij)
747 safe_allocate(this%dn_alt_ij(1:maxorbs,1:maxorbs,1:nspin,1:this%norbsets,1:this%maxneighbors))
748 this%dn_alt_ij(1:maxorbs,1:maxorbs,1:nspin,1:this%norbsets,1:this%maxneighbors) = m_zero
749 safe_deallocate_a(this%dn_alt_ii)
750 safe_allocate(this%dn_alt_ii(1:2,1:maxorbs,1:nspin,1:this%norbsets,1:this%maxneighbors))
751 this%dn_alt_ii(1:2,1:maxorbs,1:nspin,1:this%norbsets,1:this%maxneighbors) = m_zero
752 else
753 safe_deallocate_a(this%zn_ij)
754 safe_allocate(this%zn_ij(1:maxorbs,1:maxorbs,1:nspin,1:this%norbsets,1:this%maxneighbors))
755 this%zn_ij(1:maxorbs,1:maxorbs,1:nspin,1:this%norbsets,1:this%maxneighbors) = m_z0
756 safe_deallocate_a(this%zn_alt_ij)
757 safe_allocate(this%zn_alt_ij(1:maxorbs,1:maxorbs,1:nspin,1:this%norbsets,1:this%maxneighbors))
758 this%zn_alt_ij(1:maxorbs,1:maxorbs,1:nspin,1:this%norbsets,1:this%maxneighbors) = m_z0
759 safe_deallocate_a(this%zn_alt_ii)
760 safe_allocate(this%zn_alt_ii(1:2,1:maxorbs,1:nspin,1:this%norbsets,1:this%maxneighbors))
761 this%zn_alt_ii(1:2,1:maxorbs,1:nspin,1:this%norbsets,1:this%maxneighbors) = m_z0
762 end if
763 end if
764
765 ! We rebuild the phase for the orbital projection, similarly to the one of the pseudopotentials
766 ! In case of a laser field, the phase is recomputed in hamiltonian_elec_update
767 if (has_phase) then
768 call lda_u_build_phase_correction(this, space, st%d, gr%der%boundaries, namespace, kpoints)
769 else
770 if(.not. this%basisfromstates) then
771 !In case there is no phase, we perform the orthogonalization here
772 if(this%basis%orthogonalization) then
773 call dloewdin_orthogonalize(this%basis, st%d%kpt, namespace)
774 else
775 if(debug%info .and. space%is_periodic()) then
776 call dloewdin_info(this%basis, st%d%kpt, namespace)
777 end if
778 end if
779 end if
780 end if
781
782 ! Rebuild the Coulomb integrals
783 if (allocated(this%coulomb)) then
784 safe_deallocate_a(this%coulomb)
785 end if
786 if (allocated(this%zcoulomb)) then
787 safe_deallocate_a(this%zcoulomb)
788 end if
789 call lda_u_init_coulomb_integrals(this, namespace, space, gr, st, psolver)
791 pop_sub(lda_u_update_basis)
792
793 end subroutine lda_u_update_basis
794
795 ! Interface for the X(update_occ_matrices) routines
796 subroutine lda_u_update_occ_matrices(this, namespace, mesh, st, phase, energy)
797 type(lda_u_t), intent(inout) :: this
798 type(namespace_t), intent(in) :: namespace
799 class(mesh_t), intent(in) :: mesh
800 type(states_elec_t), intent(inout) :: st
801 type(phase_t), intent(in) :: phase
802 type(energy_t), intent(inout) :: energy
803
804 if (this%level == dft_u_none .or. this%freeze_occ) return
806
807 if (states_are_real(st)) then
808 call dupdate_occ_matrices(this, namespace, mesh, st, energy%dft_u)
809 else
810 if (phase%is_allocated()) then
811 call zupdate_occ_matrices(this, namespace, mesh, st, energy%dft_u, phase)
812 else
813 call zupdate_occ_matrices(this, namespace, mesh, st, energy%dft_u)
814 end if
815 end if
816
818 end subroutine lda_u_update_occ_matrices
819
820
822 subroutine lda_u_build_phase_correction(this, space, std, boundaries, namespace, kpoints, vec_pot, vec_pot_var)
823 type(lda_u_t), intent(inout) :: this
824 class(space_t), intent(in) :: space
825 type(states_elec_dim_t), intent(in) :: std
826 type(boundaries_t), intent(in) :: boundaries
827 type(namespace_t), intent(in) :: namespace
828 type(kpoints_t), intent(in) :: kpoints
829 real(real64), optional, allocatable, intent(in) :: vec_pot(:)
830 real(real64), optional, allocatable, intent(in) :: vec_pot_var(:, :)
831
832 integer :: ios
833
834 if (boundaries%spiralBC) call messages_not_implemented("DFT+U with spiral boundary conditions", &
835 namespace=namespace)
836
838
839 write(message(1), '(a)') 'Debug: Building the phase correction for DFT+U orbitals.'
840 call messages_info(1, namespace=namespace, debug_only=.true.)
841
842 do ios = 1, this%norbsets
843 call orbitalset_update_phase(this%orbsets(ios), space%dim, std%kpt, kpoints, &
844 (std%ispin==spin_polarized), vec_pot, vec_pot_var)
845 call orbitalset_update_phase_shift(this%orbsets(ios), space%dim, std%kpt, kpoints, &
846 (std%ispin==spin_polarized), vec_pot, vec_pot_var)
847 end do
848
849 if (.not. this%basisfromstates) then
850 if (this%basis%orthogonalization) then
851 call zloewdin_orthogonalize(this%basis, std%kpt, namespace)
852 else
853 if (debug%info .and. space%is_periodic()) call zloewdin_info(this%basis, std%kpt, namespace)
854 end if
855 end if
856
858
859 end subroutine lda_u_build_phase_correction
860
861 ! ---------------------------------------------------------
862 subroutine compute_acbno_u_kanamori(this, st, kanamori)
863 type(lda_u_t), intent(in) :: this
864 type(states_elec_t), intent(in) :: st
865 real(real64), intent(out) :: kanamori(:,:)
866
867 if (this%nspins == 1) then
868 if (states_are_real(st)) then
869 call dcompute_acbno_u_kanamori_restricted(this, kanamori)
870 else
871 call zcompute_acbno_u_kanamori_restricted(this, kanamori)
872 end if
873 else
874 if (states_are_real(st)) then
875 call dcompute_acbno_u_kanamori(this, kanamori)
876 else
877 call zcompute_acbno_u_kanamori(this, kanamori)
878 end if
879 end if
880
881
882 end subroutine compute_acbno_u_kanamori
883
884 ! ---------------------------------------------------------
885 subroutine lda_u_freeze_occ(this)
886 type(lda_u_t), intent(inout) :: this
887
888 this%freeze_occ = .true.
889 end subroutine lda_u_freeze_occ
890
891 ! ---------------------------------------------------------
892 subroutine lda_u_freeze_u(this)
893 type(lda_u_t), intent(inout) :: this
894
895 this%freeze_u = .true.
896 end subroutine lda_u_freeze_u
897
898 ! ---------------------------------------------------------
899 subroutine lda_u_set_effectiveu(this, Ueff)
900 type(lda_u_t), intent(inout) :: this
901 real(real64), intent(in) :: ueff(:)
902
903 integer :: ios
904
905 push_sub(lda_u_set_effectiveu)
906
907 do ios = 1,this%norbsets
908 this%orbsets(ios)%Ueff = ueff(ios)
909 end do
910
911 pop_sub(lda_u_set_effectiveu)
912 end subroutine lda_u_set_effectiveu
913
914 ! ---------------------------------------------------------
915 subroutine lda_u_get_effectiveu(this, Ueff)
916 type(lda_u_t), intent(in) :: this
917 real(real64), intent(inout) :: ueff(:)
918
919 integer :: ios
920
921 push_sub(lda_u_get_effectiveu)
922
923 do ios = 1,this%norbsets
924 ueff(ios) = this%orbsets(ios)%Ueff
925 end do
926
927 pop_sub(lda_u_get_effectiveu)
928 end subroutine lda_u_get_effectiveu
929
930 ! ---------------------------------------------------------
931 subroutine lda_u_set_effectivev(this, Veff)
932 type(lda_u_t), intent(inout) :: this
933 real(real64), intent(in) :: veff(:)
934
935 integer :: ios, ncount
936
937 push_sub(lda_u_set_effectivev)
938
939 ncount = 0
940 do ios = 1, this%norbsets
941 this%orbsets(ios)%V_ij(1:this%orbsets(ios)%nneighbors,0) = veff(ncount+1:ncount+this%orbsets(ios)%nneighbors)
942 ncount = ncount + this%orbsets(ios)%nneighbors
943 end do
944
945 pop_sub(lda_u_set_effectivev)
946 end subroutine lda_u_set_effectivev
947
948 ! ---------------------------------------------------------
949 subroutine lda_u_get_effectivev(this, Veff)
950 type(lda_u_t), intent(in) :: this
951 real(real64), intent(inout) :: veff(:)
952
953 integer :: ios, ncount
954
955 push_sub(lda_u_get_effectivev)
956
957 ncount = 0
958 do ios = 1, this%norbsets
959 veff(ncount+1:ncount+this%orbsets(ios)%nneighbors) = this%orbsets(ios)%V_ij(1:this%orbsets(ios)%nneighbors,0)
960 ncount = ncount + this%orbsets(ios)%nneighbors
961 end do
962
963 pop_sub(lda_u_get_effectivev)
964 end subroutine lda_u_get_effectivev
965
966 ! ---------------------------------------------------------
967 subroutine lda_u_write_info(this, iunit, namespace)
968 type(lda_u_t), intent(in) :: this
969 integer, optional, intent(in) :: iunit
970 type(namespace_t), optional, intent(in) :: namespace
971
972 push_sub(lda_u_write_info)
973
974 write(message(1), '(1x)')
975 call messages_info(1, iunit=iunit, namespace=namespace)
976 if (this%level == dft_u_empirical) then
977 write(message(1), '(a)') "Method:"
978 write(message(2), '(a)') " [1] Dudarev et al., Phys. Rev. B 57, 1505 (1998)"
979 call messages_info(2, iunit=iunit, namespace=namespace)
980 else
981 if (.not. this%intersite) then
982 write(message(1), '(a)') "Method:"
983 write(message(2), '(a)') " [1] Agapito et al., Phys. Rev. X 5, 011006 (2015)"
984 else
985 write(message(1), '(a)') "Method:"
986 write(message(2), '(a)') " [1] Tancogne-Dejean, and Rubio, Phys. Rev. B 102, 155117 (2020)"
987 end if
988 call messages_info(2, iunit=iunit, namespace=namespace)
989 end if
990 write(message(1), '(a)') "Implementation:"
991 write(message(2), '(a)') " [1] Tancogne-Dejean, Oliveira, and Rubio, Phys. Rev. B 69, 245133 (2017)"
992 write(message(3), '(1x)')
993 call messages_info(3, iunit=iunit, namespace=namespace)
995 pop_sub(lda_u_write_info)
996
997 end subroutine lda_u_write_info
998
999 ! ---------------------------------------------------------
1000 subroutine lda_u_loadbasis(this, namespace, space, st, mesh, mc, ierr)
1001 type(lda_u_t), intent(inout) :: this
1002 type(namespace_t), intent(in) :: namespace
1003 class(space_t), intent(in) :: space
1004 type(states_elec_t), intent(in) :: st
1005 class(mesh_t), intent(in) :: mesh
1006 type(multicomm_t), intent(in) :: mc
1007 integer, intent(out) :: ierr
1008
1009 integer :: err, wfns_file, is, ist, idim, ik, ios, iorb
1011 character(len=256) :: lines(3)
1012 character(len=256), allocatable :: restart_file(:, :)
1013 logical, allocatable :: restart_file_present(:, :)
1014 character(len=12) :: filename
1015 character(len=1) :: char
1016 character(len=50) :: str
1017 type(orbitalset_t), pointer :: os
1018 integer, allocatable :: count(:)
1019 real(real64) :: norm, center(space%dim)
1020 real(real64), allocatable :: dpsi(:,:,:)
1021 complex(real64), allocatable :: zpsi(:,:,:)
1022
1023
1024 push_sub(lda_u_loadbasis)
1025
1026 ierr = 0
1027
1028 message(1) = "Debug: Loading DFT+U basis from states."
1029 call messages_info(1, debug_only=.true.)
1030
1031 call restart_gs%init(namespace, restart_proj, restart_type_load, mc, err, mesh=mesh)
1032
1033 ! If any error occured up to this point then it is not worth continuing,
1034 ! as there something fundamentally wrong with the restart files
1035 if (err /= 0) then
1036 call restart_gs%end()
1037 message(1) = "Error loading DFT+U basis from states, cannot proceed with the calculation"
1038 call messages_fatal(1)
1039 pop_sub(lda_u_loadbasis)
1040 return
1041 end if
1042 ! only supported for obf format
1043 if (restart_gs%file_format_states /= option__restartfileformatstates__obf) then
1044 message(1) = "Error: loading DFT+U basis only supported for OBF restart format"
1045 message(2) = "Please set RestartFileFormatStates = obf and re-run the GS calculation."
1046 call messages_fatal(2)
1047 end if
1048
1049 ! open files to read
1050 wfns_file = restart_gs%open('wfns')
1051 call restart_gs%read(wfns_file, lines, 2, err)
1052 if (err /= 0) then
1053 ierr = ierr - 2**5
1054 else if (states_are_real(st)) then
1055 read(lines(2), '(a)') str
1056 if (str(2:8) == 'Complex') then
1057 message(1) = "Cannot read real states from complex wavefunctions."
1058 call messages_fatal(1, namespace=namespace)
1059 else if (str(2:5) /= 'Real') then
1060 message(1) = "Restart file 'wfns' does not specify real/complex; cannot check compatibility."
1061 call messages_warning(1, namespace=namespace)
1062 end if
1063 end if
1064 ! complex can be restarted from real, so there is no problem.
1065
1066 ! If any error occured up to this point then it is not worth continuing,
1067 ! as there something fundamentally wrong with the restart files
1068 if (err /= 0) then
1069 call restart_gs%close(wfns_file)
1070 call restart_gs%end()
1071 message(1) = "Error loading DFT+U basis from states, cannot proceed with the calculation"
1072 call messages_fatal(1)
1073 pop_sub(lda_u_loadbasis)
1074 return
1075 end if
1076
1077 safe_allocate(restart_file(1:st%d%dim, 1:st%nst))
1078 safe_allocate(restart_file_present(1:st%d%dim, 1:st%nst))
1079 restart_file_present = .false.
1080
1081 ! Next we read the list of states from the files.
1082 ! Errors in reading the information of a specific state from the files are ignored
1083 ! at this point, because later we will skip reading the wavefunction of that state.
1084 do
1085 call restart_gs%read(wfns_file, lines, 1, err)
1086 if (err == 0) then
1087 read(lines(1), '(a)') char
1088 if (char == '%') then
1089 !We reached the end of the file
1090 exit
1091 else
1092 read(lines(1), *) ik, char, ist, char, idim, char, filename
1093 end if
1094 end if
1096 if (any(this%basisstates==ist) .and. ik == 1) then
1097 restart_file(idim, ist) = trim(filename)
1098 restart_file_present(idim, ist) = .true.
1099 end if
1100 end do
1101 call restart_gs%close(wfns_file)
1102
1103 !We loop over the states we need
1104 safe_allocate(count(1:this%norbsets))
1105 count = 0
1106 do is = 1, this%maxnorbs
1107 ist = this%basisstates(is)
1108 ios = this%basisstates_os(is)
1109 count(ios) = count(ios)+1
1110 do idim = 1, st%d%dim
1111
1112 if (.not. restart_file_present(idim, ist)) then
1113 write(message(1), '(a,i3,a)') "Cannot read states ", ist, "from the projection folder"
1114 call messages_fatal(1, namespace=namespace)
1115 end if
1116
1117 if (states_are_real(st)) then
1118 call restart_gs%read_mesh_function(restart_file(idim, ist), mesh, &
1119 this%orbsets(ios)%dorb(:,idim,count(ios)), err)
1120 else
1121 call restart_gs%read_mesh_function(restart_file(idim, ist), mesh, &
1122 this%orbsets(ios)%zorb(:,idim,count(ios)), err)
1123 end if
1124 if (err /= 0) then
1125 message(1) = "Error loading mesh function "//trim(restart_file(idim, ist))
1126 call messages_fatal(1)
1127 end if
1128 end do
1129 end do
1130 safe_deallocate_a(count)
1131 safe_deallocate_a(restart_file)
1132 safe_deallocate_a(restart_file_present)
1133 call restart_gs%end()
1134
1135 ! Normalize the orbitals. This is important if we use Wannier orbitals instead of KS states
1136 if(this%basis%normalize) then
1137 do ios = 1, this%norbsets
1138 do iorb = 1, this%orbsets(ios)%norbs
1139 if (states_are_real(st)) then
1140 norm = dmf_nrm2(mesh, st%d%dim, this%orbsets(ios)%dorb(:,:,iorb))
1141 call lalg_scal(mesh%np, st%d%dim, m_one/norm, this%orbsets(ios)%dorb(:,:,iorb))
1142 else
1143 norm = zmf_nrm2(mesh, st%d%dim, this%orbsets(ios)%zorb(:,:,iorb))
1144 call lalg_scal(mesh%np, st%d%dim, m_one/norm, this%orbsets(ios)%zorb(:,:,iorb))
1145 end if
1146 end do
1147 end do
1148 end if
1149
1150 ! We rotate the orbitals in the complex plane to have them as close as possible to real functions
1151 if(states_are_complex(st) .and. st%d%dim == 1) then
1152 do ios = 1, this%norbsets
1153 do iorb = 1, this%orbsets(ios)%norbs
1154 call zmf_fix_phase(mesh, this%orbsets(ios)%zorb(:,1,iorb))
1155 end do
1156 end do
1157 end if
1158
1159 ! We determine the center of charge by computing <w|r|w>
1160 ! We could also determine the spread by \Omega = <w|r^2|w> - <w|r|w>^2
1161 do ios = 1, this%norbsets
1162 if (states_are_real(st)) then
1163 call dorbitalset_get_center_of_mass(this%orbsets(ios), space, mesh, this%latt)
1164 else
1165 call zorbitalset_get_center_of_mass(this%orbsets(ios), space, mesh, this%latt)
1166 end if
1167 end do
1168
1169 message(1) = "Debug: Converting the Wannier states to submeshes."
1170 call messages_info(1, debug_only=.true.)
1171
1172 ! We now transfer the states to a submesh centered on the center of mass of the Wannier orbitals
1173 this%max_np = 0
1174 do ios = 1, this%norbsets
1175 os => this%orbsets(ios)
1176 center = os%sphere%center
1177 safe_deallocate_a(os%sphere%center)
1178 if (states_are_real(st)) then
1179 safe_allocate(dpsi(1:mesh%np, 1:os%ndim, 1:os%norbs))
1180 dpsi(1:mesh%np, 1:os%ndim, 1:os%norbs) = os%dorb(1:mesh%np, 1:os%ndim, 1:os%norbs)
1181
1182 safe_deallocate_a(os%dorb)
1183 !We initialise the submesh corresponding to the orbital
1184 call submesh_init(os%sphere, space, mesh, this%latt, center, os%radius)
1185 safe_allocate(os%dorb(1:os%sphere%np, 1:os%ndim, 1:os%norbs))
1186 do iorb = 1, os%norbs
1187 do idim = 1, os%ndim
1188 call dsubmesh_copy_from_mesh(os%sphere, dpsi(:,idim,iorb), os%dorb(:,idim, iorb))
1189 end do
1190 end do
1191 safe_deallocate_a(dpsi)
1192 else
1193 safe_allocate(zpsi(1:mesh%np, 1:os%ndim, 1:os%norbs))
1194 zpsi(1:mesh%np, 1:os%ndim, 1:os%norbs) = os%zorb(1:mesh%np, 1:os%ndim, 1:os%norbs)
1195 safe_deallocate_a(os%zorb)
1196 !We initialise the submesh corresponding to the orbital
1197 call submesh_init(os%sphere, space, mesh, this%latt, center, os%radius)
1198 safe_allocate(os%zorb(1:os%sphere%np, 1:os%ndim, 1:os%norbs))
1199 do iorb = 1, os%norbs
1200 do idim = 1, os%ndim
1201 call zsubmesh_copy_from_mesh(os%sphere, zpsi(:,idim,iorb), os%zorb(:,idim, iorb))
1202 end do
1203 end do
1204 safe_deallocate_a(zpsi)
1205
1206 safe_allocate(os%phase(1:os%sphere%np, st%d%kpt%start:st%d%kpt%end))
1207 safe_allocate(os%eorb_submesh(1:os%sphere%np, 1:os%ndim, 1:os%norbs, st%d%kpt%start:st%d%kpt%end))
1208 end if
1209 os%use_submesh = .true. ! We are now on a submesh
1210 this%max_np = max(this%max_np, os%sphere%np)
1211 end do
1212
1213 this%basis%use_submesh = .true.
1214
1215 ! If we use GPUs, we need to transfert the orbitals on the device
1216 if (accel_is_enabled() .and. st%d%dim == 1) then
1217 do ios = 1, this%norbsets
1218 os => this%orbsets(ios)
1219
1220 os%ldorbs = max(accel_padded_size(os%sphere%np), 1)
1221 os%ldorbs_eorb = max(accel_padded_size(os%sphere%np), 1)
1222 if (states_are_real(st)) then
1223 call accel_create_buffer(os%dbuff_orb, accel_mem_read_only, type_float, os%ldorbs*os%norbs)
1224 else
1225 call accel_create_buffer(os%zbuff_orb, accel_mem_read_only, type_cmplx, os%ldorbs*os%norbs)
1226 safe_allocate(os%buff_eorb(st%d%kpt%start:st%d%kpt%end))
1227
1228 do ik= st%d%kpt%start, st%d%kpt%end
1229 call accel_create_buffer(os%buff_eorb(ik), accel_mem_read_only, type_cmplx, os%ldorbs_eorb*os%norbs)
1230 end do
1231 end if
1232
1233 call accel_create_buffer(os%sphere%buff_map, accel_mem_read_only, type_integer, max(os%sphere%np, 1))
1234 call accel_write_buffer(os%sphere%buff_map, os%sphere%np, os%sphere%map)
1235
1236 do iorb = 1, os%norbs
1237 if(states_are_complex(st)) then
1238 call accel_write_buffer(os%zbuff_orb, os%sphere%np, os%zorb(:, 1, iorb), &
1239 offset = (iorb - 1)*os%ldorbs)
1240 else
1241 call accel_write_buffer(os%dbuff_orb, os%sphere%np, os%dorb(:, 1, iorb), &
1242 offset = (iorb - 1)*os%ldorbs)
1243 end if
1244 end do
1245 end do
1246 end if
1247
1248
1249 message(1) = "Debug: Loading DFT+U basis from states done."
1250 call messages_info(1, debug_only=.true.)
1251
1252 pop_sub(lda_u_loadbasis)
1253 end subroutine lda_u_loadbasis
1254
1256 subroutine build_symmetrization_map(this, ions, gr)
1257 type(lda_u_t), intent(inout) :: this
1258 type(ions_t), intent(in) :: ions
1259 type(grid_t), intent(in) :: gr
1260
1261 integer :: nsym, iop, ios, iatom, iatom_sym, ios_sym
1262
1263 push_sub(build_symmetrization_map)
1264
1265 nsym = ions%symm%nops
1266 safe_allocate(this%inv_map_symm(1:this%norbsets, 1:nsym))
1267 this%inv_map_symm = -1
1268
1269 this%nsym = nsym
1270
1271 do ios = 1, this%norbsets
1272 iatom = this%orbsets(ios)%iatom
1273 do iop = 1, nsym
1274 iatom_sym = ions%inv_map_symm_atoms(iatom, iop)
1275
1276 do ios_sym = 1, this%norbsets
1277 if (this%orbsets(ios_sym)%iatom == iatom_sym .and. this%orbsets(ios_sym)%norbs == this%orbsets(ios)%norbs &
1278 .and. this%orbsets(ios_sym)%nn == this%orbsets(ios)%nn .and. this%orbsets(ios_sym)%ll == this%orbsets(ios)%ll &
1279 .and. is_close(this%orbsets(ios_sym)%jj, this%orbsets(ios)%jj)) then
1280 this%inv_map_symm(ios, iop) = ios_sym
1281 exit
1282 end if
1283 end do
1284 assert(this%inv_map_symm(ios, iop) > 0)
1285 end do
1286 end do
1287
1288 safe_allocate(this%symm_weight(1:this%maxnorbs, 1:this%maxnorbs, 1:this%nsym, 1:this%norbsets))
1289 this%symm_weight = m_zero
1290
1291 do ios = 1, this%norbsets
1292 ! s-orbitals
1293 if (this%orbsets(ios)%norbs == 1) then
1294 this%symm_weight(1,1, 1:this%nsym, ios) = m_one
1295 cycle
1296 end if
1297
1298 ! Not implemented yet
1299 if (this%orbsets(ios)%ndim > 1) cycle
1300
1301 call orbitals_get_symm_weight(this%orbsets(ios), ions%space, ions%latt, gr, ions%symm, this%symm_weight(:,:,:,ios))
1302 end do
1303
1305 end subroutine build_symmetrization_map
1306
1308 subroutine orbitals_get_symm_weight(os, space, latt, gr, symm, weight)
1309 type(orbitalset_t), intent(in) :: os
1310 type(space_t), intent(in) :: space
1311 type(lattice_vectors_t), intent(in) :: latt
1312 type(grid_t), intent(in) :: gr
1313 type(symmetries_t), intent(in) :: symm
1314 real(real64), intent(inout) :: weight(:,:,:)
1315
1316 integer :: im, imp, iop, mm
1317 real(real64), allocatable :: orb(:,:), orb_sym(:), ylm(:)
1318 type(submesh_t) :: sphere
1319 real(real64) :: rc, norm, origin(space%dim)
1320
1321 push_sub(orbitals_get_symm_weight)
1322
1323 assert(os%ndim == 1)
1324
1325 safe_allocate(orb_sym(1:gr%np))
1326 safe_allocate(orb(1:gr%np, 1:os%norbs))
1327
1328 assert(2*os%ll+1 == os%norbs)
1329
1330 ! We generate an artificial submesh to compute the symmetries on it
1331 ! The radius is such that we fit in 20 points
1332 rc = (50.0_real64 * m_three/m_four/m_pi*product(gr%spacing))**m_third
1333 origin = m_zero
1334 call submesh_init(sphere, space, gr, latt, origin, rc)
1335
1336 safe_allocate(ylm(1:sphere%np))
1337
1338 ! We then compute the spherical harmonics in this submesh
1339 do im = 1, os%norbs
1340 mm = im-1-os%ll
1341 call loct_ylm(sphere%np, sphere%rel_x(1,1), sphere%r(1), os%ll, mm, ylm(1))
1342 orb(:,im) = m_zero
1343 call submesh_add_to_mesh(sphere, ylm, orb(:,im))
1344 norm = dmf_nrm2(gr, orb(:,im))
1345 call lalg_scal(gr%np, m_one / norm, orb(:,im))
1346 end do
1347 safe_deallocate_a(ylm)
1348
1349 ! Then we put them on the grid, rotate them
1350 do im = 1, os%norbs
1351 do iop = 1, symm%nops
1352 call dgrid_symmetrize_single(gr, iop, orb(:,im), orb_sym)
1353
1354 do imp = 1, os%norbs
1355 weight(im, imp, iop) = dmf_dotp(gr, orb(:,imp), orb_sym, reduce=.false.)
1356 end do
1357 end do
1358 end do
1359
1360 call gr%allreduce(weight)
1361
1362 safe_deallocate_a(orb)
1363 safe_deallocate_a(orb_sym)
1364
1366 end subroutine orbitals_get_symm_weight
1367
1368#include "dft_u_noncollinear_inc.F90"
1369
1370#include "undef.F90"
1371#include "real.F90"
1372#include "lda_u_inc.F90"
1373
1374#include "undef.F90"
1375#include "complex.F90"
1376#include "lda_u_inc.F90"
1377end module lda_u_oct_m
scales a vector by a constant
Definition: lalg_basic.F90:159
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:182
subroutine, public accel_detach_buffer(this)
Clear a buffer handle without freeing device memory.
Definition: accel.F90:1049
pure logical function, public accel_is_enabled()
Definition: accel.F90:402
integer, parameter, public accel_mem_read_only
Definition: accel.F90:185
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
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:124
type(debug_t), save, public debug
Definition: debug.F90:158
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public distributed_end(this)
subroutine, public distributed_nullify(this, total)
subroutine, public distributed_init(this, total, comm, tag, scalapack_compat)
Distribute N instances across M processes of communicator comm
integer, parameter, public spinors
integer, parameter, public spin_polarized
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_four
Definition: global.F90:195
real(real64), parameter, public m_third
Definition: global.F90:198
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
complex(real64), parameter, public m_z0
Definition: global.F90:201
real(real64), parameter, public m_epsilon
Definition: global.F90:207
real(real64), parameter, public m_one
Definition: global.F90:192
real(real64), parameter, public m_three
Definition: global.F90:194
This module implements the underlying real-space grid.
Definition: grid.F90:119
subroutine, public dgrid_symmetrize_single(gr, iop, field, symm_field, suppress_warning)
Definition: grid.F90:726
integer, parameter, public dft_u_amf
Definition: lda_u.F90:210
subroutine, public zlda_u_commute_r(this, mesh, space, d, namespace, psib, gpsib)
This routine computes [r,V_lda+u] .
Definition: lda_u.F90:5088
subroutine zlda_u_allocate(this, st)
Definition: lda_u.F90:5658
subroutine, public lda_u_get_effectiveu(this, Ueff)
Definition: lda_u.F90:1011
subroutine compute_complex_coulomb_integrals(this, gr, st, psolver, namespace, space)
Definition: lda_u.F90:1485
subroutine, public dlda_u_force(this, namespace, space, mesh, st, iq, psib, grad_psib, force)
Definition: lda_u.F90:3273
subroutine, public dlda_u_apply(this, d, mesh, psib, hpsib)
Apply the +U nonlocal potential to psib and adds the result to hpsib.
Definition: lda_u.F90:1853
subroutine, public lda_u_set_effectiveu(this, Ueff)
Definition: lda_u.F90:995
subroutine build_symmetrization_map(this, ions, gr)
Builds a mapping between the orbital sets based on symmetries.
Definition: lda_u.F90:1352
subroutine dcompute_acbno_u_kanamori_restricted(this, kanamori)
This routine computes the Kanamori U, Up, and J.
Definition: lda_u.F90:2874
subroutine, public lda_u_init(this, namespace, space, level, gr, ions, st, mc, kpoints)
Definition: lda_u.F90:286
subroutine, public dlda_u_get_occupations(this, occ)
Definition: lda_u.F90:3578
subroutine, public dlda_u_rvu(this, mesh, space, d, namespace, psib, gpsib)
This routine computes .
Definition: lda_u.F90:3411
subroutine, public lda_u_update_basis(this, space, gr, ions, st, psolver, namespace, kpoints, has_phase)
Definition: lda_u.F90:791
subroutine orbitals_get_symm_weight(os, space, latt, gr, symm, weight)
Computes the weight of each rotated orbitals in the basis of the same localized subspace.
Definition: lda_u.F90:1404
subroutine dupdate_occ_matrices(this, namespace, mesh, st, lda_u_energy, phase)
This routine computes the values of the occupation matrices.
Definition: lda_u.F90:1913
subroutine, public lda_u_accel_rebuild(this, kpt)
Rebuild DFT+U accelerator buffers after intrinsic assignment.
Definition: lda_u.F90:748
subroutine, public zlda_u_apply(this, d, mesh, psib, hpsib)
Apply the +U nonlocal potential to psib and adds the result to hpsib.
Definition: lda_u.F90:3825
subroutine, public zcompute_dftu_energy(this, energy, st)
This routine computes the value of the double counting term in the DFT+U energy.
Definition: lda_u.F90:4286
integer, parameter, public dft_u_empirical
Definition: lda_u.F90:205
subroutine zcompute_acbno_u_kanamori_restricted(this, kanamori)
This routine computes the Kanamori U, Up, and J.
Definition: lda_u.F90:4846
subroutine dcompute_acbno_u_kanamori(this, kanamori)
This routine computes the Kanamori U, Up, and J.
Definition: lda_u.F90:2775
subroutine, public zlda_u_commute_r_single(this, mesh, space, d, namespace, ist, ik, psi, gpsi, has_phase)
Definition: lda_u.F90:5037
subroutine, public lda_u_freeze_occ(this)
Definition: lda_u.F90:981
integer, parameter, public dft_u_mix
Definition: lda_u.F90:210
subroutine, public lda_u_rebind_after_copy(this, ions)
Rebind non-owning pointers after intrinsic assignment of lda_u_t.
Definition: lda_u.F90:703
subroutine, public lda_u_freeze_u(this)
Definition: lda_u.F90:988
subroutine, public zlda_u_rvu(this, mesh, space, d, namespace, psib, gpsib)
This routine computes .
Definition: lda_u.F90:5413
subroutine dlda_u_allocate(this, st)
Definition: lda_u.F90:3635
subroutine, public compute_acbno_u_kanamori(this, st, kanamori)
Definition: lda_u.F90:958
subroutine, public zlda_u_update_potential(this, st)
This routine computes the potential that, once multiplied by the projector Pmm' and summed over m and...
Definition: lda_u.F90:4360
subroutine lda_u_init_coulomb_integrals(this, namespace, space, gr, st, psolver)
Definition: lda_u.F90:596
subroutine, public lda_u_write_info(this, iunit, namespace)
Definition: lda_u.F90:1063
subroutine, public dlda_u_commute_r(this, mesh, space, d, namespace, psib, gpsib)
This routine computes [r,V_lda+u] .
Definition: lda_u.F90:3116
subroutine, public zlda_u_force(this, namespace, space, mesh, st, iq, psib, grad_psib, force)
Definition: lda_u.F90:5275
subroutine dcompute_coulomb_integrals(this, namespace, space, gr, psolver)
Definition: lda_u.F90:2949
subroutine, public dlda_u_set_occupations(this, occ)
Definition: lda_u.F90:3521
subroutine, public lda_u_get_effectivev(this, Veff)
Definition: lda_u.F90:1045
subroutine lda_u_loadbasis(this, namespace, space, st, mesh, mc, ierr)
Definition: lda_u.F90:1096
subroutine, public dlda_u_commute_r_single(this, mesh, space, d, namespace, ist, ik, psi, gpsi, has_phase)
Definition: lda_u.F90:3065
subroutine, public lda_u_build_phase_correction(this, space, std, boundaries, namespace, kpoints, vec_pot, vec_pot_var)
Build the phase correction to the global phase for all orbitals.
Definition: lda_u.F90:918
subroutine, public lda_u_end(this)
Definition: lda_u.F90:657
subroutine, public lda_u_set_effectivev(this, Veff)
Definition: lda_u.F90:1027
subroutine, public lda_u_update_occ_matrices(this, namespace, mesh, st, phase, energy)
Definition: lda_u.F90:892
integer, parameter, public dft_u_acbn0
Definition: lda_u.F90:205
subroutine zcompute_acbno_u_kanamori(this, kanamori)
This routine computes the Kanamori U, Up, and J.
Definition: lda_u.F90:4747
subroutine zcompute_coulomb_integrals(this, namespace, space, gr, psolver)
Definition: lda_u.F90:4921
subroutine, public zlda_u_get_occupations(this, occ)
Definition: lda_u.F90:5601
subroutine zupdate_occ_matrices(this, namespace, mesh, st, lda_u_energy, phase)
This routine computes the values of the occupation matrices.
Definition: lda_u.F90:3885
subroutine, public dcompute_dftu_energy(this, energy, st)
This routine computes the value of the double counting term in the DFT+U energy.
Definition: lda_u.F90:2314
subroutine, public dlda_u_update_potential(this, st)
This routine computes the potential that, once multiplied by the projector Pmm' and summed over m and...
Definition: lda_u.F90:2388
subroutine, public zlda_u_set_occupations(this, occ)
Definition: lda_u.F90:5544
System information (time, memory, sysname)
Definition: loct.F90:117
subroutine, public dloewdin_orthogonalize(basis, kpt, namespace)
Definition: loewdin.F90:214
subroutine, public zloewdin_info(basis, kpt, namespace)
Definition: loewdin.F90:747
subroutine, public zloewdin_orthogonalize(basis, kpt, namespace)
Definition: loewdin.F90:539
subroutine, public dloewdin_info(basis, kpt, namespace)
Definition: loewdin.F90:415
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
This module defines various routines, operating on mesh functions.
subroutine, public zmf_fix_phase(mesh, ff)
Fix the phase of complex function.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:898
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1091
character(len=512), private msg
Definition: messages.F90:167
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1063
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
subroutine, public orbitalbasis_end(this)
subroutine, public zorbitalbasis_build(this, namespace, ions, mesh, kpt, ndim, skip_s_orb, use_all_orb, verbose)
This routine is an interface for constructing the orbital basis.
subroutine, public dorbitalbasis_build(this, namespace, ions, mesh, kpt, ndim, skip_s_orb, use_all_orb, verbose)
This routine is an interface for constructing the orbital basis.
subroutine, public dorbitalbasis_build_empty(this, mesh, ndim, norbsets, map_os, verbose)
This routine constructd an empty orbital basis.
subroutine, public zorbitalbasis_build_empty(this, mesh, ndim, norbsets, map_os, verbose)
This routine constructd an empty orbital basis.
subroutine, public orbitalbasis_init(this, namespace, periodic_dim)
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, public dorbitalset_transfer_to_device(os, kpt, use_mesh)
Allocate and transfer the orbitals to the device.
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 zorbitalset_transfer_to_device(os, kpt, use_mesh)
Allocate and transfer the orbitals to the device.
integer, parameter, public sm_poisson_psolver
integer, parameter, public sm_poisson_isf
subroutine, public zorbitalset_get_center_of_mass(os, space, mesh, latt)
subroutine, public orbitalset_init_intersite(this, namespace, space, ind, ions, der, psolver, os, nos, maxnorbs, rcut, kpt, has_phase, sm_poisson, basis_from_states, combine_j_orbitals)
integer, parameter, public sm_poisson_direct
subroutine, public dorbitalset_get_center_of_mass(os, space, mesh, latt)
integer, parameter, public sm_poisson_fft
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:615
integer, parameter, public restart_gs
Definition: restart.F90:156
integer, parameter, public restart_proj
Definition: restart.F90:156
integer, parameter, public restart_type_load
Definition: restart.F90:183
pure logical function, public states_are_complex(st)
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public zsubmesh_copy_from_mesh(this, phi, sphi, conjugate)
Definition: submesh.F90:1821
subroutine, public dsubmesh_copy_from_mesh(this, phi, sphi, conjugate)
Definition: submesh.F90:1284
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
Definition: submesh.F90:282
type(type_t), public type_float
Definition: types.F90:135
type(type_t), public type_cmplx
Definition: types.F90:136
type(type_t), public type_integer
Definition: types.F90:137
This module defines the unit system, used for input and output.
type(unit_system_t), public units_inp
the units systems for reading and writing
This class contains information about the boundary conditions.
Definition: boundaries.F90:159
Distribution of N instances over mpi_grpsize processes, for the local rank mpi_grprank....
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:171
Class to describe DFT+U parameters.
Definition: lda_u.F90:218
Describes mesh distribution to nodes.
Definition: mesh.F90:187
Stores all communicators and groups.
Definition: multicomm.F90:208
A container for the phase.
Definition: phase.F90:181
class for organizing spins and k-points
The states_elec_t class contains all electronic wave functions.
A submesh is a type of mesh, used for the projectors in the pseudopotentials It contains points on a ...
Definition: submesh.F90:177
int true(void)