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, dw
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
269 dw = (rightbound-leftbound) / (nfrequencies - 1)
272 leftbound = -search_interval
273 rightbound = search_interval
276 call modify_ot(time_steps, dt, order, ot, w, power)
277 nw_subtracted = nw_subtracted + 1
281 do i = 1, nresonances
282 do j = i + 1, nresonances
283 leftbound = wij(k) - search_interval
284 rightbound = wij(k) + search_interval
287 call modify_ot(time_steps, dt, order, ot, wij(k), power)
288 nw_subtracted = nw_subtracted + 1
293 safe_deallocate_a(wij)
294 safe_deallocate_a(c0i)
295 safe_deallocate_a(omega)
296 safe_deallocate_a(ot)
306 subroutine analyze_signal(order, omega, search_interval, final_time, nresonances, nfrequencies, &
308 integer,
intent(inout) :: order
309 real(real64),
intent(inout) :: omega
310 real(real64),
intent(inout) :: search_interval
311 real(real64),
intent(inout) :: final_time
312 integer,
intent(inout) :: nresonances
313 integer,
intent(inout) :: nfrequencies
314 real(real64),
intent(in) :: damping
317 real(real64) :: leftbound, rightbound, dw, power
318 real(real64),
allocatable :: w(:), c0I2(:)
319 integer :: nspin, i, ierr, order_in_file, nw_subtracted
320 logical :: file_exists
325 inquire(file=
"ot", exist = file_exists)
326 if (.not. file_exists)
then
327 write(
message(1),
'(a)')
"Could not find 'ot' file."
334 write(
message(1),
'(a)')
"Could not retrieve units in the 'ot' file."
344 if (search_interval >
m_zero)
then
350 leftbound = omega - search_interval
351 rightbound = omega + search_interval
353 call read_ot(namespace, nspin, order_in_file, nw_subtracted)
355 if (order_in_file /= order)
then
356 write(
message(1),
'(a)')
'Internal error in analyze_signal'
360 if (mod(order, 2) == 1)
then
361 mode = sine_transform
366 if (final_time >
m_zero)
then
368 if (total_time > dt*time_steps)
then
369 total_time = dt*time_steps
370 write(
message(1),
'(a)')
'The requested total time to process is larger than the time available in the input file.'
371 write(
message(2),
'(a,f8.4,a)')
'The time has been adjusted to ', &
375 time_steps = int(total_time / dt)
376 total_time = time_steps * dt
378 total_time = dt*time_steps
381 dw = (rightbound-leftbound) / (nfrequencies - 1)
383 safe_allocate( w(1:nresonances))
384 safe_allocate(c0i2(1:nresonances))
390 if (nw_subtracted >= nresonances)
exit
408 call modify_ot(time_steps, dt, order, ot, omega, power)
410 nw_subtracted = nw_subtracted + 1
419 safe_deallocate_a(ot)
421 safe_deallocate_a(c0i2)
435 integer,
intent(in) :: nfrequencies, nresonances
436 real(real64),
intent(in) :: dw
437 real(real64),
intent(in) :: w(nresonances), c0I2(nresonances)
438 real(real64),
intent(in) :: gamma
440 integer :: iunit, i, j
442 complex(real64) :: pol
446 iunit =
io_open(
'polarizability', namespace, action =
'write', status=
'replace', die=.false.)
447 write(iunit,
'(a)')
'# Polarizability file. Generated using the SOS formula with the following data:'
448 write(iunit,
'(a)')
'#'
450 do i = 1, nresonances
451 write(iunit,
'(a1,3e20.12)')
'#', w(i),
sqrt(abs(c0i2(i))), c0i2(i)
454 write(iunit,
'(a10,f12.6)')
'# Gamma = ', gamma
455 write(iunit,
'(a)')
'#'
457 do i = 1, nfrequencies
460 do j = 1, nresonances
461 pol = pol + c0i2(j) * (
m_one / (w(j) - e -
m_zi * gamma) +
m_one/(w(j) + e +
m_zi * gamma))
463 write(iunit,
'(3e20.12)') e, pol
475 subroutine find_resonance(omega, leftbound, rightbound, nfrequencies)
476 real(real64),
intent(inout) :: omega
477 real(real64),
intent(in) :: leftbound, rightbound
478 integer,
intent(in) :: nfrequencies
481 real(real64) :: dw, w, aw, min_aw, min_w, omega_orig, min_w1, min_w2
482 real(real64),
allocatable :: warray(:), tarray(:)
486 safe_allocate(warray(1:nfrequencies))
487 safe_allocate(tarray(1:nfrequencies))
491 dw = (rightbound-leftbound) / (nfrequencies - 1)
492 do i = 1, nfrequencies
493 w = leftbound + (i-1)*dw
501 do i = 1, nfrequencies
502 w = leftbound + (i-1)*dw
503 if (tarray(i) < min_aw)
then
511 min_w1 = min_w - 2*dw
512 min_w2 = min_w + 2*dw
516 write(
message(1),
'(a)')
'Could not find a maximum.'
518 write(
message(3),
'(a,f12.8,a,f12.8,a)')
' Search interval = [', &
520 write(
message(4),
'(a,f12.4,a)')
' Search discretization = ', &
525 safe_deallocate_a(warray)
526 safe_deallocate_a(tarray)
535 real(real64),
intent(in) :: omega
536 real(real64),
intent(out) :: power
537 integer,
intent(in) :: nw_subtracted
538 real(real64),
intent(in) :: dw, leftbound, rightbound
542 call ft(omega, power)
545 case (sine_transform)
551 write(
message(1),
'(a)')
'******************************************************************'
552 write(
message(2),
'(a,i3)')
'Resonance #', nw_subtracted + 1
559 write(
message(6),
'(a,f12.8)')
'f[O->I] = ',
m_two*omega*power
561 write(
message(8),
'(a,f12.8,a,f12.8,a)')
' Search interval = [',
units_from_atomic(
units_out%energy, leftbound),
',', &
565 write(
message(10),
'(a)')
'******************************************************************'
577 real(real64),
intent(in) :: omega
578 real(real64),
intent(out) :: power
579 integer,
intent(in) :: nw_subtracted
580 real(real64),
intent(in) :: leftbound, rightbound
581 real(real64),
intent(in) :: c01, c02
585 call ft(omega, power)
587 case (sine_transform)
594 power = power /
m_two
598 write(
message(1),
'(a)')
'******************************************************************'
599 write(
message(2),
'(a,i3)')
'Resonance #', nw_subtracted + 1
607 write(
message(1),
'(a,f12.8)')
' C(omega)/(C0i*C0j) = ', power / (c01 * c02)
612 write(
message(2),
'(a,f12.8,a,f12.8,a)')
' Search interval = [',
units_from_atomic(
units_out%energy, leftbound),
',', &
614 write(
message(3),
'(a)')
'******************************************************************'
626 integer,
intent(in) :: order
627 integer,
intent(in) :: observable(2)
629 logical :: file_exists
630 integer :: i, j, nspin, time_steps, lmax, nfiles, k, add_lm, l, m, max_add_lm
631 integer,
allocatable :: iunit(:)
632 real(real64) :: dt, lambda, dump, o0
634 real(real64),
allocatable :: q(:), mu(:), qq(:, :), c(:)
635 character(len=20) :: filename
638 real(real64),
allocatable :: multipole(:, :, :), ot(:), dipole(:, :)
647 write(filename,
'(a11,i1)')
'multipoles.', nfiles+1
648 inquire(file=trim(filename), exist = file_exists)
649 if (.not. file_exists)
exit
654 if (nfiles == 0)
then
655 write(
message(1),
'(a)')
'No multipoles.x file was found'
658 if (order > nfiles)
then
659 write(
message(1),
'(a)')
'The order that you ask for is higher than the number'
660 write(
message(2),
'(a)')
'of multipoles.x file that you supply.'
665 safe_allocate(iunit(1:nfiles))
667 write(filename,
'(a11,i1)')
'multipoles.', j
668 iunit(j) =
io_open(trim(filename), namespace, action=
'read', status=
'old')
671 safe_allocate( q(1:nfiles))
672 safe_allocate(mu(1:nfiles))
673 safe_allocate(qq(1:nfiles, 1:nfiles))
674 safe_allocate( c(1:nfiles))
680 call spectrum_mult_info(namespace, iunit(1), nspin, kick, time_steps, dt, units, lmax=lmax)
683 if (kick%n_multipoles > 0)
then
684 kick_operator%n_multipoles = kick%n_multipoles
685 safe_allocate( kick_operator%l(1:kick_operator%n_multipoles))
686 safe_allocate( kick_operator%m(1:kick_operator%n_multipoles))
687 safe_allocate(kick_operator%weight(1:kick_operator%n_multipoles))
688 do i = 1, kick_operator%n_multipoles
689 kick_operator%l(i) = kick%l(i)
690 kick_operator%m(i) = kick%m(i)
691 kick_operator%weight(i) = kick%weight(i)
694 kick_operator%n_multipoles = 3
695 safe_allocate( kick_operator%l(1:kick_operator%n_multipoles))
696 safe_allocate( kick_operator%m(1:kick_operator%n_multipoles))
697 safe_allocate(kick_operator%weight(1:kick_operator%n_multipoles))
698 kick_operator%l(1:3) = 1
699 kick_operator%m(1) = -1
700 kick_operator%m(2) = 0
701 kick_operator%m(3) = 1
709 select case (observable(1))
718 safe_allocate(obs%l(1:1))
719 safe_allocate(obs%m(1:1))
720 safe_allocate(obs%weight(1:1))
722 select case (observable(2))
736 safe_allocate(obs%l(1:1))
737 safe_allocate(obs%m(1:1))
738 safe_allocate(obs%weight(1:1))
739 obs%weight(1) =
m_one
740 obs%l(1) = observable(1)
741 obs%m(1) = observable(2)
744 lambda = abs(kick%delta_strength)
745 q(1) = kick%delta_strength / lambda
748 call spectrum_mult_info(namespace, iunit(j), nspin, kick, time_steps, dt, units, lmax=lmax)
749 q(j) = kick%delta_strength / lambda
762 if (kick%n_multipoles > 0)
then
763 lmax = maxval(kick%l(1:obs%n_multipoles))
764 max_add_lm = (lmax+1)**2-1
765 safe_allocate(multipole(1:max_add_lm, 0:time_steps, 1:nspin))
767 mp_unit = units%length**kick%l(1)
770 safe_allocate(multipole(1:3, 0:time_steps, 1:nspin))
771 mp_unit = units%length
773 safe_allocate(ot(0:time_steps))
778 safe_allocate(dipole(1:3, 1:nspin))
786 read(iunit(j), *) k, dump, dump, (multipole(add_lm, i, 1), add_lm = 1, max_add_lm)
788 read(iunit(j), *) k, dump, dump, (multipole(add_lm, i, 1), add_lm = 1, max_add_lm), &
789 dump, (multipole(add_lm, i, 2), add_lm = 1, max_add_lm)
791 read(iunit(j), *) k, dump, dump, (multipole(add_lm, i, 1), add_lm = 1, max_add_lm), &
792 dump, (multipole(add_lm, i, 2), add_lm = 1, max_add_lm), &
793 dump, (multipole(add_lm, i, 3), add_lm = 1, max_add_lm), &
794 dump, (multipole(add_lm, i, 4), add_lm = 1, max_add_lm)
796 multipole(1:max_add_lm, i, :) =
units_to_atomic(mp_unit, multipole(1:max_add_lm, i, :))
800 dipole(1:3, 1:nspin) = - multipole(1:3, i, 1:nspin)
807 do k = 1, obs%n_multipoles
812 if (l == obs%l(k) .and. m == obs%m(k))
exit lcycle
818 dump = dump + obs%weight(k) * sum(multipole(add_lm, i, 1:nspin))
821 if (i == 0) o0 = dump
822 ot(i) = ot(i) + mu(j)*(dump - o0)
827 ot = ot / lambda**order
829 call write_ot(namespace, nspin, time_steps, dt, kick, order, ot(0:time_steps), observable)
835 safe_deallocate_a(iunit)
837 safe_deallocate_a(mu)
838 safe_deallocate_a(qq)
840 safe_deallocate_a(ot)
841 safe_deallocate_a(multipole)
842 safe_deallocate_a(dipole)
850 subroutine modify_ot(time_steps, dt, order, ot, omega, power)
851 integer,
intent(in) :: time_steps
852 real(real64),
intent(in) :: dt
853 integer,
intent(in) :: order
854 real(real64),
intent(inout) :: ot(0:time_steps)
855 real(real64),
intent(in) :: omega
856 real(real64),
intent(in) :: power
862 select case (mod(order, 2))
865 ot(i) = ot(i) -
m_two * power *
sin(omega*i*dt)
870 ot(i) = ot(i) - power *
cos(omega*i*dt)
874 ot(i) = ot(i) -
m_two * power *
cos(omega*i*dt)
885 subroutine write_ot(namespace, nspin, time_steps, dt, kick, order, ot, observable)
887 integer,
intent(in) :: nspin, time_steps
888 real(real64),
intent(in) :: dt
889 type(
kick_t),
intent(in) :: kick
890 integer,
intent(in) :: order
891 real(real64),
intent(in) :: ot(0:time_steps)
892 integer,
intent(in) :: observable(2)
895 character(len=20) :: header_string
900 iunit =
io_open(
'ot', namespace, action=
'write', status=
'replace')
902 write(iunit,
'(a15,i2)')
'# nspin ', nspin
903 write(iunit,
'(a15,i2)')
'# Order ', order
904 write(iunit,
'(a28,i2)')
'# Frequencies subtracted = ', 0
905 select case (observable(1))
907 write(iunit,
'(a)')
'# Observable operator = kick operator'
908 if (kick%n_multipoles > 0)
then
914 select case (observable(2))
916 write(iunit,
'(a)')
'# O = x'
918 write(iunit,
'(a)')
'# O = y'
920 write(iunit,
'(a)')
'# O = z'
924 ot_unit =
units_out%length**observable(1)
925 write(iunit,
'(a12,i1,a1,i2,a1)')
'# (l, m) = (', observable(1),
',',observable(2),
')'
930 write(iunit,
'(a1)', advance =
'no')
'#'
931 write(iunit,
'(a20)', advance =
'no')
str_center(
't', 19)
932 write(iunit,
'(a20)', advance =
'yes')
str_center(
'<O>(t)', 20)
933 write(iunit,
'(a1)', advance =
'no')
'#'
935 write(iunit,
'(a20)', advance =
'no')
str_center(trim(header_string), 19)
937 write(iunit,
'(a20)', advance =
'yes')
str_center(trim(header_string), 20)
949 subroutine read_ot(namespace, nspin, order, nw_subtracted)
951 integer,
intent(out) :: nspin
952 integer,
intent(out) :: order
953 integer,
intent(out) :: nw_subtracted
955 integer :: iunit, i, ierr
956 character(len=100) :: line
957 character(len=12) :: dummychar
958 real(real64) :: dummy, t1, t2
963 if (
allocated(ot))
then
964 safe_deallocate_a(ot)
967 iunit =
io_open(
'ot', namespace, action=
'read', status=
'old')
969 read(iunit,
'(15x,i2)') nspin
970 read(iunit,
'(15x,i2)') order
971 read(iunit,
'(28x,i2)') nw_subtracted
972 read(iunit,
'(a)') line
974 i = index(line,
'Observable')
975 if (index(line,
'Observable') /= 0)
then
977 elseif (index(line,
'# O =') /= 0)
then
979 if (index(line,
'x') /= 0)
then
981 elseif (index(line,
'y') /= 0)
then
983 elseif (index(line,
'z') /= 0)
then
986 elseif (index(line,
'# (l, m) = ') /= 0)
then
987 read(line,
'(a12,i1,a1,i2,a1)') dummychar(1:12), observable(1), dummychar(1:1), observable(2), dummychar(1:1)
989 write(
message(1),
'(a)')
'Problem reading "ot" file: could not figure out the shape'
990 write(
message(2),
'(a)')
'of the observation operator.'
995 read(iunit,
'(a)') line
996 read(iunit,
'(a)') line
1002 write(
message(1),
'(a)')
'Could not figure out the units in file "ot".'
1006 select case (observable(1))
1008 if (kick%n_multipoles > 0)
then
1016 ot_unit =
units_out%length**observable(1)
1024 read(iunit, *,
end=100) dummy
1025 time_steps = time_steps + 1
1026 if (time_steps == 1) t1 = dummy
1027 if (time_steps == 2) t2 = dummy
1030 if (time_steps <= 0)
then
1031 write(
message(1),
'(a)')
'No data points found in file "ot".'
1035 time_steps = time_steps - 1
1036 if (time_steps > 0)
then
1044 safe_allocate(ot(0:time_steps))
1046 do i = 0, time_steps
1047 read(iunit, *) dummy, ot(i)
1058 subroutine print_omega_file(namespace, omega, search_interval, final_time, nfrequencies)
1059 type(namespace_t),
intent(in) :: namespace
1060 real(real64),
intent(inout) :: omega
1061 real(real64),
intent(inout) :: search_interval
1062 real(real64),
intent(inout) :: final_time
1063 integer,
intent(inout) :: nfrequencies
1065 integer :: iunit, i, ierr, nspin, order, nw_subtracted
1066 logical :: file_exists
1067 character(len=20) :: header_string
1068 real(real64),
allocatable :: warray(:), tarray(:)
1069 real(real64) :: leftbound, rightbound, dw, w, aw
1074 inquire(file=
"ot", exist = file_exists)
1075 if (.not. file_exists)
then
1076 write(message(1),
'(a)')
"Could not find 'ot' file."
1077 call messages_fatal(1, namespace=namespace)
1081 call unit_system_from_file(
units,
"ot", namespace, ierr)
1083 write(message(1),
'(a)')
"Could not retrieve units in the 'ot' file."
1084 call messages_fatal(1, namespace=namespace)
1087 if (omega > m_zero)
then
1088 omega = units_to_atomic(
units%energy, omega)
1093 if (search_interval > m_zero)
then
1094 search_interval = units_to_atomic(
units%energy, search_interval)
1096 search_interval = m_half
1099 leftbound = omega - search_interval
1100 rightbound = omega + search_interval
1102 safe_allocate(warray(1:nfrequencies))
1103 safe_allocate(tarray(1:nfrequencies))
1105 call read_ot(namespace, nspin, order, nw_subtracted)
1107 if (mod(order, 2) == 1)
then
1113 if (final_time > m_zero)
then
1117 write(message(1),
'(a)')
'The requested total time to process is larger than the time available in the input file.'
1118 write(message(2),
'(a,f8.4,a)')
'The time has been adjusted to ', &
1119 units_from_atomic(units_out%time,
total_time), trim(units_abbrev(units_out%time))
1120 call messages_warning(2, namespace=namespace)
1130 dw = (rightbound-leftbound) / (nfrequencies - 1)
1131 do i = 1, nfrequencies
1132 w = leftbound + (i-1)*dw
1138 iunit = io_open(
'omega', namespace, action=
'write', status=
'replace')
1139 write(iunit,
'(a15,i2)')
'# nspin ', nspin
1140 call kick_write(
kick, iunit)
1141 write(iunit,
'(a)')
'#%'
1142 write(iunit,
'(a1,a20)', advance =
'no')
'#', str_center(
"omega", 20)
1143 write(header_string,
'(a)')
'F(omega)'
1144 write(iunit,
'(a20)', advance =
'yes') str_center(trim(header_string), 20)
1145 write(iunit,
'(a1,a20)', advance =
'no')
'#', str_center(
'['//trim(units_abbrev(units_out%energy)) //
']', 20)
1147 write(iunit,
'(a)', advance =
'yes')
1149 do i = 1, nfrequencies
1150 write(iunit,
'(2e20.8)') units_from_atomic(units_out%energy, warray(i)), &
1154 call io_close(iunit)
1156 safe_deallocate_a(warray)
1157 safe_deallocate_a(tarray)
1158 safe_deallocate_a(
ot)
1180 integer :: run_mode, order, nfrequencies, ierr, nresonances
1181 real(real64) :: omega, search_interval, final_time, damping
1182 integer,
parameter :: &
1183 ANALYZE_NTHORDER_SIGNAL = 1, &
1184 generate_nthorder_signal = 2, &
1185 read_resonances_from_file = 3, &
1186 generate_omega_file = 4
1187 character(len=100) :: ffile
1188 character(kind=c_char) :: cfile(101)
1193 message(1) =
"Your Fortran compiler doesn't support command-line arguments;"
1194 message(2) =
"the oct-oscillator-strength command is not available."
1199 run_mode = analyze_nthorder_signal
1201 search_interval = -
m_one
1204 final_time = -
m_one
1209 damping = 0.1_real64/27.2114_real64
1214 order, nresonances, nfrequencies, final_time, &
1216 call string_c_to_f(cfile, ffile)
1225 select case (run_mode)
1226 case (generate_nthorder_signal)
1228 case (analyze_nthorder_signal)
1229 call analyze_signal(order, omega, search_interval, final_time, nresonances, nfrequencies, damping, &
1231 case (read_resonances_from_file)
1232 call read_resonances_file(order, ffile, global_namespace, search_interval, final_time, nfrequencies)
1233 case (generate_omega_file)
1234 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