Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_time_splitting_coeff.F90
Go to the documentation of this file.
1 #ifndef DOXYGEN_SHOULD_SKIP_THIS
2 
3 !**************************************************************
4 ! Copyright INRIA
5 ! Authors :
6 ! CALVI project team
7 !
8 ! This code SeLaLib (for Semi-Lagrangian-Library)
9 ! is a parallel library for simulating the plasma turbulence
10 ! in a tokamak.
11 !
12 ! This software is governed by the CeCILL-B license
13 ! under French law and abiding by the rules of distribution
14 ! of free software. You can use, modify and redistribute
15 ! the software under the terms of the CeCILL-B license as
16 ! circulated by CEA, CNRS and INRIA at the following URL
17 ! "http://www.cecill.info".
18 !**************************************************************
19 !
20 ! Define coefficients for time splitting with semi-lagrangian
21 !
23 !
24 
25 module sll_m_time_splitting_coeff
26 
27 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
28 #include "sll_memory.h"
29 #include "sll_working_precision.h"
30 
31  implicit none
32 
33  public :: &
34  sll_f_new_time_splitting_coeff, &
35  sll_p_lie_tv, &
36  sll_p_lie_vt, &
37  sll_p_order6_tvt, &
38  sll_p_order6_vtv, &
39  sll_p_order6vp2d_vtv, &
40  sll_p_order6vp_tvt, &
41  sll_p_order6vp_vtv, &
42  sll_p_order6vpnew1_vtv, &
43  sll_p_order6vpnew2_vtv, &
44  sll_p_order6vpnew_tvt, &
45  sll_p_order6vpot_vtv, &
46  sll_p_order6vpotnew1_vtv, &
47  sll_p_order6vpotnew2_vtv, &
48  sll_p_order6vpotnew3_vtv, &
49  sll_p_strang_tvt, &
50  sll_p_strang_vtv, &
51  sll_p_triple_jump_tvt, &
52  sll_p_triple_jump_vtv, &
53  sll_t_splitting_coeff
54 
55  private
56 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
57 
58  sll_int32, parameter :: sll_user_defined = -3
59  sll_int32, parameter :: sll_p_lie_tv = -2
60  sll_int32, parameter :: sll_p_lie_vt = -1
61  sll_int32, parameter :: sll_p_strang_tvt = 0
62  sll_int32, parameter :: sll_p_strang_vtv = 1
63  sll_int32, parameter :: sll_p_triple_jump_tvt = 2
64  sll_int32, parameter :: sll_p_triple_jump_vtv = 3
65  !sll_int32, parameter :: sll_p_order6_tvt = 4
66  sll_int32, parameter :: sll_p_order6_vtv = 5
67  sll_int32, parameter :: sll_p_order6vp_tvt = 6
68  sll_int32, parameter :: sll_p_order6vp_vtv = 7
69  sll_int32, parameter :: sll_p_order6vpnew_tvt = 8
70  sll_int32, parameter :: sll_p_order6vpnew1_vtv = 9
71  sll_int32, parameter :: sll_p_order6vpnew2_vtv = 10
72  sll_int32, parameter :: sll_p_order6vp2d_vtv = 11
73  sll_int32, parameter :: sll_p_order6vpot_vtv = 12
74  sll_int32, parameter :: sll_p_order6_tvt = 13
75  sll_int32, parameter :: sll_p_order6vpotnew1_vtv = 14
76  sll_int32, parameter :: sll_p_order6vpotnew2_vtv = 15
77  sll_int32, parameter :: sll_p_order6vpotnew3_vtv = 16
78 
79  type sll_t_splitting_coeff
80  sll_int32 :: split_case
81  sll_real64, dimension(:), pointer :: split_step
82  sll_int32 :: nb_split_step
83  logical :: split_begin_T
84  sll_int32 :: dim_split_v
85  end type sll_t_splitting_coeff
86 
87 contains
88  function sll_f_new_time_splitting_coeff( &
89  split_case, &
90  split_step, &
91  nb_split_step, &
92  split_begin_T, &
93  dt) &
94  result(split)
95  type(sll_t_splitting_coeff), pointer :: split
96  sll_int32, intent(in) :: split_case
97  sll_real64, dimension(:), intent(in), optional :: split_step
98  sll_int32, intent(in), optional :: nb_split_step
99  logical, intent(in), optional :: split_begin_T
100  sll_real64, intent(in), optional :: dt
101  sll_int32 :: ierr
102 
103  sll_allocate(split, ierr)
104  call initialize_time_splitting_coeff( &
105  split, &
106  split_case, &
107  split_step, &
108  nb_split_step, &
109  split_begin_t, &
110  dt)
111 
112  end function sll_f_new_time_splitting_coeff
113 
114  subroutine initialize_time_splitting_coeff( &
115  split, &
116  split_case, &
117  split_step, &
118  nb_split_step, &
119  split_begin_T, &
120  dt)
121 
122  type(sll_t_splitting_coeff) :: split
123  sll_int32, intent(in) :: split_case
124  sll_real64, dimension(:), intent(in), optional :: split_step
125  sll_int32, intent(in), optional :: nb_split_step
126  logical, intent(in), optional :: split_begin_T
127  sll_real64, intent(in), optional :: dt
128  sll_int32 :: ierr
129 
130  split%dim_split_V = 1
131 
132  if (split_case == sll_user_defined) then
133  if (.not. ((present(split_step)) &
134  .and. (present(nb_split_step)) &
135  .and. (present(split_begin_t)))) then
136  print *, '#provide split_step, nb_split_step,split_begin_T'
137  print *, '#in initialize_time_splitting_coeff'
138  stop
139  end if
140  else
141  if ( &
142  (present(split_step)) &
143  .or. (present(nb_split_step)) &
144  .or. (present(split_begin_t)) &
145  ) then
146  print *, '# do not provide split_step, nb_split_step, split_begin_T'
147  print *, '#in initialize_time_splitting_coeff'
148  stop
149  end if
150  end if
151 ! if(associated(split%split_step))then
152 ! print *,'#split%split_step should not be allocated'
153 ! print *,'#in initialize_time_splitting_coeff'
154 ! stop
155 ! endif
156 
157  select case (split_case)
158  case (sll_user_defined)
159  if (size(split_step) < nb_split_step) then
160  print *, '#bad size for split_step', size(split_step), nb_split_step
161  print *, '#in initialize_time_splitting_coeff'
162  stop
163  else
164  sll_allocate(split%split_step(nb_split_step), ierr)
165  split%split_step(1:nb_split_step) = split_step(1:nb_split_step)
166  split%nb_split_step = nb_split_step
167  split%split_begin_T = split_begin_t
168  end if
169  case (sll_p_lie_tv)
170  split%nb_split_step = 2
171  sll_allocate(split%split_step(split%nb_split_step), ierr)
172  split%split_begin_T = .true.
173  split%split_step(1) = 1._f64
174  split%split_step(2) = 1._f64
175  case (sll_p_lie_vt)
176  split%nb_split_step = 2
177  sll_allocate(split%split_step(split%nb_split_step), ierr)
178  split%split_begin_T = .false.
179  split%split_step(1) = 1._f64
180  split%split_step(2) = 1._f64
181  case (sll_p_strang_tvt) ! Strang splitting TVT
182  split%nb_split_step = 3
183  sll_allocate(split%split_step(split%nb_split_step), ierr)
184  split%split_begin_T = .true.
185  split%split_step(1) = 0.5_f64
186  split%split_step(2) = 1._f64
187  split%split_step(3) = split%split_step(1)
188  case (sll_p_strang_vtv) ! Strang splitting VTV
189  split%nb_split_step = 3
190  sll_allocate(split%split_step(split%nb_split_step), ierr)
191  split%split_begin_T = .false.
192  split%split_step(1) = 0.5_f64
193  split%split_step(2) = 1._f64
194  split%split_step(3) = split%split_step(1)
195  case (sll_p_triple_jump_tvt) ! triple jump TVT
196  split%nb_split_step = 7
197  sll_allocate(split%split_step(split%nb_split_step), ierr)
198  split%split_begin_T = .true.
199  split%split_step(1) = 0.675603595979829_f64
200  split%split_step(2) = 1.351207191959658_f64
201  split%split_step(3) = -0.17560359597982855_f64
202  split%split_step(4) = -1.702414383919315_f64
203  split%split_step(5) = split%split_step(3)
204  split%split_step(6) = split%split_step(2)
205  split%split_step(7) = split%split_step(1)
206  case (sll_p_triple_jump_vtv) ! triple jump VTV
207  split%nb_split_step = 7
208  sll_allocate(split%split_step(split%nb_split_step), ierr)
209  split%split_begin_T = .false.
210  split%split_step(1) = 0.675603595979829_f64
211  split%split_step(2) = 1.351207191959658_f64
212  split%split_step(3) = -0.17560359597982855_f64
213  split%split_step(4) = -1.702414383919315_f64
214  split%split_step(5) = split%split_step(3)
215  split%split_step(6) = split%split_step(2)
216  split%split_step(7) = split%split_step(1)
217  case (sll_p_order6_vtv) ! Order 6 VTV (O6-11 of Blanes)
218  split%nb_split_step = 23
219  sll_allocate(split%split_step(split%nb_split_step), ierr)
220  split%split_begin_T = .false.
221  split%split_step(1) = 0.0414649985182624_f64
222  split%split_step(2) = 0.123229775946271_f64
223  split%split_step(3) = 0.198128671918067_f64
224  split%split_step(4) = 0.290553797799558_f64
225  split%split_step(5) = -0.0400061921041533_f64
226  split%split_step(6) = -0.127049212625417_f64
227  split%split_step(7) = 0.0752539843015807_f64
228  split%split_step(8) = -0.246331761062075_f64
229  split%split_step(9) = -0.0115113874206879_f64
230  split%split_step(10) = 0.357208872795928_f64
231  split%split_step(11) = 0.23666992478693111_f64
232  split%split_step(12) = 0.20477705429147008_f64
233  split%split_step(13) = split%split_step(11)
234  split%split_step(14) = split%split_step(10)
235  split%split_step(15) = split%split_step(9)
236  split%split_step(16) = split%split_step(8)
237  split%split_step(17) = split%split_step(7)
238  split%split_step(18) = split%split_step(6)
239  split%split_step(19) = split%split_step(5)
240  split%split_step(20) = split%split_step(4)
241  split%split_step(21) = split%split_step(3)
242  split%split_step(22) = split%split_step(2)
243  split%split_step(23) = split%split_step(1)
244  case (sll_p_order6_tvt) ! Order 6 TVT (O6-14 of Blanes)
245  split%nb_split_step = 29
246  sll_allocate(split%split_step(split%nb_split_step), ierr)
247  split%split_begin_T = .true.
248  split%split_step(1) = 0.0378593198406116_f64
249  split%split_step(2) = 0.09171915262446165_f64
250  split%split_step(3) = 0.102635633102435_f64
251  split%split_step(4) = 0.183983170005006_f64
252  split%split_step(5) = -0.0258678882665587_f64
253  split%split_step(6) = -0.05653436583288827_f64
254  split%split_step(7) = 0.314241403071447_f64
255  split%split_step(8) = 0.004914688774712854_f64
256  split%split_step(9) = -0.130144459517415_f64
257  split%split_step(10) = 0.143761127168358_f64
258  split%split_step(11) = 0.106417700369543_f64
259  split%split_step(12) = 0.328567693746804_f64
260  split%split_step(13) = -0.00879424312851058_f64
261  split%split_step(14) = -0.196411466486454234_f64
262  split%split_step(15) = 0.20730506905689536_f64
263  split%split_step(16) = split%split_step(14)
264  split%split_step(17) = split%split_step(13)
265  split%split_step(18) = split%split_step(12)
266  split%split_step(19) = split%split_step(11)
267  split%split_step(20) = split%split_step(10)
268  split%split_step(21) = split%split_step(9)
269  split%split_step(22) = split%split_step(8)
270  split%split_step(23) = split%split_step(7)
271  split%split_step(24) = split%split_step(6)
272  split%split_step(25) = split%split_step(5)
273  split%split_step(26) = split%split_step(4)
274  split%split_step(27) = split%split_step(3)
275  split%split_step(28) = split%split_step(2)
276  split%split_step(29) = split%split_step(1)
277 
278  case (sll_p_order6vp_tvt) ! Order 6 for Vlasov-Poisson TVT
279  if (present(dt)) then
280  split%nb_split_step = 9
281  sll_allocate(split%split_step(split%nb_split_step), ierr)
282  split%split_begin_T = .true.
283  split%split_step(1) = 0.1095115577513980413559540_f64
284  split%split_step(2) = 0.268722208204814693684441_f64 &
285  - 2._f64*dt**2*0.000805681667096178271312_f64 &
286  + 4._f64*dt**4*0.000017695766224036466792_f64
287  split%split_step(3) = 0.4451715080955340951457244_f64
288  split%split_step(4) = 0.2312777917951853063155588_f64 &
289  - 2._f64*dt**2*0.003955911930042478239763_f64 &
290  + 4._f64*dt**4*0.000052384078562246674986_f64
291  split%split_step(5) = -0.1093661316938642730033570_f64
292  split%split_step(6) = split%split_step(4)
293  split%split_step(7) = split%split_step(3)
294  split%split_step(8) = split%split_step(2)
295  split%split_step(9) = split%split_step(1)
296  else
297  print *, '#provide dt for use of case sll_p_order6vp_tvt=', sll_p_order6vp_tvt
298  stop
299  end if
300  case (sll_p_order6vp_vtv) ! Order 6 for Vlasov-Poisson VTV
301  if (present(dt)) then
302  split%nb_split_step = 9
303  sll_allocate(split%split_step(split%nb_split_step), ierr)
304  split%split_begin_T = .false.
305  split%split_step(1) = 0.359950808794143627485664_f64 &
306  - 2._f64*dt**2*(-0.01359558332625151635_f64) &
307  + 4._f64*dt**4*(-8.562814848565929e-6_f64)
308  split%split_step(2) = 1.079852426382430882456991_f64
309  split%split_step(3) = -0.1437147273026540434771131_f64 &
310  - 2._f64*dt**2*(-0.00385637757897273261_f64) &
311  + 4._f64*dt**4*(0.0004883788785819335822_f64)
312  split%split_step(4) = -0.579852426382430882456991_f64
313  split%split_step(5) = 0.567527837017020831982899_f64 &
314  - 2._f64*dt**2*(-0.03227361602037480885_f64) &
315  + 4._f64*dt**4*0.002005141087312622342_f64
316  split%split_step(6) = split%split_step(4)
317  split%split_step(7) = split%split_step(3)
318  split%split_step(8) = split%split_step(2)
319  split%split_step(9) = split%split_step(1)
320  else
321  print *, '#provide dt for use of case sll_p_order6vp_vtv=', sll_p_order6vp_vtv
322  stop
323  end if
324  case (sll_p_order6vpnew_tvt) ! Order 6 for Vlasov-Poisson TVT (new)
325  if (present(dt)) then
326  split%nb_split_step = 9
327  sll_allocate(split%split_step(split%nb_split_step), ierr)
328  split%split_begin_T = .true.
329  split%split_step(1) = 0.1095115577513980413559540_f64
330  split%split_step(2) = 0.268722208204814693684441_f64 &
331  - 2._f64*dt**2*0.000805681667096178271312_f64 &
332  + 4._f64*dt**4*(8.643923349886021963e-6_f64) &
333  - 8._f64*dt**6*(1.4231479258353431522e-6_f64)
334  split%split_step(3) = 0.4451715080955340951457244_f64
335  split%split_step(4) = 0.2312777917951853063155588_f64 &
336  - 2._f64*dt**2*0.003955911930042478239763_f64 &
337  + 4._f64*dt**4*(0.000061435921436397119815_f64)
338  split%split_step(5) = -0.1093661316938642730033570_f64
339  split%split_step(6) = split%split_step(4)
340  split%split_step(7) = split%split_step(3)
341  split%split_step(8) = split%split_step(2)
342  split%split_step(9) = split%split_step(1)
343  else
344  print *, '#provide dt for use of case sll_p_order6vpnew_tvt=', sll_p_order6vpnew_tvt
345  stop
346  end if
347 
348  case (sll_p_order6vpnew1_vtv) ! Order 6 for Vlasov-Poisson VTV (new1)
349  if (present(dt)) then
350  split%nb_split_step = 11
351  sll_allocate(split%split_step(split%nb_split_step), ierr)
352  split%split_begin_T = .false.
353  split%split_step(1) = 0.0490864609761162454914412_f64 &
354  - 2._f64*dt**2*(0.0000697287150553050840999_f64)
355  split%split_step(2) = 0.1687359505634374224481957_f64
356  split%split_step(3) = 0.2641776098889767002001462_f64 &
357  - 2._f64*dt**2*(0.000625704827430047189169_f64) &
358  + 4._f64*dt**4*(-2.91660045768984781644e-6_f64)
359  split%split_step(4) = 0.377851589220928303880766_f64
360  split%split_step(5) = 0.1867359291349070543084126_f64 &
361  - 2._f64*dt**2*(0.00221308512404532556163_f64) &
362  + 4._f64*dt**4*(0.0000304848026170003878868_f64) &
363  - 8._f64*dt**6*(4.98554938787506812159e-7_f64)
364  split%split_step(6) = -0.0931750795687314526579244_f64
365  split%split_step(7) = split%split_step(5)
366  split%split_step(8) = split%split_step(4)
367  split%split_step(9) = split%split_step(3)
368  split%split_step(10) = split%split_step(2)
369  split%split_step(11) = split%split_step(1)
370  else
371  print *, '#provide dt for use of case sll_p_order6vpnew1_vtv=', sll_p_order6vpnew1_vtv
372  stop
373  end if
374  case (sll_p_order6vp2d_vtv) ! Order 6 for Vlasov-Poisson VTV 2D
375  if (present(dt)) then
376  split%nb_split_step = 11
377  sll_allocate(split%split_step(split%nb_split_step), ierr)
378  split%split_begin_T = .false.
379  split%split_step(1) = 0.0490864609761162454914412_f64 &
380  + 2._f64*dt**2*(0.00166171386175851683711044_f64)
381  split%split_step(2) = 0.1687359505634374224481957_f64
382  split%split_step(3) = 0.2641776098889767002001462_f64 &
383  - 2._f64*dt**2*(0.00461492847770001641230401_f64)
384  split%split_step(4) = 0.377851589220928303880766_f64
385  split%split_step(5) = 0.1867359291349070543084126_f64 &
386  + 2._f64*dt**2*(0.0000446959494108217402966857_f64)
387  split%split_step(6) = -0.0931750795687314526579244_f64
388  split%split_step(7) = split%split_step(5)
389  split%split_step(8) = split%split_step(4)
390  split%split_step(9) = split%split_step(3)
391  split%split_step(10) = split%split_step(2)
392  split%split_step(11) = split%split_step(1)
393  else
394  print *, '#provide dt for use of case sll_p_order6vp2d_vtv=', sll_p_order6vp2d_vtv
395  stop
396  end if
397  case (sll_p_order6vpot_vtv) ! Order 6 for Vlasov-Poisson VTV 2D with potential modif
398  if (present(dt)) then
399  split%nb_split_step = 11
400  split%dim_split_V = 2
401  sll_allocate(split%split_step(17), ierr)
402  split%split_begin_T = .false.
403  split%split_step(1) = 0.0490864609761162454914412_f64 &
404  + 2._f64*dt**2*(0.00166171386175851683711044_f64)
405  split%split_step(2) = dt**2*(0.00166171386175851683711044_f64)
406  split%split_step(3) = 0.1687359505634374224481957_f64
407  split%split_step(4) = 0.2641776098889767002001462_f64 &
408  - 2._f64*dt**2*(0.00461492847770001641230401_f64)
409  split%split_step(5) = -dt**2*(0.00461492847770001641230401_f64)
410  split%split_step(6) = 0.377851589220928303880766_f64
411  split%split_step(7) = 0.1867359291349070543084126_f64 &
412  + 2._f64*dt**2*(0.0000446959494108217402966857_f64)
413  split%split_step(8) = dt**2*(0.0000446959494108217402966857_f64)
414  split%split_step(9) = -0.0931750795687314526579244_f64
415  split%split_step(10) = split%split_step(7)
416  split%split_step(11) = split%split_step(8)
417  split%split_step(12) = split%split_step(6)
418  split%split_step(13) = split%split_step(4)
419  split%split_step(14) = split%split_step(5)
420  split%split_step(15) = split%split_step(3)
421  split%split_step(16) = split%split_step(1)
422  split%split_step(17) = split%split_step(2)
423  else
424  print *, '#provide dt for use of case sll_p_order6vpot_vtv=', sll_p_order6vpot_vtv
425  stop
426  end if
427 
428  case (sll_p_order6vpotnew1_vtv) ! Order 6 for Vlasov-Poisson VTV 2D with potential modif
429  !we change also sign here
430  if (present(dt)) then
431  split%nb_split_step = 9
432  split%dim_split_V = 2
433  sll_allocate(split%split_step(14), ierr)
434  split%split_begin_T = .false.
435  !b1 a1 b2 a2 b3 a2 b2 a1 b1
436  !b1
437  split%split_step(1) = 0.359950808794143627485664_f64 &
438  + 2._f64*dt**2*(0._f64)
439  split%split_step(2) = dt**2*(0._f64)
440  !a1
441  split%split_step(3) = 1.079852426382430882456991_f64
442  !b2
443  split%split_step(4) = -0.1437147273026540434771131_f64 &
444  + 2._f64*dt**2*(0.0139652542242388403673_f64)
445  split%split_step(5) = dt**2*(0.0139652542242388403673_f64)
446  !a2
447  split%split_step(6) = -0.579852426382430882456991_f64
448  !b3
449  split%split_step(7) = 0.567527837017020831982899_f64 &
450  + 2._f64*dt**2*(0.039247029382345626020_f64)
451  split%split_step(8) = dt**2*(0.039247029382345626020_f64)
452  !a2
453  split%split_step(9) = split%split_step(6)
454  !b2
455  split%split_step(10) = split%split_step(4)
456  split%split_step(11) = split%split_step(5)
457  !a1
458  split%split_step(12) = split%split_step(3)
459  !b1
460  split%split_step(13) = split%split_step(1)
461  split%split_step(14) = split%split_step(2)
462  else
463  print *, '#provide dt for use of case sll_p_order6vpotnew1_vtv=', &
464  sll_p_order6vpotnew1_vtv
465  stop
466  end if
467 
468  case (sll_p_order6vpotnew2_vtv) ! Order 6 for Vlasov-Poisson VTV 2D with potential modif
469  !warning we try with sign change of dt**2
470  !this seems to be the right choice
471  if (present(dt)) then
472  split%nb_split_step = 11
473  split%dim_split_V = 2
474  sll_allocate(split%split_step(17), ierr)
475  split%split_begin_T = .false.
476  !b1 a1 b2 a2 b3 a3 b3 a2 b2 a1 b1
477  !b1
478  split%split_step(1) = 0.086971698963920047813358_f64 &
479  + 2._f64*dt**2*(1.98364114652831655458915e-6_f64)
480  split%split_step(2) = dt**2*(1.98364114652831655458915e-6_f64)
481  !a1
482  split%split_step(3) = 0.303629319055488881944104_f64
483  !b2
484  split%split_step(4) = 0.560744966588102145251453_f64 &
485  - 2._f64*dt**2*(0.00553752115152236516667268_f64)
486  split%split_step(5) = -dt**2*(0.00553752115152236516667268_f64)
487  !a2
488  split%split_step(6) = 0.303629319055488881944104_f64
489  !b3
490  split%split_step(7) = -0.1477166655520221930648117_f64 &
491  - 2._f64*dt**2*(0.00284218110811634663914191_f64)
492  split%split_step(8) = -dt**2*(0.00284218110811634663914191_f64)
493  !a3
494  split%split_step(9) = -0.2145172762219555277764167_f64
495  !b3
496  split%split_step(10) = split%split_step(7)
497  split%split_step(11) = split%split_step(8)
498  !a2
499  split%split_step(12) = split%split_step(6)
500  !b2
501  split%split_step(13) = split%split_step(4)
502  split%split_step(14) = split%split_step(5)
503  !a1
504  split%split_step(15) = split%split_step(3)
505  !b1
506  split%split_step(16) = split%split_step(1)
507  split%split_step(17) = split%split_step(2)
508  else
509  print *, '#provide dt for use of case sll_p_order6vpot_vtv=', sll_p_order6vpot_vtv
510  stop
511  end if
512 
513  case (sll_p_order6vpotnew3_vtv) ! Order 6 for Vlasov-Poisson VTV 2D with potential modif
514  !we change also sign
515  if (present(dt)) then
516  split%nb_split_step = 13
517  split%dim_split_V = 2
518  sll_allocate(split%split_step(20), ierr)
519  split%split_begin_T = .false.
520  !b1 a1 b2 a2 b3 a3 b4 a3 b3 a2 b2 a1 b1
521  !b1
522  split%split_step(1) = 0.0482332301753032567427580_f64 &
523  + 2._f64*dt**2*(0.0002566567904012107264_f64)
524  split%split_step(2) = dt**2*(0.0002566567904012107264_f64)
525  !a1
526  split%split_step(3) = 0.2701015188126056215752542_f64
527  !b2
528  split%split_step(4) = 0.0482332301753032567427580_f64 &
529  + 2._f64*dt**2*(0.0009439771580927593579_f64)
530  split%split_step(5) = dt**2*(0.0009439771580927593579_f64)
531  !a2
532  split%split_step(6) = -0.108612186368692920020654_f64
533  !b3
534  split%split_step(7) = 0.2361392603742494444753990_f64 &
535  - 2._f64*dt**2*(0.002494619878121813220_f64)
536  split%split_step(8) = -dt**2*(0.002494619878121813220_f64)
537  !a3
538  split%split_step(9) = 0.3385106675560872984454001_f64
539  !b4
540  split%split_step(10) = 0.3347885585502880840781703_f64 &
541  - 2._f64*dt**2*(0.002670269183371982607658111_f64)
542  split%split_step(11) = -dt**2*(0.002670269183371982607658111_f64)
543  !a3
544  split%split_step(12) = split%split_step(9)
545  !b3
546  split%split_step(13) = split%split_step(7)
547  split%split_step(14) = split%split_step(8)
548  !a2
549  split%split_step(15) = split%split_step(6)
550  !b2
551  split%split_step(16) = split%split_step(4)
552  split%split_step(17) = split%split_step(5)
553  !a1
554  split%split_step(18) = split%split_step(3)
555  !b1
556  split%split_step(19) = split%split_step(1)
557  split%split_step(20) = split%split_step(2)
558  else
559  print *, '#provide dt for use of case sll_p_order6vpot_vtv=', sll_p_order6vpot_vtv
560  stop
561  end if
562 
563  case (sll_p_order6vpnew2_vtv) ! Order 6 for Vlasov-Poisson VTV (new2)
564  if (present(dt)) then
565  split%nb_split_step = 11
566  sll_allocate(split%split_step(split%nb_split_step), ierr)
567  split%split_begin_T = .false.
568  split%split_step(1) = 0.083335463273305120964507_f64 &
569  - 2._f64*dt**2*(-0.00015280483587048489661_f64) &
570  + 4._f64*dt**4*(-0.0017675734111895638156_f64) &
571  - 8._f64*dt**6*(0.00021214072262165668039_f64)
572  split%split_step(2) = 0.72431592569108212422250_f64
573  split%split_step(3) = 0.827694857845135145869413_f64 &
574  - 2._f64*dt**2*(-0.010726848627286273332_f64) &
575  + 4._f64*dt**4*(0.012324362982853212700_f64)
576  split%split_step(4) = -0.4493507217041624582458844_f64
577  split%split_step(5) = -0.4110303211184402668339201_f64 &
578  - 2._f64*dt**2*(0.014962337009932798678_f64)
579  ! +4._f64*sim%dt**4*(-0.20213028317750837901_f64)
580  split%split_step(6) = 0.4500695920261606680467717_f64
581  split%split_step(7) = split%split_step(5)
582  split%split_step(8) = split%split_step(4)
583  split%split_step(9) = split%split_step(3)
584  split%split_step(10) = split%split_step(2)
585  split%split_step(11) = split%split_step(1)
586  else
587  print *, '#provide dt for use of case sll_p_order6vpnew2_vtv=', sll_p_order6vpnew2_vtv
588  stop
589  end if
590  case default
591  print *, '#split_case not defined'
592  stop
593  end select
594  end subroutine initialize_time_splitting_coeff
595 
596 end module sll_m_time_splitting_coeff
597 
598 #endif /* DOXYGEN_SHOULD_SKIP_THIS */
    Report Typos and Errors