22 #include "sll_assert.h"
23 #include "sll_memory.h"
24 #include "sll_working_precision.h"
65 sll_real64 :: current_time = 0.0_f64
68 sll_int32 :: split_case
70 sll_int32 :: nb_split_step
72 sll_real64,
dimension(:),
pointer :: split_step
74 logical :: split_begin_t
76 sll_int32 :: dim_split_v
88 subroutine operator(this, dt)
90 sll_real64,
intent(in) :: dt
92 print *,
'This is a dummy implementation for providint the interface &
93 & in sll_operator_splitting_base.F90. &
94 & The class sll_t_operator_splitting needs to be extended &
95 & providing the data to be evolved by the operator for an actual &
99 sll_assert(storage_size(this) > 0)
100 sll_assert(dt > 0.0_f64)
134 sll_int32,
intent(in) :: split_case
135 sll_real64,
dimension(:),
intent(in),
optional :: split_step
136 sll_int32,
intent(in),
optional :: nb_split_step
137 logical,
intent(in),
optional :: split_begin_t
138 sll_real64,
intent(in),
optional :: dt
141 sll_allocate(split, ierr)
184 sll_int32,
intent(in) :: split_case
185 sll_real64,
dimension(:),
intent(in),
optional :: split_step
186 sll_int32,
intent(in),
optional :: nb_split_step
187 logical,
intent(in),
optional :: split_begin_t
188 sll_real64,
intent(in),
optional :: dt
191 split%dim_split_V = 1
192 split%split_case = split_case
195 if (.not. ((
present(split_step)) &
196 .and. (
present(nb_split_step)) &
197 .and. (
present(split_begin_t))))
then
198 print *,
'#provide split_step, nb_split_step,split_begin_T'
199 print *,
'#in initialize_operator_splitting_coeff'
204 (
present(split_step)) &
205 .or. (
present(nb_split_step)) &
206 .or. (
present(split_begin_t)) &
208 print *,
'# do not provide split_step, nb_split_step, split_begin_T'
209 print *,
'#in initialize_operator_splitting_coeff'
214 select case (split_case)
216 if (
size(split_step) < nb_split_step)
then
217 print *,
'#bad size for split_step',
size(split_step), nb_split_step
218 print *,
'#in initialize_operator_splitting_coeff'
221 sll_allocate(split%split_step(nb_split_step), ierr)
222 split%split_step(1:nb_split_step) = split_step(1:nb_split_step)
223 split%nb_split_step = nb_split_step
224 split%split_begin_T = split_begin_t
227 split%nb_split_step = 2
228 sll_allocate(split%split_step(split%nb_split_step), ierr)
229 split%split_begin_T = .true.
230 split%split_step(1) = 1._f64
231 split%split_step(2) = 1._f64
233 split%nb_split_step = 2
234 sll_allocate(split%split_step(split%nb_split_step), ierr)
235 split%split_begin_T = .false.
236 split%split_step(1) = 1._f64
237 split%split_step(2) = 1._f64
239 split%nb_split_step = 3
240 sll_allocate(split%split_step(split%nb_split_step), ierr)
241 split%split_begin_T = .true.
242 split%split_step(1) = 0.5_f64
243 split%split_step(2) = 1._f64
244 split%split_step(3) = split%split_step(1)
246 split%nb_split_step = 3
247 sll_allocate(split%split_step(split%nb_split_step), ierr)
248 split%split_begin_T = .false.
249 split%split_step(1) = 0.5_f64
250 split%split_step(2) = 1._f64
251 split%split_step(3) = split%split_step(1)
253 split%nb_split_step = 7
254 sll_allocate(split%split_step(split%nb_split_step), ierr)
255 split%split_begin_T = .true.
256 split%split_step(1) = 0.675603595979829_f64
257 split%split_step(2) = 1.351207191959658_f64
258 split%split_step(3) = -0.17560359597982855_f64
259 split%split_step(4) = -1.702414383919315_f64
260 split%split_step(5) = split%split_step(3)
261 split%split_step(6) = split%split_step(2)
262 split%split_step(7) = split%split_step(1)
264 split%nb_split_step = 7
265 sll_allocate(split%split_step(split%nb_split_step), ierr)
266 split%split_begin_T = .false.
267 split%split_step(1) = 0.675603595979829_f64
268 split%split_step(2) = 1.351207191959658_f64
269 split%split_step(3) = -0.17560359597982855_f64
270 split%split_step(4) = -1.702414383919315_f64
271 split%split_step(5) = split%split_step(3)
272 split%split_step(6) = split%split_step(2)
273 split%split_step(7) = split%split_step(1)
275 split%nb_split_step = 23
276 sll_allocate(split%split_step(split%nb_split_step), ierr)
277 split%split_begin_T = .true.
278 split%split_step(1) = 0.0414649985182624_f64
279 split%split_step(2) = 0.123229775946271_f64
280 split%split_step(3) = 0.198128671918067_f64
281 split%split_step(4) = 0.290553797799558_f64
282 split%split_step(5) = -0.0400061921041533_f64
283 split%split_step(6) = -0.127049212625417_f64
284 split%split_step(7) = 0.0752539843015807_f64
285 split%split_step(8) = -0.246331761062075_f64
286 split%split_step(9) = -0.0115113874206879_f64
287 split%split_step(10) = 0.357208872795928_f64
288 split%split_step(11) = 0.23666992478693111_f64
289 split%split_step(12) = 0.20477705429147008_f64
290 split%split_step(13) = split%split_step(11)
291 split%split_step(14) = split%split_step(10)
292 split%split_step(15) = split%split_step(9)
293 split%split_step(16) = split%split_step(8)
294 split%split_step(17) = split%split_step(7)
295 split%split_step(18) = split%split_step(6)
296 split%split_step(19) = split%split_step(5)
297 split%split_step(20) = split%split_step(4)
298 split%split_step(21) = split%split_step(3)
299 split%split_step(22) = split%split_step(2)
300 split%split_step(23) = split%split_step(1)
302 split%nb_split_step = 23
303 sll_allocate(split%split_step(split%nb_split_step), ierr)
304 split%split_begin_T = .false.
305 split%split_step(1) = 0.0414649985182624_f64
306 split%split_step(2) = 0.123229775946271_f64
307 split%split_step(3) = 0.198128671918067_f64
308 split%split_step(4) = 0.290553797799558_f64
309 split%split_step(5) = -0.0400061921041533_f64
310 split%split_step(6) = -0.127049212625417_f64
311 split%split_step(7) = 0.0752539843015807_f64
312 split%split_step(8) = -0.246331761062075_f64
313 split%split_step(9) = -0.0115113874206879_f64
314 split%split_step(10) = 0.357208872795928_f64
315 split%split_step(11) = 0.23666992478693111_f64
316 split%split_step(12) = 0.20477705429147008_f64
317 split%split_step(13) = split%split_step(11)
318 split%split_step(14) = split%split_step(10)
319 split%split_step(15) = split%split_step(9)
320 split%split_step(16) = split%split_step(8)
321 split%split_step(17) = split%split_step(7)
322 split%split_step(18) = split%split_step(6)
323 split%split_step(19) = split%split_step(5)
324 split%split_step(20) = split%split_step(4)
325 split%split_step(21) = split%split_step(3)
326 split%split_step(22) = split%split_step(2)
327 split%split_step(23) = split%split_step(1)
329 if (
present(dt))
then
330 split%nb_split_step = 9
331 sll_allocate(split%split_step(split%nb_split_step), ierr)
332 split%split_begin_T = .true.
333 split%split_step(1) = 0.1095115577513980413559540_f64
334 split%split_step(2) = 0.268722208204814693684441_f64 &
335 - 2._f64*dt**2*0.000805681667096178271312_f64 &
336 + 4._f64*dt**4*0.000017695766224036466792_f64
337 split%split_step(3) = 0.4451715080955340951457244_f64
338 split%split_step(4) = 0.2312777917951853063155588_f64 &
339 - 2._f64*dt**2*0.003955911930042478239763_f64 &
340 + 4._f64*dt**4*0.000052384078562246674986_f64
341 split%split_step(5) = -0.1093661316938642730033570_f64
342 split%split_step(6) = split%split_step(4)
343 split%split_step(7) = split%split_step(3)
344 split%split_step(8) = split%split_step(2)
345 split%split_step(9) = split%split_step(1)
351 if (
present(dt))
then
352 split%nb_split_step = 9
353 sll_allocate(split%split_step(split%nb_split_step), ierr)
354 split%split_begin_T = .false.
355 split%split_step(1) = 0.359950808794143627485664_f64 &
356 - 2._f64*dt**2*(-0.01359558332625151635_f64) &
357 + 4._f64*dt**4*(-8.562814848565929e-6_f64)
358 split%split_step(2) = 1.079852426382430882456991_f64
359 split%split_step(3) = -0.1437147273026540434771131_f64 &
360 - 2._f64*dt**2*(-0.00385637757897273261_f64) &
361 + 4._f64*dt**4*(0.0004883788785819335822_f64)
362 split%split_step(4) = -0.579852426382430882456991_f64
363 split%split_step(5) = 0.567527837017020831982899_f64 &
364 - 2._f64*dt**2*(-0.03227361602037480885_f64) &
365 + 4._f64*dt**4*0.002005141087312622342_f64
366 split%split_step(6) = split%split_step(4)
367 split%split_step(7) = split%split_step(3)
368 split%split_step(8) = split%split_step(2)
369 split%split_step(9) = split%split_step(1)
375 if (
present(dt))
then
376 split%nb_split_step = 9
377 sll_allocate(split%split_step(split%nb_split_step), ierr)
378 split%split_begin_T = .true.
379 split%split_step(1) = 0.1095115577513980413559540_f64
380 split%split_step(2) = 0.268722208204814693684441_f64 &
381 - 2._f64*dt**2*0.000805681667096178271312_f64 &
382 + 4._f64*dt**4*(8.643923349886021963e-6_f64) &
383 - 8._f64*dt**6*(1.4231479258353431522e-6_f64)
384 split%split_step(3) = 0.4451715080955340951457244_f64
385 split%split_step(4) = 0.2312777917951853063155588_f64 &
386 - 2._f64*dt**2*0.003955911930042478239763_f64 &
387 + 4._f64*dt**4*(0.000061435921436397119815_f64)
388 split%split_step(5) = -0.1093661316938642730033570_f64
389 split%split_step(6) = split%split_step(4)
390 split%split_step(7) = split%split_step(3)
391 split%split_step(8) = split%split_step(2)
392 split%split_step(9) = split%split_step(1)
399 if (
present(dt))
then
400 split%nb_split_step = 11
401 sll_allocate(split%split_step(split%nb_split_step), ierr)
402 split%split_begin_T = .false.
403 split%split_step(1) = 0.0490864609761162454914412_f64 &
404 - 2._f64*dt**2*(0.0000697287150553050840999_f64)
405 split%split_step(2) = 0.1687359505634374224481957_f64
406 split%split_step(3) = 0.2641776098889767002001462_f64 &
407 - 2._f64*dt**2*(0.000625704827430047189169_f64) &
408 + 4._f64*dt**4*(-2.91660045768984781644e-6_f64)
409 split%split_step(4) = 0.377851589220928303880766_f64
410 split%split_step(5) = 0.1867359291349070543084126_f64 &
411 - 2._f64*dt**2*(0.00221308512404532556163_f64) &
412 + 4._f64*dt**4*(0.0000304848026170003878868_f64) &
413 - 8._f64*dt**6*(4.98554938787506812159e-7_f64)
414 split%split_step(6) = -0.0931750795687314526579244_f64
415 split%split_step(7) = split%split_step(5)
416 split%split_step(8) = split%split_step(4)
417 split%split_step(9) = split%split_step(3)
418 split%split_step(10) = split%split_step(2)
419 split%split_step(11) = split%split_step(1)
425 if (
present(dt))
then
426 split%nb_split_step = 11
427 sll_allocate(split%split_step(split%nb_split_step), ierr)
428 split%split_begin_T = .false.
429 split%split_step(1) = 0.0490864609761162454914412_f64 &
430 + 2._f64*dt**2*(0.00166171386175851683711044_f64)
431 split%split_step(2) = 0.1687359505634374224481957_f64
432 split%split_step(3) = 0.2641776098889767002001462_f64 &
433 - 2._f64*dt**2*(0.00461492847770001641230401_f64)
434 split%split_step(4) = 0.377851589220928303880766_f64
435 split%split_step(5) = 0.1867359291349070543084126_f64 &
436 + 2._f64*dt**2*(0.0000446959494108217402966857_f64)
437 split%split_step(6) = -0.0931750795687314526579244_f64
438 split%split_step(7) = split%split_step(5)
439 split%split_step(8) = split%split_step(4)
440 split%split_step(9) = split%split_step(3)
441 split%split_step(10) = split%split_step(2)
442 split%split_step(11) = split%split_step(1)
448 if (
present(dt))
then
449 split%nb_split_step = 11
450 split%dim_split_V = 2
451 sll_allocate(split%split_step(17), ierr)
452 split%split_begin_T = .false.
453 split%split_step(1) = 0.0490864609761162454914412_f64 &
454 + 2._f64*dt**2*(0.00166171386175851683711044_f64)
455 split%split_step(2) = dt**2*(0.00166171386175851683711044_f64)
456 split%split_step(3) = 0.1687359505634374224481957_f64
457 split%split_step(4) = 0.2641776098889767002001462_f64 &
458 - 2._f64*dt**2*(0.00461492847770001641230401_f64)
459 split%split_step(5) = -dt**2*(0.00461492847770001641230401_f64)
460 split%split_step(6) = 0.377851589220928303880766_f64
461 split%split_step(7) = 0.1867359291349070543084126_f64 &
462 + 2._f64*dt**2*(0.0000446959494108217402966857_f64)
463 split%split_step(8) = dt**2*(0.0000446959494108217402966857_f64)
464 split%split_step(9) = -0.0931750795687314526579244_f64
465 split%split_step(10) = split%split_step(7)
466 split%split_step(11) = split%split_step(8)
467 split%split_step(12) = split%split_step(6)
468 split%split_step(13) = split%split_step(4)
469 split%split_step(14) = split%split_step(5)
470 split%split_step(15) = split%split_step(3)
471 split%split_step(16) = split%split_step(1)
472 split%split_step(17) = split%split_step(2)
478 if (
present(dt))
then
479 split%nb_split_step = 11
480 sll_allocate(split%split_step(split%nb_split_step), ierr)
481 split%split_begin_T = .false.
482 split%split_step(1) = 0.083335463273305120964507_f64 &
483 - 2._f64*dt**2*(-0.00015280483587048489661_f64) &
484 + 4._f64*dt**4*(-0.0017675734111895638156_f64) &
485 - 8._f64*dt**6*(0.00021214072262165668039_f64)
486 split%split_step(2) = 0.72431592569108212422250_f64
487 split%split_step(3) = 0.827694857845135145869413_f64 &
488 - 2._f64*dt**2*(-0.010726848627286273332_f64) &
489 + 4._f64*dt**4*(0.012324362982853212700_f64)
490 split%split_step(4) = -0.4493507217041624582458844_f64
491 split%split_step(5) = -0.4110303211184402668339201_f64 &
492 - 2._f64*dt**2*(0.014962337009932798678_f64)
494 split%split_step(6) = 0.4500695920261606680467717_f64
495 split%split_step(7) = split%split_step(5)
496 split%split_step(8) = split%split_step(4)
497 split%split_step(9) = split%split_step(3)
498 split%split_step(10) = split%split_step(2)
499 split%split_step(11) = split%split_step(1)
505 print *,
'#split_case not defined'
516 sll_real64,
intent(in) :: dt
517 sll_int32,
intent(in) :: number_time_steps
519 sll_int32 :: i, istep
523 do i = 1, number_time_steps
524 call split%operatorT(dt)
525 call split%operatorV(dt)
527 split%current_time = split%current_time + dt
530 do i = 1, number_time_steps
531 call split%operatorV(dt)
532 call split%operatorT(dt)
534 split%current_time = split%current_time + dt
540 if (split%split_begin_T)
then
544 call split%operatorT(split%split_step(istep)*dt)
545 do i = 1, number_time_steps - 1
549 do while (istep < split%nb_split_step - 2)
550 call split%operatorV(split%split_step(istep)*dt)
552 call split%operatorT(split%split_step(istep)*dt)
557 call split%operatorV(split%split_step(istep)*dt)
559 call split%operatorT(2*split%split_step(istep)*dt)
565 do while (istep < split%nb_split_step)
566 call split%operatorV(split%split_step(istep)*dt)
568 call split%operatorT(split%split_step(istep)*dt)
572 split%current_time = split%current_time + dt
577 call split%operatorV(split%split_step(istep)*dt)
578 do i = 1, number_time_steps - 1
582 do while (istep < split%nb_split_step - 2)
583 call split%operatorT(split%split_step(istep)*dt)
585 call split%operatorV(split%split_step(istep)*dt)
590 call split%operatorT(split%split_step(istep)*dt)
592 call split%operatorV(2*split%split_step(istep)*dt)
598 do while (istep < split%nb_split_step)
599 call split%operatorT(split%split_step(istep)*dt)
601 call split%operatorV(split%split_step(istep)*dt)
605 split%current_time = split%current_time + dt
Base class of operator splitting library. It is only used with particle-in-cell method.
subroutine operator(this, dt)
dummy implementation of an sll_t_operator_splitting needed to provide the interface....
integer(kind=i32), parameter, public sll_p_triple_jump_vtv
Triple jump splitting V first (order 4)
integer(kind=i32), parameter sll_order6vpnew2_vtv
Specific Vlasov-Poisson splitting V first (order 6)
integer(kind=i32), parameter sll_order6vpot_vtv
Specific Vlasov-Poisson splitting V first (order 6)
class(sll_t_operator_splitting) function, pointer new_operator_splitting(split_case, split_step, nb_split_step, split_begin_T, dt)
Returns a pointer to a heap-allocated sll_t_operator_splitting object.
integer(kind=i32), parameter sll_user_defined
user defined splitting
integer(kind=i32), parameter, public sll_p_order6_tvt
Order 6 splitting T first (order 6)
integer(kind=i32), parameter sll_order6vpnew_tvt
Specific Vlasov-Poisson splitting T first (order 6)
integer(kind=i32), parameter, public sll_p_lie_vt
Lie splitting V first x (order 1)
integer(kind=i32), parameter, public sll_p_order6_vtv
Order 6 splitting V first (order 6)
integer(kind=i32), parameter, public sll_p_lie_tv
Lie splitting T first (order 1)
subroutine, public sll_s_do_split_steps(split, dt, number_time_steps)
Apply the composition method for given number of times steps.
integer(kind=i32), parameter, public sll_p_strang_vtv
Strang splitting V first (order 2)
integer(kind=i32), parameter sll_order6vp_vtv
Specific Vlasov-Poisson splitting V first (order 6)
integer(kind=i32), parameter sll_order6vp_tvt
Specific Vlasov-Poisson splitting T first (order 6)
integer(kind=i32), parameter sll_order6vpnew1_vtv
Specific Vlasov-Poisson splitting V first (order 6)
integer(kind=i32), parameter, public sll_p_triple_jump_tvt
Triple jump splitting T first (order 4)
subroutine, public sll_s_operator_splitting_init(split, split_case, split_step, nb_split_step, split_begin_T, dt)
Initialises data for operator splitting.
integer(kind=i32), parameter sll_order6vp2d_vtv
Specific Vlasov-Poisson splitting V first (order 6)
integer(kind=i32), parameter, public sll_p_strang_tvt
Strang splitting T first (order 2)
operator splitting object