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