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