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
74 if (st%system_grp%is_root())
then
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
177 integer :: iunit, ios
181 if (st%system_grp%is_root())
then
182 iunit =
io_open(trim(dir) //
"/effectiveU", namespace, action=
'write')
185 write(iunit,
'(a,a,a,f7.3,a)')
'Hubbard U [', &
187 write(iunit,
'(a,6x,14x,a)')
' Orbital',
'U'
188 do ios = 1, this%norbsets
189 if (.not. this%basisfromstates)
then
190 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals)
then
191 if (this%orbsets(ios)%nn /= 0)
then
192 write(iunit,
'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
193 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
196 write(iunit,
'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
201 if (this%orbsets(ios)%nn /= 0)
then
202 write(iunit,
'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
203 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
204 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
207 write(iunit,
'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
209 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
219 write(iunit,
'(a,a,a,f7.3,a)')
'Hund J [', &
221 write(iunit,
'(a,6x,14x,a)')
' Orbital',
'J'
222 do ios = 1, this%norbsets
223 if (.not. this%basisfromstates)
then
224 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals)
then
225 if (this%orbsets(ios)%nn /= 0)
then
226 write(iunit,
'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
227 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
230 write(iunit,
'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
235 if (this%orbsets(ios)%nn /= 0)
then
236 write(iunit,
'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
237 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
238 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
241 write(iunit,
'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
243 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
260 type(
lda_u_t),
intent(in) :: this
262 character(len=*),
intent(in) :: dir
265 integer :: iunit, ios
266 real(real64),
allocatable :: kanamori(:,:)
270 if (st%system_grp%is_root())
then
271 safe_allocate(kanamori(1:3,1:this%norbsets))
275 iunit =
io_open(trim(dir) //
"/kanamoriU", namespace, action=
'write')
277 write(iunit,
'(a,a,a,f7.3,a)')
'Intraorbital U [', &
279 write(iunit,
'(a,6x,14x,a)')
' Orbital',
'U'
280 do ios = 1, this%norbsets
281 if (.not. this%basisfromstates)
then
282 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals)
then
283 if (this%orbsets(ios)%nn /= 0)
then
284 write(iunit,
'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
285 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
288 write(iunit,
'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
293 if (this%orbsets(ios)%nn /= 0)
then
294 write(iunit,
'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
295 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
296 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
299 write(iunit,
'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
301 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
311 write(iunit,
'(a,a,a,f7.3,a)')
'Interorbital Up [', &
313 write(iunit,
'(a,6x,14x,a)')
' Orbital',
'Up'
314 do ios = 1, this%norbsets
315 if (.not. this%basisfromstates)
then
316 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals)
then
317 if (this%orbsets(ios)%nn /= 0)
then
318 write(iunit,
'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
319 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
322 write(iunit,
'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
327 if (this%orbsets(ios)%nn /= 0)
then
328 write(iunit,
'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
329 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
330 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
333 write(iunit,
'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
335 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
344 write(iunit,
'(a,a,a,f7.3,a)')
'Hund J [', &
346 write(iunit,
'(a,6x,14x,a)')
' Orbital',
'J'
347 do ios = 1, this%norbsets
348 if (.not. this%basisfromstates)
then
349 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals)
then
350 if (this%orbsets(ios)%nn /= 0)
then
351 write(iunit,
'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
352 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
355 write(iunit,
'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
360 if (this%orbsets(ios)%nn /= 0)
then
361 write(iunit,
'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
362 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
363 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
366 write(iunit,
'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
368 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
380 safe_deallocate_a(kanamori)
391 type(
lda_u_t),
intent(in) :: this
392 character(len=*),
intent(in) :: dir
393 type(
ions_t),
intent(in) :: ions
394 class(
mesh_t),
intent(in) :: mesh
398 integer :: iunit, ia, ios, im
399 real(real64),
allocatable :: mm(:,:)
401 if (.not. st%system_grp%is_root())
return
406 iunit =
io_open(trim(dir)//
"/magnetization.xsf", namespace, action=
'write', position=
'asis')
408 if (this%nspins > 1)
then
409 safe_allocate(mm(1:mesh%box%dim, 1:ions%natoms))
412 do ios = 1, this%norbsets
413 ia = this%orbsets(ios)%iatom
414 do im = 1, this%orbsets(ios)%norbs
416 mm(3, ia) = mm(3, ia) + this%dn(im,im,1,ios) - this%dn(im,im,2,ios)
418 mm(3, ia) = mm(3, ia) + real(this%zn(im,im,1,ios) - this%zn(im,im,2,ios), real64)
420 if (this%nspins /= this%spin_channels)
then
421 mm(1, ia) = mm(1, ia) + 2*real(this%zn(im,im,3,ios), real64)
422 mm(2, ia) = mm(2, ia) - 2*aimag(this%zn(im,im,3,ios))
427 call write_xsf_geometry(iunit, ions%space, ions%latt, ions%pos, ions%atom, mesh, forces = mm)
428 safe_deallocate_a(mm)
438 type(
lda_u_t),
intent(in) :: this
439 integer,
optional,
intent(in) :: iunit
440 type(
namespace_t),
optional,
intent(in) :: namespace
446 write(
message(1),
'(a,a,a,f7.3,a)')
'Effective Hubbard U [', &
448 write(
message(2),
'(a,6x,14x,a)')
' Orbital',
'U'
451 do ios = 1, this%norbsets
452 if (.not. this%basisfromstates)
then
453 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals)
then
454 if (this%orbsets(ios)%nn /= 0)
then
455 write(
message(1),
'(i4,a10, 2x, i1, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
456 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
459 write(
message(1),
'(i4,a10, 3x, a1, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
464 if (this%orbsets(ios)%nn /= 0)
then
465 write(
message(1),
'(i4,a10, 2x, i1, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
466 this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
467 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
470 write(
message(1),
'(i4,a10, 3x, a1, i1, a2, f15.6)') ios, trim(this%orbsets(ios)%spec%get_label()), &
472 int(
m_two*(this%orbsets(ios)%jj)),
'/2', &
487 type(
lda_u_t),
intent(in) :: this
488 integer,
optional,
intent(in) :: iunit
489 type(
namespace_t),
optional,
intent(in) :: namespace
491 integer :: ios, icopies, ios2
493 if (.not. this%intersite)
return
497 write(
message(1),
'(a,a,a,f7.3,a)')
'Effective intersite V [', &
499 write(
message(2),
'(a,14x,a)')
' Orbital',
'V'
502 do ios = 1, this%norbsets
503 do icopies = 1, this%orbsets(ios)%nneighbors
504 ios2 = this%orbsets(ios)%map_os(icopies)
505 if (.not.this%basisfromstates)
then
506 if (this%orbsets(ios)%ndim == 1 .or. this%basis%combine_j_orbitals)
then
507 if (this%orbsets(ios)%nn /= 0)
then
508 write(
message(1),
'(i4,a10, 2x, i1, a1, i4, 1x, i1, a1, f7.3, f15.6)') ios, &
509 trim(this%orbsets(ios)%spec%get_label()), this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), ios2, &
510 this%orbsets(ios2)%nn,
l_notation(this%orbsets(ios2)%ll), &
515 write(
message(1),
'(i4,a10, 3x, a1, i4, 1x, a1, f7.3, f15.6)') ios, &
516 trim(this%orbsets(ios)%spec%get_label()),
l_notation(this%orbsets(ios)%ll), ios2, &
522 if (this%orbsets(ios)%nn /= 0)
then
523 write(
message(1),
'(i4,a10, 2x, i1, a1, i1, a2, i4, f7.3, f15.6)') ios, &
524 trim(this%orbsets(ios)%spec%get_label()), this%orbsets(ios)%nn,
l_notation(this%orbsets(ios)%ll), &
525 int(
m_two*(this%orbsets(ios)%jj)),
'/2', ios2, &
530 write(
message(1),
'(i4,a10, 3x, a1, i1, a2, i4, f7.3, f15.6)') ios, &
531 trim(this%orbsets(ios)%spec%get_label()),
l_notation(this%orbsets(ios)%ll), &
532 int(
m_two*(this%orbsets(ios)%jj)),
'/2', ios2, &
539 write(
message(1),
'(i4,a10, i4, f7.3, f15.6)') ios,
'states', ios2, &
553 subroutine lda_u_dump(restart, namespace, this, st, mesh, ierr)
556 type(
lda_u_t),
intent(in) :: this
558 class(
mesh_t),
intent(in) :: mesh
559 integer,
intent(out) :: ierr
561 integer :: err, occsize, ios, ncount
562 real(real64),
allocatable :: ueff(:), docc(:), veff(:)
563 complex(real64),
allocatable :: zocc(:)
569 if (restart%skip())
then
574 message(1) =
"Debug: Writing DFT+U restart."
577 occsize = this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets
580 if (this%intersite)
then
581 occsize = occsize + 2*this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets*this%maxneighbors
587 safe_allocate(docc(1:occsize))
590 call restart%write_binary(
"lda_u_occ", occsize, docc, err)
591 if (err /= 0) ierr = ierr + 1
592 safe_deallocate_a(docc)
594 safe_allocate(zocc(1:occsize))
597 call restart%write_binary(
"lda_u_occ", occsize, zocc, err)
598 if (err /= 0) ierr = ierr + 1
599 safe_deallocate_a(zocc)
604 safe_allocate(ueff(1:this%norbsets))
607 call restart%write_binary(
"lda_u_Ueff", this%norbsets, ueff, err)
608 safe_deallocate_a(ueff)
609 if (err /= 0) ierr = ierr + 1
611 if (this%intersite .and. this%maxneighbors > 0)
then
613 do ios = 1, this%norbsets
614 ncount = ncount + this%orbsets(ios)%nneighbors
616 safe_allocate(veff(1:ncount))
618 call restart%write_binary(
"lda_u_Veff", ncount, veff, err)
619 safe_deallocate_a(veff)
620 if (err /= 0) ierr = ierr + 1
626 message(1) =
"Debug: Writing DFT+U restart done."
634 subroutine lda_u_load(restart, this, st, dftu_energy, ierr, occ_only, u_only)
636 type(
lda_u_t),
intent(inout) :: this
638 real(real64),
intent(out) :: dftu_energy
639 integer,
intent(out) :: ierr
640 logical,
optional,
intent(in) :: occ_only
641 logical,
optional,
intent(in) :: u_only
643 integer :: err, occsize, ncount, ios
644 real(real64),
allocatable :: ueff(:), docc(:), veff(:)
645 complex(real64),
allocatable :: zocc(:)
651 if (restart%skip())
then
656 message(1) =
"Debug: Reading DFT+U restart."
661 safe_allocate(ueff(1:this%norbsets))
662 call restart%read_binary(
"lda_u_Ueff", this%norbsets, ueff, err)
668 safe_deallocate_a(ueff)
670 if (this%intersite .and. this%maxneighbors > 0)
then
672 do ios = 1, this%norbsets
673 ncount = ncount + this%orbsets(ios)%nneighbors
675 safe_allocate(veff(1:ncount))
676 call restart%read_binary(
"lda_u_Veff", ncount, veff, err)
682 safe_deallocate_a(veff)
688 occsize = this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets
691 if (this%intersite)
then
692 occsize = occsize + 2*this%maxnorbs*this%maxnorbs*this%nspins*this%norbsets*this%maxneighbors
698 safe_allocate(docc(1:occsize))
699 call restart%read_binary(
"lda_u_occ", occsize, docc, err)
705 safe_deallocate_a(docc)
707 safe_allocate(zocc(1:occsize))
708 call restart%read_binary(
"lda_u_occ", occsize, zocc, err)
714 safe_deallocate_a(zocc)
726 message(1) =
"Debug: Reading DFT+U restart done."
734 type(
lda_u_t),
intent(in) :: this
737 class(
mesh_t),
intent(in) :: mesh
738 integer,
intent(out) :: ierr
740 integer :: coulomb_int_file, idim1, idim2, err
741 integer :: ios, im, imp, impp, imppp
742 character(len=256) :: lines(3)
743 logical :: complex_coulomb_integrals
754 message(1) =
"Debug: Writing Coulomb integrals restart."
757 complex_coulomb_integrals = .false.
758 do ios = 1, this%norbsets
759 if (this%orbsets(ios)%ndim > 1) complex_coulomb_integrals = .
true.
762 coulomb_int_file = restart%open(
'coulomb_integrals')
767 write(lines(1),
'(a20,i21)')
"norb=", this%norbsets
768 write(lines(2),
'(a20,i21)')
"dim=", st%d%dim
769 write(lines(3),
'(a20,i21)')
"checksum=", mesh%idx%checksum
770 call restart%write(coulomb_int_file, lines, 3, err)
771 if (err /= 0) ierr = ierr - 2
773 do ios = 1, this%norbsets
774 do im = 1, this%orbsets(ios)%norbs
775 do imp = 1, this%orbsets(ios)%norbs
776 do impp = 1, this%orbsets(ios)%norbs
777 do imppp = 1, this%orbsets(ios)%norbs
778 if(.not. complex_coulomb_integrals)
then
779 write(lines(1),
'(i4,i4,i4,i4,i4,e20.12)') ios, im, imp, impp, imppp, this%coulomb(im, imp, impp, imppp, ios)
781 do idim1 = 1, st%d%dim
782 do idim2 = 1, st%d%dim
783 write(lines(1),
'(i4,i4,i4,i4,i4,i4,i4,2e20.12)') ios, im, imp, impp, imppp, idim1, idim2, &
784 real(this%zcoulomb(im, imp, impp, imppp, idim1, idim2, ios), real64) , &
785 aimag(this%zcoulomb(im, imp, impp, imppp, idim1, idim2, ios))
789 call restart%write(coulomb_int_file, lines, 1, err)
799 call restart%close(coulomb_int_file)
801 message(1) =
"Debug: Writing Coulomb integrals restart done."
character(len=1), dimension(0:3), parameter, public l_notation
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, restart, st, mesh, ierr)
subroutine, public lda_u_write_effectiveu(dir, this, st, namespace)
subroutine, public lda_u_write_magnetization(dir, this, ions, mesh, st, 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.
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
This module handles the communicators for the various parallelization strategies.
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.