Octopus
oscillator_strength.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22 use debug_oct_m
23 use global_oct_m
24 use io_oct_m
25 use kick_oct_m
26 use, intrinsic :: iso_fortran_env
27 use iso_c_binding
32 use parser_oct_m
35 use string_oct_m
36 use unit_oct_m
38
39 implicit none
40
41 integer, parameter :: &
42 SINE_TRANSFORM = 1, &
44
46 integer :: n_multipoles
47 integer, allocatable :: l(:), m(:)
48 real(real64), allocatable :: weight(:)
49 end type local_operator_t
50
51 integer :: observable(2)
52 type(unit_system_t) :: units
53 real(real64), allocatable :: ot(:)
54 type(kick_t) :: kick
55 integer :: time_steps
56 real(real64) :: total_time
57 integer :: mode
58 real(real64) :: dt
59
60contains
61
62 ! ---------------------------------------------------------
63 subroutine ft(omega, power)
64 real(real64), intent(in) :: omega
65 real(real64), intent(out) :: power
66
67 integer :: j
68 real(real64) :: x
69 power = m_zero
70
71 push_sub(ft)
72
73 select case (mode)
74
75 case (sine_transform)
76
77 do j = 0, time_steps
78 x = sin(omega*j*dt)
79 power = power + x*ot(j)
80 end do
81 if (time_steps > 0) then
82 power = power / real(time_steps, real64)
83 else
84 power = m_zero
85 end if
86
87 case (cosine_transform)
88
89 do j = 0, time_steps
90 x = cos(omega*j*dt)
91 power = power + x*ot(j)
92 end do
93 if (time_steps > 0) then
94 power = power / real(time_steps, real64)
95 else
96 power = m_zero
97 end if
98
99 end select
100
101 pop_sub(ft)
102 end subroutine ft
103
104
105 ! ---------------------------------------------------------
106 subroutine ft2(omega, power)
107 real(real64), intent(in) :: omega
108 real(real64), intent(out) :: power
109
110 integer :: j
111 real(real64) :: x
112 power = m_zero
113
114 push_sub(ft2)
115
116 select case (mode)
117
118 case (sine_transform)
119
120 do j = 0, time_steps
121 x = sin(omega*j*dt)
122 power = power + x*ot(j)
123 end do
124 ! The function should be *minus* the sine Fourier transform, since this
125 ! is the function to be minimized.
126 if (time_steps > 0) then
127 power = - (power / real(time_steps, real64))**2
128 else
129 power = m_zero
130 end if
131
132 case (cosine_transform)
133
134 do j = 0, time_steps
135 x = cos(omega*j*dt)
136 power = power + x*ot(j)
137 end do
138 if (time_steps > 0) then
139 power = - (power / real(time_steps, real64))**2
140 else
141 power = m_zero
142 end if
144 end select
145
146 pop_sub(ft2)
147 end subroutine ft2
150 ! ---------------------------------------------------------
151 subroutine local_operator_copy(o, i)
152 type(local_operator_t), intent(inout) :: o
153 type(local_operator_t), intent(inout) :: i
154
155 integer :: j
156
157 push_sub(local_operator_copy)
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))
163
164 do j = 1, o%n_multipoles
165 o%l(j) = i%l(j)
166 o%m(j) = i%m(j)
167 o%weight(j) = i%weight(j)
168 end do
169
170 pop_sub(local_operator_copy)
171 end subroutine local_operator_copy
172
173 ! ---------------------------------------------------------
174 subroutine read_resonances_file(order, ffile, namespace, search_interval, final_time, nfrequencies)
175 integer, intent(inout) :: order
176 character(len=*), intent(in) :: ffile
177 type(namespace_t), intent(in):: namespace
178 real(real64), intent(inout) :: search_interval
179 real(real64), intent(in) :: final_time
180 integer, intent(in) :: nfrequencies
181
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(:)
186
187 push_sub(read_resonances_file)
188
189 if (order /= 2) then
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.'
192 call messages_fatal(2)
193 end if
194
195 ! First, let us check that the file "ot" exists.
196 inquire(file="ot", exist = file_exists)
197 if (.not. file_exists) then
198 write(message(1),'(a)') "Could not find 'ot' file."
199 call messages_fatal(1)
200 end if
202 ! Now, we should find out which units the file "ot" has.
203 call unit_system_from_file(units, "ot", namespace, ierr)
204 if (ierr /= 0) then
205 write(message(1),'(a)') "Could not retrieve units in the 'ot' file."
206 call messages_fatal(1)
207 end if
208
209 mode = cosine_transform
210
211 iunit = io_open(trim(ffile), namespace, action='read', status='old')
212
213 call io_skip_header(iunit)
214 ! Count the number of resonances
215 nresonances = 0
216 do
217 read(iunit, *, iostat = ios) dummy, dummy
218 if (ios /= 0) exit
219 nresonances = nresonances + 1
220 end do
221
222 npairs = (nresonances*(nresonances-1))/2
223
224 safe_allocate(omega(1:nresonances))
225 safe_allocate( c0i(1:nresonances))
226 safe_allocate(wij(1:npairs))
227
228 call io_skip_header(iunit)
229 do i = 1, nresonances
230 read(iunit, *) omega(i), c0i(i)
231 end do
232
233 k = 1
234 do i = 1, nresonances
235 do j = i + 1, nresonances
236 wij(k) = omega(j) - omega(i)
237 k = k + 1
238 end do
239 end do
240
241 if (search_interval > m_zero) then
242 search_interval = units_to_atomic(units%energy, search_interval)
243 else
244 search_interval = m_half
245 end if
247 call read_ot(namespace, nspin, order_in_file, nw_subtracted)
248
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.'
251 call messages_fatal(1)
252 end if
253
254 if (final_time > m_zero) then
255 total_time = units_to_atomic(units%time, final_time)
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 ', &
260 units_from_atomic(units%time, total_time), units_abbrev(units%time)
261 call messages_warning(2, namespace=namespace)
262 end if
263 time_steps = int(total_time / dt)
264 total_time = time_steps * dt
265 else
266 total_time = dt*time_steps
267 end if
268
269 dw = (rightbound-leftbound) / (nfrequencies - 1)
270
271 ! First, subtract zero resonance...
272 leftbound = -search_interval
273 rightbound = search_interval
274 w = m_zero
275 call resonance_second_order(w, power, nw_subtracted, leftbound, rightbound, m_zero, m_zero)
276 call modify_ot(time_steps, dt, order, ot, w, power)
277 nw_subtracted = nw_subtracted + 1
278
279 ! Then, get all the others...
280 k = 1
281 do i = 1, nresonances
282 do j = i + 1, nresonances
283 leftbound = wij(k) - search_interval
284 rightbound = wij(k) + search_interval
285 call find_resonance(wij(k), leftbound, rightbound, nfrequencies)
286 call resonance_second_order(wij(k), power, nw_subtracted, leftbound, rightbound, c0i(i), c0i(j))
287 call modify_ot(time_steps, dt, order, ot, wij(k), power)
288 nw_subtracted = nw_subtracted + 1
289 k = k + 1
290 end do
291 end do
292
293 safe_deallocate_a(wij)
294 safe_deallocate_a(c0i)
295 safe_deallocate_a(omega)
296 safe_deallocate_a(ot)
297 call kick_end(kick)
298 call io_close(iunit)
299
300 pop_sub(read_resonances_file)
301 end subroutine read_resonances_file
302 ! ---------------------------------------------------------
303
304
305 ! ---------------------------------------------------------
306 subroutine analyze_signal(order, omega, search_interval, final_time, nresonances, nfrequencies, &
307 damping, namespace)
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
315 type(namespace_t), intent(in) :: namespace
316
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
321
322 push_sub(analyze_signal)
323
324 ! First, let us check that the file "ot" exists.
325 inquire(file="ot", exist = file_exists)
326 if (.not. file_exists) then
327 write(message(1),'(a)') "Could not find 'ot' file."
328 call messages_fatal(1)
329 end if
330
331 ! Now, we should find out which units the file "ot" has.
332 call unit_system_from_file(units, "ot", namespace, ierr)
333 if (ierr /= 0) then
334 write(message(1),'(a)') "Could not retrieve units in the 'ot' file."
335 call messages_fatal(1)
336 end if
337
338 if (omega > m_zero) then
339 omega = units_to_atomic(units%energy, omega)
340 else
341 omega = m_half
342 end if
343
344 if (search_interval > m_zero) then
345 search_interval = units_to_atomic(units%energy, search_interval)
346 else
347 search_interval = m_half
348 end if
349
350 leftbound = omega - search_interval
351 rightbound = omega + search_interval
352
353 call read_ot(namespace, nspin, order_in_file, nw_subtracted)
354
355 if (order_in_file /= order) then
356 write(message(1), '(a)') 'Internal error in analyze_signal'
357 call messages_fatal(1, namespace=namespace)
358 end if
359
360 if (mod(order, 2) == 1) then
361 mode = sine_transform
362 else
363 mode = cosine_transform
364 end if
365
366 if (final_time > m_zero) then
367 total_time = units_to_atomic(units%time, final_time)
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 ', &
372 units_from_atomic(units%time, total_time), units_abbrev(units%time)
373 call messages_warning(2, namespace=namespace)
374 end if
375 time_steps = int(total_time / dt)
376 total_time = time_steps * dt
377 else
378 total_time = dt*time_steps
379 end if
380
381 dw = (rightbound-leftbound) / (nfrequencies - 1)
382
383 safe_allocate( w(1:nresonances))
384 safe_allocate(c0i2(1:nresonances))
385 w = m_zero
386 c0i2 = m_zero
387
388 i = 1
389 do
390 if (nw_subtracted >= nresonances) exit
391
392 if (mode == cosine_transform .and. nw_subtracted == 0) then
393 omega = m_zero
394 else
395 call find_resonance(omega, leftbound, rightbound, nfrequencies)
396 end if
397
398 select case (order)
399 case (1)
400 call resonance_first_order(omega, power, nw_subtracted, dw, leftbound, rightbound)
401 case (2)
402 call resonance_second_order(omega, power, nw_subtracted, leftbound, rightbound, m_zero, m_zero)
403 end select
404
405 w(i) = omega
406 c0i2(i) = power
407
408 call modify_ot(time_steps, dt, order, ot, omega, power)
409
410 nw_subtracted = nw_subtracted + 1
411 i = i + 1
412 end do
413
414 select case (order)
415 case (1)
416 call write_polarizability(namespace, nfrequencies, nresonances, dw, w, c0i2, damping)
417 end select
418
419 safe_deallocate_a(ot)
420 safe_deallocate_a(w)
421 safe_deallocate_a(c0i2)
422 call kick_end(kick)
423
424 pop_sub(analyze_signal)
425 end subroutine analyze_signal
426 ! ---------------------------------------------------------
427
428
429 ! ---------------------------------------------------------
433 subroutine write_polarizability(namespace, nfrequencies, nresonances, dw, w, c0I2, gamma)
434 type(namespace_t), intent(in) :: namespace
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
439
440 integer :: iunit, i, j
441 real(real64) :: e
442 complex(real64) :: pol
443
444 push_sub(write_polarizability)
445
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)') '#'
449
450 do i = 1, nresonances
451 write(iunit, '(a1,3e20.12)') '#', w(i), sqrt(abs(c0i2(i))), c0i2(i)
452 end do
453
454 write(iunit, '(a10,f12.6)') '# Gamma = ', gamma
455 write(iunit, '(a)') '#'
456
457 do i = 1, nfrequencies
458 e = (i - 1) * dw
459 pol = m_z0
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))
462 end do
463 write(iunit, '(3e20.12)') e, pol
464 end do
465
466 call io_close(iunit)
467
468 pop_sub(write_polarizability)
469 end subroutine write_polarizability
470 ! ---------------------------------------------------------
471
472
473 ! ---------------------------------------------------------
474 ! \todo This subroutine should be simplified.
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
479
480 integer :: i, ierr
481 real(real64) :: dw, w, aw, min_aw, min_w, omega_orig, min_w1, min_w2
482 real(real64), allocatable :: warray(:), tarray(:)
483
484 push_sub(find_resonance)
485
486 safe_allocate(warray(1:nfrequencies))
487 safe_allocate(tarray(1:nfrequencies))
488
489 warray = m_zero
490 tarray = m_zero
491 dw = (rightbound-leftbound) / (nfrequencies - 1)
492 do i = 1, nfrequencies
493 w = leftbound + (i-1)*dw
494 warray(i) = w
495 call ft2(w, aw)
496 tarray(i) = aw
497 end do
498
499 min_w = omega
500 min_aw = m_zero
501 do i = 1, nfrequencies
502 w = leftbound + (i-1)*dw
503 if (tarray(i) < min_aw) then
504 min_aw = tarray(i)
505 min_w = w
506 end if
507 end do
508
509 omega_orig = omega
510 omega = min_w
511 min_w1 = min_w - 2*dw
512 min_w2 = min_w + 2*dw
513 call loct_1dminimize(min_w1, min_w2, omega, ft2, ierr)
514
515 if (ierr /= 0) then
516 write(message(1),'(a)') 'Could not find a maximum.'
517 write(message(2),'(a)')
518 write(message(3), '(a,f12.8,a,f12.8,a)') ' Search interval = [', &
519 units_from_atomic(units%energy, leftbound), ',', units_from_atomic(units%energy, rightbound), ']'
520 write(message(4), '(a,f12.4,a)') ' Search discretization = ', &
521 units_from_atomic(units%energy, dw), ' '//trim(units_abbrev(units%energy))
522 call messages_fatal(4)
523 end if
524
525 safe_deallocate_a(warray)
526 safe_deallocate_a(tarray)
527
529 end subroutine find_resonance
530 ! ---------------------------------------------------------
531
532
533 ! ---------------------------------------------------------
534 subroutine resonance_first_order(omega, power, nw_subtracted, dw, leftbound, rightbound)
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
539
540 push_sub(resonance_first_order)
541
542 call ft(omega, power)
543
544 select case (mode)
545 case (sine_transform)
546 power = power / (m_one - sin(m_two*omega*total_time)/(m_two*omega*total_time))
547 case (cosine_transform)
548 call messages_not_implemented("resonance first order cosine transform")
549 end select
550
551 write(message(1), '(a)') '******************************************************************'
552 write(message(2), '(a,i3)') 'Resonance #', nw_subtracted + 1
553 write(message(3), '(a,f12.8,a,f12.8,a)') 'omega = ', units_from_atomic(units_out%energy, omega), &
554 ' '//trim(units_abbrev(units_out%energy))//' = ', omega, ' Ha'
555 write(message(4), '(a,f12.8,a,f12.8,a)') 'C(omega) = ', units_from_atomic(units_out%length**2, power), &
556 ' '//trim(units_abbrev(units_out%length**2))//' =', power, ' b^2'
557 write(message(5), '(a,f12.8,a,f12.8,a)') '<0|P|I> = ', units_from_atomic(units_out%length, sqrt(abs(power))), &
558 ' '//trim(units_abbrev(units_out%length))//' = ', sqrt(abs(power)),' b'
559 write(message(6), '(a,f12.8)') 'f[O->I] = ', m_two*omega*power
560 write(message(7), '(a)')
561 write(message(8), '(a,f12.8,a,f12.8,a)') ' Search interval = [', units_from_atomic(units_out%energy, leftbound), ',', &
562 units_from_atomic(units_out%energy, rightbound), ']'
563 write(message(9), '(a,f12.4,a)') ' Search discretization = ', units_from_atomic(units_out%energy, dw), &
564 ' '//trim(units_abbrev(units_out%energy))
565 write(message(10), '(a)') '******************************************************************'
566 write(message(11), '(a)')
567 call messages_info(11)
568
569 pop_sub(resonance_first_order)
571 end subroutine resonance_first_order
572 ! ---------------------------------------------------------
573
574
575 ! ---------------------------------------------------------
576 subroutine resonance_second_order(omega, power, nw_subtracted, leftbound, rightbound, c01, c02)
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
582
583 push_sub(resonance_second_order)
584
585 call ft(omega, power)
586 select case (mode)
587 case (sine_transform)
588 power = power / (m_one - sin(m_two*omega*total_time)/(m_two*omega*total_time))
589 case (cosine_transform)
590 ! WARNING: there is some difference between the omega=0 case and the rest.
591 if (abs(omega) > m_epsilon) then
592 power = power / (m_one + sin(m_two*omega*total_time)/(m_two*omega*total_time))
593 else
594 power = power / m_two
595 end if
596 end select
597
598 write(message(1), '(a)') '******************************************************************'
599 write(message(2), '(a,i3)') 'Resonance #', nw_subtracted + 1
600 write(message(3), '(a,f12.8,a,f12.8,a)') 'omega = ', units_from_atomic(units_out%energy, omega), &
601 ' '//trim(units_abbrev(units_out%energy))//' = ', omega, ' Ha'
602 write(message(4), '(a,f12.8,a,f12.8,a)') 'C(omega) = ', units_from_atomic(units_out%length**3, power), &
603 ' '//trim(units_abbrev(units_out%length**3))//' = ', power, ' b^3'
604 call messages_info(4)
605
606 if (abs(c01*c02) > m_epsilon) then
607 write(message(1), '(a,f12.8)') ' C(omega)/(C0i*C0j) = ', power / (c01 * c02)
608 call messages_info(1)
609 end if
610
611 write(message(1), '(a)')
612 write(message(2), '(a,f12.8,a,f12.8,a)') ' Search interval = [', units_from_atomic(units_out%energy, leftbound), ',', &
613 units_from_atomic(units_out%energy, rightbound), ']'
614 write(message(3), '(a)') '******************************************************************'
615 write(message(4), '(a)')
616 call messages_info(4)
617
619 end subroutine resonance_second_order
620 ! ---------------------------------------------------------
621
622
623 ! ---------------------------------------------------------
624 subroutine generate_signal(namespace, order, observable)
625 type(namespace_t), intent(in) :: namespace
626 integer, intent(in) :: order
627 integer, intent(in) :: observable(2)
628
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
633 type(unit_t) :: mp_unit
634 real(real64), allocatable :: q(:), mu(:), qq(:, :), c(:)
635 character(len=20) :: filename
636 type(kick_t) :: kick
637 type(unit_system_t) :: units
638 real(real64), allocatable :: multipole(:, :, :), ot(:), dipole(:, :)
639 type(local_operator_t) :: kick_operator
640 type(local_operator_t) :: obs
641
642 push_sub(generate_signal)
643
644 ! Find out how many files do we have
645 nfiles = 0
646 do
647 write(filename,'(a11,i1)') 'multipoles.', nfiles+1
648 inquire(file=trim(filename), exist = file_exists)
649 if (.not. file_exists) exit
650 nfiles = nfiles + 1
651 end do
652
653 ! WARNING: Check that order is smaller or equal to nfiles
654 if (nfiles == 0) then
655 write(message(1),'(a)') 'No multipoles.x file was found'
656 call messages_fatal(1, namespace=namespace)
657 end if
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.'
661 call messages_fatal(2, namespace=namespace)
662 end if
663
664 ! Open the files.
665 safe_allocate(iunit(1:nfiles))
666 do j = 1, nfiles
667 write(filename,'(a11,i1)') 'multipoles.', j
668 iunit(j) = io_open(trim(filename), namespace, action='read', status='old')
669 end do
670
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))
675
676 c = m_zero
677 c(order) = m_one
678
679 ! Get the basic info from the first file
680 call spectrum_mult_info(namespace, iunit(1), nspin, kick, time_steps, dt, units, lmax=lmax)
681
682 ! Sets the kick operator...
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)
692 end do
693 else
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
702 ! WARNING: not sure if m = -1 => y, and m = 1 => x. What is the convention?
703 kick_operator%weight(1) = -sqrt((m_four*m_pi)/m_three) * kick%pol(2, kick%pol_dir)
704 kick_operator%weight(2) = sqrt((m_four*m_pi)/m_three) * kick%pol(3, kick%pol_dir)
705 kick_operator%weight(3) = -sqrt((m_four*m_pi)/m_three) * kick%pol(1, kick%pol_dir)
706 end if
707
708 ! Sets the observation operator
709 select case (observable(1))
710 case (-1)
711 ! This means that the "observation operator" should be equal
712 ! to the "perturbation operator", i.e., the kick.
713 call local_operator_copy(obs, kick_operator)
714 case (0)
715 ! This means that the observable is the dipole operator; observable(2) determines
716 ! if it is x, y or z.
717 obs%n_multipoles = 1
718 safe_allocate(obs%l(1:1))
719 safe_allocate(obs%m(1:1))
720 safe_allocate(obs%weight(1:1))
721 obs%l(1) = 1
722 select case (observable(2))
723 case (1)
724 obs%m(1) = -1
725 obs%weight(1) = -sqrt((m_four*m_pi)/m_three)
726 case (2)
727 obs%m(1) = 1
728 obs%weight(1) = sqrt((m_four*m_pi)/m_three)
729 case (3)
730 obs%m(1) = 0
731 obs%weight(1) = -sqrt((m_four*m_pi)/m_three)
732 end select
733 case default
734 ! This means that the observation operator is (l,m) = (observable(1), observable(2))
735 obs%n_multipoles = 1
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)
742 end select
743
744 lambda = abs(kick%delta_strength)
745 q(1) = kick%delta_strength / lambda
746
747 do j = 2, nfiles
748 call spectrum_mult_info(namespace, iunit(j), nspin, kick, time_steps, dt, units, lmax=lmax)
749 q(j) = kick%delta_strength / lambda
750 end do
751
752 do i = 1, nfiles
753 do j = 1, nfiles
754 qq(i, j) = q(j)**i
755 end do
756 end do
757
758 call lalg_inverse(nfiles, qq, 'dir')
759
760 mu = matmul(qq, c)
761
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))
766 ! The units have nothing to do with the perturbing kick??
767 mp_unit = units%length**kick%l(1)
768 else
769 max_add_lm = 3
770 safe_allocate(multipole(1:3, 0:time_steps, 1:nspin))
771 mp_unit = units%length
772 end if
773 safe_allocate(ot(0:time_steps))
774 multipole = m_zero
775 ot = m_zero
776 o0 = m_zero
777
778 safe_allocate(dipole(1:3, 1:nspin))
779
780 do j = 1, nfiles
781 call io_skip_header(iunit(j))
782
783 do i = 0, time_steps
784 select case (nspin)
785 case (1)
786 read(iunit(j), *) k, dump, dump, (multipole(add_lm, i, 1), add_lm = 1, max_add_lm)
787 case (2)
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)
790 case (4)
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)
795 end select
796 multipole(1:max_add_lm, i, :) = units_to_atomic(mp_unit, multipole(1:max_add_lm, i, :))
797
798 ! The dipole is treated differently in the multipoles file: first of all,
799 ! the program should have printed *minus* the dipole operator.
800 dipole(1:3, 1:nspin) = - multipole(1:3, i, 1:nspin)
801 ! And then it contains the "cartesian" dipole, opposed to the spherical dipole:
802 multipole(1, i, 1:nspin) = -sqrt(m_three/(m_four*m_pi)) * dipole(2, 1:nspin)
803 multipole(2, i, 1:nspin) = sqrt(m_three/(m_four*m_pi)) * dipole(3, 1:nspin)
804 multipole(3, i, 1:nspin) = -sqrt(m_three/(m_four*m_pi)) * dipole(1, 1:nspin)
805
806 dump = m_zero
807 do k = 1, obs%n_multipoles
808 add_lm = 1
809 l = 1
810 lcycle: do
811 do m = -l, l
812 if (l == obs%l(k) .and. m == obs%m(k)) exit lcycle
813 add_lm = add_lm + 1
814 end do
815 l = l + 1
816 end do lcycle
817 ! Warning: it should not add the nspin components?
818 dump = dump + obs%weight(k) * sum(multipole(add_lm, i, 1:nspin))
819 end do
820
821 if (i == 0) o0 = dump
822 ot(i) = ot(i) + mu(j)*(dump - o0)
823 end do
824
825 end do
826
827 ot = ot / lambda**order
828
829 call write_ot(namespace, nspin, time_steps, dt, kick, order, ot(0:time_steps), observable)
830
831 ! Close files and exit.
832 do j = 1, nfiles
833 call io_close(iunit(j))
834 end do
835 safe_deallocate_a(iunit)
836 safe_deallocate_a(q)
837 safe_deallocate_a(mu)
838 safe_deallocate_a(qq)
839 safe_deallocate_a(c)
840 safe_deallocate_a(ot)
841 safe_deallocate_a(multipole)
842 safe_deallocate_a(dipole)
843
844 pop_sub(generate_signal)
845 end subroutine generate_signal
846 ! ---------------------------------------------------------
847
848
849 ! ---------------------------------------------------------
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
857
858 integer :: i
859
860 push_sub(modify_ot)
861
862 select case (mod(order, 2))
863 case (1)
864 do i = 0, time_steps
865 ot(i) = ot(i) - m_two * power * sin(omega*i*dt)
866 end do
867 case (0)
868 if (abs(omega) <= m_epsilon) then
869 do i = 0, time_steps
870 ot(i) = ot(i) - power * cos(omega*i*dt)
871 end do
872 else
873 do i = 0, time_steps
874 ot(i) = ot(i) - m_two * power * cos(omega*i*dt)
875 end do
876 end if
877 end select
878
879 pop_sub(modify_ot)
880 end subroutine modify_ot
881 ! ---------------------------------------------------------
882
883
884 ! ---------------------------------------------------------
885 subroutine write_ot(namespace, nspin, time_steps, dt, kick, order, ot, observable)
886 type(namespace_t), intent(in) :: namespace
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)
893
894 integer :: iunit, i
895 character(len=20) :: header_string
896 type(unit_t) :: ot_unit
897
898 push_sub(write_ot)
899
900 iunit = io_open('ot', namespace, action='write', status='replace')
901
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))
906 case (-1)
907 write(iunit,'(a)') '# Observable operator = kick operator'
908 if (kick%n_multipoles > 0) then
909 ot_unit = units_out%length**kick%l(1)
910 else
911 ot_unit = units_out%length
912 end if
913 case (0)
914 select case (observable(2))
915 case (1)
916 write(iunit,'(a)') '# O = x'
917 case (2)
918 write(iunit,'(a)') '# O = y'
919 case (3)
920 write(iunit,'(a)') '# O = z'
921 end select
922 ot_unit = units_out%length
923 case default
924 ot_unit = units_out%length**observable(1)
925 write(iunit, '(a12,i1,a1,i2,a1)') '# (l, m) = (', observable(1),',',observable(2),')'
926 end select
927 call kick_write(kick, iunit)
928
929 ! Units
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') '#'
934 write(header_string, '(a)') '['//trim(units_abbrev(units_out%time))//']'
935 write(iunit,'(a20)', advance = 'no') str_center(trim(header_string), 19)
936 write(header_string, '(a)') '['//trim(units_abbrev(units_out%length))//']'
937 write(iunit,'(a20)', advance = 'yes') str_center(trim(header_string), 20)
938
939 do i = 0, time_steps
940 write(iunit, '(2e20.8)') units_from_atomic(units_out%time, i*dt), units_from_atomic(ot_unit, ot(i))
941 end do
942
943 call io_close(iunit)
944 pop_sub(write_ot)
945 end subroutine write_ot
946
947
948 ! ---------------------------------------------------------
949 subroutine read_ot(namespace, nspin, order, nw_subtracted)
950 type(namespace_t), intent(in) :: namespace
951 integer, intent(out) :: nspin
952 integer, intent(out) :: order
953 integer, intent(out) :: nw_subtracted
954
955 integer :: iunit, i, ierr
956 character(len=100) :: line
957 character(len=12) :: dummychar
958 real(real64) :: dummy, t1, t2
959 type(unit_t) :: ot_unit
960
961 push_sub(read_ot)
962
963 if (allocated(ot)) then
964 safe_deallocate_a(ot)
965 end if
966
967 iunit = io_open('ot', namespace, action='read', status='old')
968
969 read(iunit, '(15x,i2)') nspin
970 read(iunit, '(15x,i2)') order
971 read(iunit, '(28x,i2)') nw_subtracted
972 read(iunit, '(a)') line
973
974 i = index(line, 'Observable')
975 if (index(line, 'Observable') /= 0) then
976 observable(1) = -1
977 elseif (index(line, '# O =') /= 0) then
978 observable(1) = 0
979 if (index(line,'x') /= 0) then
980 observable(2) = 1
981 elseif (index(line,'y') /= 0) then
982 observable(2) = 2
983 elseif (index(line,'z') /= 0) then
984 observable(2) = 3
985 end if
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)
988 else
989 write(message(1),'(a)') 'Problem reading "ot" file: could not figure out the shape'
990 write(message(2),'(a)') 'of the observation operator.'
991 call messages_fatal(2, namespace=namespace)
992 end if
993
994 call kick_read(kick, iunit, namespace)
995 read(iunit, '(a)') line
996 read(iunit, '(a)') line
997 call io_skip_header(iunit)
998
999 ! Figure out about the units of the file
1000 call unit_system_from_file(units, "ot", namespace, ierr)
1001 if (ierr /= 0) then
1002 write(message(1), '(a)') 'Could not figure out the units in file "ot".'
1003 call messages_fatal(1, namespace=namespace)
1004 end if
1005
1006 select case (observable(1))
1007 case (-1)
1008 if (kick%n_multipoles > 0) then
1009 ot_unit = units_out%length**kick%l(1)
1010 else
1011 ot_unit = units_out%length
1012 end if
1013 case (0)
1014 ot_unit = units_out%length
1015 case default
1016 ot_unit = units_out%length**observable(1)
1017 end select
1018
1019 ! count number of time points
1020 t1 = m_zero
1021 t2 = m_zero
1022 time_steps = 0
1023 do
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
1028 end do
1029100 continue
1030 if (time_steps <= 0) then
1031 write(message(1), '(a)') 'No data points found in file "ot".'
1032 call messages_fatal(1, namespace=namespace)
1033 end if
1034
1035 time_steps = time_steps - 1
1036 if (time_steps > 0) then
1037 dt = units_to_atomic(units%time, (t2 - t1))
1038 else
1039 dt = m_one
1040 end if
1041
1042 call io_skip_header(iunit)
1043
1044 safe_allocate(ot(0:time_steps))
1045
1046 do i = 0, time_steps
1047 read(iunit, *) dummy, ot(i)
1048 ot(i) = units_to_atomic(ot_unit, ot(i))
1049 end do
1050
1051 call io_close(iunit)
1052
1053 pop_sub(read_ot)
1054 end subroutine read_ot
1055
1056
1057 ! ---------------------------------------------------------
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
1064
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
1070
1071 push_sub(print_omega_file)
1072
1073 ! First, let us check that the file "ot" exists.
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)
1078 end if
1079
1080 ! Now, we should find out which units the file "ot" has.
1081 call unit_system_from_file(units, "ot", namespace, ierr)
1082 if (ierr /= 0) then
1083 write(message(1),'(a)') "Could not retrieve units in the 'ot' file."
1084 call messages_fatal(1, namespace=namespace)
1085 end if
1086
1087 if (omega > m_zero) then
1088 omega = units_to_atomic(units%energy, omega)
1089 else
1090 omega = m_half
1091 end if
1092
1093 if (search_interval > m_zero) then
1094 search_interval = units_to_atomic(units%energy, search_interval)
1095 else
1096 search_interval = m_half
1097 end if
1098
1099 leftbound = omega - search_interval
1100 rightbound = omega + search_interval
1101
1102 safe_allocate(warray(1:nfrequencies))
1103 safe_allocate(tarray(1:nfrequencies))
1104
1105 call read_ot(namespace, nspin, order, nw_subtracted)
1106
1107 if (mod(order, 2) == 1) then
1109 else
1111 end if
1112
1113 if (final_time > m_zero) then
1114 total_time = units_to_atomic(units%time, final_time)
1115 if (total_time > dt*time_steps) 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)
1121 end if
1122 time_steps = int(total_time / dt)
1124 else
1126 end if
1127
1128 warray = m_zero
1129 tarray = m_zero
1130 dw = (rightbound-leftbound) / (nfrequencies - 1)
1131 do i = 1, nfrequencies
1132 w = leftbound + (i-1)*dw
1133 warray(i) = w
1134 call ft(w, aw)
1135 tarray(i) = aw
1136 end do
1137
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)
1146 ! Here we should print the units of the transform.
1147 write(iunit, '(a)', advance = 'yes')
1148
1149 do i = 1, nfrequencies
1150 write(iunit,'(2e20.8)') units_from_atomic(units_out%energy, warray(i)), &
1151 tarray(i)
1152 end do
1154 call io_close(iunit)
1155
1156 safe_deallocate_a(warray)
1157 safe_deallocate_a(tarray)
1158 safe_deallocate_a(ot)
1159 call kick_end(kick)
1160
1161 pop_sub(print_omega_file)
1162 end subroutine print_omega_file
1163 ! ---------------------------------------------------------
1164
1166
1167! ---------------------------------------------------------
1168! ---------------------------------------------------------
1169! ---------------------------------------------------------
1170program oscillator_strength
1172 use global_oct_m
1173 use io_oct_m
1174 use messages_oct_m
1177
1178 implicit none
1179
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)
1189
1190 ! Reads the information passed through the command line options (if available).
1191 call getopt_init(ierr)
1192 if (ierr /= 0) then
1193 message(1) = "Your Fortran compiler doesn't support command-line arguments;"
1194 message(2) = "the oct-oscillator-strength command is not available."
1195 call messages_fatal(2)
1196 end if
1197
1198 ! Set the default values
1199 run_mode = analyze_nthorder_signal
1200 omega = - m_one
1201 search_interval = - m_one
1202 order = 1
1203 nfrequencies = 1000
1204 final_time = - m_one
1205 nresonances = 1
1206 observable(1) = -1
1207 observable(2) = 0
1208 ffile = ""
1209 damping = 0.1_real64/27.2114_real64 ! This is the usual damping factor used in the literature.
1210
1211 ! Get the parameters from the command line.
1212 cfile = c_null_char
1213 call getopt_oscillator_strength(run_mode, omega, search_interval, &
1214 order, nresonances, nfrequencies, final_time, &
1215 observable(1), observable(2), damping, cfile)
1216 call string_c_to_f(cfile, ffile)
1217 call getopt_end()
1218
1219 ! Initialize stuff
1221 call parser_init()
1222 call io_init(defaults = .true.)
1223 call unit_system_init(global_namespace)
1224
1225 select case (run_mode)
1226 case (generate_nthorder_signal)
1227 call generate_signal(global_namespace, order, observable)
1228 case (analyze_nthorder_signal)
1229 call analyze_signal(order, omega, search_interval, final_time, nresonances, nfrequencies, damping, &
1230 global_namespace)
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)
1235 case default
1236 end select
1237
1238 call io_end()
1239 call parser_end()
1240 call global_end()
1241
1242end program oscillator_strength
1243! ---------------------------------------------------------
1244
1245!! Local Variables:
1246!! mode: f90
1247!! coding: utf-8
1248!! End:
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
Definition: global.F90:193
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
Definition: global.F90:458
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_four
Definition: global.F90:195
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
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...
Definition: global.F90:284
subroutine, public init_octopus_globals(comm)
Initialise Octopus-specific global constants and files. This routine performs no initialisation calls...
Definition: global.F90:396
complex(real64), parameter, public m_z0
Definition: global.F90:201
complex(real64), parameter, public m_zi
Definition: global.F90:205
real(real64), parameter, public m_epsilon
Definition: global.F90:207
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_one
Definition: global.F90:192
real(real64), parameter, public m_three
Definition: global.F90:194
Definition: io.F90:116
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...
Definition: io.F90:165
subroutine, public io_close(iunit, grp)
Definition: io.F90:467
subroutine, public io_skip_header(iunit)
Definition: io.F90:646
subroutine, public io_end()
Definition: io.F90:271
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:402
subroutine, public kick_read(kick, iunit, namespace)
Definition: kick.F90:820
subroutine, public kick_end(kick)
Definition: kick.F90:796
subroutine, public kick_write(kick, iunit, out)
Definition: kick.F90:892
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1067
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:524
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:409
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:593
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)
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 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)
Definition: spectrum.F90:2347
character(len=80) function, public str_center(s_in, l_in)
puts space around string, so that it is centered
Definition: string.F90:199
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:225
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
int true(void)