64 type(lda_u_t),
intent(in) :: this
65 character(len=*),
intent(in) :: dir
66 type(states_elec_t),
intent(in) :: st
67 type(namespace_t),
intent(in) :: namespace
69 integer :: iunit, ios, ispin, im, imp, ios2, inn
70 type(orbitalset_t),
pointer :: os, os2
75 iunit =
io_open(trim(dir) //
"/occ_matrices", namespace, action=
'write')
76 write(iunit,
'(a)')
' Occupation matrices '
78 do ios = 1, this%norbsets
79 os => this%orbsets(ios)
81 do ispin = 1,st%d%nspin
82 write(iunit,
'(a, i4, a, i4)')
' Orbital set ', ios,
' spin ', ispin
85 write(iunit,
'(1x)',advance=
'no')
88 do imp = 1, os%norbs-1
89 write(iunit,
'(f14.8)', advance=
'no') this%dn(im, imp, ispin, ios)
91 write(iunit,
'(f14.8)') this%dn(im, os%norbs, ispin, ios)
93 do imp = 1, os%norbs-1
94 write(iunit,
'(f14.8,f14.8)', advance=
'no') this%zn(im, imp, ispin, ios)
96 write(iunit,
'(f14.8,f14.8)') this%zn(im, os%norbs, ispin, ios)
104 iunit =
io_open(trim(dir) //
"/renorm_occ_matrices", namespace, action=
'write')
105 write(iunit,
'(a)')
' Renormalized occupation matrices '
107 do ios = 1, this%norbsets
108 os => this%orbsets(ios)
109 do ispin = 1, st%d%nspin
110 write(iunit,
'(a, i4, a, i4)')
' Orbital set ', ios,
' spin ', ispin
113 write(iunit,
'(1x)',advance=
'no')
116 do imp = 1, os%norbs-1
117 write(iunit,
'(f14.8)', advance=
'no') this%dn_alt(im, imp, ispin, ios)
119 write(iunit,
'(f14.8)') this%dn_alt(im, os%norbs, ispin, ios)
121 do imp = 1, os%norbs-1
122 write(iunit,
'(f14.8,f14.8)', advance=
'no') this%zn_alt(im, imp, ispin, ios)
124 write(iunit,
'(f14.8,f14.8)') this%zn_alt(im, os%norbs, ispin, ios)
132 if(this%intersite)
then
133 iunit =
io_open(trim(dir) //
"/intersite_occ_matrices", namespace, action=
'write')
134 write(iunit,
'(a)')
' Intersite occupation matrices '
135 do ios = 1, this%norbsets
136 os => this%orbsets(ios)
137 do ispin = 1, st%d%nspin
138 write(iunit,
'(a, i4, a, i4)')
' Orbital set ', ios,
' spin ', ispin
139 do inn = 1, os%nneighbors
140 write(iunit,
'(a, i4)')
' Neighbour ', inn
141 ios2 = os%map_os(inn)
142 os2 => this%orbsets(ios2)
145 write(iunit,
'(1x)',advance=
'no')
148 do imp = 1, os2%norbs - 1
149 write(iunit,
'(f14.8)', advance=
'no') this%dn_ij(im, imp, ispin, ios, inn)
151 write(iunit,
'(f14.8)') this%dn_ij(im, os2%norbs, ispin, ios, inn)
153 do imp = 1, os2%norbs - 1
154 write(iunit,
'(f14.8,f14.8)', advance=
'no') this%zn_ij(im, imp, ispin, ios, inn)
156 write(iunit,
'(f14.8,f14.8)') this%zn_ij(im, os2%norbs, ispin, ios, inn)
172 type(
lda_u_t),
intent(in) :: this
173 character(len=*),
intent(in) :: dir
176 integer :: iunit, ios
181 iunit =
io_open(trim(dir) //
"/effectiveU", namespace, action=
'write')
184 write(iunit,
'(a,a,a,f7.3,a)')
'Hubbard U [', &
186 write(iunit,
'(a,6x,14x,a)')
' Orbital',
'U'
187 do ios = 1, this%norbsets
188 if (.not. this%basisfromstates)
then
189 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals)
then
190 if (this%orbsets(ios)%nn /= 0)
then
191 write(iunit,
'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
192 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
195 write(iunit,
'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
200 if (this%orbsets(ios)%nn /= 0)
then
201 write(iunit,
'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
202 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
203 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
206 write(iunit,
'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
208 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
218 write(iunit,
'(a,a,a,f7.3,a)')
'Hund J [', &
220 write(iunit,
'(a,6x,14x,a)')
' Orbital',
'J'
221 do ios = 1, this%norbsets
222 if (.not. this%basisfromstates)
then
223 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals)
then
224 if (this%orbsets(ios)%nn /= 0)
then
225 write(iunit,
'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
226 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
229 write(iunit,
'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
234 if (this%orbsets(ios)%nn /= 0)
then
235 write(iunit,
'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
236 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
237 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
240 write(iunit,
'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
242 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
259 type(
lda_u_t),
intent(in) :: this
261 character(len=*),
intent(in) :: dir
264 integer :: iunit, ios
265 real(real64),
allocatable :: kanamori(:,:)
270 safe_allocate(kanamori(1:3,1:this%norbsets))
274 iunit =
io_open(trim(dir) //
"/kanamoriU", namespace, action=
'write')
276 write(iunit,
'(a,a,a,f7.3,a)')
'Intraorbital U [', &
278 write(iunit,
'(a,6x,14x,a)')
' Orbital',
'U'
279 do ios = 1, this%norbsets
280 if (.not. this%basisfromstates)
then
281 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals)
then
282 if (this%orbsets(ios)%nn /= 0)
then
283 write(iunit,
'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
284 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
287 write(iunit,
'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
292 if (this%orbsets(ios)%nn /= 0)
then
293 write(iunit,
'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
294 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
295 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
298 write(iunit,
'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
300 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
310 write(iunit,
'(a,a,a,f7.3,a)')
'Interorbital Up [', &
312 write(iunit,
'(a,6x,14x,a)')
' Orbital',
'Up'
313 do ios = 1, this%norbsets
314 if (.not. this%basisfromstates)
then
315 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals)
then
316 if (this%orbsets(ios)%nn /= 0)
then
317 write(iunit,
'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
318 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
321 write(iunit,
'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
326 if (this%orbsets(ios)%nn /= 0)
then
327 write(iunit,
'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
328 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
329 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
332 write(iunit,
'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
334 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
343 write(iunit,
'(a,a,a,f7.3,a)')
'Hund J [', &
345 write(iunit,
'(a,6x,14x,a)')
' Orbital',
'J'
346 do ios = 1, this%norbsets
347 if (.not. this%basisfromstates)
then
348 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals)
then
349 if (this%orbsets(ios)%nn /= 0)
then
350 write(iunit,
'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
354 write(iunit,
'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
359 if (this%orbsets(ios)%nn /= 0)
then
360 write(iunit,
'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
361 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
362 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
365 write(iunit,
'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
367 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
379 safe_deallocate_a(kanamori)
390 type(
lda_u_t),
intent(in) :: this
391 character(len=*),
intent(in) :: dir
392 type(
ions_t),
intent(in) :: ions
393 class(
mesh_t),
intent(in) :: mesh
397 integer :: iunit, ia, ios, im
398 real(real64),
allocatable :: mm(:,:)
405 iunit =
io_open(trim(dir)//
"/magnetization.xsf", namespace, action=
'write', position=
'asis')
407 if (this%nspins > 1)
then
408 safe_allocate(mm(1:mesh%box%dim, 1:ions%natoms))
411 do ios = 1, this%norbsets
412 ia = this%orbsets(ios)%iatom
413 do im = 1, this%orbsets(ios)%norbs
415 mm(3, ia) = mm(3, ia) + this%dn(im,im,1,ios) - this%dn(im,im,2,ios)
417 mm(3, ia) = mm(3, ia) + real(this%zn(im,im,1,ios) - this%zn(im,im,2,ios), real64)
419 if (this%nspins /= this%spin_channels)
then
420 mm(1, ia) = mm(1, ia) + 2*real(this%zn(im,im,3,ios), real64)
421 mm(2, ia) = mm(2, ia) - 2*aimag(this%zn(im,im,3,ios))
426 call write_xsf_geometry(iunit, ions%space, ions%latt, ions%pos, ions%atom, mesh, forces = mm)
427 safe_deallocate_a(mm)
437 type(
lda_u_t),
intent(in) :: this
438 integer,
optional,
intent(in) :: iunit
439 type(
namespace_t),
optional,
intent(in) :: namespace
445 write(
message(1),
'(a,a,a,f7.3,a)')
'Effective Hubbard U [', &
447 write(
message(2),
'(a,6x,14x,a)')
' Orbital',
'U'
450 do ios = 1, this%norbsets
451 if (.not. this%basisfromstates)
then
452 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals)
then
453 if (this%orbsets(ios)%nn /= 0)
then
454 write(
message(1),
'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
455 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
458 write(
message(1),
'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
463 if (this%orbsets(ios)%nn /= 0)
then
464 write(
message(1),
'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
465 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
466 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
469 write(
message(1),
'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
471 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
486 type(
lda_u_t),
intent(in) :: this
487 integer,
optional,
intent(in) :: iunit
488 type(
namespace_t),
optional,
intent(in) :: namespace
490 integer :: ios, icopies, ios2
492 if (.not. this%intersite)
return
496 write(
message(1),
'(a,a,a,f7.3,a)')
'Effective intersite V [', &
498 write(
message(2),
'(a,14x,a)')
' Orbital',
'V'
501 do ios = 1, this%norbsets
502 do icopies = 1, this%orbsets(ios)%nneighbors
503 ios2 = this%orbsets(ios)%map_os(icopies)
504 if (.not.this%basisfromstates)
then
505 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals)
then
506 if (this%orbsets(ios)%nn /= 0)
then
507 write(
message(1),
'(i4,a10, 2x, i1, a1, i4, 1x, i1, a1, f7.3, f15.6)') ios, &
508 trim(this%orbsets(ios)%spec%get_label()), this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), ios2, &
509 this%orbsets(ios2)%nn,
l_notation(this%orbsets(ios2)%ll), &
514 write(
message(1),
'(i4,a10, 3x, a1, i4, 1x, a1, f7.3, f15.6)') ios, &
515 trim(this%orbsets(ios)%spec%get_label()),
l_notation(this%orbsets(ios)%ll), ios2, &
521 if (this%orbsets(ios)%nn /= 0)
then
522 write(
message(1),
'(i4,a10, 2x, i1, a1, i1, a2, i4, f7.3, f15.6)') ios, &
523 trim(this%orbsets(ios)%spec%get_label()), this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
524 int(
m_two*(this%orbsets(ios)%jj)),
'/2', ios2, &
529 write(
message(1),
'(i4,a10, 3x, a1, i1, a2, i4, f7.3, f15.6)') ios, &
530 trim(this%orbsets(ios)%spec%get_label()),
l_notation(this%orbsets(ios)%ll), &
531 int(
m_two*(this%orbsets(ios)%jj)),
'/2', ios2, &
538 write(
message(1),
'(i4,a10, i4, f7.3, f15.6)') ios,
'states', ios2, &
552 subroutine lda_u_dump(restart, namespace, this, st, mesh, ierr)
555 type(
lda_u_t),
intent(in) :: this
557 class(
mesh_t),
intent(in) :: mesh
558 integer,
intent(out) :: ierr
560 integer :: err, occsize, ios, ncount
561 real(real64),
allocatable :: ueff(:), docc(:), veff(:)
562 complex(real64),
allocatable :: zocc(:)
574 message(1) =
"Debug: Writing DFT+U restart."
578 occsize = this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets
581 if (this%intersite)
then
582 occsize = occsize + 2*this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets*this%maxneighbors
588 safe_allocate(docc(1:occsize))
592 if (err /= 0) ierr = ierr + 1
593 safe_deallocate_a(docc)
595 safe_allocate(zocc(1:occsize))
599 if (err /= 0) ierr = ierr + 1
600 safe_deallocate_a(zocc)
605 safe_allocate(ueff(1:this%norbsets))
609 safe_deallocate_a(ueff)
610 if (err /= 0) ierr = ierr + 1
612 if (this%intersite .and. this%maxneighbors > 0)
then
614 do ios = 1, this%norbsets
615 ncount = ncount + this%orbsets(ios)%nneighbors
617 safe_allocate(veff(1:ncount))
620 safe_deallocate_a(veff)
621 if (err /= 0) ierr = ierr + 1
628 message(1) =
"Debug: Writing DFT+U restart done."
637 subroutine lda_u_load(restart, this, st, dftu_energy, ierr, occ_only, u_only)
639 type(
lda_u_t),
intent(inout) :: this
641 real(real64),
intent(out) :: dftu_energy
642 integer,
intent(out) :: ierr
643 logical,
optional,
intent(in) :: occ_only
644 logical,
optional,
intent(in) :: u_only
646 integer :: err, occsize, ncount, ios
647 real(real64),
allocatable :: ueff(:), docc(:), veff(:)
648 complex(real64),
allocatable :: zocc(:)
661 message(1) =
"Debug: Reading DFT+U restart."
667 safe_allocate(ueff(1:this%norbsets))
674 safe_deallocate_a(ueff)
676 if (this%intersite .and. this%maxneighbors > 0)
then
678 do ios = 1, this%norbsets
679 ncount = ncount + this%orbsets(ios)%nneighbors
681 safe_allocate(veff(1:ncount))
688 safe_deallocate_a(veff)
694 occsize = this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets
697 if (this%intersite)
then
698 occsize = occsize + 2*this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets*this%maxneighbors
704 safe_allocate(docc(1:occsize))
711 safe_deallocate_a(docc)
713 safe_allocate(zocc(1:occsize))
720 safe_deallocate_a(zocc)
733 message(1) =
"Debug: Reading DFT+U restart done."
742 type(
lda_u_t),
intent(in) :: this
746 class(
mesh_t),
intent(in) :: mesh
747 integer,
intent(out) :: ierr
749 integer :: coulomb_int_file, idim1, idim2, err
750 integer :: ios, im, imp, impp, imppp
751 character(len=256) :: lines(3)
752 logical :: complex_coulomb_integrals
765 message(1) =
"Debug: Writing Coulomb integrals restart."
769 complex_coulomb_integrals = .false.
770 do ios = 1, this%norbsets
771 if (this%orbsets(ios)%ndim > 1) complex_coulomb_integrals = .
true.
774 coulomb_int_file =
restart_open(restart,
'coulomb_integrals')
779 write(lines(1),
'(a20,i21)')
"norb=", this%norbsets
780 write(lines(2),
'(a20,i21)')
"dim=", st%d%dim
781 write(lines(3),
'(a20,i21)')
"checksum=", mesh%idx%checksum
783 if (err /= 0) ierr = ierr - 2
785 do ios = 1, this%norbsets
786 do im = 1, this%orbsets(ios)%norbs
787 do imp = 1, this%orbsets(ios)%norbs
788 do impp = 1, this%orbsets(ios)%norbs
789 do imppp = 1, this%orbsets(ios)%norbs
790 if(.not. complex_coulomb_integrals)
then
791 write(lines(1),
'(i4,i4,i4,i4,i4,e20.12)') ios, im, imp, impp, imppp, this%coulomb(im, imp, impp, imppp, ios)
793 do idim1 = 1, st%d%dim
794 do idim2 = 1, st%d%dim
795 write(lines(1),
'(i4,i4,i4,i4,i4,i4,i4,2e20.12)') ios, im, imp, impp, imppp, idim1, idim2, &
796 real(this%zcoulomb(im, imp, impp, imppp, idim1, idim2, ios), real64) , &
797 aimag(this%zcoulomb(im, imp, impp, imppp, idim1, idim2, ios))
814 message(1) =
"Debug: Writing Coulomb integrals restart done."
character(len=1), dimension(0:3), parameter, public l_notation
type(debug_t), save, public debug
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
complex(real64), parameter, public m_z0
subroutine, public write_xsf_geometry(iunit, space, latt, pos, atoms, mesh, forces, index)
for format specification see: http:
subroutine, public io_close(iunit, grp)
subroutine, public io_mkdir(fname, namespace, parents)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
subroutine, public lda_u_write_occupation_matrices(dir, this, st, namespace)
Prints the occupation matrices at the end of the scf calculation.
subroutine, public lda_u_write_kanamoriu(dir, st, this, namespace)
subroutine, public lda_u_dump(restart, namespace, this, st, mesh, ierr)
subroutine, public lda_u_write_u(this, iunit, namespace)
subroutine, public lda_u_load(restart, this, st, dftu_energy, ierr, occ_only, u_only)
subroutine lda_u_dump_coulomb_integrals(this, namespace, restart, st, mesh, ierr)
subroutine, public lda_u_write_magnetization(dir, this, ions, mesh, st, namespace)
subroutine, public lda_u_write_effectiveu(dir, this, namespace)
subroutine, public lda_u_write_v(this, iunit, namespace)
subroutine, public lda_u_get_effectiveu(this, Ueff)
subroutine, public lda_u_set_effectiveu(this, Ueff)
subroutine, public dlda_u_get_occupations(this, occ)
subroutine, public zcompute_dftu_energy(this, energy, st)
This routine computes the value of the double counting term in the DFT+U energy.
integer, parameter, public dft_u_empirical
subroutine, public compute_acbno_u_kanamori(this, st, kanamori)
subroutine, public zlda_u_update_potential(this, st)
This routine computes the potential that, once multiplied by the projector Pmm' and summed over m and...
subroutine, public dlda_u_set_occupations(this, occ)
subroutine, public lda_u_get_effectivev(this, Veff)
subroutine, public lda_u_set_effectivev(this, Veff)
integer, parameter, public dft_u_acbn0
subroutine, public zlda_u_get_occupations(this, occ)
subroutine, public dcompute_dftu_energy(this, energy, st)
This routine computes the value of the double counting term in the DFT+U energy.
subroutine, public dlda_u_update_potential(this, st)
This routine computes the potential that, once multiplied by the projector Pmm' and summed over m and...
subroutine, public zlda_u_set_occupations(this, occ)
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
type(mpi_grp_t), public mpi_world
This module handles the communicators for the various parallelization strategies.
subroutine, public restart_close(restart, iunit)
Close a file previously opened with restart_open.
subroutine, public restart_write(restart, iunit, lines, nlines, ierr)
logical pure function, public restart_skip(restart)
Returns true if the restart information should neither be read nor written. This might happen because...
integer function, public restart_open(restart, filename, status, position, silent)
Open file "filename" found inside the current restart directory. Depending on the type of restart,...
pure logical function, public states_are_real(st)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
character(len=20) pure function, public units_abbrev(this)
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
Class to describe DFT+U parameters.
Describes mesh distribution to nodes.
The states_elec_t class contains all electronic wave functions.