Multidip  1.0
Multi-photon matrix elements
multidip_smooth_rawdips.f90
Go to the documentation of this file.
1 ! Copyright 2020
2 !
3 ! For a comprehensive list of the developers that contributed to these codes
4 ! see the UK-AMOR website.
5 !
6 ! This file is part of UKRmol-out (UKRmol+ suite).
7 !
8 ! UKRmol-out is free software: you can redistribute it and/or modify
9 ! it under the terms of the GNU General Public License as published by
10 ! the Free Software Foundation, either version 3 of the License, or
11 ! (at your option) any later version.
12 !
13 ! UKRmol-out is distributed in the hope that it will be useful,
14 ! but WITHOUT ANY WARRANTY; without even the implied warranty of
15 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 ! GNU General Public License for more details.
17 !
18 ! You should have received a copy of the GNU General Public License
19 ! along with UKRmol-out (in source/COPYING). Alternatively, you can also visit
20 ! <https://www.gnu.org/licenses/>.
21 !
22 
39 
40  use iso_fortran_env, only: real64
41 
42  implicit none
43 
44  real(real64), parameter :: sigma = 1e-3 ! Gaussian filter sharpness (large values will conserve narrow features)
45  integer, parameter :: niter = 5 ! number of smoothing iterations (to compensate large out-of-trend spikes)
46 
47  real(real64), allocatable :: dips(:, :)
48  character(len=1024) :: rawname, smoothname
49  integer :: iarg, narg, dot
50 
51  narg = command_argument_count()
52 
53  if (narg == 0) then
54  print '(A)', 'Just give me a list of files to smooth...'
55  stop
56  end if
57 
58  do iarg = 1, narg
59 
60  ! get next filename from the command line
61  call get_command_argument(iarg, rawname)
62  smoothname = rawname
63 
64  ! translate 'XXX.YYY' to 'XXX-smooth.YYY'
65  dot = scan(rawname, '.', back = .true.)
66  if (dot /= 0) then
67  smoothname = smoothname(1:dot-1) // '-smooth' // smoothname(dot:)
68  end if
69 
70  print '(3A)', trim(rawname), ' -> ', trim(smoothname)
71 
72  ! perform the smoothing
73  call read_dipoles(trim(rawname), dips)
74  call smooth_dipoles(dips)
75  call write_dipoles(trim(smoothname), dips)
76 
77  end do
78 
79 contains
80 
81  subroutine read_dipoles (filename, dips)
82 
83  character(len=*), intent(in) :: filename
84  real(real64), allocatable, intent(inout) :: dips(:, :)
85 
86  integer :: i, u, ierr, n
87  real(real64) :: x, y
88 
89  open (newunit = u, file = filename, action = 'read', form = 'formatted')
90 
91  ! first, count the lines
92  n = 0
93  ierr = 0
94  do while (ierr == 0)
95  read (u, *, iostat = ierr) x, y
96  if (ierr == 0) then
97  n = n + 1
98  end if
99  end do
100 
101  ! set up dipole storage
102  if (allocated(dips)) then
103  if (size(dips, 2) /= n) then
104  deallocate (dips)
105  end if
106  end if
107  if (.not. allocated(dips)) then
108  allocate (dips(2, n))
109  end if
110 
111  ! now read from beginning
112  rewind(u)
113  do i = 1, n
114  read (u, *) dips(:, i)
115  end do
116 
117  close (u)
118 
119  end subroutine read_dipoles
120 
121 
122  subroutine write_dipoles (filename, dips)
123 
124  character(len=*), intent(in) :: filename
125  real(real64), allocatable, intent(in) :: dips(:, :)
126 
127  integer :: u
128 
129  open (newunit = u, file = filename, action = 'write', form = 'formatted')
130 
131  write (u, '(2E25.15E3)') dips
132 
133  close (u)
134 
135  end subroutine write_dipoles
136 
137 
138  subroutine smooth_dipoles (dips)
139 
140  real(real64), allocatable, intent(inout) :: dips(:, :)
141  real(real64), allocatable :: weights(:), dips0(:, :), dips1(:, :)
142  real(real64) :: weight, weight_sum, re_d, im_d
143  integer :: i, j, iter, n0, n, d
144 
145  n = size(dips, 2)
146 
147  allocate (weights(n), dips0(2, n), dips1(2, n))
148 
149  ! set up the Gaussian weights
150  do i = 1, n
151  weights(i) = exp(-sigma*(i - 1)**2)
152  dips0(:, i) = 0
153  end do
154 
155  ! count the number of leading zeros; they will not participate in smoothing
156  n0 = 1
157  do i = 1, n
158  if (any(dips(:, i) /= 0)) then
159  n0 = i
160  exit
161  end if
162  end do
163 
164  ! do the smoothing iterations
165  do iter = 1, niter
166  ! calculate entries of the smoothed dataset
167  !$omp parallel do default(none) private(i,j,weight_sum,d,re_d,im_d,weight) shared(dips,dips1,dips0,weights,n,n0)
168  do i = n0, n
169  weight_sum = 0
170  dips1(:, i) = 0
171  ! use all elements from the previous iteration
172  do j = n0, n
173  d = 1 + abs(i - j)
174  re_d = dips(1, j) - dips0(1, j)
175  im_d = dips(2, j) - dips0(2, j)
176  weight = weights(d) / sqrt(1 + re_d*re_d + im_d*im_d)
177  dips1(:, i) = dips1(:, i) + dips(:, j) * weight
178  weight_sum = weight_sum + weight
179  end do
180  dips1(:, i) = dips1(:, i) / weight_sum
181  end do
182  !$omp end parallel do
183  dips0 = dips1
184  end do
185  dips = dips0
186 
187  end subroutine smooth_dipoles
188 
189 end program multidip_smooth_rawdips
subroutine read_dipoles(filename, dips)
program multidip_smooth_rawdips
Smoothing of MULTIDIP raw transition dipoles.
subroutine smooth_dipoles(dips)
subroutine write_dipoles(filename, dips)