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