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