58 integer :: in_file, ii, jj, kk, ierr, ip_h, irow, ifreq, nrow, it
59 integer :: ik, ist, uist, istep, ikpoint, irep, out_file, iop, idim
60 integer :: time_steps, energy_steps, istart, iend, ntiter, Nreplica, Ntrans
61 real(real64) :: dt, tt, weight, kpoint(3), kpoint_sym(3), kred(3), kred_sym(3)
62 real(real64) :: xx_h_sym(3)
63 integer :: irep_h, ip_h_sym, rankmin
64 real(real64) :: start_time, dmin
65 real(real64),
allocatable :: Et(:), ftreal(:, :, :), ftimag(:, :, :), tmp(:), omega(:)
66 complex(real64),
allocatable :: Xiak(:,:,:), Yiak(:,:,:)
67 real(real64),
allocatable :: proj_r(:,:,:,:), proj_i(:,:,:,:)
68 real(real64),
allocatable :: proj_r_corr(:,:), proj_i_corr(:,:), centers(:,:)
69 complex(real64),
allocatable :: tdm(:,:), tdm_1D(:,:,:,:)
70 complex(real64),
allocatable,
target :: psi(:,:), upsi(:,:)
71 complex(real64),
allocatable :: phase(:,:,:), ftcmplx(:,:)
72 complex(real64),
pointer :: psi_sym(:,:), upsi_sym(:,:)
73 type(spectrum_t) :: spectrum
74 type(electrons_t),
pointer :: sys
75 type(batch_t) :: projb_r, projb_i, ftrealb, ftimagb
76 character(len=MAX_PATH_LEN) :: fname
77 type(states_elec_t),
pointer :: st
78 type(states_elec_t) :: gs_st
79 type(restart_t) :: restart
80 type(unit_t) :: fn_unit
81 integer :: kpt_start, kpt_end, supercell(3), nomega, ncols
83 real(real64) :: pos_h(3), norm
107 if(sys%st%d%ispin ==
spinors)
then
111 if(st%parallel_in_states)
then
115 if(sys%gr%parallel_in_domains)
then
132 safe_allocate(omega(1:nrow))
140 message(1) =
"oct-tdtdm: TDTDMFrequencies must be defined."
146 if(any(omega > spectrum%max_energy))
then
147 message(1) =
"One requested frequecy is larger than PropagationSpectrumMaxEnergy."
148 message(2) =
"Please increase the value of PropagationSpectrumMaxEnergy."
151 if(any(omega > -spectrum%min_energy))
then
152 message(1) =
"One requested frequency is larger than -PropagationSpectrumMinEnergy."
153 message(2) =
"Please decrease the value of PropagationSpectrumMinEnergy."
160 safe_deallocate_a(gs_st%node)
165 message(1) =
"oct-tdtdm: Unable to read states information."
170 safe_allocate(gs_st%occ(1:gs_st%nst, 1:gs_st%nik))
171 safe_allocate(gs_st%eigenval(1:gs_st%nst, 1:gs_st%nik))
175 safe_allocate(gs_st%node(1:gs_st%nst))
180 kpt_start = gs_st%d%kpt%start
181 kpt_end = gs_st%d%kpt%end
185 gs_st%eigenval = huge(gs_st%eigenval)
187 if(gs_st%d%ispin ==
spinors)
then
188 safe_deallocate_a(gs_st%spin)
189 safe_allocate(gs_st%spin(1:3, 1:gs_st%nst, 1:gs_st%nik))
194 if(ierr /= 0 .and. ierr /= (gs_st%st_end-gs_st%st_start+1)*(kpt_end-kpt_start+1)*gs_st%d%dim)
then
195 message(1) =
"oct-tdtdm: Unable to read wavefunctions for TDOutput."
201 in_file =
io_open(
'td.general/projections', action=
'read', status=
'old')
207 safe_allocate(tmp(1:st%nst*gs_st%nst*st%nik*2))
208 safe_allocate(proj_r(1:time_steps, 1:gs_st%nst, 1:st%nst, 1:st%nik))
209 safe_allocate(proj_i(1:time_steps, 1:gs_st%nst, 1:st%nst, 1:st%nik))
214 do ii = 1, time_steps
215 read(in_file, *) jj, tt, (tmp(kk), kk = 1, st%nst*gs_st%nst*st%nik*2)
218 do uist = 1, gs_st%nst
219 jj = (ik-1)*st%nst*gs_st%nst + (ist-1)*gs_st%nst + uist
220 proj_r(ii, uist, ist, ik) = tmp((jj-1)*2+1)
223 proj_i(ii, uist, ist, ik) = -tmp((jj-1)*2+2)
228 safe_deallocate_a(tmp)
232 write(
message(1),
'(a, i7, a)')
"oct-tdtdm: Read ", time_steps,
" steps from file '"// &
236 start_time = spectrum%start_time
243 safe_allocate(proj_r_corr(1:time_steps, 1:gs_st%nst*st%nst*(kpt_end-kpt_start+1)))
244 safe_allocate(proj_i_corr(1:time_steps, 1:gs_st%nst*st%nst*(kpt_end-kpt_start+1)))
247 do ik = kpt_start, kpt_end
249 do uist = ist+1, gs_st%nst
250 jj = (ik-kpt_start)*st%nst*gs_st%nst+(ist-1)*gs_st%nst+uist
251 do ii = 1, time_steps
252 norm =
hypot(proj_r(ii, ist, ist, ik),proj_i(ii, ist, ist, ik))
254 proj_r_corr(ii, jj) = (proj_r(ii, uist, ist, ik) * proj_r(ii, ist, ist, ik) &
255 + proj_i(ii, uist, ist, ik) * proj_i(ii, ist, ist, ik))/norm
256 proj_i_corr(ii, jj) =(-proj_r(ii, uist, ist, ik) * proj_i(ii, ist, ist, ik) &
257 + proj_i(ii, uist, ist, ik) * proj_r(ii, ist, ist, ik))/norm
259 proj_r_corr(ii, jj) =
m_zero
260 proj_i_corr(ii, jj) =
m_zero
267 safe_deallocate_a(proj_r)
268 safe_deallocate_a(proj_i)
272 istart = max(1, istart)
275 safe_allocate(ftreal(1:energy_steps, 1:st%nst*gs_st%nst*(kpt_end-kpt_start+1), 1:2))
276 safe_allocate(ftimag(1:energy_steps, 1:st%nst*gs_st%nst*(kpt_end-kpt_start+1), 1:2))
278 call batch_init(projb_r, 1, 1, st%nst*gs_st%nst*(kpt_end-kpt_start+1), proj_r_corr)
279 call batch_init(projb_i, 1, 1, st%nst*gs_st%nst*(kpt_end-kpt_start+1), proj_i_corr)
280 call batch_init(ftrealb, 1, 1, st%nst*gs_st%nst*(kpt_end-kpt_start+1), ftreal(:,:,1))
281 call batch_init(ftimagb, 1, 1, st%nst*gs_st%nst*(kpt_end-kpt_start+1), ftimag(:,:,1))
283 write(
message(1),
'(a)')
"oct-tdtdm: Fourier transforming real part of the projections"
287 istart, iend, spectrum%start_time, dt, projb_r, spectrum%min_energy, spectrum%max_energy, spectrum%energy_step, ftrealb)
290 istart, iend, spectrum%start_time, dt, projb_r, spectrum%min_energy, spectrum%max_energy, spectrum%energy_step, ftimagb)
295 safe_allocate(ftcmplx(1:energy_steps, 1:st%nst*gs_st%nst*(kpt_end-kpt_start+1)))
296 do ii = 1, st%nst*gs_st%nst*(kpt_end-kpt_start+1)
297 ftcmplx(1:energy_steps,ii) = cmplx(ftreal(1:energy_steps,ii,1), ftimag(1:energy_steps,ii,1), real64)
300 write(
message(1),
'(a)')
"oct-tdtdm: Fourier transforming imaginary part of the projections"
303 call batch_init(ftrealb, 1, 1, st%nst*gs_st%nst*(kpt_end-kpt_start+1), ftreal(:,:,2))
304 call batch_init(ftimagb, 1, 1, st%nst*gs_st%nst*(kpt_end-kpt_start+1), ftimag(:,:,2))
307 istart, iend, spectrum%start_time, dt, projb_i, spectrum%min_energy, spectrum%max_energy, spectrum%energy_step, ftrealb)
310 istart, iend, spectrum%start_time, dt, projb_i, spectrum%min_energy, spectrum%max_energy, spectrum%energy_step, ftimagb)
316 safe_deallocate_a(proj_r_corr)
317 safe_deallocate_a(proj_i_corr)
319 do ii = 1, st%nst*gs_st%nst*(kpt_end-kpt_start+1)
320 ftcmplx(1:energy_steps,ii) = ftcmplx(1:energy_steps,ii) +
m_zi*ftreal(1:energy_steps,ii,2) - ftimag(1:energy_steps,ii,2)
323 safe_deallocate_a(ftreal)
324 safe_deallocate_a(ftimag)
326 write(
message(1),
'(a)')
"oct-tdtdm: Constructing the two-particle wavefunctions."
338 if (
parse_block(sys%namespace,
'SupercellDimensions', blk) == 0)
then
340 if (ncols /= sys%space%dim)
then
341 write(
message(1),
'(a,i3,a,i3)')
'SupercellDimensions has ', ncols,
' columns but must have ', sys%space%dim
344 do ii = 1, sys%space%dim
351 supercell(1:sys%space%dim) = sys%kpoints%nik_axis(1:sys%space%dim)
354 nreplica = product(supercell(1:sys%space%dim))
357 safe_allocate(centers(1:sys%space%dim, 1:nreplica))
359 do ii = 0, supercell(1)-1
360 do jj = 0, supercell(2)-1
361 do kk = 0, supercell(3)-1
362 centers(1, irep) = -
floor((supercell(1)-1)/
m_two)+ii
363 centers(2, irep) = -
floor((supercell(2)-1)/
m_two)+jj
364 centers(3, irep) = -
floor((supercell(3)-1)/
m_two)+kk
365 centers(:, irep) = matmul(sys%ions%latt%rlattice, centers(:, irep))
373 do ik = kpt_start, kpt_end
374 ikpoint = gs_st%d%get_kpoint_index(ik)
377 safe_allocate(phase(kpt_start:kpt_end, 1:irep, 1:nreplica))
378 do irep = 1, nreplica
379 do ik = kpt_start, kpt_end
380 ikpoint = gs_st%d%get_kpoint_index(ik)
381 kpoint(1:sys%space%dim) = sys%kpoints%get_point(ikpoint)
385 if (sys%kpoints%use_symmetries)
then
392 phase(ik, ii, irep) =
exp(-
m_zi*sum(kpoint_sym(1:sys%space%dim)*centers(:, irep)))
399 if(sys%space%dim > 1)
then
405 do ist = 1, gs_st%nst
406 if(abs(gs_st%occ(ist, 1)) <
m_min_occ) cycle
408 do uist = 1, gs_st%nst
409 if(abs(gs_st%occ(uist, 1)) >
m_min_occ) cycle
410 weight = gs_st%kweights(1) * (gs_st%occ(ist, 1)-gs_st%occ(uist, 1))
416 write(
message(1),
'(a)')
"oct-tdtdm: No transition found."
417 write(
message(2),
'(a)')
"Please check that unoccupied states are included in the ground state calculation."
421 safe_allocate(xiak(1:st%nst, 1:gs_st%nst, 1:st%nik))
422 safe_allocate(yiak(1:st%nst, 1:gs_st%nst, 1:st%nik))
423 safe_allocate(et(1:ntrans*st%nik))
424 safe_allocate(psi(1:sys%gr%np, 1:gs_st%d%dim))
425 safe_allocate(upsi(1:sys%gr%np, 1:gs_st%d%dim))
427 if(sys%kpoints%use_symmetries)
then
428 safe_allocate(psi_sym(1:sys%gr%np, 1:st%d%dim))
429 safe_allocate(upsi_sym(1:sys%gr%np, 1:st%d%dim))
432 select case(sys%space%dim)
434 safe_allocate(tdm(1:sys%gr%np, 1:nreplica))
436 safe_allocate(tdm_1d(1:sys%gr%np, 1:sys%gr%np, 1:nreplica, 1:nreplica))
441 write(
message(1),
'(a, f6.4, a)')
"oct-tdtdm: Constructing the two-particle wavefunction at ", omega(ifreq),
" Ha."
444 select case(sys%space%dim)
456 it = (kpt_start-1)*ntrans + 1
458 do ik = kpt_start, kpt_end
459 ikpoint = st%d%get_kpoint_index(ik)
462 if(abs(gs_st%occ(ist, ik)) <
m_min_occ) cycle
465 if (sys%hm%phase%is_allocated())
then
466 call sys%hm%phase%apply_to_single(psi, sys%gr%np, gs_st%d%dim, ik, .false.)
469 do uist = 1, gs_st%nst
470 if(abs(gs_st%occ(uist, ik)) >
m_min_occ) cycle
475 jj = (ik-kpt_start)*st%nst*gs_st%nst+(ist-1)*gs_st%nst+uist
476 istep = int((+omega(ifreq)-spectrum%min_energy)/spectrum%energy_step)
477 xiak(ist, uist, ik) = conjg(ftcmplx(istep, jj))
478 istep = int((+omega(ifreq)-spectrum%min_energy)/spectrum%energy_step)
479 yiak(ist, uist, ik) = ftcmplx(istep, jj)
482 weight = gs_st%kweights(ik) * (gs_st%occ(ist, ik)-gs_st%occ(uist, ik)) &
487 if(sys%hm%phase%is_allocated())
then
488 call sys%hm%phase%apply_to_single(upsi, sys%gr%np, st%d%dim, ik, .false.)
494 if(sys%kpoints%use_symmetries)
then
495 do idim = 1, st%d%dim
502 xx_h_sym = sys%ions%latt%fold_into_cell(xx_h_sym)
504 assert(.not.sys%gr%parallel_in_domains)
515 select case(sys%space%dim)
517 do irep = 1, nreplica
518 call lalg_axpy(sys%gr%np, phase(ik, ii, irep) * weight &
519 * conjg(xiak(ist,uist,ik))*conjg(psi_sym(ip_h_sym,1)), upsi_sym(:, 1), tdm(:,irep))
520 call lalg_axpy(sys%gr%np, phase(ik, ii, irep) * weight &
521 * yiak(ist,uist,ik)*conjg(upsi_sym(ip_h_sym,1)), psi_sym(:, 1), tdm(:,irep))
525 do irep_h = 1, nreplica
526 do irep = 1, nreplica
527 do ip_h = 1, sys%gr%np
528 call lalg_axpy(sys%gr%np, phase(ik, ii, irep) * conjg(phase(ik, ii, irep_h)) &
529 * weight * conjg(xiak(ist,uist,ik)) * conjg(psi_sym(ip_h,1)), &
530 upsi_sym(:, 1), tdm_1d(:, ip_h, irep, irep_h))
531 call lalg_axpy(sys%gr%np, phase(ik, ii, irep) * conjg(phase(ik, ii, irep_h)) &
532 * weight * conjg(yiak(ist,uist,ik)) * conjg(upsi_sym(ip_h,1)), &
533 psi_sym(:, 1), tdm_1d(:, ip_h, irep, irep_h))
541 et(it) = gs_st%eigenval(uist, ik) - gs_st%eigenval(ist, ik)
547 if(gs_st%d%kpt%parallel)
then
548 if(sys%space%dim > 1)
then
564 safe_deallocate_a(et)
565 safe_deallocate_a(xiak)
566 safe_deallocate_a(yiak)
567 safe_deallocate_a(tdm)
568 safe_deallocate_a(tdm_1d)
570 safe_deallocate_a(psi)
571 safe_deallocate_a(upsi)
572 if(sys%kpoints%use_symmetries)
then
573 safe_deallocate_p(psi_sym)
574 safe_deallocate_p(upsi_sym)
576 safe_deallocate_a(ftcmplx)
577 safe_deallocate_a(centers)
578 safe_deallocate_a(phase)
579 safe_deallocate_a(omega)
581 safe_deallocate_p(sys)
597 real(real64),
intent(out) :: xx_h(1:sys%space%dim)
598 integer,
intent(out) :: ip_h
601 integer :: idir, rankmin
626 do idir = 1, sys%space%dim
642 do idir = 1, sys%space%dim
646 xx_h = sys%ions%latt%red_to_cart(xx_h)
648 xx_h(1:sys%space%dim) = sys%ions%pos(1:sys%space%dim, 1)
653 xx_h = sys%ions%latt%fold_into_cell(xx_h)
656 assert(.not.sys%gr%parallel_in_domains)
658 write(
message(1),
'(a, 3(1x,f7.4,a))')
"oct-tdtdm: Requesting the hole at (", xx_h(1), &
659 ",", xx_h(2),
",", xx_h(3),
")."
660 call mesh_r(sys%gr, ip_h, dmin, coords=xx_h)
661 write(
message(2),
'(a, 3(1x,f7.4,a))')
"oct-tdtdm: Setting the hole at (", xx_h(1), &
662 ",", xx_h(2),
",", xx_h(3),
")."
670 real(real64),
allocatable :: den(:,:), den_1d(:,:,:,:)
671 real(real64) :: norm, xx(3), xx_h(3)
677 select case(sys%space%dim)
679 safe_allocate(den(1:sys%gr%np, 1:nreplica))
680 do irep = 1, nreplica
682 den(ii, irep) = real(tdm(ii, irep)*conjg(tdm(ii, irep)), real64)
691 safe_allocate(den_1d(1:sys%gr%np, 1:sys%gr%np, 1:nreplica, 1:nreplica))
692 do irep_h = 1, nreplica
693 do irep = 1, nreplica
694 do ip_h = 1, sys%gr%np
696 tdm_1d(ii, ip_h, irep, irep_h) = conjg(tdm_1d(ii, ip_h, irep, irep_h))
697 den_1d(ii, ip_h, irep, irep_h) = real(tdm_1d(ii, ip_h, irep, irep_h)*conjg(tdm_1d(ii,ip_h, irep, irep_h)), real64)
704 fn_unit =
units_out%length**(-sys%space%dim)
706 select case(sys%space%dim)
708 write(fname,
'(a, f0.4)')
'tdm_density-0', omega(ifreq)
710 sys%gr, sys%space, sys%ions%latt, den, centers, supercell, fn_unit, &
711 ierr,
global_namespace, pos=sys%ions%pos, atoms=sys%ions%atom, grp = st%dom_st_kpt_mpi_grp, extra_atom=pos_h)
714 sys%gr, sys%space, sys%ions%latt, den, centers, supercell, fn_unit, &
717 safe_deallocate_a(den)
724 write(fname,
'(a, f0.4)')
'tdm_density-0', omega(ifreq)
726 sys%gr, sys%space, sys%ions%latt, &
727 den_1d(:,ip_h,:,irep_h), centers, supercell, fn_unit, ierr,
global_namespace, &
728 grp = st%dom_st_kpt_mpi_grp)
730 write(fname,
'(a, f0.4)')
'tdm_wfn-0', omega(ifreq)
732 sys%gr, sys%space, sys%ions%latt, &
733 tdm_1d(:,ip_h,:,irep_h), centers, supercell, fn_unit, ierr,
global_namespace, &
734 grp = st%dom_st_kpt_mpi_grp)
736 assert(.not.sys%gr%parallel_in_domains)
738 write(fname,
'(a, f0.4)')
'td.general/tdm_density-0', omega(ifreq)
739 iunit =
io_open(fname, action=
'write')
740 write(iunit,
'(a)', iostat=ierr)
'# r_e r_h Re(\Psi(r_e,r_h)) Im(\Psi(r_e,r_h)) |\Psi(r_e,r_h)|^2'
742 do irep_h = 1, nreplica
743 do ip_h = 1, sys%gr%np
745 + centers(1:sys%space%dim, irep_h))
747 do irep = 1, nreplica
750 + centers(1:sys%space%dim, irep))
751 write(iunit,
'(5es23.14E3)', iostat=ierr) xx(1), xx_h(1), &
752 real(units_from_atomic(fn_unit, tdm_1D(ii, ip_h, irep, irep_h)), real64) ,&
753 aimag(units_from_atomic(fn_unit, tdm_1D(ii, ip_h, irep, irep_h))), &
754 units_from_atomic(fn_unit, den_1D(ii, ip_h, irep, irep_h))
761 safe_deallocate_a(den_1d)
769 real(real64),
allocatable :: weight(:,:)
775 safe_allocate(weight(1:st%nik, 1:gs_st%nst))
780 if(abs(gs_st%occ(ist, ik)) <
m_min_occ) cycle
782 do uist = ist+1, gs_st%nst
783 if(abs(gs_st%occ(uist, ik)) >
m_min_occ) cycle
785 weight(ik, ist) = weight(ik, ist) + abs(xiak(ist, uist, ik))**2
786 weight(ik, uist) = weight(ik,uist) + abs(yiak(ist, uist, ik))**2
791 write(fname,
'(a, f0.4)')
'td.general/tdm_weights-0', omega(ifreq)
792 out_file =
io_open(fname, action=
'write')
793 write(out_file,
'(a)')
'# ik - kx - ky - kz - sum weights - eigenval and weights(ist,ik) '
795 ikpoint = st%d%get_kpoint_index(ik)
796 kpoint(1:sys%space%dim) = sys%kpoints%reduced%point1BZ(1:sys%space%dim,ikpoint)
797 write(out_file,
'(i4,4e15.6)', advance=
'no') ik, kpoint(1:3), sum(weight(ik, 1:gs_st%nst))
798 do uist = 1, gs_st%nst-1
799 write(out_file,
'(2e15.6)', advance=
'no') gs_st%eigenval(uist, ik), weight(ik, uist)
801 write(out_file,
'(e15.6)') weight(ik, uist)
805 safe_deallocate_a(weight)
initialize a batch with existing memory
constant times a vector plus a vector
scales a vector by a constant
double hypot(double __x, double __y) __attribute__((__nothrow__
double exp(double __x) __attribute__((__nothrow__
double floor(double __x) __attribute__((__nothrow__
This module implements batches of mesh functions.
This module handles the calculation mode.
type(calc_mode_par_t), public calc_mode_par
Singleton instance of parallel calculation mode.
integer, parameter, public p_strategy_states
parallelization in states
integer, parameter, public spinors
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
subroutine, public fft_all_init(namespace)
initialize the table
subroutine, public fft_all_end()
delete all plans
real(real64), parameter, public m_two
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
real(real64), parameter, public m_zero
complex(real64), parameter, public m_z0
complex(real64), parameter, public m_zi
real(real64), parameter, public m_epsilon
subroutine, public global_init(communicator)
Initialise Octopus.
real(real64), parameter, public m_one
real(real64), parameter, public m_min_occ
Minimal occupation that is considered to be non-zero.
This module implements the underlying real-space grid.
subroutine, public zgrid_symmetrize_single(gr, iop, field, symm_field)
integer(int64) function, public io_function_fill_how(where)
Use this function to quickly plot functions for debugging purposes: call dio_function_output(io_funct...
subroutine, public io_init(defaults)
If the argument defaults is present and set to true, then the routine will not try to read anything f...
subroutine, public io_close(iunit, grp)
subroutine, public io_skip_header(iunit)
subroutine, public io_end()
character(len=max_path_len) function, public io_workpath(path, namespace)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
integer pure function, public kpoints_get_num_symmetry_ops(this, ik)
integer pure function, public kpoints_get_symmetry_ops(this, ik, index)
subroutine, public kpoints_to_reduced(latt, kin, kout)
subroutine, public kpoints_to_absolute(latt, kin, kout)
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
integer function, public mesh_nearest_point(mesh, pos, dmin, rankmin)
Returns the index of the point which is nearest to a given vector position pos.
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
real(real64) function, dimension(1:mesh%box%dim), public mesh_x_global(mesh, ipg)
subroutine, public messages_end()
subroutine, public messages_not_implemented(feature, namespace)
subroutine, public messages_init(output_dir)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_input_error(namespace, var, details, row, column)
subroutine, public messages_experimental(name, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
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.
type(namespace_t), public global_namespace
logical function, public parse_is_defined(namespace, name)
subroutine, public parser_init()
Initialise the Octopus parser.
subroutine, public parser_end()
End the Octopus parser.
integer function, public parse_block(namespace, name, blk, check_varinfo_)
subroutine, public profiling_end(namespace)
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
subroutine, public restart_module_init(namespace)
subroutine, public restart_init(restart, namespace, data_type, type, mc, ierr, mesh, dir, exact)
Initializes a restart object.
integer, parameter, public restart_proj
integer, parameter, public restart_type_load
subroutine, public restart_end(restart)
subroutine, public spectrum_fix_time_limits(spectrum, time_steps, dt, istart, iend, ntiter)
subroutine, public spectrum_fourier_transform(method, transform, noise, time_start, time_end, t0, time_step, time_function, energy_start, energy_end, energy_step, energy_function)
Computes the sine, cosine, (or "exponential") Fourier transform of the real function given in the tim...
subroutine, public spectrum_init(spectrum, namespace, default_energy_step, default_max_energy)
integer, parameter, public spectrum_transform_cos
integer, parameter, public spectrum_transform_sin
subroutine, public spectrum_count_time_steps(namespace, iunit, time_steps, dt)
pure integer function, public spectrum_nenergy_steps(spectrum)
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public states_elec_distribute_nodes(st, namespace, mc)
@Brief. Distribute states over the processes for states parallelization
subroutine, public states_elec_end(st)
finalize the states_elec_t object
subroutine, public states_elec_allocate_wfns(st, mesh, wfs_type, skip, packed)
Allocates the KS wavefunctions defined within a states_elec_t structure.
subroutine, public kpoints_distribute(this, mc)
distribute k-points over the nodes in the corresponding communicator
subroutine, public states_elec_copy(stout, stin, exclude_wfns, exclude_eigenval, special)
make a (selective) copy of a states_elec_t object
subroutine, public states_elec_look(restart, nik, dim, nst, ierr)
Reads the 'states' file in the restart directory, and finds out the nik, dim, and nst contained in it...
This module handles reading and writing restart information for the states_elec_t.
subroutine, public states_elec_load(restart, namespace, space, st, mesh, kpoints, fixed_occ, ierr, iter, lr, lowest_missing, label, verbose, skip)
returns in ierr: <0 => Fatal error, or nothing read =0 => read all wavefunctions >0 => could only rea...
subroutine, public symmetries_apply_kpoint_red(this, iop, aa, bb)
type(type_t), parameter, public type_cmplx
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
subroutine, public unit_system_init(namespace)
type(unit_system_t), public units_inp
the units systems for reading and writing
Class describing the electron system.
subroutine tdtdm_excitonic_weight()
subroutine tdtdm_output_density()
subroutine tdtdm_get_hole_position(xx_h, ip_h)