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
31 use parser_oct_m
34 use string_oct_m
35 use unit_oct_m
37
38 implicit none
39
40 integer, parameter :: &
41 SINE_TRANSFORM = 1, &
43
45 integer :: n_multipoles
46 integer, allocatable :: l(:), m(:)
47 real(real64), allocatable :: weight(:)
48 end type local_operator_t
49
50 integer :: observable(2)
51 type(unit_system_t) :: units
52 real(real64), allocatable :: ot(:)
53 type(kick_t) :: kick
54 integer :: time_steps
55 real(real64) :: total_time
56 integer :: mode
57 real(real64) :: dt
58
59contains
60
61 ! ---------------------------------------------------------
62 subroutine ft(omega, power)
63 real(real64), intent(in) :: omega
64 real(real64), intent(out) :: power
65
66 integer :: j
67 real(real64) :: x
68 power = m_zero
69
70 push_sub(ft)
71
72 select case (mode)
73
74 case (sine_transform)
75
76 do j = 0, time_steps
77 x = sin(omega*j*dt)
78 power = power + x*ot(j)
79 end do
80 power = power*dt / (dt*time_steps)
81
82 case (cosine_transform)
83
84 do j = 0, time_steps
85 x = cos(omega*j*dt)
86 power = power + x*ot(j)
87 end do
88 power = power*dt / (dt*time_steps)
89
90 end select
91
92 pop_sub(ft)
93 end subroutine ft
94
95
96 ! ---------------------------------------------------------
97 subroutine ft2(omega, power)
98 real(real64), intent(in) :: omega
99 real(real64), intent(out) :: power
100
101 integer :: j
102 real(real64) :: x
103 power = m_zero
104
105 push_sub(ft2)
106
107 select case (mode)
108
109 case (sine_transform)
110
111 do j = 0, time_steps
112 x = sin(omega*j*dt)
113 power = power + x*ot(j)
114 end do
115 ! The function should be *minus* the sine Fourier transform, since this
116 ! is the function to be minimized.
117 power = - (power*dt / (dt*time_steps))**2
118
119 case (cosine_transform)
120
121 do j = 0, time_steps
122 x = cos(omega*j*dt)
123 power = power + x*ot(j)
124 end do
125 power = - (power*dt / (dt*time_steps))**2
126
127 end select
128
129 pop_sub(ft2)
130 end subroutine ft2
131
132
133 ! ---------------------------------------------------------
134 subroutine local_operator_copy(o, i)
135 type(local_operator_t), intent(inout) :: o
136 type(local_operator_t), intent(inout) :: i
138 integer :: j
141
142 o%n_multipoles = i%n_multipoles
143 safe_allocate( o%l(1:o%n_multipoles))
144 safe_allocate( o%m(1:o%n_multipoles))
145 safe_allocate(o%weight(1:o%n_multipoles))
147 do j = 1, o%n_multipoles
148 o%l(j) = i%l(j)
149 o%m(j) = i%m(j)
150 o%weight(j) = i%weight(j)
151 end do
152
153 pop_sub(local_operator_copy)
154 end subroutine local_operator_copy
156 ! ---------------------------------------------------------
157 subroutine read_resonances_file(order, ffile, namespace, search_interval, final_time, nfrequencies)
158 integer, intent(inout) :: order
159 character(len=*), intent(in) :: ffile
160 type(namespace_t), intent(in):: namespace
161 real(real64), intent(inout) :: search_interval
162 real(real64), intent(in) :: final_time
163 integer, intent(in) :: nfrequencies
164
165 real(real64) :: dummy, leftbound, rightbound, w, power
166 integer :: iunit, nresonances, ios, i, j, k, npairs, nspin, order_in_file, nw_subtracted, ierr
167 logical :: file_exists
168 real(real64), allocatable :: wij(:), omega(:), c0i(:)
169
170 push_sub(read_resonances_file)
171
172 if (order /= 2) then
173 write(message(1),'(a)') 'The run mode #3 is only compatible with the analysis of the'
174 write(message(2),'(a)') 'second-order response.'
175 call messages_fatal(2)
176 end if
177
178 ! First, let us check that the file "ot" exists.
179 inquire(file="ot", exist = file_exists)
180 if (.not. file_exists) then
181 write(message(1),'(a)') "Could not find 'ot' file."
182 call messages_fatal(1)
183 end if
184
185 ! Now, we should find out which units the file "ot" has.
186 call unit_system_from_file(units, "ot", namespace, ierr)
187 if (ierr /= 0) then
188 write(message(1),'(a)') "Could not retrieve units in the 'ot' file."
189 call messages_fatal(1)
190 end if
191
192 mode = cosine_transform
193
194 iunit = io_open(trim(ffile), namespace, action='read', status='old')
195
196 call io_skip_header(iunit)
197 ! Count the number of resonances
198 nresonances = 0
199 do
200 read(iunit, *, iostat = ios) dummy, dummy
201 if (ios /= 0) exit
202 nresonances = nresonances + 1
203 end do
204
205 npairs = (nresonances*(nresonances-1))/2
206
207 safe_allocate(omega(1:nresonances))
208 safe_allocate( c0i(1:nresonances))
209 safe_allocate(wij(1:npairs))
210
211 call io_skip_header(iunit)
212 do i = 1, nresonances
213 read(iunit, *) omega(i), c0i(i)
214 end do
215
216 k = 1
217 do i = 1, nresonances
218 do j = i + 1, nresonances
219 wij(k) = omega(j) - omega(i)
220 k = k + 1
221 end do
222 end do
223
224 if (search_interval > m_zero) then
225 search_interval = units_to_atomic(units%energy, search_interval)
226 else
227 search_interval = m_half
228 end if
229
230 call read_ot(namespace, nspin, order_in_file, nw_subtracted)
231
232 if (order_in_file /= order) then
233 write(message(1), '(a)') 'The ot file should contain the second-order response in this run mode.'
234 call messages_fatal(1)
235 end if
236
237 if (final_time > m_zero) then
238 total_time = units_to_atomic(units%time, final_time)
239 if (total_time > dt*time_steps) then
240 total_time = dt*time_steps
241 write(message(1), '(a)') 'The requested total time to process is larger than the time available in the input file.'
242 write(message(2), '(a,f8.4,a)') 'The time has been adjusted to ', &
243 units_from_atomic(units%time, total_time), units_abbrev(units%time)
244 call messages_warning(2, namespace=namespace)
245 end if
246 time_steps = int(total_time / dt)
247 total_time = time_steps * dt
248 else
249 total_time = dt*time_steps
250 end if
251
252 ! First, subtract zero resonance...
253 w = m_zero
254 call resonance_second_order(w, power, nw_subtracted, leftbound, rightbound, m_zero, m_zero)
255 call modify_ot(time_steps, dt, order, ot, w, power)
256 nw_subtracted = nw_subtracted + 1
257
258 ! Then, get all the others...
259 k = 1
260 do i = 1, nresonances
261 do j = i + 1, nresonances
262 leftbound = wij(k) - search_interval
263 rightbound = wij(k) + search_interval
264 call find_resonance(wij(k), leftbound, rightbound, nfrequencies)
265 call resonance_second_order(wij(k), power, nw_subtracted, leftbound, rightbound, c0i(i), c0i(j))
266 call modify_ot(time_steps, dt, order, ot, wij(k), power)
267 nw_subtracted = nw_subtracted + 1
268 k = k + 1
269 end do
270 end do
271
272 safe_deallocate_a(wij)
273 safe_deallocate_a(c0i)
274 safe_deallocate_a(omega)
275 call io_close(iunit)
276
277 pop_sub(read_resonances_file)
278 end subroutine read_resonances_file
279 ! ---------------------------------------------------------
280
281
282 ! ---------------------------------------------------------
283 subroutine analyze_signal(order, omega, search_interval, final_time, nresonances, nfrequencies, &
284 damping, namespace)
285 integer, intent(inout) :: order
286 real(real64), intent(inout) :: omega
287 real(real64), intent(inout) :: search_interval
288 real(real64), intent(inout) :: final_time
289 integer, intent(inout) :: nresonances
290 integer, intent(inout) :: nfrequencies
291 real(real64), intent(in) :: damping
292 type(namespace_t), intent(in) :: namespace
293
294 real(real64) :: leftbound, rightbound, dw, power
295 real(real64), allocatable :: w(:), c0I2(:)
296 integer :: nspin, i, ierr, order_in_file, nw_subtracted
297 logical :: file_exists
298
299 push_sub(analyze_signal)
300
301 ! First, let us check that the file "ot" exists.
302 inquire(file="ot", exist = file_exists)
303 if (.not. file_exists) then
304 write(message(1),'(a)') "Could not find 'ot' file."
305 call messages_fatal(1)
306 end if
307
308 ! Now, we should find out which units the file "ot" has.
309 call unit_system_from_file(units, "ot", namespace, ierr)
310 if (ierr /= 0) then
311 write(message(1),'(a)') "Could not retrieve units in the 'ot' file."
312 call messages_fatal(1)
313 end if
314
315 if (omega > m_zero) then
316 omega = units_to_atomic(units%energy, omega)
317 else
318 omega = m_half
319 end if
320
321 if (search_interval > m_zero) then
322 search_interval = units_to_atomic(units%energy, search_interval)
323 else
324 search_interval = m_half
325 end if
326
327 leftbound = omega - search_interval
328 rightbound = omega + search_interval
329
330 call read_ot(namespace, nspin, order_in_file, nw_subtracted)
331
332 if (order_in_file /= order) then
333 write(message(1), '(a)') 'Internal error in analyze_signal'
334 call messages_fatal(1, namespace=namespace)
335 end if
336
337 if (mod(order, 2) == 1) then
338 mode = sine_transform
339 else
340 mode = cosine_transform
341 end if
342
343 if (final_time > m_zero) then
344 total_time = units_to_atomic(units%time, final_time)
345 if (total_time > dt*time_steps) then
346 total_time = dt*time_steps
347 write(message(1), '(a)') 'The requested total time to process is larger than the time available in the input file.'
348 write(message(2), '(a,f8.4,a)') 'The time has been adjusted to ', &
349 units_from_atomic(units%time, total_time), units_abbrev(units%time)
350 call messages_warning(2, namespace=namespace)
351 end if
352 time_steps = int(total_time / dt)
353 total_time = time_steps * dt
354 else
355 total_time = dt*time_steps
356 end if
357
358 dw = (rightbound-leftbound) / (nfrequencies - 1)
359
360 safe_allocate( w(1:nresonances))
361 safe_allocate(c0i2(1:nresonances))
362 w = m_zero
363 c0i2 = m_zero
364
365 i = 1
366 do
367 if (nw_subtracted >= nresonances) exit
368
369 if (mode == cosine_transform .and. nw_subtracted == 0) then
370 omega = m_zero
371 else
372 call find_resonance(omega, leftbound, rightbound, nfrequencies)
373 end if
374
375 select case (order)
376 case (1)
377 call resonance_first_order(omega, power, nw_subtracted, dw, leftbound, rightbound)
378 case (2)
379 call resonance_second_order(omega, power, nw_subtracted, leftbound, rightbound, m_zero, m_zero)
380 end select
381
382 w(i) = omega
383 c0i2(i) = power
384
385 call modify_ot(time_steps, dt, order, ot, omega, power)
386
387 nw_subtracted = nw_subtracted + 1
388 i = i + 1
389 end do
390
391 select case (order)
392 case (1)
393 call write_polarizability(namespace, nfrequencies, nresonances, dw, w, c0i2, damping)
394 end select
395
396 safe_deallocate_a(ot)
397 safe_deallocate_a(w)
398 safe_deallocate_a(c0i2)
399
400 pop_sub(analyze_signal)
401 end subroutine analyze_signal
402 ! ---------------------------------------------------------
403
404
405 ! ---------------------------------------------------------
409 subroutine write_polarizability(namespace, nfrequencies, nresonances, dw, w, c0I2, gamma)
410 type(namespace_t), intent(in) :: namespace
411 integer, intent(in) :: nfrequencies, nresonances
412 real(real64), intent(in) :: dw
413 real(real64), intent(in) :: w(nresonances), c0I2(nresonances)
414 real(real64), intent(in) :: gamma
415
416 integer :: iunit, i, j
417 real(real64) :: e
418 complex(real64) :: pol
419
420 push_sub(write_polarizability)
421
422 iunit = io_open('polarizability', namespace, action = 'write', status='replace', die=.false.)
423 write(iunit, '(a)') '# Polarizability file. Generated using the SOS formula with the following data:'
424 write(iunit, '(a)') '#'
425
426 do i = 1, nresonances
427 write(iunit, '(a1,3e20.12)') '#', w(i), sqrt(abs(c0i2(i))), c0i2(i)
428 end do
429
430 write(iunit, '(a10,f12.6)') '# Gamma = ', gamma
431 write(iunit, '(a)') '#'
432
433 do i = 1, nfrequencies
434 e = (i - 1) * dw
435 pol = m_z0
436 do j = 1, nresonances
437 pol = pol + c0i2(j) * (m_one / (w(j) - e - m_zi * gamma) + m_one/(w(j) + e + m_zi * gamma))
438 end do
439 write(iunit, '(3e20.12)') e, pol
440 end do
441
442 call io_close(iunit)
443
444 pop_sub(write_polarizability)
445 end subroutine write_polarizability
446 ! ---------------------------------------------------------
447
448
449 ! ---------------------------------------------------------
450 ! \todo This subroutine should be simplified.
451 subroutine find_resonance(omega, leftbound, rightbound, nfrequencies)
452 real(real64), intent(inout) :: omega
453 real(real64), intent(in) :: leftbound, rightbound
454 integer, intent(in) :: nfrequencies
455
456 integer :: i, ierr
457 real(real64) :: dw, w, aw, min_aw, min_w, omega_orig, min_w1, min_w2
458 real(real64), allocatable :: warray(:), tarray(:)
459
460 push_sub(find_resonance)
461
462 safe_allocate(warray(1:nfrequencies))
463 safe_allocate(tarray(1:nfrequencies))
464
465 warray = m_zero
466 tarray = m_zero
467 dw = (rightbound-leftbound) / (nfrequencies - 1)
468 do i = 1, nfrequencies
469 w = leftbound + (i-1)*dw
470 warray(i) = w
471 call ft2(w, aw)
472 tarray(i) = aw
473 end do
474
475 min_w = omega
476 min_aw = m_zero
477 do i = 1, nfrequencies
478 w = leftbound + (i-1)*dw
479 if (tarray(i) < min_aw) then
480 min_aw = tarray(i)
481 min_w = w
482 end if
483 end do
484
485 omega_orig = omega
486 omega = min_w
487 min_w1 = min_w - 2*dw
488 min_w2 = min_w + 2*dw
489 call loct_1dminimize(min_w1, min_w2, omega, ft2, ierr)
490
491 if (ierr /= 0) then
492 write(message(1),'(a)') 'Could not find a maximum.'
493 write(message(2),'(a)')
494 write(message(3), '(a,f12.8,a,f12.8,a)') ' Search interval = [', &
495 units_from_atomic(units%energy, leftbound), ',', units_from_atomic(units%energy, rightbound), ']'
496 write(message(4), '(a,f12.4,a)') ' Search discretization = ', &
497 units_from_atomic(units%energy, dw), ' '//trim(units_abbrev(units%energy))
498 call messages_fatal(4)
499 end if
500
501 safe_deallocate_a(warray)
502 safe_deallocate_a(tarray)
503
504 pop_sub(find_resonance)
505 end subroutine find_resonance
506 ! ---------------------------------------------------------
507
508
509 ! ---------------------------------------------------------
510 subroutine resonance_first_order(omega, power, nw_subtracted, dw, leftbound, rightbound)
511 real(real64), intent(in) :: omega
512 real(real64), intent(out) :: power
513 integer, intent(in) :: nw_subtracted
514 real(real64), intent(in) :: dw, leftbound, rightbound
515
516 push_sub(resonance_first_order)
517
518 call ft(omega, power)
519
520 select case (mode)
521 case (sine_transform)
522 power = power / (m_one - sin(m_two*omega*total_time)/(m_two*omega*total_time))
523 case (cosine_transform)
524 call messages_not_implemented("resonance first order cosine transform")
525 end select
526
527 write(message(1), '(a)') '******************************************************************'
528 write(message(2), '(a,i3)') 'Resonance #', nw_subtracted + 1
529 write(message(3), '(a,f12.8,a,f12.8,a)') 'omega = ', units_from_atomic(units_out%energy, omega), &
530 ' '//trim(units_abbrev(units_out%energy))//' = ', omega, ' Ha'
531 write(message(4), '(a,f12.8,a,f12.8,a)') 'C(omega) = ', units_from_atomic(units_out%length**2, power), &
532 ' '//trim(units_abbrev(units_out%length**2))//' =', power, ' b^2'
533 write(message(5), '(a,f12.8,a,f12.8,a)') '<0|P|I> = ', units_from_atomic(units_out%length, sqrt(abs(power))), &
534 ' '//trim(units_abbrev(units_out%length))//' = ', sqrt(abs(power)),' b'
535 write(message(6), '(a,f12.8)') 'f[O->I] = ', m_two*omega*power
536 write(message(7), '(a)')
537 write(message(8), '(a,f12.8,a,f12.8,a)') ' Search interval = [', units_from_atomic(units_out%energy, leftbound), ',', &
538 units_from_atomic(units_out%energy, rightbound), ']'
539 write(message(9), '(a,f12.4,a)') ' Search discretization = ', units_from_atomic(units_out%energy, dw), &
540 ' '//trim(units_abbrev(units_out%energy))
541 write(message(10), '(a)') '******************************************************************'
542 write(message(11), '(a)')
543 call messages_info(11)
545 pop_sub(resonance_first_order)
546
547 end subroutine resonance_first_order
548 ! ---------------------------------------------------------
549
550
551 ! ---------------------------------------------------------
552 subroutine resonance_second_order(omega, power, nw_subtracted, leftbound, rightbound, c01, c02)
553 real(real64), intent(in) :: omega
554 real(real64), intent(out) :: power
555 integer, intent(in) :: nw_subtracted
556 real(real64), intent(in) :: leftbound, rightbound
557 real(real64), intent(in) :: c01, c02
558
559 push_sub(resonance_second_order)
560
561 call ft(omega, power)
562 select case (mode)
563 case (sine_transform)
564 power = power / (m_one - sin(m_two*omega*total_time)/(m_two*omega*total_time))
565 case (cosine_transform)
566 ! WARNING: there is some difference between the omega=0 case and the rest.
567 if (abs(omega) > m_epsilon) then
568 power = power / (m_one + sin(m_two*omega*total_time)/(m_two*omega*total_time))
569 else
570 power = power / m_two
571 end if
572 end select
573
574 write(message(1), '(a)') '******************************************************************'
575 write(message(2), '(a,i3)') 'Resonance #', nw_subtracted + 1
576 write(message(3), '(a,f12.8,a,f12.8,a)') 'omega = ', units_from_atomic(units_out%energy, omega), &
577 ' '//trim(units_abbrev(units_out%energy))//' = ', omega, ' Ha'
578 write(message(4), '(a,f12.8,a,f12.8,a)') 'C(omega) = ', units_from_atomic(units_out%length**3, power), &
579 ' '//trim(units_abbrev(units_out%length**3))//' = ', power, ' b^3'
580 call messages_info(4)
581
582 if (abs(c01*c02) > m_epsilon) then
583 write(message(1), '(a,f12.8)') ' C(omega)/(C0i*C0j) = ', power / (c01 * c02)
584 call messages_info(1)
585 end if
586
587 write(message(1), '(a)')
588 write(message(2), '(a,f12.8,a,f12.8,a)') ' Search interval = [', units_from_atomic(units_out%energy, leftbound), ',', &
589 units_from_atomic(units_out%energy, rightbound), ']'
590 write(message(3), '(a)') '******************************************************************'
591 write(message(4), '(a)')
592 call messages_info(4)
593
595 end subroutine resonance_second_order
596 ! ---------------------------------------------------------
597
598
599 ! ---------------------------------------------------------
600 subroutine generate_signal(namespace, order, observable)
601 type(namespace_t), intent(in) :: namespace
602 integer, intent(in) :: order
603 integer, intent(in) :: observable(2)
604
605 logical :: file_exists
606 integer :: i, j, nspin, time_steps, lmax, nfiles, k, add_lm, l, m, max_add_lm
607 integer, allocatable :: iunit(:)
608 real(real64) :: dt, lambda, dump, o0
609 type(unit_t) :: mp_unit
610 real(real64), allocatable :: q(:), mu(:), qq(:, :), c(:)
611 character(len=20) :: filename
612 type(kick_t) :: kick
613 type(unit_system_t) :: units
614 real(real64), allocatable :: multipole(:, :, :), ot(:), dipole(:, :)
615 type(local_operator_t) :: kick_operator
616 type(local_operator_t) :: obs
617
618 push_sub(generate_signal)
619
620 ! Find out how many files do we have
621 nfiles = 0
622 do
623 write(filename,'(a11,i1)') 'multipoles.', nfiles+1
624 inquire(file=trim(filename), exist = file_exists)
625 if (.not. file_exists) exit
626 nfiles = nfiles + 1
627 end do
628
629 ! WARNING: Check that order is smaller or equal to nfiles
630 if (nfiles == 0) then
631 write(message(1),'(a)') 'No multipoles.x file was found'
632 call messages_fatal(1, namespace=namespace)
633 end if
634 if (order > nfiles) then
635 write(message(1),'(a)') 'The order that you ask for is higher than the number'
636 write(message(2),'(a)') 'of multipoles.x file that you supply.'
637 call messages_fatal(2, namespace=namespace)
638 end if
639
640 ! Open the files.
641 safe_allocate(iunit(1:nfiles))
642 do j = 1, nfiles
643 write(filename,'(a11,i1)') 'multipoles.', j
644 iunit(j) = io_open(trim(filename), namespace, action='read', status='old')
645 end do
646
647 safe_allocate( q(1:nfiles))
648 safe_allocate(mu(1:nfiles))
649 safe_allocate(qq(1:nfiles, 1:nfiles))
650 safe_allocate( c(1:nfiles))
651
652 c = m_zero
653 c(order) = m_one
654
655 ! Get the basic info from the first file
656 call spectrum_mult_info(namespace, iunit(1), nspin, kick, time_steps, dt, units, lmax=lmax)
657
658 ! Sets the kick operator...
659 if (kick%n_multipoles > 0) then
660 kick_operator%n_multipoles = kick%n_multipoles
661 safe_allocate( kick_operator%l(1:kick_operator%n_multipoles))
662 safe_allocate( kick_operator%m(1:kick_operator%n_multipoles))
663 safe_allocate(kick_operator%weight(1:kick_operator%n_multipoles))
664 do i = 1, kick_operator%n_multipoles
665 kick_operator%l(i) = kick%l(i)
666 kick_operator%m(i) = kick%m(i)
667 kick_operator%weight(i) = kick%weight(i)
668 end do
669 else
670 kick_operator%n_multipoles = 3
671 safe_allocate( kick_operator%l(1:kick_operator%n_multipoles))
672 safe_allocate( kick_operator%m(1:kick_operator%n_multipoles))
673 safe_allocate(kick_operator%weight(1:kick_operator%n_multipoles))
674 kick_operator%l(1:3) = 1
675 kick_operator%m(1) = -1
676 kick_operator%m(2) = 0
677 kick_operator%m(3) = 1
678 ! WARNING: not sure if m = -1 => y, and m = 1 => x. What is the convention?
679 kick_operator%weight(1) = -sqrt((m_four*m_pi)/m_three) * kick%pol(2, kick%pol_dir)
680 kick_operator%weight(2) = sqrt((m_four*m_pi)/m_three) * kick%pol(3, kick%pol_dir)
681 kick_operator%weight(3) = -sqrt((m_four*m_pi)/m_three) * kick%pol(1, kick%pol_dir)
682 end if
683
684 ! Sets the observation operator
685 select case (observable(1))
686 case (-1)
687 ! This means that the "observation operator" should be equal
688 ! to the "perturbation operator", i.e., the kick.
689 call local_operator_copy(obs, kick_operator)
690 case (0)
691 ! This means that the observable is the dipole operator; observable(2) determines
692 ! if it is x, y or z.
693 obs%n_multipoles = 1
694 safe_allocate(obs%l(1:1))
695 safe_allocate(obs%m(1:1))
696 safe_allocate(obs%weight(1:1))
697 obs%l(1) = 1
698 select case (observable(2))
699 case (1)
700 obs%m(1) = -1
701 obs%weight(1) = -sqrt((m_four*m_pi)/m_three)
702 case (2)
703 obs%m(1) = 1
704 obs%weight(1) = sqrt((m_four*m_pi)/m_three)
705 case (3)
706 obs%m(1) = 0
707 obs%weight(1) = -sqrt((m_four*m_pi)/m_three)
708 end select
709 case default
710 ! This means that the observation operator is (l,m) = (observable(1), observable(2))
711 obs%n_multipoles = 1
712 safe_allocate(obs%l(1:1))
713 safe_allocate(obs%m(1:1))
714 safe_allocate(obs%weight(1:1))
715 obs%weight(1) = m_one
716 obs%l(1) = observable(1)
717 obs%m(1) = observable(2)
718 end select
719
720 lambda = abs(kick%delta_strength)
721 q(1) = kick%delta_strength / lambda
722
723 do j = 2, nfiles
724 call spectrum_mult_info(namespace, iunit(j), nspin, kick, time_steps, dt, units, lmax=lmax)
725 q(j) = kick%delta_strength / lambda
726 end do
727
728 do i = 1, nfiles
729 do j = 1, nfiles
730 qq(i, j) = q(j)**i
731 end do
732 end do
733
734 call lalg_inverse(nfiles, qq, 'dir')
735
736 mu = matmul(qq, c)
737
738 if (kick%n_multipoles > 0) then
739 lmax = maxval(kick%l(1:obs%n_multipoles))
740 max_add_lm = (lmax+1)**2-1
741 safe_allocate(multipole(1:max_add_lm, 0:time_steps, 1:nspin))
742 ! The units have nothing to do with the perturbing kick??
743 mp_unit = units%length**kick%l(1)
744 else
745 max_add_lm = 3
746 safe_allocate(multipole(1:3, 0:time_steps, 1:nspin))
747 mp_unit = units%length
748 end if
749 safe_allocate(ot(0:time_steps))
750 multipole = m_zero
751 ot = m_zero
752 o0 = m_zero
753
754 safe_allocate(dipole(1:3, 1:nspin))
755
756 do j = 1, nfiles
757 call io_skip_header(iunit(j))
758
759 do i = 0, time_steps
760 select case (nspin)
761 case (1)
762 read(iunit(j), *) k, dump, dump, (multipole(add_lm, i, 1), add_lm = 1, max_add_lm)
763 case (2)
764 read(iunit(j), *) k, dump, dump, (multipole(add_lm, i, 1), add_lm = 1, max_add_lm), &
765 dump, (multipole(add_lm, i, 2), add_lm = 1, max_add_lm)
766 case (4)
767 read(iunit(j), *) k, dump, dump, (multipole(add_lm, i, 1), add_lm = 1, max_add_lm), &
768 dump, (multipole(add_lm, i, 2), add_lm = 1, max_add_lm), &
769 dump, (multipole(add_lm, i, 3), add_lm = 1, max_add_lm), &
770 dump, (multipole(add_lm, i, 4), add_lm = 1, max_add_lm)
771 end select
772 multipole(1:max_add_lm, i, :) = units_to_atomic(mp_unit, multipole(1:max_add_lm, i, :))
773
774 ! The dipole is treated differently in the multipoles file: first of all,
775 ! the program should have printed *minus* the dipole operator.
776 dipole(1:3, 1:nspin) = - multipole(1:3, i, 1:nspin)
777 ! And then it contains the "cartesian" dipole, opposed to the spherical dipole:
778 multipole(1, i, 1:nspin) = -sqrt(m_three/(m_four*m_pi)) * dipole(2, 1:nspin)
779 multipole(2, i, 1:nspin) = sqrt(m_three/(m_four*m_pi)) * dipole(3, 1:nspin)
780 multipole(3, i, 1:nspin) = -sqrt(m_three/(m_four*m_pi)) * dipole(1, 1:nspin)
781
782 dump = m_zero
783 do k = 1, obs%n_multipoles
784 add_lm = 1
785 l = 1
786 lcycle: do
787 do m = -l, l
788 if (l == obs%l(k) .and. m == obs%m(k)) exit lcycle
789 add_lm = add_lm + 1
790 end do
791 l = l + 1
792 end do lcycle
793 ! Warning: it should not add the nspin components?
794 dump = dump + obs%weight(k) * sum(multipole(add_lm, i, 1:nspin))
795 end do
796
797 if (i == 0) o0 = dump
798 ot(i) = ot(i) + mu(j)*(dump - o0)
799 end do
800
801 end do
802
803 ot = ot / lambda**order
804
805 call write_ot(namespace, nspin, time_steps, dt, kick, order, ot(0:time_steps), observable)
806
807 ! Close files and exit.
808 do j = 1, nfiles
809 call io_close(iunit(j))
810 end do
811 safe_deallocate_a(iunit)
812 safe_deallocate_a(q)
813 safe_deallocate_a(mu)
814 safe_deallocate_a(qq)
815 safe_deallocate_a(c)
816 safe_deallocate_a(ot)
817 safe_deallocate_a(multipole)
818 safe_deallocate_a(dipole)
819
820 pop_sub(generate_signal)
821 end subroutine generate_signal
822 ! ---------------------------------------------------------
823
824
825 ! ---------------------------------------------------------
826 subroutine modify_ot(time_steps, dt, order, ot, omega, power)
827 integer, intent(in) :: time_steps
828 real(real64), intent(in) :: dt
829 integer, intent(in) :: order
830 real(real64), intent(inout) :: ot(0:time_steps)
831 real(real64), intent(in) :: omega
832 real(real64), intent(in) :: power
833
834 integer :: i
835
836 push_sub(modify_ot)
837
838 select case (mod(order, 2))
839 case (1)
840 do i = 0, time_steps
841 ot(i) = ot(i) - m_two * power * sin(omega*i*dt)
842 end do
843 case (0)
844 if (abs(omega) <= m_epsilon) then
845 do i = 0, time_steps
846 ot(i) = ot(i) - power * cos(omega*i*dt)
847 end do
848 else
849 do i = 0, time_steps
850 ot(i) = ot(i) - m_two * power * cos(omega*i*dt)
851 end do
852 end if
853 end select
854
855 pop_sub(modify_ot)
856 end subroutine modify_ot
857 ! ---------------------------------------------------------
858
859
860 ! ---------------------------------------------------------
861 subroutine write_ot(namespace, nspin, time_steps, dt, kick, order, ot, observable)
862 type(namespace_t), intent(in) :: namespace
863 integer, intent(in) :: nspin, time_steps
864 real(real64), intent(in) :: dt
865 type(kick_t), intent(in) :: kick
866 integer, intent(in) :: order
867 real(real64), intent(in) :: ot(0:time_steps)
868 integer, intent(in) :: observable(2)
869
870 integer :: iunit, i
871 character(len=20) :: header_string
872 type(unit_t) :: ot_unit
873
874 push_sub(write_ot)
875
876 iunit = io_open('ot', namespace, action='write', status='replace')
877
878 write(iunit, '(a15,i2)') '# nspin ', nspin
879 write(iunit, '(a15,i2)') '# Order ', order
880 write(iunit, '(a28,i2)') '# Frequencies subtracted = ', 0
881 select case (observable(1))
882 case (-1)
883 write(iunit,'(a)') '# Observable operator = kick operator'
884 if (kick%n_multipoles > 0) then
885 ot_unit = units_out%length**kick%l(1)
886 else
887 ot_unit = units_out%length
888 end if
889 case (0)
890 select case (observable(2))
891 case (1)
892 write(iunit,'(a)') '# O = x'
893 case (2)
894 write(iunit,'(a)') '# O = y'
895 case (3)
896 write(iunit,'(a)') '# O = z'
897 end select
898 ot_unit = units_out%length
899 case default
900 ot_unit = units_out%length**observable(1)
901 write(iunit, '(a12,i1,a1,i2,a1)') '# (l, m) = (', observable(1),',',observable(2),')'
902 end select
903 call kick_write(kick, iunit)
904
905 ! Units
906 write(iunit,'(a1)', advance = 'no') '#'
907 write(iunit,'(a20)', advance = 'no') str_center('t', 19)
908 write(iunit,'(a20)', advance = 'yes') str_center('<O>(t)', 20)
909 write(iunit,'(a1)', advance = 'no') '#'
910 write(header_string, '(a)') '['//trim(units_abbrev(units_out%time))//']'
911 write(iunit,'(a20)', advance = 'no') str_center(trim(header_string), 19)
912 write(header_string, '(a)') '['//trim(units_abbrev(units_out%length))//']'
913 write(iunit,'(a20)', advance = 'yes') str_center(trim(header_string), 20)
914
915 do i = 0, time_steps
916 write(iunit, '(2e20.8)') units_from_atomic(units_out%time, i*dt), units_from_atomic(ot_unit, ot(i))
917 end do
918
919 call io_close(iunit)
920 pop_sub(write_ot)
921 end subroutine write_ot
922
923
924 ! ---------------------------------------------------------
925 subroutine read_ot(namespace, nspin, order, nw_subtracted)
926 type(namespace_t), intent(in) :: namespace
927 integer, intent(out) :: nspin
928 integer, intent(out) :: order
929 integer, intent(out) :: nw_subtracted
930
931 integer :: iunit, i, ierr
932 character(len=100) :: line
933 character(len=12) :: dummychar
934 real(real64) :: dummy, t1, t2
935 type(unit_t) :: ot_unit
936
937 push_sub(read_ot)
938
939 iunit = io_open('ot', namespace, action='read', status='old')
940
941 read(iunit, '(15x,i2)') nspin
942 read(iunit, '(15x,i2)') order
943 read(iunit, '(28x,i2)') nw_subtracted
944 read(iunit, '(a)') line
945
946 i = index(line, 'Observable')
947 if (index(line, 'Observable') /= 0) then
948 observable(1) = -1
949 elseif (index(line, '# O =') /= 0) then
950 observable(1) = 0
951 if (index(line,'x') /= 0) then
952 observable(2) = 1
953 elseif (index(line,'y') /= 0) then
954 observable(2) = 2
955 elseif (index(line,'z') /= 0) then
956 observable(2) = 3
957 end if
958 elseif (index(line, '# (l, m) = ') /= 0) then
959 read(line,'(a12,i1,a1,i2,a1)') dummychar(1:12), observable(1), dummychar(1:1), observable(2), dummychar(1:1)
960 else
961 write(message(1),'(a)') 'Problem reading "ot" file: could not figure out the shape'
962 write(message(2),'(a)') 'of the observation operator.'
963 call messages_fatal(2, namespace=namespace)
964 end if
965
966 call kick_read(kick, iunit, namespace)
967 read(iunit, '(a)') line
968 read(iunit, '(a)') line
969 call io_skip_header(iunit)
970
971 ! Figure out about the units of the file
972 call unit_system_from_file(units, "ot", namespace, ierr)
973 if (ierr /= 0) then
974 write(message(1), '(a)') 'Could not figure out the units in file "ot".'
975 call messages_fatal(1, namespace=namespace)
976 end if
977
978 select case (observable(1))
979 case (-1)
980 if (kick%n_multipoles > 0) then
981 ot_unit = units_out%length**kick%l(1)
982 else
983 ot_unit = units_out%length
984 end if
985 case (0)
986 ot_unit = units_out%length
987 case default
988 ot_unit = units_out%length**observable(1)
989 end select
990
991 ! count number of time_steps
992 time_steps = 0
993 do
994 read(iunit, *, end=100) dummy
995 time_steps = time_steps + 1
996 if (time_steps == 1) t1 = dummy
997 if (time_steps == 2) t2 = dummy
998 end do
999100 continue
1000 dt = units_to_atomic(units%time, (t2 - t1)) ! units_out is OK
1001
1002 call io_skip_header(iunit)
1003
1004 safe_allocate(ot(0:time_steps))
1005
1006 do i = 0, time_steps
1007 read(iunit, *) dummy, ot(i)
1008 ot(i) = units_to_atomic(ot_unit, ot(i))
1009 end do
1010
1011 pop_sub(read_ot)
1012 end subroutine read_ot
1013
1014
1015 ! ---------------------------------------------------------
1016 subroutine print_omega_file(namespace, omega, search_interval, final_time, nfrequencies)
1017 type(namespace_t), intent(in) :: namespace
1018 real(real64), intent(inout) :: omega
1019 real(real64), intent(inout) :: search_interval
1020 real(real64), intent(inout) :: final_time
1021 integer, intent(inout) :: nfrequencies
1022
1023 integer :: iunit, i, ierr, nspin, order, nw_subtracted
1024 logical :: file_exists
1025 character(len=20) :: header_string
1026 real(real64), allocatable :: warray(:), tarray(:)
1027 real(real64) :: leftbound, rightbound, dw, w, aw
1028
1029 push_sub(print_omega_file)
1030
1031 ! First, let us check that the file "ot" exists.
1032 inquire(file="ot", exist = file_exists)
1033 if (.not. file_exists) then
1034 write(message(1),'(a)') "Could not find 'ot' file."
1035 call messages_fatal(1, namespace=namespace)
1036 end if
1037
1038 ! Now, we should find out which units the file "ot" has.
1039 call unit_system_from_file(units, "ot", namespace, ierr)
1040 if (ierr /= 0) then
1041 write(message(1),'(a)') "Could not retrieve units in the 'ot' file."
1042 call messages_fatal(1, namespace=namespace)
1043 end if
1044
1045 if (omega > m_zero) then
1046 omega = units_to_atomic(units%energy, omega)
1047 else
1048 omega = m_half
1049 end if
1050
1051 if (search_interval > m_zero) then
1052 search_interval = units_to_atomic(units%energy, search_interval)
1053 else
1054 search_interval = m_half
1055 end if
1056
1057 leftbound = omega - search_interval
1058 rightbound = omega + search_interval
1059
1060 safe_allocate(warray(1:nfrequencies))
1061 safe_allocate(tarray(1:nfrequencies))
1062
1063 call read_ot(namespace, nspin, order, nw_subtracted)
1064
1065 if (mod(order, 2) == 1) then
1067 else
1069 end if
1070
1071 if (final_time > m_zero) then
1072 total_time = units_to_atomic(units%time, final_time)
1073 if (total_time > dt*time_steps) then
1075 write(message(1), '(a)') 'The requested total time to process is larger than the time available in the input file.'
1076 write(message(2), '(a,f8.4,a)') 'The time has been adjusted to ', &
1077 units_from_atomic(units_out%time, total_time), trim(units_abbrev(units_out%time))
1078 call messages_warning(2, namespace=namespace)
1079 end if
1080 time_steps = int(total_time / dt)
1082 else
1084 end if
1085
1086 warray = m_zero
1087 tarray = m_zero
1088 dw = (rightbound-leftbound) / (nfrequencies - 1)
1089 do i = 1, nfrequencies
1090 w = leftbound + (i-1)*dw
1091 warray(i) = w
1092 call ft(w, aw)
1093 tarray(i) = aw
1094 end do
1095
1096 iunit = io_open('omega', namespace, action='write', status='replace')
1097 write(iunit, '(a15,i2)') '# nspin ', nspin
1098 call kick_write(kick, iunit)
1099 write(iunit, '(a)') '#%'
1100 write(iunit, '(a1,a20)', advance = 'no') '#', str_center("omega", 20)
1101 write(header_string,'(a)') 'F(omega)'
1102 write(iunit, '(a20)', advance = 'yes') str_center(trim(header_string), 20)
1103 write(iunit, '(a1,a20)', advance = 'no') '#', str_center('['//trim(units_abbrev(units_out%energy)) // ']', 20)
1104 ! Here we should print the units of the transform.
1105 write(iunit, '(a)', advance = 'yes')
1106
1107 do i = 1, nfrequencies
1108 write(iunit,'(2e20.8)') units_from_atomic(units_out%energy, warray(i)), &
1109 tarray(i)
1110 end do
1111
1112 safe_deallocate_a(warray)
1113 safe_deallocate_a(tarray)
1114
1115 pop_sub(print_omega_file)
1116 end subroutine print_omega_file
1117 ! ---------------------------------------------------------
1118
1120
1121! ---------------------------------------------------------
1122! ---------------------------------------------------------
1123! ---------------------------------------------------------
1124program oscillator_strength
1126 use global_oct_m
1127 use io_oct_m
1128 use messages_oct_m
1130
1131 implicit none
1132
1133 integer :: run_mode, order, nfrequencies, ierr, nresonances
1134 real(real64) :: omega, search_interval, final_time, damping
1135 integer, parameter :: &
1136 ANALYZE_NTHORDER_SIGNAL = 1, &
1137 generate_nthorder_signal = 2, &
1138 read_resonances_from_file = 3, &
1139 generate_omega_file = 4
1140 character(len=100) :: ffile
1141
1142 ! Reads the information passed through the command line options (if available).
1143 call getopt_init(ierr)
1144 if (ierr /= 0) then
1145 message(1) = "Your Fortran compiler doesn't support command-line arguments;"
1146 message(2) = "the oct-oscillator-strength command is not available."
1147 call messages_fatal(2)
1148 end if
1149
1150 ! Set the default values
1151 run_mode = analyze_nthorder_signal
1152 omega = - m_one
1153 search_interval = - m_one
1154 order = 1
1155 nfrequencies = 1000
1156 final_time = - m_one
1157 nresonances = 1
1158 observable(1) = -1
1159 observable(2) = 0
1160 ffile = ""
1161 damping = 0.1_real64/27.2114_real64 ! This is the usual damping factor used in the literature.
1162
1163 ! Get the parameters from the command line.
1164 call getopt_oscillator_strength(run_mode, omega, search_interval, &
1165 order, nresonances, nfrequencies, final_time, &
1166 observable(1), observable(2), damping, ffile)
1167 call getopt_end()
1168
1169 ! Initialize stuff
1171 call parser_init()
1172 call io_init(defaults = .true.)
1173
1174 select case (run_mode)
1175 case (generate_nthorder_signal)
1176 call generate_signal(global_namespace, order, observable)
1177 case (analyze_nthorder_signal)
1178 call analyze_signal(order, omega, search_interval, final_time, nresonances, nfrequencies, damping, &
1179 global_namespace)
1180 case (read_resonances_from_file)
1181 call read_resonances_file(order, ffile, global_namespace, search_interval, final_time, nfrequencies)
1182 case (generate_omega_file)
1183 call print_omega_file(global_namespace, omega, search_interval, final_time, nfrequencies)
1184 case default
1185 end select
1186
1187 call io_end()
1188 call parser_end()
1189 call global_end()
1190
1191end program oscillator_strength
1192! ---------------------------------------------------------
1193
1194!! Local Variables:
1195!! mode: f90
1196!! coding: utf-8
1197!! 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:190
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
Definition: global.F90:382
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_four
Definition: global.F90:192
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:186
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:265
subroutine, public init_octopus_globals(comm)
Initialise Octopus-specific global constants and files. This routine performs no initialisation calls...
Definition: global.F90:354
complex(real64), parameter, public m_z0
Definition: global.F90:198
complex(real64), parameter, public m_zi
Definition: global.F90:202
real(real64), parameter, public m_epsilon
Definition: global.F90:204
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
real(real64), parameter, public m_three
Definition: global.F90:191
Definition: io.F90:114
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:161
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
subroutine, public io_skip_header(iunit)
Definition: io.F90:597
subroutine, public io_end()
Definition: io.F90:258
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
subroutine, public kick_read(kick, iunit, namespace)
Definition: kick.F90:819
subroutine, public kick_write(kick, iunit, out)
Definition: kick.F90:891
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1114
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:538
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:415
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
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)
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:2339
character(len=80) function, public str_center(s_in, l_in)
puts space around string, so that it is centered
Definition: string.F90:174
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
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...
program oscillator_strength
int true(void)