Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_reduction.F90
Go to the documentation of this file.
2 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 #include "sll_memory.h"
4 #include "sll_working_precision.h"
5 
6  implicit none
7 
8  public :: &
14 
15  private
16 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
17 
18  abstract interface
19  function sll_integration_discrete_1d(data, Npts, delta, params) result(res)
21  sll_real64 :: res
22  sll_real64, dimension(:), intent(in) :: data
23  sll_int32, intent(in) :: npts
24  sll_real64, intent(in) :: delta
25  sll_real64, dimension(:), intent(in), optional :: params
26  end function sll_integration_discrete_1d
27  end interface
28 
29 contains
30 
31  !--------------------------------------------------
32  ! Generic function which can be used for computing charge density
33  ! by default we suppose a uniform mesh node based
34  ! and trapezoid quadrature formula is used
35  ! as an option we can use a way to integrate
36  ! by providing a function whose signature is sll_integration_discrete_1d
37  !---------------------------------------------------
38 
40  data_4d, &
41  data_3d, &
42  Npts1, &
43  Npts2, &
44  Npts3, &
45  Npts4, &
46  delta4, &
47  integration_func, &
48  integration_func_params)
49 
50  sll_real64, dimension(:, :, :, :), intent(in) :: data_4d
51  sll_real64, dimension(:, :, :), intent(out) :: data_3d
52  sll_int32, intent(in) :: npts1
53  sll_int32, intent(in) :: npts2
54  sll_int32, intent(in) :: npts3
55  sll_int32, intent(in) :: npts4
56  sll_real64, intent(in) :: delta4
57  procedure(sll_integration_discrete_1d), optional :: integration_func
58  sll_real64, dimension(:), optional :: integration_func_params
59  sll_int32 :: i1, i2, i3, i4
60  sll_real64 :: tmp
61 
62  if (npts1 > size(data_4d, 1)) then
63  print *, '#Problem for size1 in sll_s_compute_reduction_4d_to_3d_direction4'
64  print *, '#Npts=', (/npts1, npts2, npts3, npts4/)
65  print *, '#size(data_4d)=', (/ &
66  size(data_4d, 1), &
67  size(data_4d, 2), &
68  size(data_4d, 3), &
69  size(data_4d, 4)/)
70  stop
71  end if
72  if (npts2 > size(data_4d, 2)) then
73  print *, '#Problem for size2 in sll_s_compute_reduction_4d_to_3d_direction4'
74  stop
75  end if
76  if (npts3 > size(data_4d, 3)) then
77  print *, '#Problem for size3 in sll_s_compute_reduction_4d_to_3d_direction4'
78  stop
79  end if
80  if (npts4 > size(data_4d, 4)) then
81  print *, '#Problem for size3 in sll_s_compute_reduction_4d_to_3d_direction4'
82  stop
83  end if
84  if (.not. (present(integration_func))) then
85  do i3 = 1, npts3
86  do i2 = 1, npts2
87  do i1 = 1, npts1
88  tmp = 0.5_f64*(data_4d(i1, i2, i3, 1) &
89  + data_4d(i1, i2, i3, npts4))
90  do i4 = 2, npts4 - 1
91  tmp = tmp + data_4d(i1, i2, i3, i4)
92  end do
93  data_3d(i1, i2, i3) = tmp*delta4
94  !data_3d(i1,i2,i3) = delta4*sum(data_4d(i1,i2,i3,1:Npts4-1))
95  end do
96  end do
97  end do
98  else
99  do i3 = 1, npts3
100  do i2 = 1, npts2
101  do i1 = 1, npts1
102  data_3d(i1, i2, i3) = integration_func( &
103  data_4d(i1, i2, i3, :), &
104  npts4, &
105  delta4, &
106  integration_func_params)
107  end do
108  end do
109  end do
110  end if
112 
114  data_4d, &
115  data_3d, &
116  Npts1, &
117  Npts2, &
118  Npts3, &
119  Npts4, &
120  delta4, &
121  integration_func, &
122  integration_func_params)
123 
124  sll_real64, dimension(:, :, :, :), intent(in) :: data_4d
125  sll_real64, dimension(:, :, :), intent(inout) :: data_3d
126  sll_int32, intent(in) :: npts1
127  sll_int32, intent(in) :: npts2
128  sll_int32, intent(in) :: npts3
129  sll_int32, intent(in) :: npts4
130  sll_real64, intent(in) :: delta4
131  procedure(sll_integration_discrete_1d), optional :: integration_func
132  sll_real64, dimension(:), optional :: integration_func_params
133  sll_int32 :: i1, i2, i3, i4
134  sll_real64 :: tmp
135 
136  if (npts1 > size(data_4d, 1)) then
137  print *, '#Problem for size1 in sll_s_compute_reduction_4d_to_3d_direction4'
138  print *, '#Npts=', (/npts1, npts2, npts3, npts4/)
139  print *, '#size(data_4d)=', (/ &
140  size(data_4d, 1), &
141  size(data_4d, 2), &
142  size(data_4d, 3), &
143  size(data_4d, 4)/)
144  stop
145  end if
146  if (npts2 > size(data_4d, 2)) then
147  print *, '#Problem for size2 in sll_s_compute_reduction_4d_to_3d_direction4'
148  stop
149  end if
150  if (npts3 > size(data_4d, 3)) then
151  print *, '#Problem for size3 in sll_s_compute_reduction_4d_to_3d_direction4'
152  stop
153  end if
154  if (npts4 > size(data_4d, 4)) then
155  print *, '#Problem for size3 in sll_s_compute_reduction_4d_to_3d_direction4'
156  stop
157  end if
158  if (.not. (present(integration_func))) then
159  do i3 = 1, npts3
160  do i2 = 1, npts2
161  do i1 = 1, npts1
162  tmp = 0.5_f64*(data_4d(i1, i2, i3, 1) &
163  + data_4d(i1, i2, i3, npts4))
164  do i4 = 2, npts4 - 1
165  tmp = tmp + data_4d(i1, i2, i3, i4)
166  end do
167  data_3d(i1, i2, i3) = data_3d(i1, i2, i3) + tmp*delta4
168  !data_3d(i1,i2,i3) = delta4*sum(data_4d(i1,i2,i3,1:Npts4-1))
169  end do
170  end do
171  end do
172  else
173  do i3 = 1, npts3
174  do i2 = 1, npts2
175  do i1 = 1, npts1
176  data_3d(i1, i2, i3) = data_3d(i1, i2, i3) + integration_func( &
177  data_4d(i1, i2, i3, :), &
178  npts4, &
179  delta4, &
180  integration_func_params)
181  end do
182  end do
183  end do
184  end if
186 
188  data_4d, &
189  data_2d, &
190  Npts1, &
191  Npts2, &
192  Npts3, &
193  Npts4, &
194  delta3, &
195  delta4, &
196  integration_func, &
197  integration_func_params)
198 
199  sll_real64, dimension(:, :, :, :), intent(in) :: data_4d
200  sll_real64, dimension(:, :), intent(out) :: data_2d
201  sll_int32, intent(in) :: npts1
202  sll_int32, intent(in) :: npts2
203  sll_int32, intent(in) :: npts3
204  sll_int32, intent(in) :: npts4
205  sll_real64, intent(in) :: delta3
206  sll_real64, intent(in) :: delta4
207  procedure(sll_integration_discrete_1d), optional :: integration_func
208  sll_real64, dimension(:), optional :: integration_func_params
209  sll_int32 :: i1, i2, i3, i4
210  sll_real64 :: tmp
211 
212  if (npts1 > size(data_4d, 1)) then
213  print *, '#Problem for size1 in sll_s_compute_reduction_4d_to_2d_direction34'
214  print *, 'Npts1=', npts1, 'size(data_4d,1)', size(data_4d, 1)
215  stop
216  end if
217  if (npts2 > size(data_4d, 2)) then
218  print *, '#Problem for size2 in sll_s_compute_reduction_4d_to_2d_direction34'
219  stop
220  end if
221  if (npts3 > size(data_4d, 3)) then
222  print *, '#Problem for size3 in sll_s_compute_reduction_4d_to_2d_direction34'
223  stop
224  end if
225  if (npts4 > size(data_4d, 4)) then
226  print *, '#Problem for size3 in sll_s_compute_reduction_4d_to_2d_direction34'
227  stop
228  end if
229  if (.not. (present(integration_func))) then
230  do i2 = 1, npts2
231  do i1 = 1, npts1
232  data_2d(i1, i2) = 0._f64
233  i3 = 1
234  tmp = 0.5_f64*(data_4d(i1, i2, i3, 1) &
235  + data_4d(i1, i2, i3, npts4))
236  do i4 = 2, npts4 - 1
237  tmp = tmp + data_4d(i1, i2, i3, i4)
238  end do
239  tmp = tmp*delta4
240  data_2d(i1, i2) = data_2d(i1, i2) + 0.5_f64*tmp
241 
242  do i3 = 2, npts3 - 1
243  tmp = 0.5_f64*(data_4d(i1, i2, i3, 1) &
244  + data_4d(i1, i2, i3, npts4))
245  do i4 = 2, npts4 - 1
246  tmp = tmp + data_4d(i1, i2, i3, i4)
247  end do
248  tmp = tmp*delta4
249  data_2d(i1, i2) = data_2d(i1, i2) + tmp
250  end do
251 
252  i3 = npts3
253  tmp = 0.5_f64*(data_4d(i1, i2, i3, 1) &
254  + data_4d(i1, i2, i3, npts4))
255  do i4 = 2, npts4 - 1
256  tmp = tmp + data_4d(i1, i2, i3, i4)
257  end do
258  tmp = tmp*delta4
259  data_2d(i1, i2) = data_2d(i1, i2) + 0.5_f64*tmp
260 
261  data_2d(i1, i2) = data_2d(i1, i2)*delta3
262  end do
263  end do
264  else
265  if (present(integration_func_params)) then
266  print *, 'integration_func_params is present'
267  end if
268  print *, '#not implemented yet'
269  print *, '#in sll_s_compute_reduction_4d_to_2d_direction34'
270  stop
271  end if
273 
275  data_4d, &
276  data_2d, &
277  Npts1, &
278  Npts2, &
279  Npts3, &
280  Npts4, &
281  delta1, &
282  delta2, &
283  integration_func, &
284  integration_func_params)
285 
286  sll_real64, dimension(:, :, :, :), intent(in) :: data_4d
287  sll_real64, dimension(:, :), intent(out) :: data_2d
288  sll_int32, intent(in) :: npts1
289  sll_int32, intent(in) :: npts2
290  sll_int32, intent(in) :: npts3
291  sll_int32, intent(in) :: npts4
292  sll_real64, intent(in) :: delta1
293  sll_real64, intent(in) :: delta2
294  procedure(sll_integration_discrete_1d), optional :: integration_func
295  sll_real64, dimension(:), optional :: integration_func_params
296  sll_int32 :: i1, i2, i3, i4
297  sll_real64 :: tmp
298 
299  if (npts1 > size(data_4d, 1)) then
300  print *, '#Problem for size1 in compute_reduction_4d_to_2d_direction12'
301  print *, 'Npts1=', npts1, 'size(data_4d,1)', size(data_4d, 1)
302  stop
303  end if
304  if (npts2 > size(data_4d, 2)) then
305  print *, '#Problem for size2 in compute_reduction_4d_to_2d_direction12'
306  print *, 'Npts2=', npts2, 'size(data_4d,2)', size(data_4d, 2)
307  stop
308  end if
309  if (npts3 > size(data_4d, 3)) then
310  print *, '#Problem for size3 in compute_reduction_4d_to_2d_direction12'
311  print *, 'Npts3=', npts3, 'size(data_4d,3)', size(data_4d, 3)
312  stop
313  end if
314  if (npts4 > size(data_4d, 4)) then
315  print *, '#Problem for size3 in compute_reduction_4d_to_2d_direction12'
316  print *, 'Npts4=', npts4, 'size(data_4d,4)', size(data_4d, 4)
317  stop
318  end if
319  if (.not. (present(integration_func))) then
320  do i4 = 1, npts4
321  do i3 = 1, npts3
322  data_2d(i3, i4) = 0._f64
323  i1 = 1
324  tmp = 0.5_f64*(data_4d(i1, 1, i3, i4) &
325  + data_4d(i1, npts2, i3, i4))
326  do i2 = 2, npts2 - 1
327  tmp = tmp + data_4d(i1, i2, i3, i4)
328  end do
329  tmp = tmp*delta2
330  data_2d(i3, i4) = data_2d(i3, i4) + 0.5_f64*tmp
331 
332  do i1 = 2, npts1 - 1
333  tmp = 0.5_f64*(data_4d(i1, 1, i3, i4) &
334  + data_4d(i1, npts2, i3, i4))
335  do i2 = 2, npts2 - 1
336  tmp = tmp + data_4d(i1, i2, i3, i4)
337  end do
338  tmp = tmp*delta2
339  data_2d(i3, i4) = data_2d(i3, i4) + tmp
340  end do
341 
342  i1 = npts1
343  tmp = 0.5_f64*(data_4d(i1, 1, i3, i4) &
344  + data_4d(i1, npts2, i3, i4))
345  do i2 = 2, npts2 - 1
346  tmp = tmp + data_4d(i1, i2, i3, i4)
347  end do
348  tmp = tmp*delta2
349  data_2d(i3, i4) = data_2d(i3, i4) + 0.5_f64*tmp
350 
351  data_2d(i3, i4) = data_2d(i3, i4)*delta1
352  end do
353  end do
354  else
355  if (present(integration_func_params)) then
356  print *, 'integration_func_params is present'
357  end if
358  print *, '#not implemented yet'
359  print *, '#in compute_reduction_4d_to_2d_direction12'
360  stop
361  end if
363 
365  data_4d, &
366  data_diag_2d, &
367  Npts1, &
368  Npts2, &
369  Npts3, &
370  Npts4, &
371  delta1, &
372  delta2, &
373  integration_func, &
374  integration_func_params)
375 
376  sll_real64, dimension(:, :, :, :), intent(in) :: data_4d
377  sll_real64, dimension(:, :, :), intent(out) :: data_diag_2d
378  sll_int32, intent(in) :: npts1
379  sll_int32, intent(in) :: npts2
380  sll_int32, intent(in) :: npts3
381  sll_int32, intent(in) :: npts4
382  sll_real64, intent(in) :: delta1
383  sll_real64, intent(in) :: delta2
384  procedure(sll_integration_discrete_1d), optional :: integration_func
385  sll_real64, dimension(:), optional :: integration_func_params
386  sll_int32 :: i1, i2, i3, i4
387  sll_real64 :: tmp
388  sll_real64 :: tmp_l1
389  sll_real64 :: tmp_l2
390 
391  if (npts1 > size(data_4d, 1)) then
392  print *, '#Problem for size1 in compute_reduction_4d_to_2d_direction12'
393  print *, 'Npts1=', npts1, 'size(data_4d,1)', size(data_4d, 1)
394  stop
395  end if
396  if (npts2 > size(data_4d, 2)) then
397  print *, '#Problem for size2 in compute_reduction_4d_to_2d_direction12'
398  print *, 'Npts2=', npts2, 'size(data_4d,2)', size(data_4d, 2)
399  stop
400  end if
401  if (npts3 > size(data_4d, 3)) then
402  print *, '#Problem for size3 in compute_reduction_4d_to_2d_direction12'
403  print *, 'Npts3=', npts3, 'size(data_4d,3)', size(data_4d, 3)
404  stop
405  end if
406  if (npts4 > size(data_4d, 4)) then
407  print *, '#Problem for size3 in compute_reduction_4d_to_2d_direction12'
408  print *, 'Npts4=', npts4, 'size(data_4d,4)', size(data_4d, 4)
409  stop
410  end if
411  if (.not. (present(integration_func))) then
412  do i4 = 1, npts4
413  do i3 = 1, npts3
414  data_diag_2d(i3, i4, 1:3) = 0._f64
415  i1 = 1
416  tmp = 0.5_f64*(data_4d(i1, 1, i3, i4) &
417  + data_4d(i1, npts2, i3, i4))
418  tmp_l1 = 0.5_f64*(abs(data_4d(i1, 1, i3, i4)) &
419  + abs(data_4d(i1, npts2, i3, i4)))
420  tmp_l2 = 0.5_f64*((data_4d(i1, 1, i3, i4))**2 &
421  + (data_4d(i1, npts2, i3, i4))**2)
422  do i2 = 2, npts2 - 1
423  tmp = tmp + data_4d(i1, i2, i3, i4)
424  tmp_l1 = tmp_l1 + abs(data_4d(i1, i2, i3, i4))
425  tmp_l2 = tmp_l2 + data_4d(i1, i2, i3, i4)**2
426  end do
427  tmp = tmp*delta2
428  tmp_l1 = tmp_l1*delta2
429  tmp_l2 = tmp_l2*delta2
430  data_diag_2d(i3, i4, 1) = data_diag_2d(i3, i4, 1) + 0.5_f64*tmp
431  data_diag_2d(i3, i4, 2) = data_diag_2d(i3, i4, 2) + 0.5_f64*tmp_l1
432  data_diag_2d(i3, i4, 3) = data_diag_2d(i3, i4, 3) + 0.5_f64*tmp_l2
433 
434  do i1 = 2, npts1 - 1
435  tmp = 0.5_f64*(data_4d(i1, 1, i3, i4) &
436  + data_4d(i1, npts2, i3, i4))
437  tmp_l1 = 0.5_f64*(abs(data_4d(i1, 1, i3, i4)) &
438  + abs(data_4d(i1, npts2, i3, i4)))
439  tmp_l2 = 0.5_f64*((data_4d(i1, 1, i3, i4))**2 &
440  + (data_4d(i1, npts2, i3, i4))**2)
441  do i2 = 2, npts2 - 1
442  tmp = tmp + data_4d(i1, i2, i3, i4)
443  tmp_l1 = tmp_l1 + abs(data_4d(i1, i2, i3, i4))
444  tmp_l2 = tmp_l2 + data_4d(i1, i2, i3, i4)**2
445  end do
446  tmp = tmp*delta2
447  tmp_l1 = tmp_l1*delta2
448  tmp_l2 = tmp_l2*delta2
449  data_diag_2d(i3, i4, 1) = data_diag_2d(i3, i4, 1) + tmp
450  data_diag_2d(i3, i4, 2) = data_diag_2d(i3, i4, 2) + tmp_l1
451  data_diag_2d(i3, i4, 3) = data_diag_2d(i3, i4, 3) + tmp_l2
452  end do
453 
454  i1 = npts1
455  tmp = 0.5_f64*(data_4d(i1, 1, i3, i4) &
456  + data_4d(i1, npts2, i3, i4))
457  tmp_l1 = 0.5_f64*(abs(data_4d(i1, 1, i3, i4)) &
458  + abs(data_4d(i1, npts2, i3, i4)))
459  tmp_l2 = 0.5_f64*((data_4d(i1, 1, i3, i4))**2 &
460  + (data_4d(i1, npts2, i3, i4))**2)
461  do i2 = 2, npts2 - 1
462  tmp = tmp + data_4d(i1, i2, i3, i4)
463  tmp_l1 = tmp_l1 + abs(data_4d(i1, i2, i3, i4))
464  tmp_l2 = tmp_l2 + data_4d(i1, i2, i3, i4)**2
465  end do
466  tmp = tmp*delta2
467  tmp_l1 = tmp_l1*delta2
468  tmp_l2 = tmp_l2*delta2
469  data_diag_2d(i3, i4, 1) = data_diag_2d(i3, i4, 1) + 0.5_f64*tmp
470  data_diag_2d(i3, i4, 2) = data_diag_2d(i3, i4, 2) + 0.5_f64*tmp_l1
471  data_diag_2d(i3, i4, 3) = data_diag_2d(i3, i4, 3) + 0.5_f64*tmp_l2
472 
473  data_diag_2d(i3, i4, 1:3) = data_diag_2d(i3, i4, 1:3)*delta1
474  end do
475  end do
476  else
477  if (present(integration_func_params)) then
478  print *, 'integration_func_params is present'
479  end if
480  print *, '#not implemented yet'
481  print *, '#in compute_reduction_4d_to_2d_direction12'
482  stop
483  end if
485 
487  data_2d, &
488  res, &
489  Npts1, &
490  Npts2, &
491  delta1, &
492  delta2, &
493  integration_func, &
494  integration_func_params)
495 
496  sll_real64, dimension(:, :), intent(in) :: data_2d
497  sll_real64, intent(out) :: res
498  sll_int32, intent(in) :: npts1
499  sll_int32, intent(in) :: npts2
500  sll_real64, intent(in) :: delta1
501  sll_real64, intent(in) :: delta2
502  procedure(sll_integration_discrete_1d), optional :: integration_func
503  sll_real64, dimension(:), optional :: integration_func_params
504  sll_int32 :: i1, i2
505  sll_real64 :: tmp
506 
507  if (npts1 > size(data_2d, 1)) then
508  print *, '#Problem for size1 in sll_s_compute_reduction_2d_to_0d'
509  print *, 'Npts1=', npts1, 'size(data_2d,1)', size(data_2d, 1)
510  stop
511  end if
512  if (npts2 > size(data_2d, 2)) then
513  print *, '#Problem for size2 in sll_s_compute_reduction_2d_to_0d'
514  print *, 'Npts2=', npts2, 'size(data_2d,2)', size(data_2d, 2)
515  stop
516  end if
517  if (.not. (present(integration_func))) then
518  res = 0._f64
519  i1 = 1
520  tmp = 0.5_f64*(data_2d(i1, 1) &
521  + data_2d(i1, npts2))
522  do i2 = 2, npts2 - 1
523  tmp = tmp + data_2d(i1, i2)
524  end do
525  tmp = tmp*delta2
526  res = res + 0.5_f64*tmp
527 
528  do i1 = 2, npts1 - 1
529  tmp = 0.5_f64*(data_2d(i1, 1) &
530  + data_2d(i1, npts2))
531  do i2 = 2, npts2 - 1
532  tmp = tmp + data_2d(i1, i2)
533  end do
534  tmp = tmp*delta2
535  res = res + tmp
536  end do
537 
538  i1 = npts1
539  tmp = 0.5_f64*(data_2d(i1, 1) &
540  + data_2d(i1, npts2))
541  do i2 = 2, npts2 - 1
542  tmp = tmp + data_2d(i1, i2)
543  end do
544  tmp = tmp*delta2
545  res = res + 0.5_f64*tmp
546 
547  res = res*delta1
548  else
549  if (present(integration_func_params)) then
550  print *, 'integration_func_params is present'
551  end if
552  print *, '#not implemented yet'
553  print *, '#in sll_s_compute_reduction_2d_to_0d'
554  stop
555  end if
556  end subroutine sll_s_compute_reduction_2d_to_0d
557 
558  function sll_f_compute_integral_trapezoid_1d(data, Npts, delta, func_params) result(res)
559  sll_real64, dimension(:), intent(in) :: data
560  sll_int32, intent(in) :: npts
561  sll_real64, intent(in) :: delta
562  sll_real64, dimension(:), intent(in), optional :: func_params
563  sll_real64 :: res
564  sll_int32 :: i
565 
566  if (present(func_params)) then
567  if (size(func_params) > 100) then
568  print *, '#size of func_params is >100'
569  end if
570  end if
571 
572  res = 0.5*(data(1) + data(npts))
573  do i = 2, npts - 1
574  res = res + data(i)
575  end do
576  res = res*delta
578 
579  function compute_integral_conservative_1d(data, Npts, node_positions) result(res)
580  sll_real64, dimension(:), intent(in) :: data
581  sll_int32, intent(in) :: npts
582  sll_real64, dimension(:), intent(in) :: node_positions
583  sll_real64 :: res
584  sll_int32 :: i
585 
586  res = 0._f64
587  do i = 1, npts - 1
588  res = res + data(i)*(node_positions(i + 1) - node_positions(i))
589  end do
591 
592 end module sll_m_reduction
real(kind=f64) function compute_integral_conservative_1d(data, Npts, node_positions)
subroutine, public sll_s_compute_reduction_2d_to_0d(data_2d, res, Npts1, Npts2, delta1, delta2, integration_func, integration_func_params)
real(kind=f64) function, public sll_f_compute_integral_trapezoid_1d(data, Npts, delta, func_params)
subroutine, public sll_s_compute_reduction_diag_4d_to_2d_direction12(data_4d, data_diag_2d, Npts1, Npts2, Npts3, Npts4, delta1, delta2, integration_func, integration_func_params)
subroutine, public sll_s_compute_reduction_4d_to_2d_direction34(data_4d, data_2d, Npts1, Npts2, Npts3, Npts4, delta3, delta4, integration_func, integration_func_params)
subroutine compute_reduction_4d_to_2d_direction12(data_4d, data_2d, Npts1, Npts2, Npts3, Npts4, delta1, delta2, integration_func, integration_func_params)
subroutine, public sll_s_compute_reduction_4d_to_3d_direction4(data_4d, data_3d, Npts1, Npts2, Npts3, Npts4, delta4, integration_func, integration_func_params)
subroutine compute_reduction_4d_to_3d_direction4_accumulate(data_4d, data_3d, Npts1, Npts2, Npts3, Npts4, delta4, integration_func, integration_func_params)
Module to select the kind parameter.
    Report Typos and Errors