Octopus
maxwell_function.F90
Go to the documentation of this file.
1!! Copyright (C) 2019 R. Jestaedt, H. Appel, F. Bonafe, M. Oliveira, N. Tancogne-Dejean
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
26 use iso_c_binding
27 use debug_oct_m
28 use fft_oct_m
29 use global_oct_m
30 use io_oct_m
32 use math_oct_m
33 use mpi_oct_m
34 use parser_oct_m
37 use string_oct_m
38 use unit_oct_m
41
42 implicit none
43
44 public :: &
45 mxf_t, &
51 mxf, &
53
54 integer, public, parameter :: &
55 MXF_EMPTY = 10001, &
56 mxf_const_wave = 10002, &
57 mxf_const_phase = 10004, &
58 mxf_gaussian_wave = 10005, &
59 mxf_cosinoidal_wave = 10006, &
60 mxf_logistic_wave = 10007, &
61 mxf_trapezoidal_wave = 10008, &
62 mxf_from_file = 10009, &
63 mxf_numerical = 10010, &
64 mxf_from_expr = 10011, &
65 mxf_fourier_series = 10012, &
66 mxf_zero_fourier = 10013
67
68 type mxf_t
69 integer :: mode = mxf_empty
70 real(real64) :: k_vector(3) = m_zero
71 real(real64) :: r0(3) = m_zero
72 real(real64) :: width = m_zero
73 real(real64) :: a0 = m_zero
74 real(real64) :: dx = m_zero
75 real(real64) :: init_x = m_zero
76 real(real64) :: final_x = m_zero
77 real(real64) :: growth = m_zero
78 integer :: niter = 0
79 integer :: nfreqs = 0
80
81 type(spline_t) :: amplitude
82 character(len=200) :: expression
83 type(fft_t) :: fft_handler
84 end type mxf_t
85
86 interface mxf
87 module procedure mxf_eval
88 end interface mxf
89
90contains
91
92 !------------------------------------------------------------
94 subroutine mxf_read(f, namespace, function_name, ierr)
95 type(mxf_t), intent(inout) :: f
96 type(namespace_t), intent(in) :: namespace
97 character(len=*), intent(in) :: function_name
98 integer, intent(out) :: ierr
99
100 type(block_t) :: blk
101 integer :: nrows, ncols, i, function_type, idim
102 character(len=100) :: row_name, function_expression
103 real(real64) :: a0, r0(3), growth, width, k_vector(3)
104
105 push_sub(mxf_read)
106
107 !%Variable MaxwellFunctions
108 !%Type block
109 !%Section Maxwell
110 !%Description
111 !% This block specifies the shape of a "spatial-dependent function", such as the
112 !% envelope needed when using the <tt>MaxwellFunctions</tt> block. Each line in the block
113 !% specifies one function. The first element of each line will be a string
114 !% that defines the name of the function. The second element specifies which type
115 !% of function we are using; in the following we provide an example for each of the
116 !% possible types. In the following, we will use</br>
117 !% <math>
118 !% \mathbf{r} = \begin{pmatrix} {\rm x} \\\\ {\rm y} \\\\ {\rm z} \end{pmatrix},
119 !% \mathbf{r}_0 = \begin{pmatrix} {\rm x0} \\\\ {\rm y0} \\\\ {\rm z0} \end{pmatrix},
120 !% \mathbf{k} = \begin{pmatrix} {\rm kx} \\\\ {\rm ky} \\\\ {\rm kz} \end{pmatrix},
121 !% w = {\rm width},
122 !% g = {\rm growth}
123 !% </math>
124 !%Option mxf_const_wave 10002
125 !%
126 !% <tt>%MaxwellFunctions
127 !% <br>&nbsp;&nbsp; "function-name" | mxf_const_wave | kx | ky | kz | x0 | y0 | z0
128 !% <br>%</tt>
129 !%
130 !% The function is constant plane wave <math> f(\mathbf{r}) = \cos( \mathbf{k} (\mathbf{r} - \mathbf{r}_0) ) </math>
131 !%
132 !%Option mxf_const_phase 10004
133 !%
134 !% <tt>%MaxwellFunctions
135 !% <br>&nbsp;&nbsp; "function-name" | mxf_const_phase | kx | ky | kz | x0 | y0 | z0
136 !% <br>%</tt>
137 !%
138 !% The function is a constant phase of <math> f(\mathbf{r}) = (\mathbf{k} \mathbf{r}_0) </math>
139 !%
140 !%Option mxf_gaussian_wave 10005
141 !%
142 !% <tt>%MaxwellFunctions
143 !% <br>&nbsp;&nbsp; "function-name" | mxf_gaussian_wave | kx | ky | kz | x0 | y0 | z0 | width
144 !% <br>%</tt>
145 !%
146 !% The function is a Gaussian,
147 !% <math> f(\mathbf{r}) = \exp( -( \mathbf{k} (\mathbf{r}-\mathbf{r}_0) )^2 / (2 {w}^2) ) </math>
148 !%
149 !%Option mxf_cosinoidal_wave 10006
150 !%
151 !% <tt>%MaxwellFunctions
152 !% <br>&nbsp;&nbsp; "function-name" | mxf_cosinoidal_wave | kx | ky | kz | x0 | y0 | z0 | width
153 !% <br>%</tt>
154 !%
155 !% <math>
156 !% f(\mathbf{r}) =
157 !% \begin{cases}
158 !% \cos( \frac{\pi}{2} \frac{ \mathbf{k}(\mathbf{r}-\mathbf{r}_0) - 2 w}{w} + \pi )
159 !% & \text{if } | \mathbf{k}(\mathbf{r}-\mathbf{r}_0) | \le w \\\\
160 !% 0 & \text{otherwise}
161 !% \end{cases}
162 !% </math>
163 !%
164 !%Option mxf_logistic_wave 10007
165 !%
166 !% <tt>%MaxwellFunctions
167 !% <br>&nbsp;&nbsp; "function-name" | mxf_logistic_wave | kx | ky | kz | x0 | y0 | z0 | growth | width
168 !% <br>%</tt>
169 !%
170 !% The function is a logistic function, <br>
171 !% <math>
172 !% f(\mathbf{r}) =
173 !% \left( \frac{1}{1 + e^{g \left( \frac{\mathbf{k} \cdot (\mathbf{r-r_0})}{|\mathbf{k}|} + \frac{w}{2} \right)}} \right)
174 !% \cdot
175 !% \left( \frac{1}{1 + e^{-g \left( \frac{\mathbf{k} \cdot (\mathbf{r-r_0})}{|\mathbf{k}|} - \frac{w}{2} \right)}} \right)
176 !% \cdot
177 !% {\left(1 + e^{\frac{g \cdot w}{2}}\right)}^2
178 !% </math></br>
179 !%
180 !%Option mxf_trapezoidal_wave 10008
181 !%
182 !% <tt>%MaxwellFunctions
183 !% <br>&nbsp;&nbsp; "function-name" | mxf_trapezoidal_wave | kx | ky | kz | x0 | y0 | z0 | growth | width
184 !% <br>%
185 !% </tt>
186 !%
187 !% The function is a logistic function,
188 !% <br>
189 !% <math>
190 !% f(\mathbf{r}) =
191 !% \begin{cases}
192 !% 0 & \text{if }
193 !% \mathbf{k} \cdot (\mathbf{r} - \mathbf{r}_0) / |\mathbf{k}| \leq -\frac{w}{2} \\\\
194 !% 1 & \text{if }
195 !% -\frac{w}{2} + \frac{1}{g} < \mathbf{k} \cdot (\mathbf{r} - \mathbf{r}_0) / |\mathbf{k}| < \frac{w}{2} - \frac{1}{g} \\\\
196 !% 1 - g\left[\mathbf{k} \cdot (\mathbf{r} - \mathbf{r}_0) / |\mathbf{k}| - \left(\frac{w}{2} - \frac{1}{g}\right)\right] &
197 !% \text{if } \frac{w}{2} - \frac{1}{g} \leq \mathbf{k} \cdot (\mathbf{r} - \mathbf{r}_0) / |\mathbf{k}| < \frac{w}{2} \\\\
198 !% 0 &
199 !% \text{if } \mathbf{k} \cdot (\mathbf{r} - \mathbf{r}_0) / |\mathbf{k}| \geq \frac{w}{2}
200 !% \end{cases}
201 !% </math></br>
202 !%Option mxf_from_expr 10011
203 !%
204 !% <tt>%MaxwellFunctions
205 !% <br>&nbsp;&nbsp; "function-name" | mxf_from_expr | kx | ky | kz | "expression"
206 !% <br>%</tt>
207 !%
208 !% The temporal shape of the field is given as an expression (e.g., <tt>cos(2.0*x-3*y+4*z)</tt>. The
209 !% letter <i>x</i>, <i>y</i>, <i>z</i> means spatial coordinates, obviously.
210 !% The expression is used to construct the function <i>f</i>
211 !% that defines the field.
212 !%End
213 ierr = -3
214 if (parse_block(namespace, 'MaxwellFunctions', blk) /= 0) then
215 ierr = -1
216 pop_sub(mxf_read)
217 return
218 end if
219
220 nrows = parse_block_n(blk)
221 row_loop: do i = 1, nrows
222 call parse_block_string(blk, i-1, 0, row_name)
223 if (trim(row_name) == trim(function_name)) then
224
225 ncols = parse_block_cols(blk, i-1)
226 call parse_block_integer(blk, i-1, 1, function_type)
227
228 a0 = m_one
229 r0 = m_zero
230 width = m_zero
231 k_vector = m_zero
232 select case (function_type)
233 case (mxf_const_wave)
234 do idim = 1, 3
235 call parse_block_float(blk, i-1, 1+idim, k_vector(idim), unit_one/units_inp%length)
236 end do
237 do idim = 1, 3
238 call parse_block_float(blk, i-1, 4+idim, r0(idim), units_inp%length)
239 end do
240 call mxf_init_const_wave(f, a0, k_vector, r0)
241 case (mxf_const_phase)
242 do idim = 1, 3
243 call parse_block_float(blk, i-1, 1+idim, k_vector(idim), unit_one/units_inp%length)
244 end do
245 do idim = 1, 3
246 call parse_block_float(blk, i-1, 4+idim, r0(idim), units_inp%length)
247 end do
248 call mxf_init_const_phase(f, a0, k_vector, r0)
249 case (mxf_gaussian_wave)
250 do idim = 1, 3
251 call parse_block_float(blk, i-1, 1+idim, k_vector(idim), unit_one/units_inp%length)
252 end do
253 do idim = 1, 3
254 call parse_block_float(blk, i-1, 4+idim, r0(idim), units_inp%length)
255 end do
256 call parse_block_float(blk, i-1, 8, width, units_inp%length)
257 call mxf_init_gaussian_wave(f, a0, k_vector, r0, width)
259 do idim = 1, 3
260 call parse_block_float(blk, i-1, 1+idim, k_vector(idim), unit_one/units_inp%length)
261 end do
262 do idim = 1, 3
263 call parse_block_float(blk, i-1, 4+idim, r0(idim), units_inp%length)
264 end do
265 call parse_block_float(blk, i-1, 8, width, units_inp%length)
266 call mxf_init_cosinoidal_wave(f, a0, k_vector, r0, width)
267 case (mxf_logistic_wave)
268 do idim = 1, 3
269 call parse_block_float(blk, i-1, 1+idim, k_vector(idim), unit_one/units_inp%length)
270 end do
271 do idim = 1, 3
272 call parse_block_float(blk, i-1, 4+idim, r0(idim), units_inp%length)
273 end do
274 call parse_block_float(blk, i-1, 8, growth, units_inp%length)
275 call parse_block_float(blk, i-1, 9, width, units_inp%length)
276 call mxf_init_logistic_wave(f, a0, k_vector, r0, growth, width)
278 do idim = 1, 3
279 call parse_block_float(blk, i-1, 1+idim, k_vector(idim), unit_one/units_inp%length)
280 end do
281 do idim = 1, 3
282 call parse_block_float(blk, i-1, 4+idim, r0(idim), units_inp%length)
283 end do
284 call parse_block_float(blk, i-1, 8, growth, units_inp%length)
285 call parse_block_float(blk, i-1, 9, width, units_inp%length)
286 call mxf_init_trapezoidal_wave(f, a0, k_vector, r0, growth, width)
287 case (mxf_from_expr)
288 do idim = 1, 3
289 call parse_block_float(blk, i-1, 1+idim, k_vector(idim), unit_one/units_inp%length)
290 end do
291 call parse_block_string(blk, i-1, 5, function_expression, convert_to_c=.true.)
292 call mxf_init_fromexpr(f, k_vector, trim(function_expression))
293 case default
294 ierr = -2
295 call parse_block_end(blk)
296 pop_sub(mxf_read)
297 return
298 end select
299
300 ierr = 0
301 exit row_loop
302 end if
303 end do row_loop
304
305 call parse_block_end(blk)
306
307 pop_sub(mxf_read)
308 end subroutine mxf_read
309 !------------------------------------------------------------
310
311
312 !------------------------------------------------------------
313 subroutine mxf_init_const_wave(f, a0, k_vector, r0)
314 type(mxf_t), intent(inout) :: f
315 real(real64), intent(in) :: a0, k_vector(3), r0(3)
316
317 push_sub(mxf_init_const_wave)
318
319 f%mode = mxf_const_wave
320 f%a0 = a0
321 f%k_vector = k_vector
322 f%r0 = r0
323
324 pop_sub(mxf_init_const_wave)
325 end subroutine mxf_init_const_wave
326 !------------------------------------------------------------
327
328
329 !------------------------------------------------------------
330 subroutine mxf_init_const_phase(f, a0, k_vector, r0)
331 type(mxf_t), intent(inout) :: f
332 real(real64), intent(in) :: a0, k_vector(3), r0(3)
333
334 push_sub(mxf_init_const_phase)
335
336 f%mode = mxf_const_phase
337 f%a0 = a0
338 f%k_vector = k_vector
339 f%r0 = r0
340
341 pop_sub(mxf_init_const_phase)
342 end subroutine mxf_init_const_phase
343 !------------------------------------------------------------
344
345
346 !------------------------------------------------------------
347 subroutine mxf_init_gaussian_wave(f, a0, k_vector, r0, width)
348 type(mxf_t), intent(inout) :: f
349 real(real64), intent(in) :: a0, k_vector(3), r0(3), width
350
351 push_sub(mxf_init_gaussian_wave)
352
353 f%mode = mxf_gaussian_wave
354 f%a0 = a0
355 f%k_vector = k_vector
356 f%r0 = r0
357 f%width = width
358
360 end subroutine mxf_init_gaussian_wave
361 !------------------------------------------------------------
362
363
364 !------------------------------------------------------------
365 subroutine mxf_init_cosinoidal_wave(f, a0, k_vector, r0, width)
366 type(mxf_t), intent(inout) :: f
367 real(real64), intent(in) :: a0, k_vector(3), r0(3), width
368
370
372 f%a0 = a0
373 f%k_vector = k_vector
374 f%r0 = r0
375 f%width = width
376
378 end subroutine mxf_init_cosinoidal_wave
379 !------------------------------------------------------------
380
381
382 !------------------------------------------------------------
383 subroutine mxf_init_logistic_wave(f, a0, k_vector, r0, growth, width)
384 type(mxf_t), intent(inout) :: f
385 real(real64), intent(in) :: a0, k_vector(3), r0(3), growth, width
386
387 push_sub(mxf_init_logistic_wave)
388
389 f%mode = mxf_logistic_wave
390 f%a0 = a0
391 f%k_vector = k_vector
392 f%r0 = r0
393 f%growth = growth
394 f%width = width
395
397 end subroutine
398 !------------------------------------------------------------
399
400
401 !------------------------------------------------------------
402 subroutine mxf_init_trapezoidal_wave(f, a0, k_vector, r0, growth, width)
403 type(mxf_t), intent(inout) :: f
404 real(real64), intent(in) :: a0, k_vector(3), r0(3), growth, width
405
407
409 f%a0 = a0
410 f%k_vector = k_vector
411 f%r0 = r0
412 f%growth = growth
413 f%width = width
414
416 end subroutine
417 !------------------------------------------------------------
418
419
420 !------------------------------------------------------------
421 subroutine mxf_init_fromexpr(f, k_vector, expression)
422 type(mxf_t), intent(inout) :: f
423 real(real64), intent(in) :: k_vector(3)
424 character(len=*), intent(in) :: expression
426 push_sub(mxf_init_fromexpr)
427
428 f%mode = mxf_from_expr
429 f%k_vector = k_vector
430 f%expression = trim(expression)
431
432 pop_sub(mxf_init_fromexpr)
433 end subroutine mxf_init_fromexpr
434 !------------------------------------------------------------
435
440 !------------------------------------------------------------
441 complex(real64) function mxf_envelope_eval(f, x) result(env)
442 type(mxf_t), intent(in) :: f
443 real(real64), intent(in) :: x(:)
444 real(real64) :: xx, limit_1, limit_2, limit_3, limit_4, r_re, r_im, rr
445 complex(real64) :: r
446 integer :: xdim
447
448 xdim = size(x)
449
450 select case (f%mode)
451 case (mxf_const_wave)
452
453 env = f%a0
454
455 case (mxf_const_phase)
456
457 env = f%a0 * sum(f%k_vector(1:xdim)*(x(:) - f%r0(1:xdim)))
458
459 case (mxf_gaussian_wave)
461 r = exp(-((sum(f%k_vector(1:xdim)*(x(:) - f%r0(1:xdim))) &
462 / norm2(f%k_vector(1:xdim)) )**2 / (m_two*f%width**2)))
463 env = f%a0 * r
464
466
467 r = m_zero
468 if (abs( sum(f%k_vector(1:xdim)*(x(:) - f%r0(1:xdim))/norm2(f%k_vector(1:xdim)))) <= f%width) then
469 r = - cos((m_pi/m_two) * ((sum(f%k_vector(1:xdim)*(x(:) - f%r0(1:xdim))) &
470 / norm2(f%k_vector(1:xdim)) - m_two*f%width) / f%width))
471 end if
472 env = f%a0 * r
473
474 case (mxf_logistic_wave)
475
476 r = m_one/(m_one + exp(f%growth*(sum(f%k_vector(1:xdim)*(x(:) - f%r0(1:xdim))) &
477 / norm2(f%k_vector(1:xdim)) - f%width/m_two))) &
478 + m_one/(m_one + exp(-f%growth*(sum(f%k_vector(1:xdim)*(x(:) - f%r0(1:xdim))) &
479 / norm2(f%k_vector(1:xdim)) + f%width/m_two))) - m_one
480 env = f%a0 * r
481
483
484 xx = sum(f%k_vector(1:xdim)*(x(:) - f%r0(1:xdim))) / norm2(f%k_vector(1:xdim))
485 limit_1 = - f%width/m_two - m_one/f%growth
486 limit_2 = - f%width/m_two
487 limit_3 = f%width/m_two
488 limit_4 = f%width/m_two + m_one/f%growth
489 if ((xx > limit_1) .and. (xx <= limit_2)) then
490 r = m_one + f%growth * (sum(f%k_vector(1:xdim)*(x(:) - f%r0(1:xdim)))/norm2(f%k_vector(1:xdim)) + f%width/m_two)
491 else if ((xx > limit_2) .and. (xx <= limit_3)) then
492 r = m_one
493 else if ((xx > limit_3) .and. (xx <= limit_4)) then
494 r = m_one - f%growth * (sum(f%k_vector(1:xdim)*(x(:) - f%r0(1:xdim)))/norm2(f%k_vector(1:xdim)) - f%width/m_two)
495 else
496 r = m_zero
497 end if
498 env = f%a0 * r
499
500 case (mxf_from_expr)
501 ! here x(:) is actually x(:) - v(:)*t
502 rr = norm2(x)
503 call parse_expression(r_re, r_im, 3, x(:), rr, m_zero, trim(f%expression))
504 env = cmplx(r_re, r_im, real64)
505
506 case default
507
508 env = m_zero
509
510 end select
511
512 end function mxf_envelope_eval
513 !------------------------------------------------------------
514
521 !------------------------------------------------------------
522 complex(real64) function mxf_eval(f, x, phi) result(y)
523 type(mxf_t), intent(in) :: f
524 real(real64), intent(in) :: x(:)
525 real(real64), optional, intent(in) :: phi
526
527 real(real64) :: phi_
528 integer :: xdim
529
530 ! no push_sub because it is called too frequently
531
532 phi_ = optional_default(phi, m_zero)
533
534 xdim = size(x)
535
537 if (f%mode /= mxf_const_phase) then
538 y = y * exp(m_zi * (sum(f%k_vector(1:xdim)*(x(:) - f%r0(1:xdim))) + phi_))
539 endif
540
541 end function mxf_eval
542 !------------------------------------------------------------
543
544end module maxwell_function_oct_m
545
546!! Local Variables:
547!! mode: f90
548!! coding: utf-8
549!! End:
double exp(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
Definition: fft.F90:120
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
complex(real64), parameter, public m_zi
Definition: global.F90:205
real(real64), parameter, public m_one
Definition: global.F90:192
Definition: io.F90:116
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
integer, parameter, public mxf_fourier_series
integer, parameter, public mxf_logistic_wave
complex(real64) function mxf_envelope_eval(f, x)
Evaluation of envelope itself.
complex(real64) function mxf_eval(f, x, phi)
Evaluation of spatial envelope Functions.
integer, parameter, public mxf_trapezoidal_wave
subroutine, public mxf_read(f, namespace, function_name, ierr)
This function initializes "f" from the MXFunctions block.
integer, parameter, public mxf_from_file
integer, parameter, public mxf_cosinoidal_wave
integer, parameter, public mxf_zero_fourier
subroutine, public mxf_init_const_wave(f, a0, k_vector, r0)
integer, parameter, public mxf_numerical
subroutine mxf_init_trapezoidal_wave(f, a0, k_vector, r0, growth, width)
subroutine, public mxf_init_const_phase(f, a0, k_vector, r0)
integer, parameter, public mxf_gaussian_wave
integer, parameter, public mxf_const_wave
subroutine, public mxf_init_gaussian_wave(f, a0, k_vector, r0, width)
subroutine, public mxf_init_fromexpr(f, k_vector, expression)
integer, parameter, public mxf_const_phase
subroutine, public mxf_init_cosinoidal_wave(f, a0, k_vector, r0, width)
integer, parameter, public mxf_from_expr
subroutine mxf_init_logistic_wave(f, a0, k_vector, r0, growth, width)
subroutine, public parse_block_string(blk, l, c, res, convert_to_c)
Definition: parser.F90:810
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:615
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
This module defines the unit system, used for input and output.
type(unit_system_t), public units_inp
the units systems for reading and writing
type(unit_t), public unit_one
some special units required for particular quantities
static double f(double w, void *p)
int true(void)