39 use,
intrinsic :: iso_fortran_env
89 integer,
parameter :: &
90 INITRHO_PARAMAGNETIC = 1, &
98 logical :: complex_ylms
99 logical :: initialized
102 integer,
allocatable :: atom(:)
103 integer,
allocatable :: level(:)
104 integer,
allocatable :: ddim(:)
105 logical :: alternative
106 logical :: derivative
107 integer,
allocatable :: cst(:, :)
108 integer,
allocatable :: ck(:, :)
109 real(4),
allocatable :: dbuff_single(:, :, :, :)
110 complex(4),
allocatable :: zbuff_single(:, :, :, :)
111 real(real64),
allocatable :: dbuff(:, :, :, :)
112 complex(real64),
allocatable :: zbuff(:, :, :, :)
113 logical :: save_memory
114 logical :: initialized_orbitals
115 real(real64) :: orbital_scale_factor
118 logical,
allocatable :: is_empty(:)
124 real(real64),
allocatable :: radius(:)
125 real(real64) :: lapdist
129 integer,
allocatable :: basis_atom(:)
130 integer,
allocatable :: basis_orb(:)
131 integer,
allocatable :: atom_orb_basis(:, :)
132 integer,
allocatable :: norb_atom(:)
134 integer :: lsize(1:2)
135 integer :: nproc(1:2)
136 integer :: myroc(1:2)
137 integer :: desc(1:BLACS_DLEN)
138 logical,
allocatable :: calc_atom(:)
139 real(real64) :: diag_tol
140 type(submesh_t),
allocatable :: sphere(:)
141 type(batch_t),
allocatable :: orbitals(:)
142 logical,
allocatable :: is_orbital_initialized(:)
149 real(real64),
parameter,
public :: default_eigenval = 1.0e10_real64
154 subroutine lcao_init(this, namespace, space, gr, ions, st, st_start)
155 type(lcao_t),
intent(out) :: this
156 type(namespace_t),
intent(in) :: namespace
157 type(electron_space_t),
intent(in) :: space
158 type(grid_t),
intent(in) :: gr
159 type(ions_t),
intent(in) :: ions
160 type(states_elec_t),
intent(in) :: st
161 integer,
intent(in) :: st_start
163 integer :: ia, n, iorb, jj, maxj, idim
164 integer :: ii, ll, mm, norbs, ii_tmp
165 integer :: mode_default
166 real(real64) :: max_orb_radius, maxradius
171 this%initialized = .
true.
176 mode_default = option__lcaostart__lcao_states
177 if (ions%only_user_def) mode_default = option__lcaostart__lcao_none
234 if (this%mode == option__lcaostart__lcao_states_batch .and. st_start > 1)
then
235 message(1) =
"LCAOStart = lcao_states_batch not compatible with this run."
236 message(2) =
"Please use LCAOStart = lcao_states instead"
241 this%alternative = this%mode == option__lcaostart__lcao_states_batch
243 if (this%mode == option__lcaostart__lcao_none)
then
253 if (st%d%ispin ==
spinors .and. this%alternative)
then
254 message(1) =
"LCAOStart = lcao_states_batch is not working for spinors."
257 if (space%is_periodic() .and. this%alternative)
then
275 call parse_variable(namespace,
'LCAOComplexYlms', .false., this%complex_ylms)
277 this%complex_ylms = .false.
288 call parse_variable(namespace,
'LCAOSaveMemory', .false., this%save_memory)
291 if (
debug%info .and. st%system_grp%is_root())
then
292 call io_mkdir(
'debug/lcao', namespace)
293 iunit_o =
io_open(
'debug/lcao/orbitals', namespace, action=
'write')
294 write(iunit_o,
'(7a6)')
'iorb',
'atom',
'level',
'i',
'l',
'm',
'spin'
297 if (.not. this%alternative)
then
318 call parse_variable(namespace,
'LCAOMaximumOrbitalRadius', 20.0_real64, max_orb_radius, unit =
units_inp%length)
324 do ia = 1, ions%natoms
325 maxj = max(maxj, ions%atom(ia)%species%get_niwfs())
326 this%maxorbs = this%maxorbs + ions%atom(ia)%species%get_niwfs()
329 this%maxorbs = this%maxorbs*st%d%dim
331 if (this%maxorbs == 0)
then
332 call messages_write(
'The are no atomic orbitals available, cannot do LCAO.')
334 this%mode = option__lcaostart__lcao_none
341 safe_allocate( this%atom(1:this%maxorbs))
342 safe_allocate(this%level(1:this%maxorbs))
343 safe_allocate( this%ddim(1:this%maxorbs))
345 safe_allocate(this%is_empty(1:this%maxorbs))
346 this%is_empty = .false.
349 this%lapdist = maxval(abs(gr%idx%enlarge)*gr%spacing)
352 safe_allocate(this%radius(1:ions%natoms))
354 do ia = 1, ions%natoms
355 norbs = ions%atom(ia)%species%get_niwfs()
359 call ions%atom(ia)%species%get_iwf_ilm(iorb, 1, ii, ll, mm)
361 if(ions%atom(ia)%species%is_full())
call ions%atom(ia)%species%get_iwf_n( iorb, 1, ii)
362 maxradius = max(maxradius, ions%atom(ia)%species%get_iwf_radius( ii, is = 1))
365 maxradius = min(maxradius,
m_two*maxval(gr%box%bounding_box_l(1:space%dim)))
367 this%radius(ia) = maxradius
399 do ia = 1, ions%natoms
401 if (jj > ions%atom(ia)%species%get_niwfs()) cycle
402 call ions%atom(ia)%species%get_iwf_ilm(jj, idim, ii, ll, mm)
404 if(ions%atom(ia)%species%is_full())
then
406 call ions%atom(ia)%species%get_iwf_n( ii_tmp, 1, ii)
408 if (this%orbital_scale_factor*ions%atom(ia)%species%get_iwf_radius( ii, is = 1) >= max_orb_radius) cycle
411 this%level(iorb) = jj
412 this%ddim(iorb) = idim
414 if (
debug%info .and. st%system_grp%is_root())
then
415 write(iunit_o,
'(7i6)') iorb, this%atom(iorb), this%level(iorb), ii, ll, mm, this%ddim(iorb)
423 if (
debug%info .and. st%system_grp%is_root())
then
428 if (this%maxorbs /= iorb - 1)
then
433 call messages_write(
' orbitals cannot be used for the LCAO calculation,')
435 call messages_write(
' their radii exceeds LCAOMaximumOrbitalRadius (')
440 this%maxorbs = iorb - 1
443 if (this%maxorbs < st%nst)
then
444 call messages_write(
'Cannot do LCAO for all states because there are not enough atomic orbitals.')
484 if (n > 0 .and. n <= st%nst .and. st%nst <= this%maxorbs)
then
487 else if (n > st%nst .and. n <= this%maxorbs)
then
490 else if (n == 0)
then
492 this%norbs = min(this%maxorbs, 2*st%nst)
495 this%norbs = this%maxorbs
498 assert(this%norbs <= this%maxorbs)
500 safe_allocate(this%cst(1:this%norbs, 1:st%d%spin_channels))
501 safe_allocate(this%ck(1:this%norbs, 1:st%d%spin_channels))
502 this%initialized_orbitals = .false.
512 integer :: iatom, iorb, norbs
513 real(real64) :: maxradius
516 integer :: jatom, jorb, jbasis, ilbasis, jlbasis, proc(1:2),
info, nbl
557 call parse_variable(namespace,
'LCAOExtraOrbitals', .false., this%derivative)
564 if (this%derivative)
then
567 if (st%nst * st%smear%el_per_state > st%qtot)
then
568 message(1) =
"Lower-lying empty states may be missed with LCAOExtraOrbitals."
581 call parse_variable(namespace,
'LCAODiagTol', 1e-10_real64, this%diag_tol)
583 if (this%derivative)
then
589 safe_allocate(this%sphere(1:ions%natoms))
590 safe_allocate(this%orbitals(1:ions%natoms))
591 safe_allocate(this%is_orbital_initialized(1:ions%natoms))
592 this%is_orbital_initialized = .false.
594 safe_allocate(this%norb_atom(1:ions%natoms))
598 do iatom = 1, ions%natoms
599 this%norb_atom(iatom) = this%mult*ions%atom(iatom)%species%get_niwfs()
600 this%maxorb = max(this%maxorb, ions%atom(iatom)%species%get_niwfs())
601 this%norbs = this%norbs + ions%atom(iatom)%species%get_niwfs()
604 this%maxorb = this%maxorb*this%mult
605 this%norbs = this%norbs*this%mult
607 safe_allocate(this%basis_atom(1:this%norbs))
608 safe_allocate(this%basis_orb(1:this%norbs))
609 safe_allocate(this%atom_orb_basis(1:ions%natoms, 1:this%maxorb))
614 do iatom = 1, ions%natoms
615 norbs = ions%atom(iatom)%species%get_niwfs()
617 do iorb = 1, this%mult*norbs
619 this%atom_orb_basis(iatom, iorb) = ibasis
620 this%basis_atom(ibasis) = iatom
621 this%basis_orb(ibasis) = iorb
624 if (
debug%info .and. st%system_grp%is_root())
then
625 call ions%atom(iatom)%species%get_iwf_ilm(iorb, 1, ii, ll, mm)
626 write(iunit_o,
'(7i6)') ibasis, iatom, iorb, ii, ll, mm, 1
631 if (
debug%info .and. st%system_grp%is_root())
then
636 this%lapdist = maxval(abs(gr%idx%enlarge)*gr%spacing)
639 safe_allocate(this%radius(1:ions%natoms))
641 do iatom = 1, ions%natoms
642 norbs = ions%atom(iatom)%species%get_niwfs()
646 call ions%atom(iatom)%species%get_iwf_ilm(iorb, 1, ii, ll, mm)
648 if(ions%atom(iatom)%species%is_full())
call ions%atom(iatom)%species%get_iwf_n( iorb, 1, ii)
649 maxradius = max(maxradius, ions%atom(iatom)%species%get_iwf_radius( ii, is = 1))
652 if (this%derivative) maxradius = maxradius + this%lapdist
654 maxradius = min(maxradius,
m_two*maxval(gr%box%bounding_box_l(1:space%dim)))
656 this%radius(iatom) = maxradius
659 safe_allocate(this%calc_atom(1:ions%natoms))
660 this%calc_atom = .
true.
663#ifndef HAVE_SCALAPACK
664 this%parallel = .false.
666 this%parallel = (st%parallel_in_states .or. gr%parallel_in_domains) &
669 if (this%parallel)
then
670 nbl = min(16, this%norbs)
673 this%lsize(1) = max(1,
numroc(this%norbs, nbl, st%dom_st_proc_grid%myrow, 0, st%dom_st_proc_grid%nprow))
674 this%lsize(2) = max(1,
numroc(this%norbs, nbl, st%dom_st_proc_grid%mycol, 0, st%dom_st_proc_grid%npcol))
676 this%nproc(1) = st%dom_st_proc_grid%nprow
677 this%nproc(2) = st%dom_st_proc_grid%npcol
678 this%myroc(1) = st%dom_st_proc_grid%myrow
679 this%myroc(2) = st%dom_st_proc_grid%mycol
681 call descinit(this%desc(1), this%norbs, this%norbs, nbl, nbl, 0, 0, &
682 st%dom_st_proc_grid%context, this%lsize(1),
info)
685 write(
message(1),
'(a,i6)')
'descinit for BLACS failed with error code ',
info
689 this%calc_atom = .false.
690 do iatom = 1, ions%natoms
691 ibasis = this%atom_orb_basis(iatom, 1)
693 do jatom = 1, ions%natoms
694 jbasis = this%atom_orb_basis(jatom, 1)
696 do iorb = 1, this%norb_atom(iatom)
697 do jorb = 1, this%norb_atom(jatom)
699 ilbasis, jlbasis, proc(1), proc(2))
701 this%calc_atom(this%basis_atom(jbasis)) = &
702 this%calc_atom(this%basis_atom(jbasis)) .or. proc(2) == this%myroc(2)
720 subroutine lcao_run(namespace, space, gr, ions, ext_partners, st, ks, hm, st_start, lmm_r, known_lower_bound)
723 type(
grid_t),
intent(in) :: gr
724 type(
ions_t),
intent(in) :: ions
727 type(
v_ks_t),
intent(inout) :: ks
729 integer,
optional,
intent(in) :: st_start
730 real(real64),
optional,
intent(in) :: lmm_r
731 logical,
optional,
intent(out) :: known_lower_bound
734 integer :: st_start_random, required_min_nst
736 logical :: is_orbital_dependent
740 if (
present(known_lower_bound)) known_lower_bound = .false.
742 if (
present(st_start))
then
745 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
746 calc_eigenval=.not.
present(st_start), calc_current=.false.)
748 assert(st_start >= 1)
749 if (st_start > st%nst)
then
750 if (
present(known_lower_bound)) known_lower_bound = .
true.
770 if (.not.
present(st_start))
then
771 call lcao_guess_density(lcao, namespace, st, gr, hm, ions, st%qtot, st%d%ispin, st%rho)
774 assert(
present(lmm_r))
780 message(1) =
'Info: Setting up Hamiltonian.'
784 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval=.false., &
785 calc_current=.false., calc_energy=.false., force_semilocal=is_orbital_dependent)
788 st%eigenval = default_eigenval
797 if (
present(st_start))
then
798 write(
message(1),
'(a,i8,a)')
'Performing LCAO for states ', st_start,
' and above'
802 call lcao_wf(lcao, st, gr, ions, hm, namespace, start = st_start)
808 select case (st%d%ispin)
810 required_min_nst = int(st%qtot/2)
812 required_min_nst = int(st%qtot/2)
814 required_min_nst = int(st%qtot)
817 if (lcao%norbs <= required_min_nst .and. lcao%norbs < st%nst)
then
818 st%eigenval(lcao%norbs+1:,:) =
m_zero
822 if (.not.
present(st_start))
then
827 if (lcao%mode == option__lcaostart__lcao_full)
then
828 call v_ks_h_setup(namespace, space, gr, ions, ext_partners, st, ks, hm, &
829 calc_eigenval = .false., calc_current=.false.)
831 assert(
present(lmm_r))
837 if (
present(known_lower_bound) .and. lcao%norbs >= st%nst) known_lower_bound = .
true.
840 if (.not. lcao_done .or. lcao%norbs < st%nst)
then
843 st_start_random = lcao%norbs + 1
847 if (
present(st_start)) st_start_random = max(st_start, st_start_random)
849 if (st_start_random > 1)
then
850 write(
message(1),
'(a,i8,a)')
'Generating random wavefunctions for states ', st_start_random,
' and above'
859 if (.not. st%fixed_spins)
then
867 if (.not. lcao_done)
then
870 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
871 calc_eigenval=.not.
present(st_start), calc_current=.false., force_semilocal=is_orbital_dependent)
872 if (.not.
present(st_start))
then
878 else if (
present(st_start))
then
880 if (st_start > 1)
then
897 type(
lcao_t),
intent(inout) :: this
901 safe_deallocate_a(this%calc_atom)
902 safe_deallocate_a(this%norb_atom)
903 safe_deallocate_a(this%basis_atom)
904 safe_deallocate_a(this%basis_orb)
905 safe_deallocate_a(this%atom_orb_basis)
906 safe_deallocate_a(this%radius)
907 safe_deallocate_a(this%sphere)
908 safe_deallocate_a(this%orbitals)
910 safe_deallocate_a(this%atom)
911 safe_deallocate_a(this%level)
912 safe_deallocate_a(this%ddim)
913 safe_deallocate_a(this%cst)
914 safe_deallocate_a(this%ck)
915 safe_deallocate_a(this%dbuff_single)
916 safe_deallocate_a(this%zbuff_single)
917 safe_deallocate_a(this%dbuff)
918 safe_deallocate_a(this%zbuff)
920 this%initialized = .false.
926 subroutine lcao_wf(this, st, gr, ions, hm, namespace, start)
927 type(
lcao_t),
intent(inout) :: this
929 type(
grid_t),
intent(in) :: gr
930 type(
ions_t),
intent(in) :: ions
933 integer,
optional,
intent(in) :: start
937 assert(this%initialized)
943 if (
present(start)) start_ = start
945 if (this%alternative)
then
947 call dlcao_alt_wf(this, st, gr, ions, hm, namespace, start_)
949 call zlcao_alt_wf(this, st, gr, ions, hm, namespace, start_)
953 call dlcao_wf(this, st, gr, ions, hm, namespace, start_)
955 call zlcao_wf(this, st, gr, ions, hm, namespace, start_)
966 type(
lcao_t),
intent(in) :: this
970 available = this%initialized .and. this%mode /= option__lcaostart__lcao_none &
979 type(
lcao_t),
intent(in) :: this
980 integer,
intent(in) :: ig
981 integer,
intent(in) :: jg
982 integer,
intent(out) :: il
983 integer,
intent(out) :: jl
984 integer,
intent(out) :: prow
985 integer,
intent(out) :: pcol
989 call infog2l(ig, jg, this%desc(1), this%nproc(1), this%nproc(2), this%myroc(1), this%myroc(2), &
1006 type(
lcao_t),
intent(inout) :: this
1007 integer,
intent(in) :: iatom
1011 if (this%is_orbital_initialized(iatom))
then
1012 call this%orbitals(iatom)%end()
1013 this%is_orbital_initialized(iatom) = .false.
1023 type(
lcao_t),
intent(inout) :: this
1025 class(
mesh_t),
intent(in) :: mesh
1027 integer,
intent(in) :: iatom
1028 integer,
intent(in) :: spin_channels
1029 real(real64),
intent(inout) :: rho(:, :)
1031 real(real64),
allocatable :: dorbital(:, :)
1032 complex(real64),
allocatable :: zorbital(:, :)
1033 real(real64),
allocatable :: factors(:)
1034 real(real64) :: factor, aa
1035 integer :: iorb, ip, ii, ll, mm, ispin
1036 type(
ps_t),
pointer :: ps
1037 logical :: use_stored_orbitals
1043 use_stored_orbitals = ions%atom(iatom)%species%is_ps() &
1045 .and. st%d%dim == 1 .and. .not. ions%space%is_periodic()
1049 if (use_stored_orbitals)
then
1051 assert(.not. ions%space%is_periodic())
1053 select type(spec=>ions%atom(iatom)%species)
1060 if (.not. this%alternative)
then
1063 safe_allocate(dorbital(1:mesh%np, 1:st%d%dim))
1065 safe_allocate(zorbital(1:mesh%np, 1:st%d%dim))
1068 do iorb = 1, this%norbs
1069 if (iatom /= this%atom(iorb)) cycle
1071 call ions%atom(iatom)%species%get_iwf_ilm(this%level(iorb), 1, ii, ll, mm)
1072 factor = ps%conf%occ(ii, 1)/(
m_two*ll +
m_one)
1075 call dget_ao(this, st, mesh, ions, iorb, 1, dorbital, use_psi = .
true.)
1078 rho(ip, 1) = rho(ip, 1) + factor*dorbital(ip, 1)**2
1081 call zget_ao(this, st, mesh, ions, iorb, 1, zorbital, use_psi = .
true.)
1084 rho(ip, 1) = rho(ip, 1) + factor*abs(zorbital(ip, 1))**2
1090 safe_deallocate_a(dorbital)
1091 safe_deallocate_a(zorbital)
1099 safe_allocate(factors(1:this%norb_atom(iatom)/this%mult))
1101 do iorb = 1, this%norb_atom(iatom)/this%mult
1102 call ions%atom(iatom)%species%get_iwf_ilm(iorb, 1, ii, ll, mm)
1103 factors(iorb) = ps%conf%occ(ii, 1)/(
m_two*ll +
m_one)
1107 do ip = 1, this%sphere(iatom)%np
1109 do iorb = 1, this%norb_atom(iatom)/this%mult
1110 aa = aa + factors(iorb)*this%orbitals(iatom)%dff_linear(ip, iorb)**2
1112 rho(this%sphere(iatom)%map(ip), 1) = rho(this%sphere(iatom)%map(ip), 1) + aa
1115 safe_deallocate_a(factors)
1121 ions%pos(:, iatom), mesh, spin_channels, rho)
1126 do ispin = 1, spin_channels
1129 rho(ip, ispin) = max(rho(ip, ispin),
m_zero)
1139 type(
lcao_t),
intent(inout) :: this
1142 type(
grid_t),
intent(in) :: gr
1144 type(
ions_t),
intent(in) :: ions
1145 real(real64),
intent(in) :: qtot
1146 integer,
intent(in) :: ispin
1147 real(real64),
contiguous,
intent(out) :: rho(:, :)
1149 integer :: ia, is, idir, ip, m_dim
1152 real(real64) :: rr, rnd, phi, theta, lmag, n1, n2
1153 real(real64),
allocatable :: atom_rho(:,:), mag(:,:)
1154 real(real64),
parameter :: tol_min_mag = 1.0e-20_real64
1158 if (st%d%spin_channels == 1)
then
1159 this%gmd_opt = initrho_paramagnetic
1192 message(1) =
"GuessMagnetDensity cannot be used for Hartree-Fock and generalized Kohn-Sham calculation."
1193 message(2) =
"Please perform a LDA or GGA calculation first and restart from this calculation."
1217 if (
parse_block(namespace,
'AtomsMagnetDirection', blk) < 0)
then
1218 message(1) =
"AtomsMagnetDirection block is not defined."
1223 message(1) =
"The number of rows in the AtomsMagnetDirection block does not equal the number of atoms."
1233 safe_allocate(mag(1:m_dim, 1:ions%natoms))
1234 do ia = 1, ions%natoms
1238 if (abs(mag(idir, ia)) < tol_min_mag) mag(idir, ia) =
m_zero
1246 safe_allocate(atom_rho(1:gr%np, 1:st%d%spin_channels))
1247 select case (this%gmd_opt)
1248 case (initrho_paramagnetic)
1249 do ia = ions%atoms_dist%start, ions%atoms_dist%end
1254 if (st%d%spin_channels == 2)
then
1257 rho(ip, 1) =
m_half*(rho(ip, 1) + rho(ip, 2))
1258 rho(ip, 2) = rho(ip, 1)
1263 do ia = ions%atoms_dist%start, ions%atoms_dist%end
1265 rho(1:gr%np, 1:2) = rho(1:gr%np, 1:2) + atom_rho(1:gr%np, 1:2)
1269 do ia = ions%atoms_dist%start, ions%atoms_dist%end
1284 call lalg_axpy(gr%np, (lmag - n1 + n2)/
m_two/n2, atom_rho(:, 2), atom_rho(:, 1))
1293 rho(1:gr%np, 1:2) = rho(1:gr%np, 1:2) + atom_rho(1:gr%np, 1:2)
1295 rho(1:gr%np, 1) = rho(1:gr%np, 1) + atom_rho(1:gr%np, 2)
1296 rho(1:gr%np, 2) = rho(1:gr%np, 2) + atom_rho(1:gr%np, 1)
1298 elseif (ispin ==
spinors)
then
1309 do ia = ions%atoms_dist%start, ions%atoms_dist%end
1317 lmag = norm2(mag(:, ia))
1318 if (lmag > n1 + n2)
then
1319 mag = mag*(n1 + n2)/lmag
1327 atom_rho(ip, 1) =
m_half*(atom_rho(ip, 1) + atom_rho(ip, 2))
1328 rho(ip, 1) = rho(ip, 1) + atom_rho(ip, 1)
1329 rho(ip, 2) = rho(ip, 2) + atom_rho(ip, 1)
1336 if (n1 - n2 < lmag)
then
1337 call lalg_axpy(gr%np, (lmag - n1 + n2)/
m_two/n2, atom_rho(:, 2), atom_rho(:, 1))
1339 elseif (n1 - n2 > lmag)
then
1340 call lalg_axpy(gr%np, (n1 - n2 - lmag)/
m_two/n1, atom_rho(:, 1), atom_rho(:, 2))
1348 if (mag(1, ia) >
m_zero)
then
1355 elseif (ispin ==
spinors)
then
1365 if (ions%atoms_dist%parallel)
then
1366 do is = 1, st%d%nspin
1367 call lalg_copy(gr%np, rho(:,is), atom_rho(:,1))
1368 call ions%atoms_dist%mpi_grp%allreduce(atom_rho(1, 1), rho(1, is), gr%np, mpi_double_precision, mpi_sum)
1374 write(
message(1),
'(a,f13.6)')
'Info: Unnormalized total charge = ', rr
1379 call lalg_scal(gr%np, st%d%nspin, qtot/rr, rho)
1381 write(
message(1),
'(a,f13.6)')
'Info: Renormalized total charge = ', rr
1386 if (st%symmetrize_density)
then
1387 do is = 1, st%d%nspin
1392 safe_deallocate_a(atom_rho)
1393 safe_deallocate_a(mag)
1399 real(real64) function integrated_charge_density(gr, st, rho) result(rr)
1400 type(
grid_t),
intent(in) :: gr
1402 real(real64),
intent(in) :: rho(:,:)
1407 do is = 1, st%d%spin_channels
1411 call gr%allreduce(rr)
1416 class(mesh_t),
intent(in) :: mesh
1417 real(real64),
intent(inout) :: rho(:,:)
1418 real(real64),
intent(in) :: atom_rho(:,:)
1419 real(real64),
intent(in) :: theta, phi
1422 real(real64) :: half_theta
1426 half_theta = m_half * theta
1430 rho(ip, 1) = rho(ip, 1) +
cos(half_theta)**2*atom_rho(ip, 1) +
sin(half_theta)**2*atom_rho(ip, 2)
1431 rho(ip, 2) = rho(ip, 2) +
sin(half_theta)**2*atom_rho(ip, 1) +
cos(half_theta)**2*atom_rho(ip, 2)
1432 rho(ip, 3) = rho(ip, 3) +
cos(half_theta)*
sin(half_theta)*
cos(phi)*(atom_rho(ip, 1)-atom_rho(ip, 2))
1433 rho(ip, 4) = rho(ip, 4) -
cos(half_theta)*
sin(half_theta)*
sin(phi)*(atom_rho(ip, 1)-atom_rho(ip, 2))
1442 type(
lcao_t),
intent(inout) :: this
1443 type(namespace_t),
intent(in) :: namespace
1444 type(states_elec_t),
intent(inout) :: st
1445 type(grid_t),
intent(in) :: gr
1446 type(ions_t),
intent(in) :: ions
1447 integer,
optional,
intent(in) :: start
1453 if (.not. this%alternative)
then
1454 if (states_are_real(st))
then
1460 if (states_are_real(st))
then
1473 real(real64),
intent(in) :: mag(3)
1474 real(real64),
intent(in) :: lmag
1475 real(real64),
intent(out) :: theta
1476 real(real64),
intent(out) :: phi
1478 assert(lmag > m_zero)
1480 theta =
acos(max(-1.0_real64, min(1.0_real64, mag(3)/lmag)))
1482 if (abs(mag(1)) <= m_epsilon)
then
1483 if (abs(mag(2)) <= m_epsilon)
then
1485 elseif (mag(2) < m_zero)
then
1486 phi = m_pi*m_three*m_half
1492 phi =
atan2(mag(2), mag(1))
1507 type(states_elec_t),
intent(inout) :: st
1508 type(grid_t),
intent(in) :: gr
1509 integer,
intent(in) :: ist_start
1510 integer,
intent(in) :: gmd_opt
1512 integer :: ik, ist, ip
1513 real(real64),
allocatable :: md(:,:), up(:)
1514 complex(real64),
allocatable :: down(:), psi(:,:)
1515 real(real64) :: lmag, theta, phi, norm, mm(3)
1517 if (st%d%ispin /= spinors)
return
1522 safe_allocate(psi(1:gr%np, 1:st%d%dim))
1524 select case (gmd_opt)
1527 do ik = st%d%kpt%start, st%d%kpt%end
1528 do ist = ist_start, st%st_end
1529 call states_elec_get_state(st, gr, ist, ik, psi)
1532 psi(ip, 1) = psi(ip, 1) + psi(ip, 2)
1533 psi(ip, 2) = psi(ip, 1)
1535 call zmf_normalize(gr, st%d%dim, psi)
1536 call states_elec_set_state(st, gr, ist, ik, psi)
1542 do ik = st%d%kpt%start, st%d%kpt%end
1543 do ist = ist_start, st%st_end
1544 call states_elec_get_state(st, gr, ist, ik, psi)
1547 psi(ip, 1) = psi(ip, 1) + psi(ip, 2)
1550 call zmf_normalize(gr, st%d%dim, psi)
1551 call states_elec_set_state(st, gr, ist, ik, psi)
1557 safe_allocate(md(1:gr%np, 1:3))
1558 safe_allocate(up(1:gr%np))
1559 safe_allocate(down(1:gr%np))
1560 call magnetic_density(gr, st%d, st%rho, md)
1562 mm(1) = dmf_integrate(gr, md(:, 1), reduce = .false.)
1563 mm(2) = dmf_integrate(gr, md(:, 2), reduce = .false.)
1564 mm(3) = dmf_integrate(gr, md(:, 3), reduce = .false.)
1566 if (gr%parallel_in_domains)
then
1567 call gr%allreduce(mm)
1573 if (lmag > 1.0e-2_real64)
then
1575 up(:) =
cos(theta/m_two)
1576 down(:) =
exp(m_zi * phi) *
sin(theta/m_two)
1583 lmag = norm2(md(ip, 1:3))
1585 up(ip) =
cos(theta/m_two)
1586 down(ip) =
exp(m_zi * phi) *
sin(theta/m_two)
1588 safe_deallocate_a(md)
1592 do ik = st%d%kpt%start, st%d%kpt%end
1593 do ist = ist_start, st%st_end
1594 call states_elec_get_state(st, gr, ist, ik, psi)
1595 call zmf_normalize(gr, st%d%dim, psi)
1598 norm =
sqrt(real(psi(ip, 1)*conjg(psi(ip, 1)) + psi(ip, 2)*conjg(psi(ip, 2))))
1599 psi(ip, 1) = norm * up(ip)
1600 psi(ip, 2) = norm * down(ip)
1602 call states_elec_set_state(st, gr, ist, ik, psi)
1605 safe_deallocate_a(up)
1606 safe_deallocate_a(down)
1609 safe_deallocate_a(psi)
1616#include "lcao_inc.F90"
1619#include "complex.F90"
1620#include "lcao_inc.F90"
constant times a vector plus a vector
Copies a vector x, to a vector y.
scales a vector by a constant
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
double acos(double __x) __attribute__((__nothrow__
double exp(double __x) __attribute__((__nothrow__
double sin(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
double atan2(double __y, double __x) __attribute__((__nothrow__
This module implements batches of mesh functions.
This module contains interfaces for BLACS routines Interfaces are from http:
This module provides the BLACS processor grid.
logical pure function, public blacs_proc_grid_null(this)
Module implementing boundary conditions in Octopus.
type(debug_t), save, public debug
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
integer, parameter, public spin_polarized
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_pi
some mathematical constants
integer, parameter, public hartree_fock
integer, parameter, public generalized_kohn_sham_dft
integer, parameter, public kohn_sham_dft
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public m_one
integer, parameter, public hartree
This module implements the underlying real-space grid.
subroutine, public dgrid_symmetrize_scalar_field(gr, field, suppress_warning)
This module defines classes and functions for interaction partners.
subroutine, public io_close(iunit, grp)
subroutine, public io_mkdir(fname, namespace, parents)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
A module to handle KS potential, without the external potential.
This module contains interfaces for LAPACK routines.
subroutine, public lcao_init_orbitals(this, namespace, st, gr, ions, start)
subroutine, public lcao_wf(this, st, gr, ions, hm, namespace, start)
subroutine zlcao_alt_init_orbitals(this, namespace, st, gr, ions, start)
subroutine zget_ao(this, st, mesh, ions, iorb, ispin, ao, use_psi)
subroutine, public lcao_end(this)
subroutine zlcao_alt_wf(this, st, gr, ions, hm, namespace, start)
The alternative implementation.
subroutine dlcao_wf(this, st, gr, ions, hm, namespace, start)
subroutine dlcao_alt_init_orbitals(this, namespace, st, gr, ions, start)
subroutine zlcao_wf(this, st, gr, ions, hm, namespace, start)
subroutine rotate_random_states_to_local_frame(st, gr, ist_start, gmd_opt)
Rotate the spinors with band index >= ist_start to the local frame of the magnetization.
subroutine, public lcao_run(namespace, space, gr, ions, ext_partners, st, ks, hm, st_start, lmm_r, known_lower_bound)
subroutine dlcao_alt_wf(this, st, gr, ions, hm, namespace, start)
The alternative implementation.
integer, parameter initrho_userdef
subroutine get_angles_from_magnetization(mag, lmag, theta, phi)
Given a magnetization vector, and its norm, this returns the corresponding inclination and azimuthal ...
integer, parameter initrho_random
subroutine lcao_alt_end_orbital(this, iatom)
This function deallocates a set of an atomic orbitals for an atom. It can be called when the batch is...
subroutine lcao_local_index(this, ig, jg, il, jl, prow, pcol)
subroutine lcao_guess_density(this, namespace, st, gr, hm, ions, qtot, ispin, rho)
builds a density which is the sum of the atomic densities
subroutine dinit_orbitals(this, namespace, st, gr, ions, start)
real(real64) function integrated_charge_density(gr, st, rho)
Computes the integral of rho, summed over spin channels.
integer, parameter initrho_paramagnetic
subroutine dget_ao(this, st, mesh, ions, iorb, ispin, ao, use_psi)
subroutine dlcao_alt_get_orbital(this, sphere, ions, ispin, iatom, norbs)
This function generates the set of an atomic orbitals for an atom and stores it in the batch orbitalb...
subroutine zinit_orbitals(this, namespace, st, gr, ions, start)
subroutine lcao_atom_density(this, st, mesh, ions, iatom, spin_channels, rho)
subroutine, public lcao_init(this, namespace, space, gr, ions, st, st_start)
subroutine accumulate_rotated_density(mesh, rho, atom_rho, theta, phi)
integer, parameter initrho_ferromagnetic
logical function, public lcao_is_available(this)
Returns true if LCAO can be done.
System information (time, memory, sysname)
subroutine, public compute_and_write_magnetic_moments(gr, st, phase, ep, ions, lmm_r, calc_orb_moments, iunit, namespace)
Computes and prints the global and local magnetic moments.
This module is intended to contain "only mathematical" functions and procedures.
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public messages_obsolete_variable(namespace, name, rep)
subroutine, public messages_new_line()
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_input_error(namespace, var, details, row, column)
subroutine, public messages_experimental(name, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
logical function, public parse_is_defined(namespace, name)
integer function, public parse_block(namespace, name, blk, check_varinfo_)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
integer(int64), parameter, public splitmix64_321
This module contains interfaces for ScaLAPACK routines Interfaces are from http:
integer, parameter, public smear_semiconductor
integer, parameter, public smear_fixed_occ
subroutine, public species_atom_density(species, namespace, space, latt, pos, mesh, spin_channels, rho)
pure logical function, public states_are_complex(st)
pure logical function, public states_are_real(st)
subroutine, public states_elec_orthogonalize(st, namespace, mesh)
Orthonormalizes nst orbitals in mesh (honours state parallelization).
This module defines routines to write information about states.
subroutine, public states_elec_write_eigenvalues(nst, st, space, kpoints, error, st_start, compact, iunit, namespace)
write the eigenvalues for some states to a file.
subroutine, public states_elec_fermi(st, namespace, mesh, compute_spin)
calculate the Fermi level for the states in this object
subroutine, public states_elec_generate_random(st, mesh, kpoints, ist_start_, ist_end_, ikpt_start_, ikpt_end_, normalized)
randomize states
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
type(unit_system_t), public units_inp
the units systems for reading and writing
subroutine, public v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval, time, calc_energy, calc_current, force_semilocal)
subroutine, public v_ks_h_setup(namespace, space, gr, ions, ext_partners, st, ks, hm, calc_eigenval, calc_current)
logical pure function, public xc_is_orbital_dependent(xcs)
Is the xc family orbital dependent.
integer, parameter, public sic_pz_oep
Perdew-Zunger SIC (OEP way)
Extension of space that contains the knowledge of the spin dimension.
Description of the grid, containing information on derivatives, stencil, and symmetries.
Describes mesh distribution to nodes.
A type storing the information and data about a pseudopotential.
The states_elec_t class contains all electronic wave functions.