26 use,
intrinsic :: iso_fortran_env
41 integer,
parameter :: &
46 integer :: n_multipoles
47 integer,
allocatable :: l(:), m(:)
48 real(real64),
allocatable :: weight(:)
51 integer :: observable(2)
52 type(unit_system_t) :: units
53 real(real64),
allocatable :: ot(:)
56 real(real64) :: total_time
63 subroutine ft(omega, power)
64 real(real64),
intent(in) :: omega
65 real(real64),
intent(out) :: power
79 power = power + x*ot(j)
81 if (time_steps > 0)
then
82 power = power / real(time_steps, real64)
91 power = power + x*ot(j)
93 if (time_steps > 0)
then
94 power = power / real(time_steps, real64)
106 subroutine ft2(omega, power)
107 real(real64),
intent(in) :: omega
108 real(real64),
intent(out) :: power
118 case (sine_transform)
122 power = power + x*ot(j)
126 if (time_steps > 0)
then
127 power = - (power / real(time_steps, real64))**2
136 power = power + x*ot(j)
138 if (time_steps > 0)
then
139 power = - (power / real(time_steps, real64))**2
159 o%n_multipoles = i%n_multipoles
160 safe_allocate( o%l(1:o%n_multipoles))
161 safe_allocate( o%m(1:o%n_multipoles))
162 safe_allocate(o%weight(1:o%n_multipoles))
164 do j = 1, o%n_multipoles
167 o%weight(j) = i%weight(j)
174 subroutine read_resonances_file(order, ffile, namespace, search_interval, final_time, nfrequencies)
175 integer,
intent(inout) :: order
176 character(len=*),
intent(in) :: ffile
178 real(real64),
intent(inout) :: search_interval
179 real(real64),
intent(in) :: final_time
180 integer,
intent(in) :: nfrequencies
182 real(real64) :: dummy, leftbound, rightbound, w, power
183 integer :: iunit, nresonances, ios, i, j, k, npairs, nspin, order_in_file, nw_subtracted, ierr
184 logical :: file_exists
185 real(real64),
allocatable :: wij(:), omega(:), c0i(:)
190 write(
message(1),
'(a)')
'The run mode #3 is only compatible with the analysis of the'
191 write(
message(2),
'(a)')
'second-order response.'
196 inquire(file=
"ot", exist = file_exists)
197 if (.not. file_exists)
then
198 write(
message(1),
'(a)')
"Could not find 'ot' file."
205 write(
message(1),
'(a)')
"Could not retrieve units in the 'ot' file."
211 iunit =
io_open(trim(ffile), namespace, action=
'read', status=
'old')
217 read(iunit, *, iostat = ios) dummy, dummy
219 nresonances = nresonances + 1
222 npairs = (nresonances*(nresonances-1))/2
224 safe_allocate(omega(1:nresonances))
225 safe_allocate( c0i(1:nresonances))
226 safe_allocate(wij(1:npairs))
229 do i = 1, nresonances
230 read(iunit, *) omega(i), c0i(i)
234 do i = 1, nresonances
235 do j = i + 1, nresonances
236 wij(k) = omega(j) - omega(i)
241 if (search_interval >
m_zero)
then
247 call read_ot(namespace, nspin, order_in_file, nw_subtracted)
249 if (order_in_file /= order)
then
250 write(
message(1),
'(a)')
'The ot file should contain the second-order response in this run mode.'
254 if (final_time >
m_zero)
then
256 if (total_time > dt*time_steps)
then
257 total_time = dt*time_steps
258 write(
message(1),
'(a)')
'The requested total time to process is larger than the time available in the input file.'
259 write(
message(2),
'(a,f8.4,a)')
'The time has been adjusted to ', &
263 time_steps = int(total_time / dt)
264 total_time = time_steps * dt
266 total_time = dt*time_steps
270 leftbound = -search_interval
271 rightbound = search_interval
274 call modify_ot(time_steps, dt, order, ot, w, power)
275 nw_subtracted = nw_subtracted + 1
279 do i = 1, nresonances
280 do j = i + 1, nresonances
281 leftbound = wij(k) - search_interval
282 rightbound = wij(k) + search_interval
285 call modify_ot(time_steps, dt, order, ot, wij(k), power)
286 nw_subtracted = nw_subtracted + 1
291 safe_deallocate_a(wij)
292 safe_deallocate_a(c0i)
293 safe_deallocate_a(omega)
294 safe_deallocate_a(ot)
304 subroutine analyze_signal(order, omega, search_interval, final_time, nresonances, nfrequencies, &
306 integer,
intent(inout) :: order
307 real(real64),
intent(inout) :: omega
308 real(real64),
intent(inout) :: search_interval
309 real(real64),
intent(inout) :: final_time
310 integer,
intent(inout) :: nresonances
311 integer,
intent(inout) :: nfrequencies
312 real(real64),
intent(in) :: damping
315 real(real64) :: leftbound, rightbound, dw, power
316 real(real64),
allocatable :: w(:), c0I2(:)
317 integer :: nspin, i, ierr, order_in_file, nw_subtracted
318 logical :: file_exists
323 inquire(file=
"ot", exist = file_exists)
324 if (.not. file_exists)
then
325 write(
message(1),
'(a)')
"Could not find 'ot' file."
332 write(
message(1),
'(a)')
"Could not retrieve units in the 'ot' file."
342 if (search_interval >
m_zero)
then
348 leftbound = omega - search_interval
349 rightbound = omega + search_interval
351 call read_ot(namespace, nspin, order_in_file, nw_subtracted)
353 if (order_in_file /= order)
then
354 write(
message(1),
'(a)')
'Internal error in analyze_signal'
358 if (mod(order, 2) == 1)
then
359 mode = sine_transform
364 if (final_time >
m_zero)
then
366 if (total_time > dt*time_steps)
then
367 total_time = dt*time_steps
368 write(
message(1),
'(a)')
'The requested total time to process is larger than the time available in the input file.'
369 write(
message(2),
'(a,f8.4,a)')
'The time has been adjusted to ', &
373 time_steps = int(total_time / dt)
374 total_time = time_steps * dt
376 total_time = dt*time_steps
379 dw = (rightbound-leftbound) / (nfrequencies - 1)
381 safe_allocate( w(1:nresonances))
382 safe_allocate(c0i2(1:nresonances))
388 if (nw_subtracted >= nresonances)
exit
406 call modify_ot(time_steps, dt, order, ot, omega, power)
408 nw_subtracted = nw_subtracted + 1
417 safe_deallocate_a(ot)
419 safe_deallocate_a(c0i2)
433 integer,
intent(in) :: nfrequencies, nresonances
434 real(real64),
intent(in) :: dw
435 real(real64),
intent(in) :: w(nresonances), c0I2(nresonances)
436 real(real64),
intent(in) :: gamma
438 integer :: iunit, i, j
440 complex(real64) :: pol
444 iunit =
io_open(
'polarizability', namespace, action =
'write', status=
'replace', die=.false.)
445 write(iunit,
'(a)')
'# Polarizability file. Generated using the SOS formula with the following data:'
446 write(iunit,
'(a)')
'#'
448 do i = 1, nresonances
449 write(iunit,
'(a1,3e20.12)')
'#', w(i),
sqrt(abs(c0i2(i))), c0i2(i)
452 write(iunit,
'(a10,f12.6)')
'# Gamma = ', gamma
453 write(iunit,
'(a)')
'#'
455 do i = 1, nfrequencies
458 do j = 1, nresonances
459 pol = pol + c0i2(j) * (
m_one / (w(j) - e -
m_zi * gamma) +
m_one/(w(j) + e +
m_zi * gamma))
461 write(iunit,
'(3e20.12)') e, pol
473 subroutine find_resonance(omega, leftbound, rightbound, nfrequencies)
474 real(real64),
intent(inout) :: omega
475 real(real64),
intent(in) :: leftbound, rightbound
476 integer,
intent(in) :: nfrequencies
479 real(real64) :: dw, w, aw, min_aw, min_w, omega_orig, min_w1, min_w2
480 real(real64),
allocatable :: warray(:), tarray(:)
484 safe_allocate(warray(1:nfrequencies))
485 safe_allocate(tarray(1:nfrequencies))
489 dw = (rightbound-leftbound) / (nfrequencies - 1)
490 do i = 1, nfrequencies
491 w = leftbound + (i-1)*dw
499 do i = 1, nfrequencies
500 w = leftbound + (i-1)*dw
501 if (tarray(i) < min_aw)
then
509 min_w1 = min_w - 2*dw
510 min_w2 = min_w + 2*dw
514 write(
message(1),
'(a)')
'Could not find a maximum.'
516 write(
message(3),
'(a,f12.8,a,f12.8,a)')
' Search interval = [', &
518 write(
message(4),
'(a,f12.4,a)')
' Search discretization = ', &
523 safe_deallocate_a(warray)
524 safe_deallocate_a(tarray)
533 real(real64),
intent(in) :: omega
534 real(real64),
intent(out) :: power
535 integer,
intent(in) :: nw_subtracted
536 real(real64),
intent(in) :: dw, leftbound, rightbound
540 call ft(omega, power)
543 case (sine_transform)
549 write(
message(1),
'(a)')
'******************************************************************'
550 write(
message(2),
'(a,i3)')
'Resonance #', nw_subtracted + 1
557 write(
message(6),
'(a,f12.8)')
'f[O->I] = ',
m_two*omega*power
559 write(
message(8),
'(a,f12.8,a,f12.8,a)')
' Search interval = [',
units_from_atomic(
units_out%energy, leftbound),
',', &
563 write(
message(10),
'(a)')
'******************************************************************'
575 real(real64),
intent(in) :: omega
576 real(real64),
intent(out) :: power
577 integer,
intent(in) :: nw_subtracted
578 real(real64),
intent(in) :: leftbound, rightbound
579 real(real64),
intent(in) :: c01, c02
583 call ft(omega, power)
585 case (sine_transform)
592 power = power /
m_two
596 write(
message(1),
'(a)')
'******************************************************************'
597 write(
message(2),
'(a,i3)')
'Resonance #', nw_subtracted + 1
605 write(
message(1),
'(a,f12.8)')
' C(omega)/(C0i*C0j) = ', power / (c01 * c02)
610 write(
message(2),
'(a,f12.8,a,f12.8,a)')
' Search interval = [',
units_from_atomic(
units_out%energy, leftbound),
',', &
612 write(
message(3),
'(a)')
'******************************************************************'
624 integer,
intent(in) :: order
625 integer,
intent(in) :: observable(2)
627 logical :: file_exists
628 integer :: i, j, nspin, time_steps, lmax, nfiles, k, add_lm, l, m, max_add_lm
629 integer,
allocatable :: iunit(:)
630 real(real64) :: dt, lambda, dump, o0
632 real(real64),
allocatable :: q(:), mu(:), qq(:, :), c(:)
633 character(len=20) :: filename
636 real(real64),
allocatable :: multipole(:, :, :), ot(:), dipole(:, :)
645 write(filename,
'(a11,i1)')
'multipoles.', nfiles+1
646 inquire(file=trim(filename), exist = file_exists)
647 if (.not. file_exists)
exit
652 if (nfiles == 0)
then
653 write(
message(1),
'(a)')
'No multipoles.x file was found'
656 if (order > nfiles)
then
657 write(
message(1),
'(a)')
'The order that you ask for is higher than the number'
658 write(
message(2),
'(a)')
'of multipoles.x file that you supply.'
663 safe_allocate(iunit(1:nfiles))
665 write(filename,
'(a11,i1)')
'multipoles.', j
666 iunit(j) =
io_open(trim(filename), namespace, action=
'read', status=
'old')
669 safe_allocate( q(1:nfiles))
670 safe_allocate(mu(1:nfiles))
671 safe_allocate(qq(1:nfiles, 1:nfiles))
672 safe_allocate( c(1:nfiles))
678 call spectrum_mult_info(namespace, iunit(1), nspin, kick, time_steps, dt, units, lmax=lmax)
681 if (kick%n_multipoles > 0)
then
682 kick_operator%n_multipoles = kick%n_multipoles
683 safe_allocate( kick_operator%l(1:kick_operator%n_multipoles))
684 safe_allocate( kick_operator%m(1:kick_operator%n_multipoles))
685 safe_allocate(kick_operator%weight(1:kick_operator%n_multipoles))
686 do i = 1, kick_operator%n_multipoles
687 kick_operator%l(i) = kick%l(i)
688 kick_operator%m(i) = kick%m(i)
689 kick_operator%weight(i) = kick%weight(i)
692 kick_operator%n_multipoles = 3
693 safe_allocate( kick_operator%l(1:kick_operator%n_multipoles))
694 safe_allocate( kick_operator%m(1:kick_operator%n_multipoles))
695 safe_allocate(kick_operator%weight(1:kick_operator%n_multipoles))
696 kick_operator%l(1:3) = 1
697 kick_operator%m(1) = -1
698 kick_operator%m(2) = 0
699 kick_operator%m(3) = 1
707 select case (observable(1))
716 safe_allocate(obs%l(1:1))
717 safe_allocate(obs%m(1:1))
718 safe_allocate(obs%weight(1:1))
720 select case (observable(2))
734 safe_allocate(obs%l(1:1))
735 safe_allocate(obs%m(1:1))
736 safe_allocate(obs%weight(1:1))
737 obs%weight(1) =
m_one
738 obs%l(1) = observable(1)
739 obs%m(1) = observable(2)
742 lambda = abs(kick%delta_strength)
743 q(1) = kick%delta_strength / lambda
746 call spectrum_mult_info(namespace, iunit(j), nspin, kick, time_steps, dt, units, lmax=lmax)
747 q(j) = kick%delta_strength / lambda
760 if (kick%n_multipoles > 0)
then
761 lmax = maxval(kick%l(1:obs%n_multipoles))
762 max_add_lm = (lmax+1)**2-1
763 safe_allocate(multipole(1:max_add_lm, 0:time_steps, 1:nspin))
765 mp_unit = units%length**kick%l(1)
768 safe_allocate(multipole(1:3, 0:time_steps, 1:nspin))
769 mp_unit = units%length
771 safe_allocate(ot(0:time_steps))
776 safe_allocate(dipole(1:3, 1:nspin))
784 read(iunit(j), *) k, dump, dump, (multipole(add_lm, i, 1), add_lm = 1, max_add_lm)
786 read(iunit(j), *) k, dump, dump, (multipole(add_lm, i, 1), add_lm = 1, max_add_lm), &
787 dump, (multipole(add_lm, i, 2), add_lm = 1, max_add_lm)
789 read(iunit(j), *) k, dump, dump, (multipole(add_lm, i, 1), add_lm = 1, max_add_lm), &
790 dump, (multipole(add_lm, i, 2), add_lm = 1, max_add_lm), &
791 dump, (multipole(add_lm, i, 3), add_lm = 1, max_add_lm), &
792 dump, (multipole(add_lm, i, 4), add_lm = 1, max_add_lm)
794 multipole(1:max_add_lm, i, :) =
units_to_atomic(mp_unit, multipole(1:max_add_lm, i, :))
798 dipole(1:3, 1:nspin) = - multipole(1:3, i, 1:nspin)
805 do k = 1, obs%n_multipoles
810 if (l == obs%l(k) .and. m == obs%m(k))
exit lcycle
816 dump = dump + obs%weight(k) * sum(multipole(add_lm, i, 1:nspin))
819 if (i == 0) o0 = dump
820 ot(i) = ot(i) + mu(j)*(dump - o0)
825 ot = ot / lambda**order
827 call write_ot(namespace, nspin, time_steps, dt, kick, order, ot(0:time_steps), observable)
833 safe_deallocate_a(iunit)
835 safe_deallocate_a(mu)
836 safe_deallocate_a(qq)
838 safe_deallocate_a(ot)
839 safe_deallocate_a(multipole)
840 safe_deallocate_a(dipole)
848 subroutine modify_ot(time_steps, dt, order, ot, omega, power)
849 integer,
intent(in) :: time_steps
850 real(real64),
intent(in) :: dt
851 integer,
intent(in) :: order
852 real(real64),
intent(inout) :: ot(0:time_steps)
853 real(real64),
intent(in) :: omega
854 real(real64),
intent(in) :: power
860 select case (mod(order, 2))
863 ot(i) = ot(i) -
m_two * power *
sin(omega*i*dt)
868 ot(i) = ot(i) - power *
cos(omega*i*dt)
872 ot(i) = ot(i) -
m_two * power *
cos(omega*i*dt)
883 subroutine write_ot(namespace, nspin, time_steps, dt, kick, order, ot, observable)
885 integer,
intent(in) :: nspin, time_steps
886 real(real64),
intent(in) :: dt
887 type(
kick_t),
intent(in) :: kick
888 integer,
intent(in) :: order
889 real(real64),
intent(in) :: ot(0:time_steps)
890 integer,
intent(in) :: observable(2)
893 character(len=20) :: header_string
898 iunit =
io_open(
'ot', namespace, action=
'write', status=
'replace')
900 write(iunit,
'(a15,i2)')
'# nspin ', nspin
901 write(iunit,
'(a15,i2)')
'# Order ', order
902 write(iunit,
'(a28,i2)')
'# Frequencies subtracted = ', 0
903 select case (observable(1))
905 write(iunit,
'(a)')
'# Observable operator = kick operator'
906 if (kick%n_multipoles > 0)
then
912 select case (observable(2))
914 write(iunit,
'(a)')
'# O = x'
916 write(iunit,
'(a)')
'# O = y'
918 write(iunit,
'(a)')
'# O = z'
922 ot_unit =
units_out%length**observable(1)
923 write(iunit,
'(a12,i1,a1,i2,a1)')
'# (l, m) = (', observable(1),
',',observable(2),
')'
928 write(iunit,
'(a1)', advance =
'no')
'#'
929 write(iunit,
'(a20)', advance =
'no')
str_center(
't', 19)
930 write(iunit,
'(a20)', advance =
'yes')
str_center(
'<O>(t)', 20)
931 write(iunit,
'(a1)', advance =
'no')
'#'
933 write(iunit,
'(a20)', advance =
'no')
str_center(trim(header_string), 19)
935 write(iunit,
'(a20)', advance =
'yes')
str_center(trim(header_string), 20)
947 subroutine read_ot(namespace, nspin, order, nw_subtracted)
949 integer,
intent(out) :: nspin
950 integer,
intent(out) :: order
951 integer,
intent(out) :: nw_subtracted
953 integer :: iunit, i, ierr
954 character(len=100) :: line
955 character(len=12) :: dummychar
956 real(real64) :: dummy, t1, t2
961 if (
allocated(ot))
then
962 safe_deallocate_a(ot)
965 iunit =
io_open(
'ot', namespace, action=
'read', status=
'old')
967 read(iunit,
'(15x,i2)') nspin
968 read(iunit,
'(15x,i2)') order
969 read(iunit,
'(28x,i2)') nw_subtracted
970 read(iunit,
'(a)') line
972 i = index(line,
'Observable')
973 if (index(line,
'Observable') /= 0)
then
975 elseif (index(line,
'# O =') /= 0)
then
977 if (index(line,
'x') /= 0)
then
979 elseif (index(line,
'y') /= 0)
then
981 elseif (index(line,
'z') /= 0)
then
984 elseif (index(line,
'# (l, m) = ') /= 0)
then
985 read(line,
'(a12,i1,a1,i2,a1)') dummychar(1:12), observable(1), dummychar(1:1), observable(2), dummychar(1:1)
987 write(
message(1),
'(a)')
'Problem reading "ot" file: could not figure out the shape'
988 write(
message(2),
'(a)')
'of the observation operator.'
993 read(iunit,
'(a)') line
994 read(iunit,
'(a)') line
1000 write(
message(1),
'(a)')
'Could not figure out the units in file "ot".'
1004 select case (observable(1))
1006 if (kick%n_multipoles > 0)
then
1014 ot_unit =
units_out%length**observable(1)
1022 read(iunit, *,
end=100) dummy
1023 time_steps = time_steps + 1
1024 if (time_steps == 1) t1 = dummy
1025 if (time_steps == 2) t2 = dummy
1028 if (time_steps <= 0)
then
1029 write(
message(1),
'(a)')
'No data points found in file "ot".'
1033 time_steps = time_steps - 1
1034 if (time_steps > 0)
then
1042 safe_allocate(ot(0:time_steps))
1044 do i = 0, time_steps
1045 read(iunit, *) dummy, ot(i)
1056 subroutine print_omega_file(namespace, omega, search_interval, final_time, nfrequencies)
1057 type(namespace_t),
intent(in) :: namespace
1058 real(real64),
intent(inout) :: omega
1059 real(real64),
intent(inout) :: search_interval
1060 real(real64),
intent(inout) :: final_time
1061 integer,
intent(inout) :: nfrequencies
1063 integer :: iunit, i, ierr, nspin, order, nw_subtracted
1064 logical :: file_exists
1065 character(len=20) :: header_string
1066 real(real64),
allocatable :: warray(:), tarray(:)
1067 real(real64) :: leftbound, rightbound, dw, w, aw
1072 inquire(file=
"ot", exist = file_exists)
1073 if (.not. file_exists)
then
1074 write(message(1),
'(a)')
"Could not find 'ot' file."
1075 call messages_fatal(1, namespace=namespace)
1079 call unit_system_from_file(
units,
"ot", namespace, ierr)
1081 write(message(1),
'(a)')
"Could not retrieve units in the 'ot' file."
1082 call messages_fatal(1, namespace=namespace)
1085 if (omega > m_zero)
then
1086 omega = units_to_atomic(
units%energy, omega)
1091 if (search_interval > m_zero)
then
1092 search_interval = units_to_atomic(
units%energy, search_interval)
1094 search_interval = m_half
1097 leftbound = omega - search_interval
1098 rightbound = omega + search_interval
1100 safe_allocate(warray(1:nfrequencies))
1101 safe_allocate(tarray(1:nfrequencies))
1103 call read_ot(namespace, nspin, order, nw_subtracted)
1105 if (mod(order, 2) == 1)
then
1111 if (final_time > m_zero)
then
1115 write(message(1),
'(a)')
'The requested total time to process is larger than the time available in the input file.'
1116 write(message(2),
'(a,f8.4,a)')
'The time has been adjusted to ', &
1117 units_from_atomic(units_out%time,
total_time), trim(units_abbrev(units_out%time))
1118 call messages_warning(2, namespace=namespace)
1128 dw = (rightbound-leftbound) / (nfrequencies - 1)
1129 do i = 1, nfrequencies
1130 w = leftbound + (i-1)*dw
1136 iunit = io_open(
'omega', namespace, action=
'write', status=
'replace')
1137 write(iunit,
'(a15,i2)')
'# nspin ', nspin
1138 call kick_write(
kick, iunit)
1139 write(iunit,
'(a)')
'#%'
1140 write(iunit,
'(a1,a20)', advance =
'no')
'#', str_center(
"omega", 20)
1141 write(header_string,
'(a)')
'F(omega)'
1142 write(iunit,
'(a20)', advance =
'yes') str_center(trim(header_string), 20)
1143 write(iunit,
'(a1,a20)', advance =
'no')
'#', str_center(
'['//trim(units_abbrev(units_out%energy)) //
']', 20)
1145 write(iunit,
'(a)', advance =
'yes')
1147 do i = 1, nfrequencies
1148 write(iunit,
'(2e20.8)') units_from_atomic(units_out%energy, warray(i)), &
1152 call io_close(iunit)
1154 safe_deallocate_a(warray)
1155 safe_deallocate_a(tarray)
1156 safe_deallocate_a(
ot)
1178 integer :: run_mode, order, nfrequencies, ierr, nresonances
1179 real(real64) :: omega, search_interval, final_time, damping
1180 integer,
parameter :: &
1181 ANALYZE_NTHORDER_SIGNAL = 1, &
1182 generate_nthorder_signal = 2, &
1183 read_resonances_from_file = 3, &
1184 generate_omega_file = 4
1185 character(len=100) :: ffile
1186 character(kind=c_char) :: cfile(101)
1191 message(1) =
"Your Fortran compiler doesn't support command-line arguments;"
1192 message(2) =
"the oct-oscillator-strength command is not available."
1197 run_mode = analyze_nthorder_signal
1199 search_interval = -
m_one
1202 final_time = -
m_one
1207 damping = 0.1_real64/27.2114_real64
1212 order, nresonances, nfrequencies, final_time, &
1214 call string_c_to_f(cfile, ffile)
1223 select case (run_mode)
1224 case (generate_nthorder_signal)
1226 case (analyze_nthorder_signal)
1227 call analyze_signal(order, omega, search_interval, final_time, nresonances, nfrequencies, damping, &
1229 case (read_resonances_from_file)
1230 call read_resonances_file(order, ffile, global_namespace, search_interval, final_time, nfrequencies)
1231 case (generate_omega_file)
1232 call print_omega_file(global_namespace, omega, search_interval, final_time, nfrequencies)
double sin(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
subroutine, public getopt_init(ierr)
Initializes the getopt machinery. Must be called before attempting to parse the options....
subroutine, public getopt_end
real(real64), parameter, public m_two
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
real(real64), parameter, public m_pi
some mathematical constants
type(mpi_comm), parameter, public serial_dummy_comm
Alias MPI_COMM_UNDEFINED for the specific use case of initialising Octopus utilities with no MPI supp...
subroutine, public init_octopus_globals(comm)
Initialise Octopus-specific global constants and files. This routine performs no initialisation calls...
complex(real64), parameter, public m_z0
complex(real64), parameter, public m_zi
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public m_one
real(real64), parameter, public m_three
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()
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
subroutine, public kick_read(kick, iunit, namespace)
subroutine, public kick_end(kick)
subroutine, public kick_write(kick, iunit, out)
subroutine, public messages_not_implemented(feature, namespace)
subroutine, public messages_warning(no_lines, all_nodes, namespace)
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_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
subroutine ft(omega, power)
integer, dimension(2) observable
subroutine write_ot(namespace, nspin, time_steps, dt, kick, order, ot, observable)
subroutine generate_signal(namespace, order, observable)
subroutine ft2(omega, power)
real(real64), dimension(:), allocatable ot
subroutine read_ot(namespace, nspin, order, nw_subtracted)
subroutine resonance_first_order(omega, power, nw_subtracted, dw, leftbound, rightbound)
type(unit_system_t) units
integer, parameter sine_transform
integer, parameter cosine_transform
subroutine resonance_second_order(omega, power, nw_subtracted, leftbound, rightbound, c01, c02)
subroutine analyze_signal(order, omega, search_interval, final_time, nresonances, nfrequencies, damping, namespace)
subroutine modify_ot(time_steps, dt, order, ot, omega, power)
subroutine print_omega_file(namespace, omega, search_interval, final_time, nfrequencies)
subroutine find_resonance(omega, leftbound, rightbound, nfrequencies)
subroutine local_operator_copy(o, i)
subroutine read_resonances_file(order, ffile, namespace, search_interval, final_time, nfrequencies)
subroutine write_polarizability(namespace, nfrequencies, nresonances, dw, w, c0I2, gamma)
Implements the SOS formula of the polarizability, and writes down to the "polarizability" file the re...
subroutine, public spectrum_mult_info(namespace, iunit, nspin, kick, time_steps, dt, file_units, lmax)
character(len=80) function, public str_center(s_in, l_in)
puts space around string, so that it is centered
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
subroutine, public unit_system_from_file(uu, fname, namespace, ierr)
This is a very primitive procedure that attempts to find out which units were used to write an octopu...
subroutine, public unit_system_init(namespace)
program oscillator_strength