Octopus
filter.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
21module filter_oct_m
22 use debug_oct_m
23 use global_oct_m
24 use io_oct_m
27 use parser_oct_m
29 use string_oct_m
31
32 implicit none
33
34 private
35
36 public :: &
37 filter_t, &
43
44 type filter_t
45 private
46 integer :: no_filters
47 type(tdf_t), allocatable :: f(:)
48 character(len=1024), allocatable :: expression(:)
49 integer, allocatable :: domain(:)
50 end type filter_t
51
52 integer, parameter, public :: &
53 filter_freq = 1, &
54 filter_time = 2
55
56contains
57
58
59 ! ---------------------------------------------------------
60 subroutine filter_init(steps, namespace, dt, filter)
61 integer, intent(in) :: steps
62 type(namespace_t), intent(in) :: namespace
63 real(real64), intent(in) :: dt
64 type(filter_t), intent(inout) :: filter
65
66 type(block_t) :: blk
67 integer :: i, no_f
68
69 push_sub(filter_init)
70
71 filter%no_filters = 0
72
73 !%Variable OCTFilter
74 !%Type block
75 !%Section Calculation Modes::Optimal Control
76 !%Description
77 !% The block <tt>OCTFilter</tt> describes the type and shape of the filter function
78 !% that are applied to the optimized laser field in each iteration.
79 !% The filter forces the laser field to obtain the given form in frequency space.
80 !% Each line of the block describes a filter; this way you can actually have more
81 !% than one filter function (<i>e.g.</i> a filter in time and two in frequency space).
82 !% The filters are applied in the given order, <i>i.e.</i>, first the filter specified
83 !% by the first line is applied, then second line.
84 !% The syntax of each line is, then:
85 !%
86 !% <tt>%OCTFilter
87 !% <br>&nbsp;&nbsp;domain | function
88 !% <br>%</tt>
89 !%
90 !%
91 !% Possible arguments for domain are:
92 !%
93 !% (i) <tt>frequency_filter</tt>: Specifies a spectral filter.
94 !%
95 !% (ii) <tt>time_filter</tt>: DISABLED IN THIS VERSION.
96 !%
97 !% Example:
98 !%
99 !% <tt>%OCTFilter
100 !% <br>&nbsp;&nbsp;time | "exp(-80*( w + 0.1567 )^2 ) + exp(-80*( w - 0.1567 )^2 )"
101 !% <br>%</tt>
102 !%
103 !% Be careful that also the negative-frequency component is filtered since the resulting
104 !% field has to be real-valued.
105 !%
106 !%Option frequency_filter 1
107 !% The filter is applied in the frequency domain.
108 !%End
109 if (parse_block(namespace, 'OCTFilter', blk) == 0) then
110 no_f = parse_block_n(blk)
111
112 if (no_f <= 0) then
113 pop_sub(filter_init)
114 return
115 end if
116
117 filter%no_filters = no_f
118 safe_allocate(filter%f(1:no_f))
119 safe_allocate(filter%expression(1:no_f))
120 safe_allocate(filter%domain(1:no_f))
121
122 do i = 1, no_f
123 call parse_block_integer(blk, i-1, 0, filter%domain(i))
124 call parse_block_string(blk, i-1, 1, filter%expression(i))
125 call conv_to_c_string(filter%expression(i))
126 call tdf_init_numerical(filter%f(i), steps, dt, -m_one)
127 call tdf_numerical_to_fourier(filter%f(i))
128 end do
129 call build_filter(filter)
130
131 call parse_block_end(blk)
132 end if
133
134 pop_sub(filter_init)
135 end subroutine filter_init
136 ! ---------------------------------------------------------
138
139 ! ---------------------------------------------------------
140 subroutine filter_apply(f, filter)
141 type(tdf_t), intent(inout) :: f
142 type(filter_t), intent(in) :: filter
143
144 integer :: no_f, i, j, nfreqs
146 push_sub(filter_apply)
147
148 no_f = filter%no_filters
149 if (no_f <= 0) then
150 pop_sub(filter_apply)
151 return
152 end if
154 do i = 1, no_f
155 write(message(1),'(a,i2)') "Info: Applying filter "
156 call messages_info(1)
157
158 nfreqs = tdf_nfreqs(f)
159
160 select case (filter%domain(i))
161 case (filter_freq)
163 call tdf_set_numerical(f, 1, tdf(f, 1)*tdf(filter%f(i), 1))
164 do j = 2, nfreqs !+ 1
165 call tdf_set_numerical(f, j, tdf(f, j)*tdf(filter%f(i), j) / sqrt(m_two))
166 end do
167 do j = nfreqs + 1, 2*nfreqs - 1
168 call tdf_set_numerical(f, j, tdf(f, j)*tdf(filter%f(i), j-nfreqs+1) / sqrt(m_two))
169 end do
171
172 case default
173 message(1) = "...I don't know this filter type..."
174 call messages_fatal(1)
175 end select
176
177 end do
178
179 pop_sub(filter_apply)
180 end subroutine filter_apply
181 ! ---------------------------------------------------------
182
183
184 !------------------------------------------------
185 subroutine build_filter(filter)
186 type(filter_t), intent(inout) :: filter
187
188 integer :: i, ip, j, nfreqs
189 real(real64) :: f_re, f_im
190 complex(real64), allocatable :: ff(:)
191 real(real64), allocatable :: grid(:)
192
193 push_sub(build_filter)
194
195 do i = 1, filter%no_filters
196
197 nfreqs = tdf_nfreqs(filter%f(i))
198
199 safe_allocate(ff(1:tdf_nfreqs(filter%f(i))))
200 safe_allocate(grid(1:tdf_nfreqs(filter%f(i))))
201 ff = m_z0
202 grid = m_zero
203
204 select case (filter%domain(i))
205 case (filter_freq)
206 call tdf_fourier_grid(filter%f(i), grid)
207 ff = m_z1
208 do ip = 1, tdf_nfreqs(filter%f(i))
209 call parse_expression(f_re, f_im, "w", real(grid(ip), real64), filter%expression(i))
210 ff(ip) = f_re + m_zi*f_im
211 end do
212
213 case default
214 write(message(1),'(a)') "Unknown choice of domain for filter."
215 write(message(2),'(a)') "Choose: time or freq ."
216 call messages_fatal(2)
217 end select
218
219 call tdf_set_numerical(filter%f(i), 1, real(ff(1), real64))
220 do j = 2, nfreqs !+ 1
221 call tdf_set_numerical(filter%f(i), j, sqrt(m_two)*real(ff(j), real64))
222 end do
223 ! WARNING: all the sine coefficients (imaginary part of ff) should be null.
224 do j = nfreqs + 1, 2*nfreqs - 1
225 call tdf_set_numerical(filter%f(i), j, -sqrt(m_two)*aimag(ff(j-nfreqs+1)))
226 end do
227
228 safe_deallocate_a(ff)
229 safe_deallocate_a(grid)
230
231 end do
232
234 end subroutine build_filter
235 ! ---------------------------------------------------------
236
237
238 ! ---------------------------------------------------------
239 subroutine filter_write(filter, namespace)
240 type(filter_t), intent(in) :: filter
241 type(namespace_t), intent(in) :: namespace
242
243 integer :: kk, iunit, i, max_iter
244 character(len=80) :: filename
245 real(real64), allocatable :: wgrid(:)
246
247 push_sub(filter_write)
248
249 if (filter%no_filters <= 0) then
250 pop_sub(filter_write)
251 return
252 end if
253
254 do kk = 1, filter%no_filters
255 write(filename,'(a,i2.2)') oct_dir//'filter', kk
256 max_iter = tdf_niter(filter%f(kk))
257 iunit = io_open(filename, namespace, action='write')
258 safe_allocate(wgrid(1:max_iter/2+1))
259 call tdf_fourier_grid(filter%f(kk), wgrid)
260 do i = 1, max_iter/2+1
261 write(iunit, '(3es30.16e4)') wgrid(i), tdf(filter%f(kk), i)
262 end do
263 safe_deallocate_a(wgrid)
264 call io_close(iunit)
265 end do
266
267 pop_sub(filter_write)
268 end subroutine filter_write
269 ! ---------------------------------------------------------
270
271
272 ! ----------------------------------------------------------------
273 subroutine filter_end(filter)
274 type(filter_t), intent(inout) :: filter
275 integer :: i
276
277 push_sub(filter_end)
279 if (filter%no_filters <= 0) then
280 pop_sub(filter_end)
281 return
282 end if
283
284 do i = 1, filter%no_filters
285 call tdf_end(filter%f(i))
286 end do
287 safe_deallocate_a(filter%f)
288 filter%no_filters = 0
289 safe_deallocate_a(filter%expression)
290 safe_deallocate_a(filter%domain)
291
292 pop_sub(filter_end)
293 end subroutine filter_end
294 ! ---------------------------------------------------------
295
296
297 ! ---------------------------------------------------------
298 integer function filter_number(filter)
299 type(filter_t), intent(in) :: filter
300
301 push_sub(filter_number)
302 filter_number = filter%no_filters
303
304 pop_sub(filter_number)
305 end function filter_number
306 ! ---------------------------------------------------------
307
308
309end module filter_oct_m
310
311!! Local Variables:
312!! mode: f90
313!! coding: utf-8
314!! End:
double sqrt(double __x) __attribute__((__nothrow__
integer, parameter, public filter_time
Definition: filter.F90:145
integer function, public filter_number(filter)
Definition: filter.F90:392
subroutine, public filter_init(steps, namespace, dt, filter)
Definition: filter.F90:154
subroutine, public filter_apply(f, filter)
Definition: filter.F90:234
subroutine, public filter_write(filter, namespace)
Definition: filter.F90:333
subroutine, public filter_end(filter)
Definition: filter.F90:367
subroutine build_filter(filter)
Definition: filter.F90:279
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_zero
Definition: global.F90:187
complex(real64), parameter, public m_z0
Definition: global.F90:197
complex(real64), parameter, public m_zi
Definition: global.F90:201
complex(real64), parameter, public m_z1
Definition: global.F90:198
real(real64), parameter, public m_one
Definition: global.F90:188
character(len= *), parameter, public oct_dir
Definition: global.F90:258
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:468
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:395
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
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:420
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public conv_to_c_string(str)
converts to c string
Definition: string.F90:252
subroutine, public tdf_numerical_to_fourier(f)
Definition: tdfunction.F90:634
integer pure function, public tdf_nfreqs(f)
Definition: tdfunction.F90:372
integer pure function, public tdf_niter(f)
Definition: tdfunction.F90:364
subroutine, public tdf_end(f)
Definition: tdfunction.F90:974
subroutine, public tdf_fourier_grid(f, wgrid)
Definition: tdfunction.F90:607
subroutine, public tdf_init_numerical(f, niter, dt, omegamax, initval, rep)
Definition: tdfunction.F90:552
subroutine, public tdf_fourier_to_numerical(f)
Definition: tdfunction.F90:669