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
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 ! First, subtract zero resonance...
270 leftbound = -search_interval
271 rightbound = search_interval
272 w = m_zero
273 call resonance_second_order(w, power, nw_subtracted, leftbound, rightbound, m_zero, m_zero)
274 call modify_ot(time_steps, dt, order, ot, w, power)
275 nw_subtracted = nw_subtracted + 1
276
277 ! Then, get all the others...
278 k = 1
279 do i = 1, nresonances
280 do j = i + 1, nresonances
281 leftbound = wij(k) - search_interval
282 rightbound = wij(k) + search_interval
283 call find_resonance(wij(k), leftbound, rightbound, nfrequencies)
284 call resonance_second_order(wij(k), power, nw_subtracted, leftbound, rightbound, c0i(i), c0i(j))
285 call modify_ot(time_steps, dt, order, ot, wij(k), power)
286 nw_subtracted = nw_subtracted + 1
287 k = k + 1
288 end do
289 end do
290
291 safe_deallocate_a(wij)
292 safe_deallocate_a(c0i)
293 safe_deallocate_a(omega)
294 safe_deallocate_a(ot)
295 call kick_end(kick)
296 call io_close(iunit)
297
298 pop_sub(read_resonances_file)
299 end subroutine read_resonances_file
300 ! ---------------------------------------------------------
301
302
303 ! ---------------------------------------------------------
304 subroutine analyze_signal(order, omega, search_interval, final_time, nresonances, nfrequencies, &
305 damping, namespace)
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
313 type(namespace_t), intent(in) :: namespace
314
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
319
320 push_sub(analyze_signal)
321
322 ! First, let us check that the file "ot" exists.
323 inquire(file="ot", exist = file_exists)
324 if (.not. file_exists) then
325 write(message(1),'(a)') "Could not find 'ot' file."
326 call messages_fatal(1)
327 end if
328
329 ! Now, we should find out which units the file "ot" has.
330 call unit_system_from_file(units, "ot", namespace, ierr)
331 if (ierr /= 0) then
332 write(message(1),'(a)') "Could not retrieve units in the 'ot' file."
333 call messages_fatal(1)
334 end if
335
336 if (omega > m_zero) then
337 omega = units_to_atomic(units%energy, omega)
338 else
339 omega = m_half
340 end if
341
342 if (search_interval > m_zero) then
343 search_interval = units_to_atomic(units%energy, search_interval)
344 else
345 search_interval = m_half
346 end if
347
348 leftbound = omega - search_interval
349 rightbound = omega + search_interval
350
351 call read_ot(namespace, nspin, order_in_file, nw_subtracted)
352
353 if (order_in_file /= order) then
354 write(message(1), '(a)') 'Internal error in analyze_signal'
355 call messages_fatal(1, namespace=namespace)
356 end if
357
358 if (mod(order, 2) == 1) then
359 mode = sine_transform
360 else
361 mode = cosine_transform
362 end if
363
364 if (final_time > m_zero) then
365 total_time = units_to_atomic(units%time, final_time)
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 ', &
370 units_from_atomic(units%time, total_time), units_abbrev(units%time)
371 call messages_warning(2, namespace=namespace)
372 end if
373 time_steps = int(total_time / dt)
374 total_time = time_steps * dt
375 else
376 total_time = dt*time_steps
377 end if
378
379 dw = (rightbound-leftbound) / (nfrequencies - 1)
380
381 safe_allocate( w(1:nresonances))
382 safe_allocate(c0i2(1:nresonances))
383 w = m_zero
384 c0i2 = m_zero
385
386 i = 1
387 do
388 if (nw_subtracted >= nresonances) exit
389
390 if (mode == cosine_transform .and. nw_subtracted == 0) then
391 omega = m_zero
392 else
393 call find_resonance(omega, leftbound, rightbound, nfrequencies)
394 end if
395
396 select case (order)
397 case (1)
398 call resonance_first_order(omega, power, nw_subtracted, dw, leftbound, rightbound)
399 case (2)
400 call resonance_second_order(omega, power, nw_subtracted, leftbound, rightbound, m_zero, m_zero)
401 end select
402
403 w(i) = omega
404 c0i2(i) = power
405
406 call modify_ot(time_steps, dt, order, ot, omega, power)
407
408 nw_subtracted = nw_subtracted + 1
409 i = i + 1
410 end do
411
412 select case (order)
413 case (1)
414 call write_polarizability(namespace, nfrequencies, nresonances, dw, w, c0i2, damping)
415 end select
416
417 safe_deallocate_a(ot)
418 safe_deallocate_a(w)
419 safe_deallocate_a(c0i2)
420 call kick_end(kick)
421
422 pop_sub(analyze_signal)
423 end subroutine analyze_signal
424 ! ---------------------------------------------------------
425
426
427 ! ---------------------------------------------------------
431 subroutine write_polarizability(namespace, nfrequencies, nresonances, dw, w, c0I2, gamma)
432 type(namespace_t), intent(in) :: namespace
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
437
438 integer :: iunit, i, j
439 real(real64) :: e
440 complex(real64) :: pol
441
442 push_sub(write_polarizability)
443
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)') '#'
447
448 do i = 1, nresonances
449 write(iunit, '(a1,3e20.12)') '#', w(i), sqrt(abs(c0i2(i))), c0i2(i)
450 end do
451
452 write(iunit, '(a10,f12.6)') '# Gamma = ', gamma
453 write(iunit, '(a)') '#'
454
455 do i = 1, nfrequencies
456 e = (i - 1) * dw
457 pol = m_z0
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))
460 end do
461 write(iunit, '(3e20.12)') e, pol
462 end do
463
464 call io_close(iunit)
465
466 pop_sub(write_polarizability)
467 end subroutine write_polarizability
468 ! ---------------------------------------------------------
469
470
471 ! ---------------------------------------------------------
472 ! \todo This subroutine should be simplified.
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
477
478 integer :: i, ierr
479 real(real64) :: dw, w, aw, min_aw, min_w, omega_orig, min_w1, min_w2
480 real(real64), allocatable :: warray(:), tarray(:)
481
482 push_sub(find_resonance)
483
484 safe_allocate(warray(1:nfrequencies))
485 safe_allocate(tarray(1:nfrequencies))
486
487 warray = m_zero
488 tarray = m_zero
489 dw = (rightbound-leftbound) / (nfrequencies - 1)
490 do i = 1, nfrequencies
491 w = leftbound + (i-1)*dw
492 warray(i) = w
493 call ft2(w, aw)
494 tarray(i) = aw
495 end do
496
497 min_w = omega
498 min_aw = m_zero
499 do i = 1, nfrequencies
500 w = leftbound + (i-1)*dw
501 if (tarray(i) < min_aw) then
502 min_aw = tarray(i)
503 min_w = w
504 end if
505 end do
506
507 omega_orig = omega
508 omega = min_w
509 min_w1 = min_w - 2*dw
510 min_w2 = min_w + 2*dw
511 call loct_1dminimize(min_w1, min_w2, omega, ft2, ierr)
512
513 if (ierr /= 0) then
514 write(message(1),'(a)') 'Could not find a maximum.'
515 write(message(2),'(a)')
516 write(message(3), '(a,f12.8,a,f12.8,a)') ' Search interval = [', &
517 units_from_atomic(units%energy, leftbound), ',', units_from_atomic(units%energy, rightbound), ']'
518 write(message(4), '(a,f12.4,a)') ' Search discretization = ', &
519 units_from_atomic(units%energy, dw), ' '//trim(units_abbrev(units%energy))
520 call messages_fatal(4)
521 end if
522
523 safe_deallocate_a(warray)
524 safe_deallocate_a(tarray)
525
527 end subroutine find_resonance
528 ! ---------------------------------------------------------
529
530
531 ! ---------------------------------------------------------
532 subroutine resonance_first_order(omega, power, nw_subtracted, dw, leftbound, rightbound)
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
537
538 push_sub(resonance_first_order)
539
540 call ft(omega, power)
541
542 select case (mode)
543 case (sine_transform)
544 power = power / (m_one - sin(m_two*omega*total_time)/(m_two*omega*total_time))
545 case (cosine_transform)
546 call messages_not_implemented("resonance first order cosine transform")
547 end select
548
549 write(message(1), '(a)') '******************************************************************'
550 write(message(2), '(a,i3)') 'Resonance #', nw_subtracted + 1
551 write(message(3), '(a,f12.8,a,f12.8,a)') 'omega = ', units_from_atomic(units_out%energy, omega), &
552 ' '//trim(units_abbrev(units_out%energy))//' = ', omega, ' Ha'
553 write(message(4), '(a,f12.8,a,f12.8,a)') 'C(omega) = ', units_from_atomic(units_out%length**2, power), &
554 ' '//trim(units_abbrev(units_out%length**2))//' =', power, ' b^2'
555 write(message(5), '(a,f12.8,a,f12.8,a)') '<0|P|I> = ', units_from_atomic(units_out%length, sqrt(abs(power))), &
556 ' '//trim(units_abbrev(units_out%length))//' = ', sqrt(abs(power)),' b'
557 write(message(6), '(a,f12.8)') 'f[O->I] = ', m_two*omega*power
558 write(message(7), '(a)')
559 write(message(8), '(a,f12.8,a,f12.8,a)') ' Search interval = [', units_from_atomic(units_out%energy, leftbound), ',', &
560 units_from_atomic(units_out%energy, rightbound), ']'
561 write(message(9), '(a,f12.4,a)') ' Search discretization = ', units_from_atomic(units_out%energy, dw), &
562 ' '//trim(units_abbrev(units_out%energy))
563 write(message(10), '(a)') '******************************************************************'
564 write(message(11), '(a)')
565 call messages_info(11)
566
567 pop_sub(resonance_first_order)
569 end subroutine resonance_first_order
570 ! ---------------------------------------------------------
571
572
573 ! ---------------------------------------------------------
574 subroutine resonance_second_order(omega, power, nw_subtracted, leftbound, rightbound, c01, c02)
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
580
581 push_sub(resonance_second_order)
582
583 call ft(omega, power)
584 select case (mode)
585 case (sine_transform)
586 power = power / (m_one - sin(m_two*omega*total_time)/(m_two*omega*total_time))
587 case (cosine_transform)
588 ! WARNING: there is some difference between the omega=0 case and the rest.
589 if (abs(omega) > m_epsilon) then
590 power = power / (m_one + sin(m_two*omega*total_time)/(m_two*omega*total_time))
591 else
592 power = power / m_two
593 end if
594 end select
595
596 write(message(1), '(a)') '******************************************************************'
597 write(message(2), '(a,i3)') 'Resonance #', nw_subtracted + 1
598 write(message(3), '(a,f12.8,a,f12.8,a)') 'omega = ', units_from_atomic(units_out%energy, omega), &
599 ' '//trim(units_abbrev(units_out%energy))//' = ', omega, ' Ha'
600 write(message(4), '(a,f12.8,a,f12.8,a)') 'C(omega) = ', units_from_atomic(units_out%length**3, power), &
601 ' '//trim(units_abbrev(units_out%length**3))//' = ', power, ' b^3'
602 call messages_info(4)
603
604 if (abs(c01*c02) > m_epsilon) then
605 write(message(1), '(a,f12.8)') ' C(omega)/(C0i*C0j) = ', power / (c01 * c02)
606 call messages_info(1)
607 end if
608
609 write(message(1), '(a)')
610 write(message(2), '(a,f12.8,a,f12.8,a)') ' Search interval = [', units_from_atomic(units_out%energy, leftbound), ',', &
611 units_from_atomic(units_out%energy, rightbound), ']'
612 write(message(3), '(a)') '******************************************************************'
613 write(message(4), '(a)')
614 call messages_info(4)
615
617 end subroutine resonance_second_order
618 ! ---------------------------------------------------------
619
620
621 ! ---------------------------------------------------------
622 subroutine generate_signal(namespace, order, observable)
623 type(namespace_t), intent(in) :: namespace
624 integer, intent(in) :: order
625 integer, intent(in) :: observable(2)
626
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
631 type(unit_t) :: mp_unit
632 real(real64), allocatable :: q(:), mu(:), qq(:, :), c(:)
633 character(len=20) :: filename
634 type(kick_t) :: kick
635 type(unit_system_t) :: units
636 real(real64), allocatable :: multipole(:, :, :), ot(:), dipole(:, :)
637 type(local_operator_t) :: kick_operator
638 type(local_operator_t) :: obs
639
640 push_sub(generate_signal)
641
642 ! Find out how many files do we have
643 nfiles = 0
644 do
645 write(filename,'(a11,i1)') 'multipoles.', nfiles+1
646 inquire(file=trim(filename), exist = file_exists)
647 if (.not. file_exists) exit
648 nfiles = nfiles + 1
649 end do
650
651 ! WARNING: Check that order is smaller or equal to nfiles
652 if (nfiles == 0) then
653 write(message(1),'(a)') 'No multipoles.x file was found'
654 call messages_fatal(1, namespace=namespace)
655 end if
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.'
659 call messages_fatal(2, namespace=namespace)
660 end if
661
662 ! Open the files.
663 safe_allocate(iunit(1:nfiles))
664 do j = 1, nfiles
665 write(filename,'(a11,i1)') 'multipoles.', j
666 iunit(j) = io_open(trim(filename), namespace, action='read', status='old')
667 end do
668
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))
673
674 c = m_zero
675 c(order) = m_one
676
677 ! Get the basic info from the first file
678 call spectrum_mult_info(namespace, iunit(1), nspin, kick, time_steps, dt, units, lmax=lmax)
679
680 ! Sets the kick operator...
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)
690 end do
691 else
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
700 ! WARNING: not sure if m = -1 => y, and m = 1 => x. What is the convention?
701 kick_operator%weight(1) = -sqrt((m_four*m_pi)/m_three) * kick%pol(2, kick%pol_dir)
702 kick_operator%weight(2) = sqrt((m_four*m_pi)/m_three) * kick%pol(3, kick%pol_dir)
703 kick_operator%weight(3) = -sqrt((m_four*m_pi)/m_three) * kick%pol(1, kick%pol_dir)
704 end if
705
706 ! Sets the observation operator
707 select case (observable(1))
708 case (-1)
709 ! This means that the "observation operator" should be equal
710 ! to the "perturbation operator", i.e., the kick.
711 call local_operator_copy(obs, kick_operator)
712 case (0)
713 ! This means that the observable is the dipole operator; observable(2) determines
714 ! if it is x, y or z.
715 obs%n_multipoles = 1
716 safe_allocate(obs%l(1:1))
717 safe_allocate(obs%m(1:1))
718 safe_allocate(obs%weight(1:1))
719 obs%l(1) = 1
720 select case (observable(2))
721 case (1)
722 obs%m(1) = -1
723 obs%weight(1) = -sqrt((m_four*m_pi)/m_three)
724 case (2)
725 obs%m(1) = 1
726 obs%weight(1) = sqrt((m_four*m_pi)/m_three)
727 case (3)
728 obs%m(1) = 0
729 obs%weight(1) = -sqrt((m_four*m_pi)/m_three)
730 end select
731 case default
732 ! This means that the observation operator is (l,m) = (observable(1), observable(2))
733 obs%n_multipoles = 1
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)
740 end select
741
742 lambda = abs(kick%delta_strength)
743 q(1) = kick%delta_strength / lambda
744
745 do j = 2, nfiles
746 call spectrum_mult_info(namespace, iunit(j), nspin, kick, time_steps, dt, units, lmax=lmax)
747 q(j) = kick%delta_strength / lambda
748 end do
749
750 do i = 1, nfiles
751 do j = 1, nfiles
752 qq(i, j) = q(j)**i
753 end do
754 end do
755
756 call lalg_inverse(nfiles, qq, 'dir')
757
758 mu = matmul(qq, c)
759
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))
764 ! The units have nothing to do with the perturbing kick??
765 mp_unit = units%length**kick%l(1)
766 else
767 max_add_lm = 3
768 safe_allocate(multipole(1:3, 0:time_steps, 1:nspin))
769 mp_unit = units%length
770 end if
771 safe_allocate(ot(0:time_steps))
772 multipole = m_zero
773 ot = m_zero
774 o0 = m_zero
775
776 safe_allocate(dipole(1:3, 1:nspin))
777
778 do j = 1, nfiles
779 call io_skip_header(iunit(j))
780
781 do i = 0, time_steps
782 select case (nspin)
783 case (1)
784 read(iunit(j), *) k, dump, dump, (multipole(add_lm, i, 1), add_lm = 1, max_add_lm)
785 case (2)
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)
788 case (4)
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)
793 end select
794 multipole(1:max_add_lm, i, :) = units_to_atomic(mp_unit, multipole(1:max_add_lm, i, :))
795
796 ! The dipole is treated differently in the multipoles file: first of all,
797 ! the program should have printed *minus* the dipole operator.
798 dipole(1:3, 1:nspin) = - multipole(1:3, i, 1:nspin)
799 ! And then it contains the "cartesian" dipole, opposed to the spherical dipole:
800 multipole(1, i, 1:nspin) = -sqrt(m_three/(m_four*m_pi)) * dipole(2, 1:nspin)
801 multipole(2, i, 1:nspin) = sqrt(m_three/(m_four*m_pi)) * dipole(3, 1:nspin)
802 multipole(3, i, 1:nspin) = -sqrt(m_three/(m_four*m_pi)) * dipole(1, 1:nspin)
803
804 dump = m_zero
805 do k = 1, obs%n_multipoles
806 add_lm = 1
807 l = 1
808 lcycle: do
809 do m = -l, l
810 if (l == obs%l(k) .and. m == obs%m(k)) exit lcycle
811 add_lm = add_lm + 1
812 end do
813 l = l + 1
814 end do lcycle
815 ! Warning: it should not add the nspin components?
816 dump = dump + obs%weight(k) * sum(multipole(add_lm, i, 1:nspin))
817 end do
818
819 if (i == 0) o0 = dump
820 ot(i) = ot(i) + mu(j)*(dump - o0)
821 end do
822
823 end do
824
825 ot = ot / lambda**order
826
827 call write_ot(namespace, nspin, time_steps, dt, kick, order, ot(0:time_steps), observable)
828
829 ! Close files and exit.
830 do j = 1, nfiles
831 call io_close(iunit(j))
832 end do
833 safe_deallocate_a(iunit)
834 safe_deallocate_a(q)
835 safe_deallocate_a(mu)
836 safe_deallocate_a(qq)
837 safe_deallocate_a(c)
838 safe_deallocate_a(ot)
839 safe_deallocate_a(multipole)
840 safe_deallocate_a(dipole)
841
842 pop_sub(generate_signal)
843 end subroutine generate_signal
844 ! ---------------------------------------------------------
845
846
847 ! ---------------------------------------------------------
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
855
856 integer :: i
857
858 push_sub(modify_ot)
859
860 select case (mod(order, 2))
861 case (1)
862 do i = 0, time_steps
863 ot(i) = ot(i) - m_two * power * sin(omega*i*dt)
864 end do
865 case (0)
866 if (abs(omega) <= m_epsilon) then
867 do i = 0, time_steps
868 ot(i) = ot(i) - power * cos(omega*i*dt)
869 end do
870 else
871 do i = 0, time_steps
872 ot(i) = ot(i) - m_two * power * cos(omega*i*dt)
873 end do
874 end if
875 end select
876
877 pop_sub(modify_ot)
878 end subroutine modify_ot
879 ! ---------------------------------------------------------
880
881
882 ! ---------------------------------------------------------
883 subroutine write_ot(namespace, nspin, time_steps, dt, kick, order, ot, observable)
884 type(namespace_t), intent(in) :: namespace
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)
891
892 integer :: iunit, i
893 character(len=20) :: header_string
894 type(unit_t) :: ot_unit
895
896 push_sub(write_ot)
897
898 iunit = io_open('ot', namespace, action='write', status='replace')
899
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))
904 case (-1)
905 write(iunit,'(a)') '# Observable operator = kick operator'
906 if (kick%n_multipoles > 0) then
907 ot_unit = units_out%length**kick%l(1)
908 else
909 ot_unit = units_out%length
910 end if
911 case (0)
912 select case (observable(2))
913 case (1)
914 write(iunit,'(a)') '# O = x'
915 case (2)
916 write(iunit,'(a)') '# O = y'
917 case (3)
918 write(iunit,'(a)') '# O = z'
919 end select
920 ot_unit = units_out%length
921 case default
922 ot_unit = units_out%length**observable(1)
923 write(iunit, '(a12,i1,a1,i2,a1)') '# (l, m) = (', observable(1),',',observable(2),')'
924 end select
925 call kick_write(kick, iunit)
926
927 ! Units
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') '#'
932 write(header_string, '(a)') '['//trim(units_abbrev(units_out%time))//']'
933 write(iunit,'(a20)', advance = 'no') str_center(trim(header_string), 19)
934 write(header_string, '(a)') '['//trim(units_abbrev(units_out%length))//']'
935 write(iunit,'(a20)', advance = 'yes') str_center(trim(header_string), 20)
936
937 do i = 0, time_steps
938 write(iunit, '(2e20.8)') units_from_atomic(units_out%time, i*dt), units_from_atomic(ot_unit, ot(i))
939 end do
940
941 call io_close(iunit)
942 pop_sub(write_ot)
943 end subroutine write_ot
944
945
946 ! ---------------------------------------------------------
947 subroutine read_ot(namespace, nspin, order, nw_subtracted)
948 type(namespace_t), intent(in) :: namespace
949 integer, intent(out) :: nspin
950 integer, intent(out) :: order
951 integer, intent(out) :: nw_subtracted
952
953 integer :: iunit, i, ierr
954 character(len=100) :: line
955 character(len=12) :: dummychar
956 real(real64) :: dummy, t1, t2
957 type(unit_t) :: ot_unit
958
959 push_sub(read_ot)
960
961 if (allocated(ot)) then
962 safe_deallocate_a(ot)
963 end if
964
965 iunit = io_open('ot', namespace, action='read', status='old')
966
967 read(iunit, '(15x,i2)') nspin
968 read(iunit, '(15x,i2)') order
969 read(iunit, '(28x,i2)') nw_subtracted
970 read(iunit, '(a)') line
971
972 i = index(line, 'Observable')
973 if (index(line, 'Observable') /= 0) then
974 observable(1) = -1
975 elseif (index(line, '# O =') /= 0) then
976 observable(1) = 0
977 if (index(line,'x') /= 0) then
978 observable(2) = 1
979 elseif (index(line,'y') /= 0) then
980 observable(2) = 2
981 elseif (index(line,'z') /= 0) then
982 observable(2) = 3
983 end if
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)
986 else
987 write(message(1),'(a)') 'Problem reading "ot" file: could not figure out the shape'
988 write(message(2),'(a)') 'of the observation operator.'
989 call messages_fatal(2, namespace=namespace)
990 end if
991
992 call kick_read(kick, iunit, namespace)
993 read(iunit, '(a)') line
994 read(iunit, '(a)') line
995 call io_skip_header(iunit)
996
997 ! Figure out about the units of the file
998 call unit_system_from_file(units, "ot", namespace, ierr)
999 if (ierr /= 0) then
1000 write(message(1), '(a)') 'Could not figure out the units in file "ot".'
1001 call messages_fatal(1, namespace=namespace)
1002 end if
1003
1004 select case (observable(1))
1005 case (-1)
1006 if (kick%n_multipoles > 0) then
1007 ot_unit = units_out%length**kick%l(1)
1008 else
1009 ot_unit = units_out%length
1010 end if
1011 case (0)
1012 ot_unit = units_out%length
1013 case default
1014 ot_unit = units_out%length**observable(1)
1015 end select
1016
1017 ! count number of time points
1018 t1 = m_zero
1019 t2 = m_zero
1020 time_steps = 0
1021 do
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
1026 end do
1027100 continue
1028 if (time_steps <= 0) then
1029 write(message(1), '(a)') 'No data points found in file "ot".'
1030 call messages_fatal(1, namespace=namespace)
1031 end if
1032
1033 time_steps = time_steps - 1
1034 if (time_steps > 0) then
1035 dt = units_to_atomic(units%time, (t2 - t1))
1036 else
1037 dt = m_one
1038 end if
1039
1040 call io_skip_header(iunit)
1041
1042 safe_allocate(ot(0:time_steps))
1043
1044 do i = 0, time_steps
1045 read(iunit, *) dummy, ot(i)
1046 ot(i) = units_to_atomic(ot_unit, ot(i))
1047 end do
1048
1049 call io_close(iunit)
1050
1051 pop_sub(read_ot)
1052 end subroutine read_ot
1053
1054
1055 ! ---------------------------------------------------------
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
1062
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
1068
1069 push_sub(print_omega_file)
1070
1071 ! First, let us check that the file "ot" exists.
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)
1076 end if
1077
1078 ! Now, we should find out which units the file "ot" has.
1079 call unit_system_from_file(units, "ot", namespace, ierr)
1080 if (ierr /= 0) then
1081 write(message(1),'(a)') "Could not retrieve units in the 'ot' file."
1082 call messages_fatal(1, namespace=namespace)
1083 end if
1084
1085 if (omega > m_zero) then
1086 omega = units_to_atomic(units%energy, omega)
1087 else
1088 omega = m_half
1089 end if
1090
1091 if (search_interval > m_zero) then
1092 search_interval = units_to_atomic(units%energy, search_interval)
1093 else
1094 search_interval = m_half
1095 end if
1096
1097 leftbound = omega - search_interval
1098 rightbound = omega + search_interval
1099
1100 safe_allocate(warray(1:nfrequencies))
1101 safe_allocate(tarray(1:nfrequencies))
1102
1103 call read_ot(namespace, nspin, order, nw_subtracted)
1104
1105 if (mod(order, 2) == 1) then
1107 else
1109 end if
1110
1111 if (final_time > m_zero) then
1112 total_time = units_to_atomic(units%time, final_time)
1113 if (total_time > dt*time_steps) 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)
1119 end if
1120 time_steps = int(total_time / dt)
1122 else
1124 end if
1125
1126 warray = m_zero
1127 tarray = m_zero
1128 dw = (rightbound-leftbound) / (nfrequencies - 1)
1129 do i = 1, nfrequencies
1130 w = leftbound + (i-1)*dw
1131 warray(i) = w
1132 call ft(w, aw)
1133 tarray(i) = aw
1134 end do
1135
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)
1144 ! Here we should print the units of the transform.
1145 write(iunit, '(a)', advance = 'yes')
1146
1147 do i = 1, nfrequencies
1148 write(iunit,'(2e20.8)') units_from_atomic(units_out%energy, warray(i)), &
1149 tarray(i)
1150 end do
1152 call io_close(iunit)
1153
1154 safe_deallocate_a(warray)
1155 safe_deallocate_a(tarray)
1156 safe_deallocate_a(ot)
1157 call kick_end(kick)
1158
1159 pop_sub(print_omega_file)
1160 end subroutine print_omega_file
1161 ! ---------------------------------------------------------
1162
1164
1165! ---------------------------------------------------------
1166! ---------------------------------------------------------
1167! ---------------------------------------------------------
1168program oscillator_strength
1170 use global_oct_m
1171 use io_oct_m
1172 use messages_oct_m
1175
1176 implicit none
1177
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)
1187
1188 ! Reads the information passed through the command line options (if available).
1189 call getopt_init(ierr)
1190 if (ierr /= 0) then
1191 message(1) = "Your Fortran compiler doesn't support command-line arguments;"
1192 message(2) = "the oct-oscillator-strength command is not available."
1193 call messages_fatal(2)
1194 end if
1195
1196 ! Set the default values
1197 run_mode = analyze_nthorder_signal
1198 omega = - m_one
1199 search_interval = - m_one
1200 order = 1
1201 nfrequencies = 1000
1202 final_time = - m_one
1203 nresonances = 1
1204 observable(1) = -1
1205 observable(2) = 0
1206 ffile = ""
1207 damping = 0.1_real64/27.2114_real64 ! This is the usual damping factor used in the literature.
1208
1209 ! Get the parameters from the command line.
1210 cfile = c_null_char
1211 call getopt_oscillator_strength(run_mode, omega, search_interval, &
1212 order, nresonances, nfrequencies, final_time, &
1213 observable(1), observable(2), damping, cfile)
1214 call string_c_to_f(cfile, ffile)
1215 call getopt_end()
1216
1217 ! Initialize stuff
1219 call parser_init()
1220 call io_init(defaults = .true.)
1221 call unit_system_init(global_namespace)
1222
1223 select case (run_mode)
1224 case (generate_nthorder_signal)
1225 call generate_signal(global_namespace, order, observable)
1226 case (analyze_nthorder_signal)
1227 call analyze_signal(order, omega, search_interval, final_time, nresonances, nfrequencies, damping, &
1228 global_namespace)
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)
1233 case default
1234 end select
1235
1236 call io_end()
1237 call parser_end()
1238 call global_end()
1239
1240end program oscillator_strength
1241! ---------------------------------------------------------
1242
1243!! Local Variables:
1244!! mode: f90
1245!! coding: utf-8
1246!! 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:202
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
Definition: global.F90:494
real(real64), parameter, public m_zero
Definition: global.F90:200
real(real64), parameter, public m_four
Definition: global.F90:204
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:198
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:294
subroutine, public init_octopus_globals(comm)
Initialise Octopus-specific global constants and files. This routine performs no initialisation calls...
Definition: global.F90:432
complex(real64), parameter, public m_z0
Definition: global.F90:210
complex(real64), parameter, public m_zi
Definition: global.F90:214
real(real64), parameter, public m_epsilon
Definition: global.F90:216
real(real64), parameter, public m_half
Definition: global.F90:206
real(real64), parameter, public m_one
Definition: global.F90:201
real(real64), parameter, public m_three
Definition: global.F90:203
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:1068
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
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:176
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)