Octopus
lcao.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
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 lcao_oct_m
22 use atom_oct_m
24 use batch_oct_m
25 use blacs_oct_m
28 use comm_oct_m
29 use debug_oct_m
31 use global_oct_m
32 use grid_oct_m
35 use io_oct_m
37 use ions_oct_m
38 use, intrinsic :: iso_fortran_env
41 use lapack_oct_m
42 use loct_oct_m
44 use math_oct_m
45 use mesh_oct_m
48 use mpi_oct_m
51 use parser_oct_m
53 use ps_oct_m
57 use smear_oct_m
58 use space_oct_m
67 use unit_oct_m
69 use v_ks_oct_m
72 use xc_oct_m
73 use xc_sic_oct_m
74
75 implicit none
76
77 private
78 public :: &
79 lcao_t, &
80 lcao_init, &
82 lcao_wf, &
83 lcao_run, &
84 lcao_end, &
87
88 type lcao_t
89 private
90 integer :: mode
91 logical :: complex_ylms
92 logical :: initialized
93 integer :: norbs
94 integer :: maxorbs
95 integer, allocatable :: atom(:)
96 integer, allocatable :: level(:)
97 integer, allocatable :: ddim(:)
98 logical :: alternative
99 logical :: derivative
100 integer, allocatable :: cst(:, :)
101 integer, allocatable :: ck(:, :)
102 real(4), allocatable :: dbuff_single(:, :, :, :)
103 complex(4), allocatable :: zbuff_single(:, :, :, :)
104 real(real64), allocatable :: dbuff(:, :, :, :)
105 complex(real64), allocatable :: zbuff(:, :, :, :)
106 logical :: save_memory
107 logical :: initialized_orbitals
108 real(real64) :: orbital_scale_factor
109
111 logical, allocatable :: is_empty(:)
112 ! !< This occurs in domain parallelization if the atom is not
113 ! !< in the local domain
116 logical :: keep_orb
117 real(real64), allocatable :: radius(:)
118 real(real64) :: lapdist
119 integer :: mult
120 integer :: maxorb
122 integer, allocatable :: basis_atom(:)
123 integer, allocatable :: basis_orb(:)
124 integer, allocatable :: atom_orb_basis(:, :)
125 integer, allocatable :: norb_atom(:)
126 logical :: parallel
127 integer :: lsize(1:2)
128 integer :: nproc(1:2)
129 integer :: myroc(1:2)
130 integer :: desc(1:BLACS_DLEN)
131 logical, allocatable :: calc_atom(:)
132 real(real64) :: diag_tol
133 type(submesh_t), allocatable :: sphere(:)
134 type(batch_t), allocatable :: orbitals(:)
135 logical, allocatable :: is_orbital_initialized(:)
136 end type lcao_t
137
138
139 integer, parameter :: &
140 INITRHO_PARAMAGNETIC = 1, &
142 initrho_random = 3, &
143 initrho_userdef = 77
144
145contains
146
147 ! ---------------------------------------------------------
148 subroutine lcao_init(this, namespace, space, gr, ions, st)
149 type(lcao_t), intent(out) :: this
150 type(namespace_t), intent(in) :: namespace
151 type(electron_space_t), intent(in) :: space
152 type(grid_t), intent(in) :: gr
153 type(ions_t), intent(in) :: ions
154 type(states_elec_t), intent(in) :: st
155
156 integer :: ia, n, iorb, jj, maxj, idim
157 integer :: ii, ll, mm, norbs, ii_tmp
158 integer :: mode_default
159 real(real64) :: max_orb_radius, maxradius
160 integer :: iunit_o
161
162 push_sub(lcao_init)
163
164 this%initialized = .true.
165
166 ! The initial LCAO calculation is done by default if we have species representing atoms.
167 ! Otherwise, it is not the default value and has to be enforced in the input file.
168 mode_default = option__lcaostart__lcao_states
169 if (ions%only_user_def) mode_default = option__lcaostart__lcao_none
170
171 !%Variable LCAOStart
172 !%Type integer
173 !%Section SCF::LCAO
174 !%Description
175 !% Before starting a SCF calculation, <tt>Octopus</tt> can perform
176 !% a linear combination of atomic orbitals (LCAO) calculation.
177 !% These can provide <tt>Octopus</tt> with a good set
178 !% of initial wavefunctions and with a new guess for the density.
179 !% (Up to the current version, only a minimal basis set is used.)
180 !% The default is <tt>lcao_states</tt> if at least one species representing an atom is present.
181 !% The default is <tt>lcao_none</tt> if all species are <tt>species_user_defined</tt>,
182 !% <tt>species_charge_density</tt>, <tt>species_from_file</tt>, or <tt>species_jellium_slab</tt>.
183 !%
184 !% The initial guess densities for LCAO are taken from the atomic orbitals for pseudopotential species;
185 !% from the natural charge density for <tt>species_charge_density</tt>, <tt>species_point</tt>,
186 !% <tt>species_jellium</tt>, and <tt>species_jellium_slab</tt>;
187 !% or uniform for <tt>species_full_delta</tt>, <tt>species_full_gaussian</tt>,
188 !% <tt>species_user_defined</tt>, or <tt>species_from_file</tt>.
189 !% Pseudopotential species use the pseudo-wavefunctions as orbitals, full-potential atomic species
190 !% (<tt>species_full_delta</tt> and <tt>species_full_gaussian</tt>) use hydrogenic wavefunctions, and
191 !% others use harmonic-oscillator wavefunctions.
192 !%
193 !% Note: Some pseudopotential files (CPI, FHI for example) do not
194 !% contain full information about the orbitals. In this case,
195 !% Octopus generates the starting density from the normalized
196 !% square root of the local potential. If no orbitals are
197 !% available at all from the pseudopotential files, Octopus will
198 !% not be able to perform an LCAO and the initial states will be
199 !% randomized.
200 !%
201 !%Option lcao_none 0
202 !% Do not perform a LCAO calculation before the SCF cycle. Instead use random wavefunctions.
203 !%Option lcao_states 2
204 !% Do a LCAO calculation before the SCF cycle and use the resulting wavefunctions as
205 !% initial wavefunctions without changing the guess density.
206 !% This will speed up the convergence of the eigensolver during the first SCF iterations.
207 !%Option lcao_full 3
208 !% Do a LCAO calculation before the SCF cycle and use the LCAO wavefunctions to build a new
209 !% guess density and a new KS potential.
210 !% Using the LCAO density as a new guess density may improve the convergence, but can
211 !% also slow it down or yield wrong results (especially for spin-polarized calculations).
212 !%Option lcao_states_batch 4
213 !% (Experimental) Same as lcao_states, but faster for GPUs, large systems and parallel in states.
214 !% It is not working for spinors, however.
215 !%End
216 call parse_variable(namespace, 'LCAOStart', mode_default, this%mode)
217 if (.not. varinfo_valid_option('LCAOStart', this%mode)) call messages_input_error(namespace, 'LCAOStart')
219 call messages_print_var_option('LCAOStart', this%mode, namespace=namespace)
221 !TODO: Remove the alternative flag
222 this%alternative = this%mode == option__lcaostart__lcao_states_batch
224 if (this%mode == option__lcaostart__lcao_none) then
225 pop_sub(lcao_init)
226 return
227 end if
229 call messages_obsolete_variable(namespace, 'LCAOAlternative')
230
231
232 ! DAS: For spinors, you will always get magnetization in (1, 0, 0) direction, and the
233 ! eigenvalues will be incorrect. This is due to confusion between spins and spinors in the code.
234 if (st%d%ispin == spinors .and. this%alternative) then
235 message(1) = "LCAOStart = lcao_states_batch is not working for spinors."
236 call messages_fatal(1, namespace=namespace)
237 end if
238 if (space%is_periodic() .and. this%alternative) then
239 call messages_experimental("LCAOStart = lcao_states_batch in periodic systems", namespace=namespace)
240 ! specifically, if you get the message about submesh radius > box size, results will probably be totally wrong.
241 end if
242
243 !%Variable LCAOComplexYlms
244 !%Type logical
245 !%Default false
246 !%Section SCF::LCAO
247 !%Description
248 !% If set to true, and using complex states, complex spherical harmonics will be used, <i>i.e.</i>
249 !% with <math>e^{\pm i m \phi}</math>.
250 !% If false, real spherical harmonics with <math>\sin(m \phi)</math> or <math>\cos(m \phi)</math> are used.
251 !% This variable will make it more likely to get states that are eigenvectors of the <math>L_z</math>
252 !% operator, with a definite angular momentum.
253 !%End
254
255 if (states_are_complex(st)) then
256 call parse_variable(namespace, 'LCAOComplexYlms', .false., this%complex_ylms)
257 else
258 this%complex_ylms = .false.
259 end if
260
261 !%Variable LCAOSaveMemory
262 !%Type logical
263 !%Default false
264 !%Section SCF::LCAO
265 !%Description
266 !% If set to true, the LCAO will allocate extra memory needed in single precision instead of
267 !% double precision.
268 !%End
269 call parse_variable(namespace, 'LCAOSaveMemory', .false., this%save_memory)
270
271
272 if (debug%info .and. mpi_grp_is_root(mpi_world)) then
273 call io_mkdir('debug/lcao', namespace)
274 iunit_o = io_open('debug/lcao/orbitals', namespace, action='write')
275 write(iunit_o,'(7a6)') 'iorb', 'atom', 'level', 'i', 'l', 'm', 'spin'
276 end if
277
278 if (.not. this%alternative) then
279
280 !%Variable LCAOScaleFactor
281 !%Type float
282 !%Default 1.0
283 !%Section SCF::LCAO
284 !%Description
285 !% The coordinates of the atomic orbitals used by the LCAO
286 !% procedure will be rescaled by the value of this variable. 1.0 means no rescaling.
287 !%End
288 call parse_variable(namespace, 'LCAOScaleFactor', m_one, this%orbital_scale_factor)
289 call messages_print_var_value('LCAOScaleFactor', this%orbital_scale_factor, namespace=namespace)
290
291 !%Variable LCAOMaximumOrbitalRadius
292 !%Type float
293 !%Default 20.0 a.u.
294 !%Section SCF::LCAO
295 !%Description
296 !% The LCAO procedure will ignore orbitals that have an
297 !% extent greater that this value.
298 !%End
299 call parse_variable(namespace, 'LCAOMaximumOrbitalRadius', 20.0_real64, max_orb_radius, unit = units_inp%length)
300 call messages_print_var_value('LCAOMaximumOrbitalRadius', max_orb_radius, units_out%length, namespace=namespace)
301
302 ! count the number of orbitals available
303 maxj = 0
304 this%maxorbs = 0
305 do ia = 1, ions%natoms
306 maxj = max(maxj, ions%atom(ia)%species%get_niwfs())
307 this%maxorbs = this%maxorbs + ions%atom(ia)%species%get_niwfs()
308 end do
309
310 this%maxorbs = this%maxorbs*st%d%dim
311
312 if (this%maxorbs == 0) then
313 call messages_write('The are no atomic orbitals available, cannot do LCAO.')
314 call messages_warning(namespace=namespace)
315 this%mode = option__lcaostart__lcao_none
316 pop_sub(lcao_init)
317 return
318 end if
319
320 ! generate tables to know which indices each atomic orbital has
321
322 safe_allocate( this%atom(1:this%maxorbs))
323 safe_allocate(this%level(1:this%maxorbs))
324 safe_allocate( this%ddim(1:this%maxorbs))
325
326 safe_allocate(this%is_empty(1:this%maxorbs))
327 this%is_empty = .false.
328
329 ! this is determined by the stencil we are using and the spacing
330 this%lapdist = maxval(abs(gr%idx%enlarge)*gr%spacing)
331
332 ! calculate the radius of each orbital
333 safe_allocate(this%radius(1:ions%natoms))
334
335 do ia = 1, ions%natoms
336 norbs = ions%atom(ia)%species%get_niwfs()
337
338 maxradius = m_zero
339 do iorb = 1, norbs
340 call ions%atom(ia)%species%get_iwf_ilm(iorb, 1, ii, ll, mm)
341 ! For all-electron species, we need to use the principal quantum number
342 if(ions%atom(ia)%species%is_full()) call ions%atom(ia)%species%get_iwf_n( iorb, 1, ii)
343 maxradius = max(maxradius, ions%atom(ia)%species%get_iwf_radius( ii, is = 1))
344 end do
345
346 maxradius = min(maxradius, m_two*maxval(gr%box%bounding_box_l(1:space%dim)))
347
348 this%radius(ia) = maxradius
349 end do
350
351
352 ! Each atom provides niwfs pseudo-orbitals (this number is given in
353 ! ions%atom(ia)%species%niwfs for atom number ia). This number is
354 ! actually multiplied by two in case of spin-unrestricted or spinors
355 ! calculations.
356 !
357 ! The pseudo-orbitals are placed in order in the following way (Natoms
358 ! is the total number of atoms).
359 !
360 ! n = 1 => first orbital of atom 1,
361 ! n = 2 => first orbital of atom 2.
362 ! n = 3 => first orbital of atom 3.
363 ! ....
364 ! n = Natoms => first orbital of atom Natoms
365 ! n = Natoms + 1 = > second orbital of atom 1
366 ! ....
367 !
368 ! If at some point in this loop an atom pseudo cannot provide the corresponding
369 ! orbital (because the niws orbitals have been exhausted), it moves on to the following
370 ! atom.
371 !
372 ! In the spinors case, it changes a bit:
373 !
374 ! n = 1 => first spin-up orbital of atom 1, assigned to the spin-up component of the spinor.
375 ! n = 2 => first spin-down orbital of atom 1, assigned to the spin-down component of the spinor.
376 ! n = 3 => first spin-up orbital of atom 2, assigned to the spin-up component of the spinor.
377
378 iorb = 1
379 do jj = 1, maxj
380 do ia = 1, ions%natoms
381 do idim = 1,st%d%dim
382 if (jj > ions%atom(ia)%species%get_niwfs()) cycle
383 call ions%atom(ia)%species%get_iwf_ilm(jj, idim, ii, ll, mm)
384 ! For all-electron species, we need to use the principal quantum number
385 if(ions%atom(ia)%species%is_full()) then
386 ii_tmp = ii
387 call ions%atom(ia)%species%get_iwf_n( ii_tmp, 1, ii)
388 end if
389 if (this%orbital_scale_factor*ions%atom(ia)%species%get_iwf_radius( ii, is = 1) >= max_orb_radius) cycle
390
391 this%atom(iorb) = ia
392 this%level(iorb) = jj
393 this%ddim(iorb) = idim
394
395 if (debug%info .and. mpi_grp_is_root(mpi_world)) then
396 write(iunit_o,'(7i6)') iorb, this%atom(iorb), this%level(iorb), ii, ll, mm, this%ddim(iorb)
397 end if
398
399 iorb = iorb + 1
400 end do
401 end do
402 end do
403
404 if (debug%info .and. mpi_grp_is_root(mpi_world)) then
405 call io_close(iunit_o)
406 end if
407
408 ! some orbitals might have been removed because of their radii
409 if (this%maxorbs /= iorb - 1) then
410 call messages_write('Info: ')
411 call messages_write(this%maxorbs - iorb + 1)
412 call messages_write(' of ')
413 call messages_write(this%maxorbs)
414 call messages_write(' orbitals cannot be used for the LCAO calculation,')
415 call messages_new_line()
416 call messages_write(' their radii exceeds LCAOMaximumOrbitalRadius (')
417 call messages_write(max_orb_radius, units = units_out%length, fmt = '(f6.1)')
418 call messages_write(').')
419 call messages_warning(namespace=namespace)
420
421 this%maxorbs = iorb - 1
422 end if
423
424 if (this%maxorbs < st%nst) then
425 call messages_write('Cannot do LCAO for all states because there are not enough atomic orbitals.')
426 call messages_new_line()
427
428 call messages_write('Required: ')
429 call messages_write(st%nst)
430 call messages_write('. Available: ')
431 call messages_write(this%maxorbs)
432 call messages_write('. ')
433 call messages_write(st%nst - this%maxorbs)
434 call messages_write(' orbitals will be randomized.')
435 call messages_warning(namespace=namespace)
436 end if
437
438 !%Variable LCAODimension
439 !%Type integer
440 !%Section SCF::LCAO
441 !%Description
442 !% (Only applies if <tt>LCAOStart /= lcao_states_batch</tt>.)
443 !% Before starting the SCF cycle, an initial LCAO calculation can be performed
444 !% in order to obtain reasonable initial guesses for spin-orbitals and densities.
445 !% For this purpose, the code calculates a number of atomic orbitals.
446 !% The number available for a species described by a pseudopotential is all the
447 !% orbitals up the maximum angular momentum in the pseudopotential, minus any orbitals that
448 !% are found to be unbound. For non-pseudopotential species, the number is equal to
449 !% twice the valence charge.
450 !% The default dimension for the LCAO basis
451 !% set will be the sum of all these numbers, or twice the number of required orbitals
452 !% for the full calculation, whichever is less.
453 !%
454 !% This dimension however can be changed by making use of this
455 !% variable. Note that <tt>LCAODimension</tt> cannot be smaller than the
456 !% number of orbitals needed in the full calculation -- if
457 !% <tt>LCAODimension</tt> is smaller, it will be silently increased to meet
458 !% this requirement. In the same way, if <tt>LCAODimension</tt> is larger
459 !% than the available number of atomic orbitals, it will be
460 !% reduced. If you want to use the largest possible number, set
461 !% <tt>LCAODimension</tt> to a negative number.
462 !%End
463 call parse_variable(namespace, 'LCAODimension', 0, n)
464
465 if (n > 0 .and. n <= st%nst .and. st%nst <= this%maxorbs) then
466 ! n was made too small
467 this%norbs = st%nst
468 else if (n > st%nst .and. n <= this%maxorbs) then
469 ! n is a reasonable value
470 this%norbs = n
471 else if (n == 0) then
472 ! using the default
473 this%norbs = min(this%maxorbs, 2*st%nst)
474 else
475 ! n was negative, or greater than maxorbs
476 this%norbs = this%maxorbs
477 end if
478
479 assert(this%norbs <= this%maxorbs)
480
481 safe_allocate(this%cst(1:this%norbs, 1:st%d%spin_channels))
482 safe_allocate(this%ck(1:this%norbs, 1:st%d%spin_channels))
483 this%initialized_orbitals = .false.
484 else
485 call lcao2_init()
486 end if
487
488 pop_sub(lcao_init)
489
490 contains
491
492 subroutine lcao2_init()
493 integer :: iatom, iorb, norbs
494 real(real64) :: maxradius
495 integer :: ibasis
496#ifdef HAVE_SCALAPACK
497 integer :: jatom, jorb, jbasis, ilbasis, jlbasis, proc(1:2), info, nbl
498#endif
499 push_sub(lcao_init.lcao2_init)
500
501 call messages_write('Info: Using LCAO alternative implementation.')
502 call messages_info(namespace=namespace)
503
504 call messages_experimental('LCAO alternative implementation', namespace=namespace)
505
506 !%Variable LCAOKeepOrbitals
507 !%Type logical
508 !%Default yes
509 !%Section SCF::LCAO
510 !%Description
511 !% Only applies if <tt>LCAOStart = lcao_states_batch</tt>.
512 !% If set to yes (the default) Octopus keeps atomic orbitals in
513 !% memory during the LCAO procedure. If set to no, the orbitals
514 !% are generated each time that they are needed, increasing
515 !% computational time but saving memory.
516 !%
517 !% When set to yes, Octopus prints the amount of memory per node
518 !% that is required to store the orbitals.
519 !%
520 !%End
521 call parse_variable(namespace, 'LCAOKeepOrbitals', .true., this%keep_orb)
522
523 !%Variable LCAOExtraOrbitals
524 !%Type logical
525 !%Default false
526 !%Section SCF::LCAO
527 !%Description
528 !% Only applies if <tt>LCAOStart = lcao_states_batch</tt>, and all species are pseudopotentials.
529 !% (experimental) If this variable is set to yes, the LCAO
530 !% procedure will add an extra set of numerical orbitals (by
531 !% using the derivative of the radial part of the original
532 !% orbitals). Note that this corresponds roughly to adding orbitals
533 !% with higher principal quantum numbers, but the same angular momentum.
534 !% This option may cause problems for unoccupied states since you may miss
535 !% some lower-lying states which correspond to higher angular momenta instead
536 !% of higher principal quantum number.
537 !%End
538 call parse_variable(namespace, 'LCAOExtraOrbitals', .false., this%derivative)
539
540 ! DAS: if you calculate the Na atom this way, spin-polarized, with just one unoccupied state,
541 ! you will obtain states (up and down) which are actually the 10th states if you start with
542 ! random wavefunctions! We really need to implement taking the derivative of the angular part
543 ! instead to be sure of getting decent results!
544
545 if (this%derivative) then
546 call messages_experimental('LCAO extra orbitals', namespace=namespace)
547
548 if (st%nst * st%smear%el_per_state > st%qtot) then
549 message(1) = "Lower-lying empty states may be missed with LCAOExtraOrbitals."
550 call messages_warning(1, namespace=namespace)
551 end if
552 end if
553
554 !%Variable LCAODiagTol
555 !%Type float
556 !%Default 1e-10
557 !%Section SCF::LCAO
558 !%Description
559 !% Only applies if <tt>LCAOStart = lcao_states_batch</tt>.
560 !% The tolerance for the diagonalization of the LCAO Hamiltonian.
561 !%End
562 call parse_variable(namespace, 'LCAODiagTol', 1e-10_real64, this%diag_tol)
563
564 if (this%derivative) then
565 this%mult = 2
566 else
567 this%mult = 1
568 end if
569
570 safe_allocate(this%sphere(1:ions%natoms))
571 safe_allocate(this%orbitals(1:ions%natoms))
572 safe_allocate(this%is_orbital_initialized(1:ions%natoms))
573 this%is_orbital_initialized = .false.
574
575 safe_allocate(this%norb_atom(1:ions%natoms))
576
577 this%maxorb = 0
578 this%norbs = 0
579 do iatom = 1, ions%natoms
580 this%norb_atom(iatom) = this%mult*ions%atom(iatom)%species%get_niwfs()
581 this%maxorb = max(this%maxorb, ions%atom(iatom)%species%get_niwfs())
582 this%norbs = this%norbs + ions%atom(iatom)%species%get_niwfs()
583 end do
584
585 this%maxorb = this%maxorb*this%mult
586 this%norbs = this%norbs*this%mult
587
588 safe_allocate(this%basis_atom(1:this%norbs))
589 safe_allocate(this%basis_orb(1:this%norbs))
590 safe_allocate(this%atom_orb_basis(1:ions%natoms, 1:this%maxorb))
591
592 ! Initialize the mapping between indices
593
594 ibasis = 0
595 do iatom = 1, ions%natoms
596 norbs = ions%atom(iatom)%species%get_niwfs()
597
598 do iorb = 1, this%mult*norbs
599 ibasis = ibasis + 1
600 this%atom_orb_basis(iatom, iorb) = ibasis
601 this%basis_atom(ibasis) = iatom
602 this%basis_orb(ibasis) = iorb
603
604 ! no stored spin index in alternative mode
605 if (debug%info .and. mpi_grp_is_root(mpi_world)) then
606 call ions%atom(iatom)%species%get_iwf_ilm(iorb, 1, ii, ll, mm)
607 write(iunit_o,'(7i6)') ibasis, iatom, iorb, ii, ll, mm, 1
608 end if
609 end do
610 end do
611
612 if (debug%info .and. mpi_grp_is_root(mpi_world)) then
613 call io_close(iunit_o)
614 end if
615
616 ! this is determined by the stencil we are using and the spacing
617 this%lapdist = maxval(abs(gr%idx%enlarge)*gr%spacing)
618
619 ! calculate the radius of each orbital
620 safe_allocate(this%radius(1:ions%natoms))
621
622 do iatom = 1, ions%natoms
623 norbs = ions%atom(iatom)%species%get_niwfs()
624
625 maxradius = m_zero
626 do iorb = 1, norbs
627 call ions%atom(iatom)%species%get_iwf_ilm(iorb, 1, ii, ll, mm)
628 ! For all-electron species, we need to use the principal quantum number
629 if(ions%atom(iatom)%species%is_full()) call ions%atom(iatom)%species%get_iwf_n( iorb, 1, ii)
630 maxradius = max(maxradius, ions%atom(iatom)%species%get_iwf_radius( ii, is = 1))
631 end do
632
633 if (this%derivative) maxradius = maxradius + this%lapdist
634
635 maxradius = min(maxradius, m_two*maxval(gr%box%bounding_box_l(1:space%dim)))
636
637 this%radius(iatom) = maxradius
638 end do
639
640 safe_allocate(this%calc_atom(1:ions%natoms))
641 this%calc_atom = .true.
642
643 ! initialize parallel data
644#ifndef HAVE_SCALAPACK
645 this%parallel = .false.
646#else
647 this%parallel = (st%parallel_in_states .or. gr%parallel_in_domains) &
648 .and. .not. blacs_proc_grid_null(st%dom_st_proc_grid)
649
650 if (this%parallel) then
651 nbl = min(16, this%norbs)
652
653 ! The size of the distributed matrix in each node
654 this%lsize(1) = max(1, numroc(this%norbs, nbl, st%dom_st_proc_grid%myrow, 0, st%dom_st_proc_grid%nprow))
655 this%lsize(2) = max(1, numroc(this%norbs, nbl, st%dom_st_proc_grid%mycol, 0, st%dom_st_proc_grid%npcol))
656
657 this%nproc(1) = st%dom_st_proc_grid%nprow
658 this%nproc(2) = st%dom_st_proc_grid%npcol
659 this%myroc(1) = st%dom_st_proc_grid%myrow
660 this%myroc(2) = st%dom_st_proc_grid%mycol
661
662 call descinit(this%desc(1), this%norbs, this%norbs, nbl, nbl, 0, 0, &
663 st%dom_st_proc_grid%context, this%lsize(1), info)
664
665 if (info /= 0) then
666 write(message(1), '(a,i6)') 'descinit for BLACS failed with error code ', info
667 call messages_fatal(1, namespace=namespace)
668 end if
669
670 this%calc_atom = .false.
671 do iatom = 1, ions%natoms
672 ibasis = this%atom_orb_basis(iatom, 1)
673
674 do jatom = 1, ions%natoms
675 jbasis = this%atom_orb_basis(jatom, 1)
676
677 do iorb = 1, this%norb_atom(iatom)
678 do jorb = 1, this%norb_atom(jatom)
679 call lcao_local_index(this, ibasis - 1 + iorb, jbasis - 1 + jorb, &
680 ilbasis, jlbasis, proc(1), proc(2))
681
682 this%calc_atom(this%basis_atom(jbasis)) = &
683 this%calc_atom(this%basis_atom(jbasis)) .or. proc(2) == this%myroc(2)
684
685 end do
686 end do
687
688 end do
689 end do
690
691 end if
692#endif
693
694 pop_sub(lcao_init.lcao2_init)
695 end subroutine lcao2_init
696
697 end subroutine lcao_init
698
699
700 ! ---------------------------------------------------------
701 subroutine lcao_run(namespace, space, gr, ions, ext_partners, st, ks, hm, st_start, lmm_r)
702 type(namespace_t), intent(in) :: namespace
703 type(electron_space_t), intent(in) :: space
704 type(grid_t), intent(in) :: gr
705 type(ions_t), intent(in) :: ions
706 type(partner_list_t), intent(in) :: ext_partners
707 type(states_elec_t), intent(inout) :: st
708 type(v_ks_t), intent(inout) :: ks
709 type(hamiltonian_elec_t), intent(inout) :: hm
710 integer, optional, intent(in) :: st_start
711 real(real64), optional, intent(in) :: lmm_r
712
713 type(lcao_t) :: lcao
714 integer :: st_start_random, required_min_nst
715 logical :: lcao_done
716 logical :: is_orbital_dependent
717
718 push_sub(lcao_run)
719
720 if (present(st_start)) then
721 ! If we are doing unocc calculation, do not mess with the correct eigenvalues
722 ! of the occupied states.
723 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
724 calc_eigenval=.not. present(st_start), calc_current=.false.)
725
726 assert(st_start >= 1)
727 if (st_start > st%nst) then ! nothing to be done in LCAO
728 pop_sub(lcao_run)
729 return
730 end if
731 end if
732
733 call profiling_in('LCAO_RUN')
734
735 call lcao_init(lcao, namespace, space, gr, ions, st)
736
737 call lcao_init_orbitals(lcao, namespace, st, gr, ions, start = st_start)
738
739 ! By default, we want to use vxc for the LCAO.
740 ! However, we can only do this if vxc depends on the density only
741 ! For cases like MGGA or hybrids, OEP, we cannot do this
742 is_orbital_dependent = (ks%theory_level == hartree .or. ks%theory_level == hartree_fock &
743 .or. (ks%theory_level == kohn_sham_dft .and. xc_is_orbital_dependent(ks%xc)) &
744 .or. (ks%theory_level == generalized_kohn_sham_dft .and. xc_is_orbital_dependent(ks%xc)) &
745 .or. ks%sic%level == sic_pz_oep)
746
747 if (.not. present(st_start)) then
748 call lcao_guess_density(lcao, namespace, st, gr, hm, ions, st%qtot, st%d%ispin, st%d%spin_channels, st%rho)
749
750 if (st%d%ispin > unpolarized) then
751 assert(present(lmm_r))
752 call write_magnetic_moments(gr, st, ions, gr%der%boundaries, lmm_r, namespace=namespace)
753 end if
754
755 ! set up Hamiltonian (we do not call v_ks_h_setup here because we do not want to
756 ! overwrite the guess density)
757 message(1) = 'Info: Setting up Hamiltonian.'
758 call messages_info(1, namespace=namespace)
759
760 ! get the effective potential (we don`t need the eigenvalues yet)
761 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval=.false., &
762 calc_current=.false., calc_energy=.false., force_semilocal=is_orbital_dependent)
763 ! eigenvalues have nevertheless to be initialized to something
764 ! This value must be larger that the highest eigenvalue from LCAO to get the correct occupations
765 st%eigenval = 1e10_real64
766 end if
767
768 lcao_done = .false.
769
770 ! after initialized, can check that LCAO is possible
771 if (lcao_is_available(lcao)) then
772 lcao_done = .true.
773
774 if (present(st_start)) then
775 write(message(1),'(a,i8,a)') 'Performing LCAO for states ', st_start, ' and above'
776 call messages_info(1, namespace=namespace)
777 end if
778
779 call lcao_wf(lcao, st, gr, ions, hm, namespace, start = st_start)
780
781 ! In some rare cases, like bad pseudopotentials with unbound orbitals,
782 ! we might not have enough orbitals to go up to the Fermi energy
783 ! In this case, we cannot set the eigenvales to a huge value, as
784 ! the smearing will not converge. We set them to zero then
785 select case (st%d%ispin)
786 case (unpolarized)
787 required_min_nst = int(st%qtot/2)
788 case (spin_polarized)
789 required_min_nst = int(st%qtot/2)
790 case (spinors)
791 required_min_nst = int(st%qtot)
792 end select
793 if (st%smear%method /= smear_fixed_occ .and. st%smear%method /= smear_semiconductor) then
794 if (lcao%norbs <= required_min_nst .and. lcao%norbs < st%nst) then
795 st%eigenval(lcao%norbs+1:,:) = m_zero
796 end if
797 end if
798
799 if (.not. present(st_start)) then
800 call states_elec_fermi(st, namespace, gr)
801 call states_elec_write_eigenvalues(min(st%nst, lcao%norbs), st, space, hm%kpoints, namespace=namespace)
802
803 ! Update the density and the Hamiltonian
804 if (lcao%mode == option__lcaostart__lcao_full) then
805 call v_ks_h_setup(namespace, space, gr, ions, ext_partners, st, ks, hm, &
806 calc_eigenval = .false., calc_current=.false.)
807 if (st%d%ispin > unpolarized) then
808 assert(present(lmm_r))
809 call write_magnetic_moments(gr, st, ions, gr%der%boundaries, lmm_r, namespace=namespace)
810 end if
811 end if
812 end if
813 end if
814
815 if (.not. lcao_done .or. lcao%norbs < st%nst) then
816
817 if (lcao_done) then
818 st_start_random = lcao%norbs + 1
819 else
820 st_start_random = 1
821 end if
822 if (present(st_start)) st_start_random = max(st_start, st_start_random)
823
824 if (st_start_random > 1) then
825 write(message(1),'(a,i8,a)') 'Generating random wavefunctions for states ', st_start_random, ' and above'
826 call messages_info(1, namespace=namespace)
827 end if
828
829 ! Randomly generate the initial wavefunctions.
830 call states_elec_generate_random(st, gr, hm%kpoints, ist_start_ = st_start_random, normalized = .false.)
831
832 call messages_write('Orthogonalizing wavefunctions.')
833 call messages_info(namespace=namespace)
834 call states_elec_orthogonalize(st, namespace, gr)
835
836 if (.not. lcao_done) then
837 ! If we are doing unocc calculation, do not mess with the correct eigenvalues and occupations
838 ! of the occupied states.
839 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
840 calc_eigenval=.not. present(st_start), calc_current=.false.) ! get potentials
841 if (.not. present(st_start)) then
842 call states_elec_fermi(st, namespace, gr) ! occupations
843 end if
844
845 end if
846
847 else if (present(st_start)) then
848
849 if (st_start > 1) then
850 call messages_write('Orthogonalizing wavefunctions.')
851 call messages_info(namespace=namespace)
852 call states_elec_orthogonalize(st, namespace, gr)
853 end if
854
855 end if
856
857 call lcao_end(lcao)
858
859
860 call profiling_out('LCAO_RUN')
861 pop_sub(lcao_run)
862 end subroutine lcao_run
863
864 ! ---------------------------------------------------------
865 subroutine lcao_end(this)
866 type(lcao_t), intent(inout) :: this
867
868 push_sub(lcao_end)
869
870 safe_deallocate_a(this%calc_atom)
871 safe_deallocate_a(this%norb_atom)
872 safe_deallocate_a(this%basis_atom)
873 safe_deallocate_a(this%basis_orb)
874 safe_deallocate_a(this%atom_orb_basis)
875 safe_deallocate_a(this%radius)
876 safe_deallocate_a(this%sphere)
877 safe_deallocate_a(this%orbitals)
878
879 safe_deallocate_a(this%atom)
880 safe_deallocate_a(this%level)
881 safe_deallocate_a(this%ddim)
882 safe_deallocate_a(this%cst)
883 safe_deallocate_a(this%ck)
884 safe_deallocate_a(this%dbuff_single)
885 safe_deallocate_a(this%zbuff_single)
886 safe_deallocate_a(this%dbuff)
887 safe_deallocate_a(this%zbuff)
888
889 this%initialized = .false.
890 pop_sub(lcao_end)
891 end subroutine lcao_end
892
893
894 ! ---------------------------------------------------------
895 subroutine lcao_wf(this, st, gr, ions, hm, namespace, start)
896 type(lcao_t), intent(inout) :: this
897 type(states_elec_t), intent(inout) :: st
898 type(grid_t), intent(in) :: gr
899 type(ions_t), intent(in) :: ions
900 type(hamiltonian_elec_t), intent(in) :: hm
901 type(namespace_t), intent(in) :: namespace
902 integer, optional, intent(in) :: start
903
904 integer :: start_
905
906 assert(this%initialized)
907
908 call profiling_in("LCAO")
909 push_sub(lcao_wf)
910
911 start_ = 1
912 if (present(start)) start_ = start
913
914 if (this%alternative) then
915 if (states_are_real(st)) then
916 call dlcao_alt_wf(this, st, gr, ions, hm, namespace, start_)
917 else
918 call zlcao_alt_wf(this, st, gr, ions, hm, namespace, start_)
919 end if
920 else
921 if (states_are_real(st)) then
922 call dlcao_wf(this, st, gr, ions, hm, namespace, start_)
923 else
924 call zlcao_wf(this, st, gr, ions, hm, namespace, start_)
925 end if
926 end if
927 pop_sub(lcao_wf)
928 call profiling_out("LCAO")
929 end subroutine lcao_wf
930
931
932 ! ---------------------------------------------------------
933 logical function lcao_is_available(this) result(available)
934 type(lcao_t), intent(in) :: this
935
936 push_sub(lcao_is_available)
937
938 available = this%initialized .and. this%mode /= option__lcaostart__lcao_none
939
940 pop_sub(lcao_is_available)
941 end function lcao_is_available
942
943
944 ! ---------------------------------------------------------
945 integer function lcao_num_orbitals(this) result(norbs)
946 type(lcao_t), intent(in) :: this
947
948 push_sub(lcao_num_orbitals)
949 norbs = this%norbs
950
951 pop_sub(lcao_num_orbitals)
952 end function lcao_num_orbitals
953
954 ! ---------------------------------------------------------
955
956 subroutine lcao_local_index(this, ig, jg, il, jl, prow, pcol)
957 type(lcao_t), intent(in) :: this
958 integer, intent(in) :: ig
959 integer, intent(in) :: jg
960 integer, intent(out) :: il
961 integer, intent(out) :: jl
962 integer, intent(out) :: prow
963 integer, intent(out) :: pcol
964
965 ! no PUSH_SUB, called too often
966#ifdef HAVE_SCALAPACK
967 call infog2l(ig, jg, this%desc(1), this%nproc(1), this%nproc(2), this%myroc(1), this%myroc(2), &
968 il, jl, prow, pcol)
969#else
970 il = ig
971 jl = jg
972 prow = 0
973 pcol = 0
974#endif
975
976 end subroutine lcao_local_index
977
978 ! ---------------------------------------------------------
979
983 subroutine lcao_alt_end_orbital(this, iatom)
984 type(lcao_t), intent(inout) :: this
985 integer, intent(in) :: iatom
986
987 push_sub(lcao_alt_end_orbital)
989 if (this%is_orbital_initialized(iatom)) then
990 call this%orbitals(iatom)%end()
991 this%is_orbital_initialized(iatom) = .false.
992 end if
993
994 pop_sub(lcao_alt_end_orbital)
995
996 end subroutine lcao_alt_end_orbital
997
998 ! ---------------------------------------------------------
999
1000 subroutine lcao_atom_density(this, st, mesh, ions, iatom, spin_channels, rho)
1001 type(lcao_t), intent(inout) :: this
1002 type(states_elec_t), intent(in) :: st
1003 class(mesh_t), intent(in) :: mesh
1004 type(ions_t), target, intent(in) :: ions
1005 integer, intent(in) :: iatom
1006 integer, intent(in) :: spin_channels
1007 real(real64), intent(inout) :: rho(:, :)
1008
1009 real(real64), allocatable :: dorbital(:, :)
1010 complex(real64), allocatable :: zorbital(:, :)
1011 real(real64), allocatable :: factors(:)
1012 real(real64) :: factor, aa
1013 integer :: iorb, ip, ii, ll, mm, ispin
1014 type(ps_t), pointer :: ps
1015 logical :: use_stored_orbitals
1016
1017 push_sub(lcao_atom_density)
1018
1019 rho = m_zero
1020
1021 use_stored_orbitals = ions%atom(iatom)%species%is_ps() &
1022 .and. states_are_real(st) .and. spin_channels == 1 .and. lcao_is_available(this) &
1023 .and. st%d%dim == 1 .and. .not. ions%space%is_periodic()
1024
1025
1026 ! we can use the orbitals we already have calculated
1027 if (use_stored_orbitals) then
1028 !There is no periodic copies here, so this will not work for periodic systems
1029 assert(.not. ions%space%is_periodic())
1030
1031 select type(spec=>ions%atom(iatom)%species)
1032 class is(pseudopotential_t)
1033 ps => spec%ps
1034 class default
1035 assert(.false.)
1036 end select
1037
1038 if (.not. this%alternative) then
1039
1040 if (states_are_real(st)) then
1041 safe_allocate(dorbital(1:mesh%np, 1:st%d%dim))
1042 else
1043 safe_allocate(zorbital(1:mesh%np, 1:st%d%dim))
1044 end if
1045
1046 do iorb = 1, this%norbs
1047 if (iatom /= this%atom(iorb)) cycle
1048
1049 call ions%atom(iatom)%species%get_iwf_ilm(this%level(iorb), 1, ii, ll, mm)
1050 factor = ps%conf%occ(ii, 1)/(m_two*ll + m_one)
1051
1052 if (states_are_real(st)) then
1053 call dget_ao(this, st, mesh, ions, iorb, 1, dorbital, use_psi = .true.)
1054 !$omp parallel do
1055 do ip = 1, mesh%np
1056 rho(ip, 1) = rho(ip, 1) + factor*dorbital(ip, 1)**2
1057 end do
1058 else
1059 call zget_ao(this, st, mesh, ions, iorb, 1, zorbital, use_psi = .true.)
1060 !$omp parallel do
1061 do ip = 1, mesh%np
1062 rho(ip, 1) = rho(ip, 1) + factor*abs(zorbital(ip, 1))**2
1063 end do
1064 end if
1065
1066 end do
1067
1068 safe_deallocate_a(dorbital)
1069 safe_deallocate_a(zorbital)
1070
1071 else
1072
1073 ! for simplicity, always use real ones here.
1074 call dlcao_alt_get_orbital(this, this%sphere(iatom), ions, 1, iatom, this%norb_atom(iatom))
1075
1076 ! the extra orbitals with the derivative are not relevant here, hence we divide by this%mult
1077 safe_allocate(factors(1:this%norb_atom(iatom)/this%mult))
1078
1079 do iorb = 1, this%norb_atom(iatom)/this%mult
1080 call ions%atom(iatom)%species%get_iwf_ilm(iorb, 1, ii, ll, mm)
1081 factors(iorb) = ps%conf%occ(ii, 1)/(m_two*ll + m_one)
1082 end do
1083
1084 !$omp parallel do private(ip, aa, iorb) if(.not. this%sphere(iatom)%overlap)
1085 do ip = 1, this%sphere(iatom)%np
1086 aa = m_zero
1087 do iorb = 1, this%norb_atom(iatom)/this%mult
1088 aa = aa + factors(iorb)*this%orbitals(iatom)%dff_linear(ip, iorb)**2
1089 end do
1090 rho(this%sphere(iatom)%map(ip), 1) = rho(this%sphere(iatom)%map(ip), 1) + aa
1091 end do
1092
1093 safe_deallocate_a(factors)
1094
1095 end if
1096
1097 else
1098 call species_atom_density(ions%atom(iatom)%species, ions%namespace, ions%space, ions%latt, &
1099 ions%pos(:, iatom), mesh, spin_channels, rho)
1100 end if
1101
1102 ! The above code can sometimes return negative values of the density. Here we avoid introducing
1103 ! them in the calculation of v_s, mixing, ...
1104 do ispin = 1, spin_channels
1105 !$omp parallel do simd
1106 do ip = 1, mesh%np
1107 rho(ip, ispin) = max(rho(ip, ispin), m_zero)
1108 end do
1109 end do
1110
1111 pop_sub(lcao_atom_density)
1112 end subroutine lcao_atom_density
1113
1114 ! ---------------------------------------------------------
1116 subroutine lcao_guess_density(this, namespace, st, gr, hm, ions, qtot, ispin, spin_channels, rho)
1117 type(lcao_t), intent(inout) :: this
1118 type(namespace_t), intent(in) :: namespace
1119 type(states_elec_t), intent(in) :: st
1120 type(grid_t), intent(in) :: gr
1121 type(hamiltonian_elec_t), intent(in) :: hm
1122 type(ions_t), intent(in) :: ions
1123 real(real64), intent(in) :: qtot
1124 integer, intent(in) :: ispin, spin_channels
1125 real(real64), contiguous, intent(out) :: rho(:, :)
1126
1127 integer :: ia, is, idir, gmd_opt, ip
1128 integer(int64), save :: iseed = splitmix64_321
1129 type(block_t) :: blk
1130 real(real64) :: rr, rnd, phi, theta, mag(1:3), lmag, n1, n2,arg
1131 real(real64), allocatable :: atom_rho(:,:)
1132 logical :: parallelized_in_atoms
1133
1134 push_sub(lcao_guess_density)
1135
1136 parallelized_in_atoms = .false.
1137
1138 if (spin_channels == 1) then
1139 gmd_opt = initrho_paramagnetic
1140 else
1141 !%Variable GuessMagnetDensity
1142 !%Type integer
1143 !%Default ferromagnetic
1144 !%Section SCF::LCAO
1145 !%Description
1146 !% The guess density for the SCF cycle is just the sum of all the atomic densities.
1147 !% When performing spin-polarized or non-collinear-spin calculations this option sets
1148 !% the guess magnetization density.
1149 !%
1150 !% For anti-ferromagnetic configurations, the <tt>user_defined</tt> option should be used.
1151 !%
1152 !% Note that if the <tt>paramagnetic</tt> option is used, the final ground state will also be
1153 !% paramagnetic, but the same is not true for the other options.
1154 !%Option paramagnetic 1
1155 !% Magnetization density is zero.
1156 !%Option ferromagnetic 2
1157 !% Magnetization density is the sum of the atomic magnetization densities.
1158 !%Option random 3
1159 !% Each atomic magnetization density is randomly rotated.
1160 !%Option user_defined 77
1161 !% The atomic magnetization densities are rotated so that the magnetization
1162 !% vector has the same direction as a vector provided by the user. In this case,
1163 !% the <tt>AtomsMagnetDirection</tt> block has to be set.
1164 !%End
1165 call parse_variable(namespace, 'GuessMagnetDensity', initrho_ferromagnetic, gmd_opt)
1166 if (.not. varinfo_valid_option('GuessMagnetDensity', gmd_opt)) call messages_input_error(namespace, 'GuessMagnetDensity')
1167 call messages_print_var_option('GuessMagnetDensity', gmd_opt, namespace=namespace)
1168 end if
1169
1170 if (parse_is_defined(namespace, 'GuessMagnetDensity') .and. (hm%theory_level == hartree_fock &
1171 .or. hm%theory_level == generalized_kohn_sham_dft)) then
1172 message(1) = "GuessMagnetDensity cannot be used for Hartree-Fock and generalized Kohn-Sham calculation."
1173 message(2) = "Please perform a LDA or GGA calculation first and restart from this calculation."
1174 call messages_fatal(2, namespace=namespace)
1175 end if
1176
1177 rho = m_zero
1178 select case (gmd_opt)
1179 case (initrho_paramagnetic)
1180 safe_allocate(atom_rho(1:gr%np, 1:spin_channels))
1181
1182 parallelized_in_atoms = .true.
1183
1184 do ia = ions%atoms_dist%start, ions%atoms_dist%end
1185 call lcao_atom_density(this, st, gr, ions, ia, spin_channels, atom_rho)
1186 call lalg_axpy(gr%np, spin_channels, m_one, atom_rho, rho)
1187 end do
1188
1189 if (spin_channels == 2) then
1190 !$omp parallel do
1191 do ip = 1, gr%np
1192 rho(ip, 1) = m_half*(rho(ip, 1) + rho(ip, 2))
1193 rho(ip, 2) = rho(ip, 1)
1194 end do
1195 end if
1196
1198 safe_allocate(atom_rho(1:gr%np, 1:2))
1199
1200 parallelized_in_atoms = .true.
1201
1202 rho = m_zero
1203 do ia = ions%atoms_dist%start, ions%atoms_dist%end
1204 call lcao_atom_density(this, st, gr, ions, ia, 2, atom_rho(1:gr%np, 1:2))
1205 rho(1:gr%np, 1:2) = rho(1:gr%np, 1:2) + atom_rho(1:gr%np, 1:2)
1206 end do
1207
1208 case (initrho_random) ! Randomly oriented spins
1209 safe_allocate(atom_rho(1:gr%np, 1:2))
1210 do ia = 1, ions%natoms
1211 call lcao_atom_density(this, st, gr, ions, ia, 2, atom_rho)
1212
1213 if (ispin == spin_polarized) then
1214 call quickrnd(iseed, rnd)
1215 rnd = rnd - m_half
1216 if (rnd > m_zero) then
1217 rho(1:gr%np, 1:2) = rho(1:gr%np, 1:2) + atom_rho(1:gr%np, 1:2)
1218 else
1219 rho(1:gr%np, 1) = rho(1:gr%np, 1) + atom_rho(1:gr%np, 2)
1220 rho(1:gr%np, 2) = rho(1:gr%np, 2) + atom_rho(1:gr%np, 1)
1221 end if
1222 elseif (ispin == spinors) then
1223 call quickrnd(iseed, phi)
1224 call quickrnd(iseed, theta)
1225 phi = phi*m_two*m_pi
1226 theta = theta*m_pi*m_half
1227
1228 call accumulate_rotated_density(gr, rho, atom_rho, theta, phi)
1229 end if
1230 end do
1231
1232 case (initrho_userdef) ! User-defined
1233
1234 !%Variable AtomsMagnetDirection
1235 !%Type block
1236 !%Section SCF::LCAO
1237 !%Description
1238 !% This option is only used when <tt>GuessMagnetDensity</tt> is
1239 !% set to <tt>user_defined</tt>. It provides a direction for the
1240 !% magnetization vector of each atom when building the guess
1241 !% density. In order to do that, the user should specify the
1242 !% coordinates of a vector that has the desired direction and
1243 !% norm. Note that it is necessary to maintain the ordering in
1244 !% which the species were defined in the coordinates
1245 !% specifications.
1246 !%
1247 !% For spin-polarized calculations, the vectors should have only
1248 !% one component; for non-collinear-spin calculations, they
1249 !% should have three components. If the norm of the vector is greater
1250 !% than the number of valence electrons in the atom, it will be rescaled
1251 !% to this number, which is the maximum possible magnetization.
1252 !%End
1253 if (parse_block(namespace, 'AtomsMagnetDirection', blk) < 0) then
1254 message(1) = "AtomsMagnetDirection block is not defined."
1255 call messages_fatal(1, namespace=namespace)
1256 end if
1257
1258 if (parse_block_n(blk) /= ions%natoms) then
1259 message(1) = "AtomsMagnetDirection block has the wrong number of rows."
1260 call messages_fatal(1, namespace=namespace)
1261 end if
1262
1263 safe_allocate(atom_rho(1:gr%np, 1:2))
1264 do ia = 1, ions%natoms
1265 !Read from AtomsMagnetDirection block
1266 if (ispin == spin_polarized) then
1267 call parse_block_float(blk, ia-1, 0, mag(1))
1268 mag(2:3) = m_zero !Else, this is unitialized and lead to a FPE in the case (lmag > n1+n2)
1269 lmag = abs(mag(1))
1270 elseif (ispin == spinors) then
1271 do idir = 1, 3
1272 call parse_block_float(blk, ia-1, idir-1, mag(idir))
1273 if (abs(mag(idir)) < 1.0e-20_real64) mag(idir) = m_zero
1274 end do
1275 lmag = norm2(mag(1:3))
1276 end if
1277
1278 !Get atomic density
1279 call lcao_atom_density(this, st, gr, ions, ia, 2, atom_rho)
1280
1281 !Scale magnetization density
1282 n1 = dmf_integrate(gr, atom_rho(:, 1))
1283 n2 = dmf_integrate(gr, atom_rho(:, 2))
1284
1285 if (lmag > n1 + n2) then
1286 mag = mag*(n1 + n2)/lmag
1287 lmag = n1 + n2
1288 elseif (abs(lmag) <= m_epsilon) then
1289 if (abs(n1 - n2) <= m_epsilon) then
1290 call lalg_axpy(gr%np, 2, m_one, atom_rho, rho)
1291 else
1292 !$omp parallel do simd
1293 do ip = 1, gr%np
1294 atom_rho(ip, 1) = m_half*(atom_rho(ip, 1) + atom_rho(ip, 2))
1295 rho(ip, 1) = rho(ip, 1) + atom_rho(ip, 1)
1296 rho(ip, 2) = rho(ip, 2) + atom_rho(ip, 1)
1297 end do
1298 end if
1299 cycle
1300 end if
1301
1302 if (.not. is_close(n1 - n2, lmag) .and. abs(n2) > m_epsilon) then
1303 if (n1 - n2 < lmag) then
1304 call lalg_axpy(gr%np, (lmag - n1 + n2)/m_two/n2, atom_rho(:, 2), atom_rho(:, 1))
1305 call lalg_scal(gr%np, (n1 + n2 - lmag)/m_two/n2, atom_rho(:, 2))
1306 elseif (n1 - n2 > lmag) then
1307 call lalg_axpy(gr%np, (n1 - n2 - lmag)/m_two/n1, atom_rho(:, 1), atom_rho(:, 2))
1308 call lalg_scal(gr%np, (n1 + n2 + lmag)/m_two/n1, atom_rho(:, 1))
1309 end if
1310 end if
1311
1312 !Rotate magnetization density
1313 if (ispin == spin_polarized) then
1314
1315 if (mag(1) > m_zero) then
1316 call lalg_axpy(gr%np, 2, m_one, atom_rho, rho)
1317 else
1318 call lalg_axpy(gr%np, m_one, atom_rho(:,2), rho(:,1))
1319 call lalg_axpy(gr%np, m_one, atom_rho(:,1), rho(:,2))
1320 end if
1321
1322 elseif (ispin == spinors) then
1323 assert(lmag > m_zero)
1324 theta = acos(mag(3)/lmag)
1325 if (abs(mag(1)) <= m_epsilon) then
1326 if (abs(mag(2)) <= m_epsilon) then
1327 phi = m_zero
1328 elseif (mag(2) < m_zero) then
1329 phi = m_pi*m_three*m_half
1330 elseif (mag(2) > m_zero) then
1331 phi = m_pi*m_half
1332 end if
1333 else
1334 ! In some rare cases this can be larger than one
1335 arg = mag(1)/sin(theta)/lmag
1336 if (abs(arg) > m_one) arg = sign(m_one, arg)
1337 phi = acos(arg)
1338 if (mag(2) < m_zero) then
1339 phi = m_two*m_pi - phi
1340 end if
1341 end if
1342 theta = m_half*theta
1343 call accumulate_rotated_density(gr, rho, atom_rho, theta, phi)
1344 end if
1345 end do
1346
1347 call parse_block_end(blk)
1348
1349 end select
1350
1351
1352 if (ions%atoms_dist%parallel .and. parallelized_in_atoms) then
1353 ! NOTE: if random or user_defined are made parallelized in atoms, below should be st%d%nspin instead of spin_channels
1354 do is = 1, spin_channels
1355 call lalg_copy(gr%np, rho(:,is), atom_rho(:,1))
1356 call ions%atoms_dist%mpi_grp%allreduce(atom_rho(1, 1), rho(1, is), gr%np, mpi_double_precision, mpi_sum)
1357 end do
1358 end if
1359
1360 ! we now renormalize the density (necessary if we have a charged system)
1361 rr = m_zero
1362 do is = 1, spin_channels
1363 rr = rr + dmf_integrate(gr, rho(:, is), reduce = .false.)
1364 end do
1365 if (gr%parallel_in_domains) then
1366 call gr%allreduce(rr)
1367 end if
1368
1369 write(message(1),'(a,f13.6)')'Info: Unnormalized total charge = ', rr
1370 call messages_info(1, namespace=namespace)
1371
1372 if (abs(rr) > m_epsilon) then ! We only renormalize if the density is not zero
1373 rr = qtot / rr
1374 rho = rr * rho
1375 end if
1376 rr = m_zero
1377 do is = 1, spin_channels
1378 rr = rr + dmf_integrate(gr, rho(:, is), reduce = .false.)
1379 end do
1380 if (gr%parallel_in_domains) then
1381 call gr%allreduce(rr)
1382 end if
1383
1384 write(message(1),'(a,f13.6)')'Info: Renormalized total charge = ', rr
1385 call messages_info(1, namespace=namespace)
1386
1387 if (st%symmetrize_density) then
1388 do is = 1, st%d%nspin
1389 call dgrid_symmetrize_scalar_field(gr, rho(:, is))
1390 end do
1391 end if
1392
1393 safe_deallocate_a(atom_rho)
1394 pop_sub(lcao_guess_density)
1395 end subroutine lcao_guess_density
1396
1397 ! ---------------------------------------------------------
1398 subroutine accumulate_rotated_density(mesh, rho, atom_rho, theta, phi)
1399 class(mesh_t), intent(in) :: mesh
1400 real(real64), intent(inout) :: rho(:,:)
1401 real(real64), intent(in) :: atom_rho(:,:)
1402 real(real64), intent(in) :: theta, phi
1403
1404 integer :: ip
1405
1407
1408 !$omp parallel do simd
1409 do ip = 1, mesh%np
1410 rho(ip, 1) = rho(ip, 1) + cos(theta)**2*atom_rho(ip, 1) + sin(theta)**2*atom_rho(ip, 2)
1411 rho(ip, 2) = rho(ip, 2) + sin(theta)**2*atom_rho(ip, 1) + cos(theta)**2*atom_rho(ip, 2)
1412 rho(ip, 3) = rho(ip, 3) + cos(theta)*sin(theta)*cos(phi)*(atom_rho(ip, 1)-atom_rho(ip, 2))
1413 rho(ip, 4) = rho(ip, 4) - cos(theta)*sin(theta)*sin(phi)*(atom_rho(ip, 1)-atom_rho(ip, 2))
1414 end do
1415
1417 end subroutine accumulate_rotated_density
1418
1419 ! ---------------------------------------------------------
1420
1421 subroutine lcao_init_orbitals(this, namespace, st, gr, ions, start)
1422 type(lcao_t), intent(inout) :: this
1423 type(namespace_t), intent(in) :: namespace
1424 type(states_elec_t), intent(inout) :: st
1425 type(grid_t), intent(in) :: gr
1426 type(ions_t), intent(in) :: ions
1427 integer, optional, intent(in) :: start
1428
1429 if (.not. lcao_is_available(this)) return
1430
1431 push_sub(lcao_init_orbitals)
1432
1433 if (.not. this%alternative) then
1434 if (states_are_real(st)) then
1435 call dinit_orbitals(this, namespace, st, gr, ions, start)
1436 else
1437 call zinit_orbitals(this, namespace, st, gr, ions, start)
1438 end if
1439 else
1440 if (states_are_real(st)) then
1441 call dlcao_alt_init_orbitals(this, namespace, st, gr, ions, start)
1442 else
1443 call zlcao_alt_init_orbitals(this, namespace, st, gr, ions, start)
1444 end if
1445
1446 end if
1447
1448 pop_sub(lcao_init_orbitals)
1449 end subroutine lcao_init_orbitals
1450
1451#include "undef.F90"
1452#include "real.F90"
1453#include "lcao_inc.F90"
1454
1455#include "undef.F90"
1456#include "complex.F90"
1457#include "lcao_inc.F90"
1458
1459
1460end module lcao_oct_m
1461
1462!! Local Variables:
1463!! mode: f90
1464!! coding: utf-8
1465!! End:
subroutine info()
Definition: em_resp.F90:1096
constant times a vector plus a vector
Definition: lalg_basic.F90:170
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:185
scales a vector by a constant
Definition: lalg_basic.F90:156
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:180
double acos(double __x) __attribute__((__nothrow__
double sin(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
subroutine lcao2_init()
Definition: lcao.F90:586
This module implements batches of mesh functions.
Definition: batch.F90:133
This module contains interfaces for BLACS routines Interfaces are from http:
Definition: blacs.F90:27
This module provides the BLACS processor grid.
logical pure function, public blacs_proc_grid_null(this)
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:122
type(debug_t), save, public debug
Definition: debug.F90:154
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
integer, parameter, public spin_polarized
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
real(real64), parameter, public m_epsilon
Definition: global.F90:203
real(real64), parameter, public m_half
Definition: global.F90:193
real(real64), parameter, public m_one
Definition: global.F90:188
real(real64), parameter, public m_three
Definition: global.F90:190
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public dgrid_symmetrize_scalar_field(gr, field, suppress_warning)
Definition: grid.F90:661
integer, parameter, public generalized_kohn_sham_dft
integer, parameter, public hartree
integer, parameter, public hartree_fock
integer, parameter, public kohn_sham_dft
This module defines classes and functions for interaction partners.
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:468
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:354
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:395
This module contains interfaces for LAPACK routines.
Definition: lapack.F90:118
subroutine, public lcao_init_orbitals(this, namespace, st, gr, ions, start)
Definition: lcao.F90:1515
subroutine, public lcao_wf(this, st, gr, ions, hm, namespace, start)
Definition: lcao.F90:989
subroutine zlcao_alt_init_orbitals(this, namespace, st, gr, ions, start)
Definition: lcao.F90:3358
subroutine zget_ao(this, st, mesh, ions, iorb, ispin, ao, use_psi)
Definition: lcao.F90:3326
subroutine, public lcao_end(this)
Definition: lcao.F90:959
subroutine, public lcao_run(namespace, space, gr, ions, ext_partners, st, ks, hm, st_start, lmm_r)
Definition: lcao.F90:795
subroutine zlcao_alt_wf(this, st, gr, ions, hm, namespace, start)
The alternative implementation.
Definition: lcao.F90:3411
subroutine dlcao_wf(this, st, gr, ions, hm, namespace, start)
Definition: lcao.F90:1685
subroutine dlcao_alt_init_orbitals(this, namespace, st, gr, ions, start)
Definition: lcao.F90:2082
subroutine, public lcao_init(this, namespace, space, gr, ions, st)
Definition: lcao.F90:242
subroutine zlcao_wf(this, st, gr, ions, hm, namespace, start)
Definition: lcao.F90:2953
subroutine dlcao_alt_wf(this, st, gr, ions, hm, namespace, start)
The alternative implementation.
Definition: lcao.F90:2135
subroutine lcao_guess_density(this, namespace, st, gr, hm, ions, qtot, ispin, spin_channels, rho)
builds a density which is the sum of the atomic densities
Definition: lcao.F90:1210
integer, parameter initrho_userdef
Definition: lcao.F90:232
integer, parameter initrho_random
Definition: lcao.F90:232
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...
Definition: lcao.F90:1077
subroutine lcao_local_index(this, ig, jg, il, jl, prow, pcol)
Definition: lcao.F90:1050
integer function, public lcao_num_orbitals(this)
Definition: lcao.F90:1039
subroutine dinit_orbitals(this, namespace, st, gr, ions, start)
Definition: lcao.F90:1941
subroutine dget_ao(this, st, mesh, ions, iorb, ispin, ao, use_psi)
Definition: lcao.F90:2050
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...
Definition: lcao.F90:2757
subroutine zinit_orbitals(this, namespace, st, gr, ions, start)
Definition: lcao.F90:3217
subroutine lcao_atom_density(this, st, mesh, ions, iatom, spin_channels, rho)
Definition: lcao.F90:1094
subroutine accumulate_rotated_density(mesh, rho, atom_rho, theta, phi)
Definition: lcao.F90:1492
integer, parameter initrho_ferromagnetic
Definition: lcao.F90:232
logical function, public lcao_is_available(this)
Definition: lcao.F90:1027
subroutine, public write_magnetic_moments(mesh, st, ions, boundaries, lmm_r, iunit, namespace)
Definition: magnetic.F90:211
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines various routines, operating on mesh functions.
real(real64) function, public dmf_integrate(mesh, ff, mask, reduce)
Integrate a function on the mesh.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1057
subroutine, public messages_new_line()
Definition: messages.F90:1146
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1097
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:624
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
Definition: ps.F90:114
integer(int64), parameter, public splitmix64_321
Definition: quickrnd.F90:134
This module contains interfaces for ScaLAPACK routines Interfaces are from http:
Definition: scalapack.F90:131
integer, parameter, public smear_semiconductor
Definition: smear.F90:171
integer, parameter, public smear_fixed_occ
Definition: smear.F90:171
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.
Definition: unit.F90:132
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)
Definition: v_ks.F90:732
subroutine, public v_ks_h_setup(namespace, space, gr, ions, ext_partners, st, ks, hm, calc_eigenval, calc_current)
Definition: v_ks.F90:679
Definition: xc.F90:114
logical pure function, public xc_is_orbital_dependent(xcs)
Is the xc family orbital dependent.
Definition: xc.F90:541
integer, parameter, public sic_pz_oep
Perdew-Zunger SIC (OEP way)
Definition: xc_sic.F90:148
Extension of space that contains the knowledge of the spin dimension.
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
Describes mesh distribution to nodes.
Definition: mesh.F90:186
The states_elec_t class contains all electronic wave functions.
int true(void)