Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_parallel_array_initializer.F90
Go to the documentation of this file.
2 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 #include "sll_working_precision.h"
4 
5  use sll_m_cartesian_meshes, only: &
9 
13 
16 
17  use sll_m_remapper, only: &
18  sll_o_compute_local_sizes, &
19  sll_t_layout_2d, &
20  sll_t_layout_4d, &
21  sll_o_local_to_global
22 
23  implicit none
24 
25  public :: &
30 
31 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
32 
33  ! The parallel array initializer module helps initialize multidimensional
34  ! arrays that are distributed among multiple processors. Their distribution
35  ! should be described by a layout object (see the documentation for the
36  ! sll_remap module. The initializer assumes that the array is connected
37  ! to one or more logical meshes and possibly, their corresponding
38  ! transformations. Given the various ways in which the coordinates can
39  ! be grouped in terms of their associated logical meshes and transformations,
40  ! there are multiple versions of the initialization subroutines. The common
41  ! parameter is a user-provided function with the signature described in
42  ! the following abstract type.
43 
47  !module procedure sll_4d_parallel_array_initializer_cartesian_logical_2d_2d
49  end interface
50 
53  end interface
54 
58  end interface
59 
60 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
61 
62 contains
63 
64 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
65 
67  layout, &
68  mesh2d, &
69  array, &
70  func, &
71  func_params)
72 
73  type(sll_t_layout_2d), pointer :: layout
74  type(sll_t_cartesian_mesh_2d), pointer :: mesh2d
75  sll_real64, dimension(:, :), intent(out) :: array
76  procedure(sll_i_scalar_initializer_2d) :: func
77  sll_real64, dimension(:), intent(in) :: func_params
78 
79  sll_int32 :: i
80  sll_int32 :: j
81  sll_int32 :: loc_size_x1
82  sll_int32 :: loc_size_x2
83  sll_real64 :: delta1
84  sll_real64 :: delta2
85  sll_real64 :: eta1_min
86  sll_real64 :: eta2_min
87  sll_real64 :: eta1
88  sll_real64 :: eta2
89  !sll_real64 :: x1
90  !sll_real64 :: x2
91  sll_int32, dimension(1:2) :: gi ! global indices in the distributed array
92 
93  if (.not. associated(layout)) then
94  print *, '#sll_o_2d_parallel_array_initializer_cartesian error: ', &
95  '#passed layout is uninitialized.'
96  end if
97 
98  if (.not. associated(mesh2d)) then
99  print *, '#sll_o_2d_parallel_array_initializer_cartesian error: ', &
100  '#passed mesh2d_eta1_eta2 argument is uninitialized.'
101  end if
102 
103  call sll_o_compute_local_sizes(layout, loc_size_x1, loc_size_x2)
104 
105  if (size(array, 1) .lt. loc_size_x1) then
106  print *, '#sll_o_2d_parallel_array_initializer_cartesian error: ', &
107  '#first dimension of passed array is inconsistent with ', &
108  '#the size contained in the passed layout.'
109  end if
110 
111  if (size(array, 2) .lt. loc_size_x2) then
112  print *, '#sll_o_2d_parallel_array_initializer_cartesian error: ', &
113  '#second dimension of passed array is inconsistent with ', &
114  '#the size contained in the passed layout.'
115  end if
116 
117  eta1_min = mesh2d%eta1_min
118  eta2_min = mesh2d%eta2_min
119  delta1 = mesh2d%delta_eta1
120  delta2 = mesh2d%delta_eta2
121 
122  ! This initializes a node-centered array. The loop should be repeated
123  ! below if cell-centered or if arbitrary positions are specified.
124 
125  do j = 1, loc_size_x2
126  do i = 1, loc_size_x1
127  gi(:) = sll_o_local_to_global(layout, (/i, j/))
128  eta1 = eta1_min + real(gi(1) - 1, f64)*delta1
129  eta2 = eta2_min + real(gi(2) - 1, f64)*delta2
130  array(i, j) = func(eta1, eta2, func_params)
131  end do
132  end do
133 
135 
137  layout, &
138  x1_array, &
139  x2_array, &
140  array, &
141  func, &
142  func_params)
143 
144  type(sll_t_layout_2d), pointer :: layout
145  sll_real64, dimension(:), intent(in) :: x1_array
146  sll_real64, dimension(:), intent(in) :: x2_array
147  sll_real64, dimension(:, :), intent(out) :: array
148  procedure(sll_i_scalar_initializer_2d), pointer :: func
149  sll_real64, dimension(:), intent(in) :: func_params
150 
151  sll_int32 :: i
152  sll_int32 :: j
153  sll_int32 :: loc_size_x1
154  sll_int32 :: loc_size_x2
155  sll_real64 :: eta1
156  sll_real64 :: eta2
157  sll_int32, dimension(1:2) :: gi ! global indices in the distributed array
158 
159  if (.not. associated(layout)) then
160  print *, '#sll_o_2d_parallel_array_initializer_cartesian error: ', &
161  '#passed layout is uninitialized.'
162  end if
163 
164 ! if( .not. associated(mesh2d) ) then
165 ! print *, '#sll_o_2d_parallel_array_initializer_cartesian error: ', &
166 ! '#passed mesh2d_eta1_eta2 argument is uninitialized.'
167 ! end if
168 
169  call sll_o_compute_local_sizes(layout, loc_size_x1, loc_size_x2)
170 
171  if (size(array, 1) .lt. loc_size_x1) then
172  print *, '#sll_o_2d_parallel_array_initializer_cartesian error: ', &
173  '#first dimension of passed array is inconsistent with ', &
174  '#the size contained in the passed layout.'
175  end if
176 
177  if (size(array, 2) .lt. loc_size_x2) then
178  print *, '#sll_o_2d_parallel_array_initializer_cartesian error: ', &
179  '#second dimension of passed array is inconsistent with ', &
180  '#the size contained in the passed layout.'
181  end if
182 
183  do j = 1, loc_size_x2
184  do i = 1, loc_size_x1
185  gi(:) = sll_o_local_to_global(layout, (/i, j/))
186  eta1 = x1_array(gi(1))
187  eta2 = x2_array(gi(2))
188  array(i, j) = func(eta1, eta2, func_params)
189  end do
190  end do
191 
193 
195  layout, &
196  eta1_min, &
197  eta2_min, &
198  eta3_min, &
199  eta4_min, &
200  delta1, &
201  delta2, &
202  delta3, &
203  delta4, &
204  array, &
205  func, &
206  func_params)
207  type(sll_t_layout_4d), pointer :: layout
208  sll_real64, intent(in) :: eta1_min
209  sll_real64, intent(in) :: eta2_min
210  sll_real64, intent(in) :: eta3_min
211  sll_real64, intent(in) :: eta4_min
212  sll_real64, intent(in) :: delta1
213  sll_real64, intent(in) :: delta2
214  sll_real64, intent(in) :: delta3
215  sll_real64, intent(in) :: delta4
216  sll_real64, dimension(:, :, :, :), intent(out) :: array
217  procedure(sll_i_scalar_initializer_4d), pointer :: func
218  sll_real64, dimension(:), intent(in) :: func_params
219  sll_int32 :: i
220  sll_int32 :: j
221  sll_int32 :: k
222  sll_int32 :: l
223  sll_int32 :: loc_size_x1
224  sll_int32 :: loc_size_x2
225  sll_int32 :: loc_size_x3
226  sll_int32 :: loc_size_x4
227  sll_real64 :: eta1
228  sll_real64 :: eta2
229  sll_real64 :: eta3
230  sll_real64 :: eta4
231  !sll_real64 :: x1
232  !sll_real64 :: x2
233  !sll_real64 :: x3
234  !sll_real64 :: x4
235  sll_int32, dimension(1:4) :: gi ! global indices in the distributed array
236 
237  if (.not. associated(layout)) then
238  print *, 'sll_o_4d_parallel_array_initializer error: ', &
239  'passed layout is uninitialized.'
240  end if
241 
242  call sll_o_compute_local_sizes(layout, loc_size_x1, loc_size_x2, loc_size_x3, &
243  loc_size_x4)
244 
245  if (size(array, 1) .lt. loc_size_x1) then
246  print *, 'sll_o_4d_parallel_array_initializer error: ', &
247  'first dimension of passed array is inconsistent with ', &
248  'the size contained in the passed layout.'
249  end if
250 
251  if (size(array, 2) .lt. loc_size_x2) then
252  print *, 'sll_o_4d_parallel_array_initializer error: ', &
253  'second dimension of passed array is inconsistent with ', &
254  'the size contained in the passed layout.'
255  end if
256 
257  if (size(array, 3) .lt. loc_size_x3) then
258  print *, 'sll_o_4d_parallel_array_initializer error: ', &
259  'third dimension of passed array is inconsistent with ', &
260  'the size contained in the passed layout.'
261  end if
262 
263  if (size(array, 4) .lt. loc_size_x4) then
264  print *, 'sll_o_4d_parallel_array_initializer error: ', &
265  'fourth dimension of passed array is inconsistent with ', &
266  'the size contained in the passed layout.'
267  end if
268 
269  do l = 1, loc_size_x4
270  do k = 1, loc_size_x3
271  do j = 1, loc_size_x2
272  do i = 1, loc_size_x1
273  gi(:) = sll_o_local_to_global(layout, (/i, j, k, l/))
274  eta1 = eta1_min + real(gi(1) - 1, f64)*delta1
275  eta2 = eta2_min + real(gi(2) - 1, f64)*delta2
276  eta3 = eta3_min + real(gi(3) - 1, f64)*delta3
277  eta4 = eta4_min + real(gi(4) - 1, f64)*delta4
278  array(i, j, k, l) = func(eta1, eta2, eta3, eta4, func_params)
279  end do
280  end do
281  end do
282  end do
283 
285 
287  layout, &
288  mesh1d_eta1, &
289  mesh1d_eta2, &
290  mesh1d_eta3, &
291  mesh1d_eta4, &
292  array, &
293  func, &
294  func_params)
295  type(sll_t_layout_4d), pointer :: layout
296  type(sll_t_cartesian_mesh_1d), pointer :: mesh1d_eta1
297  type(sll_t_cartesian_mesh_1d), pointer :: mesh1d_eta2
298  type(sll_t_cartesian_mesh_1d), pointer :: mesh1d_eta3
299  type(sll_t_cartesian_mesh_1d), pointer :: mesh1d_eta4
300  sll_real64, dimension(:, :, :, :), intent(out) :: array
301  procedure(sll_i_scalar_initializer_4d), pointer :: func
302  sll_real64, dimension(:), intent(in) :: func_params
303  sll_real64 :: eta1_min
304  sll_real64 :: eta2_min
305  sll_real64 :: eta3_min
306  sll_real64 :: eta4_min
307  sll_real64 :: delta1
308  sll_real64 :: delta2
309  sll_real64 :: delta3
310  sll_real64 :: delta4
311 
312  eta1_min = mesh1d_eta1%eta_min
313  delta1 = mesh1d_eta1%delta_eta
314  eta2_min = mesh1d_eta2%eta_min
315  delta2 = mesh1d_eta2%delta_eta
316  eta3_min = mesh1d_eta3%eta_min
317  delta3 = mesh1d_eta3%delta_eta
318  eta4_min = mesh1d_eta4%eta_min
319  delta4 = mesh1d_eta4%delta_eta
320 
322  layout, &
323  eta1_min, &
324  eta2_min, &
325  eta3_min, &
326  eta4_min, &
327  delta1, &
328  delta2, &
329  delta3, &
330  delta4, &
331  array, &
332  func, &
333  func_params)
334 
336 
338  layout, &
339  mesh4d, &
340  array, &
341  func, &
342  func_params)
343 
344  type(sll_t_layout_4d), pointer :: layout
345  type(sll_t_cartesian_mesh_4d), pointer :: mesh4d
346  sll_real64, dimension(:, :, :, :), intent(out) :: array
347  procedure(sll_i_scalar_initializer_4d), pointer :: func
348  sll_real64, dimension(:), intent(in) :: func_params
349 
350  sll_real64 :: delta1
351  sll_real64 :: delta2
352  sll_real64 :: delta3
353  sll_real64 :: delta4
354  sll_real64 :: eta1_min
355  sll_real64 :: eta2_min
356  sll_real64 :: eta3_min
357  sll_real64 :: eta4_min
358 
359  eta1_min = mesh4d%eta1_min
360  eta2_min = mesh4d%eta2_min
361  eta3_min = mesh4d%eta3_min
362  eta4_min = mesh4d%eta4_min
363  delta1 = mesh4d%delta_eta1
364  delta2 = mesh4d%delta_eta2
365  delta3 = mesh4d%delta_eta3
366  delta4 = mesh4d%delta_eta4
367 
369  layout, &
370  eta1_min, &
371  eta2_min, &
372  eta3_min, &
373  eta4_min, &
374  delta1, &
375  delta2, &
376  delta3, &
377  delta4, &
378  array, &
379  func, &
380  func_params)
381 
383 
384  ! it is convenient for example, to separate a 4D mesh into two parts, one can remain
385  ! cartesian while the other may be transformed.
386 
388  layout, &
389  mesh2d_eta1_eta2, &
390  mesh2d_eta3_eta4, &
391  array, &
392  func, &
393  func_params, &
394  transf_x1_x2, &
395  transf_x3_x4)
396 
397  type(sll_t_layout_4d), pointer :: layout
398  type(sll_t_cartesian_mesh_2d), pointer :: mesh2d_eta1_eta2
399  type(sll_t_cartesian_mesh_2d), pointer :: mesh2d_eta3_eta4
400  sll_real64, dimension(:, :, :, :), intent(out) :: array
401  procedure(sll_i_scalar_initializer_4d), pointer :: func
402  sll_real64, dimension(:), intent(in) :: func_params
403 
404  class(sll_c_coordinate_transformation_2d_base), pointer, optional :: &
405  transf_x1_x2
406 
407  class(sll_c_coordinate_transformation_2d_base), pointer, optional :: &
408  transf_x3_x4
409 
410  sll_int32 :: i
411  sll_int32 :: j
412  sll_int32 :: k
413  sll_int32 :: l
414  sll_int32 :: loc_size_x1
415  sll_int32 :: loc_size_x2
416  sll_int32 :: loc_size_x3
417  sll_int32 :: loc_size_x4
418  sll_real64 :: delta1
419  sll_real64 :: delta2
420  sll_real64 :: delta3
421  sll_real64 :: delta4
422  sll_real64 :: eta1_min
423  sll_real64 :: eta2_min
424  sll_real64 :: eta3_min
425  sll_real64 :: eta4_min
426  sll_real64 :: eta1
427  sll_real64 :: eta2
428  sll_real64 :: eta3
429  sll_real64 :: eta4
430  sll_real64 :: x1
431  sll_real64 :: x2
432  sll_real64 :: x3
433  sll_real64 :: x4
434  sll_int32 :: case_selector
435  sll_int32, dimension(1:4) :: gi ! global indices in the distributed array
436 
437  if (.not. associated(layout)) then
438  print *, 'sll_o_4d_parallel_array_initializer error: ', &
439  'passed layout is uninitialized.'
440  end if
441 
442  if (.not. associated(mesh2d_eta1_eta2)) then
443  print *, 'sll_o_4d_parallel_array_initializer error: ', &
444  'passed mesh2d_eta1_eta2 argument is uninitialized.'
445  end if
446 
447  if (.not. associated(mesh2d_eta3_eta4)) then
448  print *, 'sll_o_4d_parallel_array_initializer error: ', &
449  'passed mesh2d_eta3_eta4 argument is uninitialized.'
450  end if
451 
452  call sll_o_compute_local_sizes(layout, loc_size_x1, loc_size_x2, loc_size_x3, &
453  loc_size_x4)
454 
455  if (size(array, 1) .lt. loc_size_x1) then
456  print *, 'sll_o_4d_parallel_array_initializer error: ', &
457  'first dimension of passed array is inconsistent with ', &
458  'the size contained in the passed layout.'
459  end if
460 
461  if (size(array, 2) .lt. loc_size_x2) then
462  print *, 'sll_o_4d_parallel_array_initializer error: ', &
463  'second dimension of passed array is inconsistent with ', &
464  'the size contained in the passed layout.'
465  end if
466 
467  if (size(array, 3) .lt. loc_size_x3) then
468  print *, 'sll_o_4d_parallel_array_initializer error: ', &
469  'third dimension of passed array is inconsistent with ', &
470  'the size contained in the passed layout.'
471  end if
472 
473  if (size(array, 4) .lt. loc_size_x4) then
474  print *, 'sll_o_4d_parallel_array_initializer error: ', &
475  'fourth dimension of passed array is inconsistent with ', &
476  'the size contained in the passed layout.'
477  end if
478 
479  ! if(present( transf_x1_x2 ) ) then
480  ! if(.not. associated(transf_x1_x2%get_cartesian_mesh(),mesh2d_eta1_eta2) )&
481  ! then
482 
483  ! print *, 'sll_o_4d_parallel_array_initializer warning: ', &
484  ! 'the mesh associated to the transf_x1_x2 transformation ', &
485  ! 'is not the same as the mesh2d_eta1_eta2 logical mesh. ', &
486  ! 'Unless the parameters of these meshes are the same, ', &
487  ! 'bad things will happen.'
488  ! end if
489  ! end if
490 
491  ! if(present( transf_x3_x4 ) ) then
492  ! if(.not. associated(transf_x3_x4%get_cartesian_mesh(),mesh2d_eta3_eta4) )&
493  ! then
494 
495  ! print *, 'sll_o_4d_parallel_array_initializer warning: ', &
496  ! 'the mesh associated to the transf_x3_x4 transformation ', &
497  ! 'is not the same as the mesh2d_eta3_eta4 logical mesh. ', &
498  ! 'Unless the parameters of these meshes are the same, ', &
499  ! 'bad things will happen.'
500  ! end if
501  ! end if
502 
503  case_selector = 0
504 
505  if (present(transf_x1_x2)) then
506  case_selector = case_selector + 1
507  end if
508 
509  if (present(transf_x3_x4)) then
510  case_selector = case_selector + 2
511  end if
512  eta1_min = mesh2d_eta1_eta2%eta1_min
513  eta2_min = mesh2d_eta1_eta2%eta2_min
514  eta3_min = mesh2d_eta3_eta4%eta1_min
515  eta4_min = mesh2d_eta3_eta4%eta2_min
516  delta1 = mesh2d_eta1_eta2%delta_eta1
517  delta2 = mesh2d_eta1_eta2%delta_eta2
518  delta3 = mesh2d_eta3_eta4%delta_eta1
519  delta4 = mesh2d_eta3_eta4%delta_eta2
520 
521  ! This initializes a node-centered array. The loop should be repeated
522  ! below if cell-centered or if arbitrary positions are specified.
523 
524  select case (case_selector)
525 
526  case (0) ! none of the transformations was provided
527  do l = 1, loc_size_x4
528  do k = 1, loc_size_x3
529  do j = 1, loc_size_x2
530  do i = 1, loc_size_x1
531  gi(:) = sll_o_local_to_global(layout, (/i, j, k, l/))
532  eta1 = eta1_min + real(gi(1) - 1, f64)*delta1
533  eta2 = eta2_min + real(gi(2) - 1, f64)*delta2
534  eta3 = eta3_min + real(gi(3) - 1, f64)*delta3
535  eta4 = eta4_min + real(gi(4) - 1, f64)*delta4
536  array(i, j, k, l) = func(eta1, eta2, eta3, eta4, func_params)
537  end do
538  end do
539  end do
540  end do
541  case (1) ! Only the x1,x2 transfomation was provided.
542  do l = 1, loc_size_x4
543  do k = 1, loc_size_x3
544  do j = 1, loc_size_x2
545  do i = 1, loc_size_x1
546  gi(:) = sll_o_local_to_global(layout, (/i, j, k, l/))
547  eta1 = eta1_min + real(gi(1) - 1, f64)*delta1
548  eta2 = eta2_min + real(gi(2) - 1, f64)*delta2
549  eta3 = eta3_min + real(gi(3) - 1, f64)*delta3
550  eta4 = eta4_min + real(gi(4) - 1, f64)*delta4
551  x1 = transf_x1_x2%x1(eta1, eta2)
552  x2 = transf_x1_x2%x2(eta1, eta2)
553  array(i, j, k, l) = func(x1, x2, eta3, eta4, func_params)
554  end do
555  end do
556  end do
557  end do
558 
559  case (2) ! Only the x3,x4 transformation was provided
560  do l = 1, loc_size_x4
561  do k = 1, loc_size_x3
562  do j = 1, loc_size_x2
563  do i = 1, loc_size_x1
564  gi(:) = sll_o_local_to_global(layout, (/i, j, k, l/))
565  eta1 = eta1_min + real(gi(1) - 1, f64)*delta1
566  eta2 = eta2_min + real(gi(2) - 1, f64)*delta2
567  eta3 = eta3_min + real(gi(3) - 1, f64)*delta3
568  eta4 = eta4_min + real(gi(4) - 1, f64)*delta4
569  x3 = transf_x3_x4%x1(eta3, eta4)
570  x4 = transf_x3_x4%x2(eta3, eta4)
571  array(i, j, k, l) = func(eta1, eta2, x3, x4, func_params)
572  end do
573  end do
574  end do
575  end do
576  case (3) ! Both transformations were provided.
577  do l = 1, loc_size_x4
578  do k = 1, loc_size_x3
579  do j = 1, loc_size_x2
580  do i = 1, loc_size_x1
581  gi(:) = sll_o_local_to_global(layout, (/i, j, k, l/))
582  eta1 = eta1_min + real(gi(1) - 1, f64)*delta1
583  eta2 = eta2_min + real(gi(2) - 1, f64)*delta2
584  eta3 = eta3_min + real(gi(3) - 1, f64)*delta3
585  eta4 = eta4_min + real(gi(4) - 1, f64)*delta4
586  x1 = transf_x1_x2%x1(eta1, eta2)
587  x2 = transf_x1_x2%x2(eta1, eta2)
588  x3 = transf_x3_x4%x1(eta3, eta4)
589  x4 = transf_x3_x4%x2(eta3, eta4)
590  array(i, j, k, l) = func(x1, x2, x3, x4, func_params)
591  end do
592  end do
593  end do
594  end do
595  end select
596 
598 
599 ! subroutine sll_o_4d_parallel_array_initializer( &
600 ! layout, &
601 ! mesh2d_eta1_eta2, &
602 ! mesh2d_eta3_eta4, &
603 ! array, &
604 ! func, &
605 ! func_params, &
606 ! transf_x1_x2, &
607 ! transf_x3_x4 )
608 !
609 ! type(sll_t_layout_4d), pointer :: layout
610 ! type(sll_t_cartesian_mesh_2d), pointer :: mesh2d_eta1_eta2
611 ! type(sll_t_cartesian_mesh_2d), pointer :: mesh2d_eta3_eta4
612 ! sll_real64, dimension(:,:,:,:), intent(out) :: array
613 ! procedure(sll_i_scalar_initializer_4d) :: func
614 ! sll_real64, dimension(:), optional :: func_params
615 !
616 ! class(sll_c_coordinate_transformation_2d_base), pointer, optional :: &
617 ! transf_x1_x2
618 !
619 ! class(sll_c_coordinate_transformation_2d_base), pointer, optional :: &
620 ! transf_x3_x4
621 !
622 ! sll_int32 :: i
623 ! sll_int32 :: j
624 ! sll_int32 :: k
625 ! sll_int32 :: l
626 ! sll_int32 :: loc_size_x1
627 ! sll_int32 :: loc_size_x2
628 ! sll_int32 :: loc_size_x3
629 ! sll_int32 :: loc_size_x4
630 ! sll_real64 :: delta1
631 ! sll_real64 :: delta2
632 ! sll_real64 :: delta3
633 ! sll_real64 :: delta4
634 ! sll_real64 :: eta1_min
635 ! sll_real64 :: eta2_min
636 ! sll_real64 :: eta3_min
637 ! sll_real64 :: eta4_min
638 ! sll_real64 :: eta1
639 ! sll_real64 :: eta2
640 ! sll_real64 :: eta3
641 ! sll_real64 :: eta4
642 ! sll_real64 :: x1
643 ! sll_real64 :: x2
644 ! sll_real64 :: x3
645 ! sll_real64 :: x4
646 ! sll_int32 :: case_selector
647 ! sll_int32, dimension(1:4) :: gi ! global indices in the distributed array
648 !
649 ! if( .not. associated(layout) ) then
650 ! print *, 'sll_o_4d_parallel_array_initializer error: ', &
651 ! 'passed layout is uninitialized.'
652 ! end if
653 !
654 ! if( .not. associated(mesh2d_eta1_eta2) ) then
655 ! print *, 'sll_m_parallel_array_initializer.F90_array_initializer error: ', &
656 ! 'passed mesh2d_eta1_eta2 argument is uninitialized.'
657 ! end if
658 !
659 ! if( .not. associated(mesh2d_eta3_eta4) ) then
660 ! print *, 'sll_o_4d_parallel_array_initializer error: ', &
661 ! 'passed mesh2d_eta3_eta4 argument is uninitialized.'
662 ! end if
663 !
664 ! call sll_o_compute_local_sizes( layout, loc_size_x1, loc_size_x2, loc_size_x3, &
665 ! loc_size_x4 )
666 !
667 ! if( size(array,1) .lt. loc_size_x1 ) then
668 ! print *, 'sll_o_4d_parallel_array_initializer error: ', &
669 ! 'first dimension of passed array is inconsistent with ', &
670 ! 'the size contained in the passed layout.'
671 ! end if
672 !
673 ! if( size(array,2) .lt. loc_size_x2 ) then
674 ! print *, 'sll_o_4d_parallel_array_initializer error: ', &
675 ! 'second dimension of passed array is inconsistent with ', &
676 ! 'the size contained in the passed layout.'
677 ! end if
678 !
679 ! if( size(array,3) .lt. loc_size_x3 ) then
680 ! print *, 'sll_o_4d_parallel_array_initializer error: ', &
681 ! 'third dimension of passed array is inconsistent with ', &
682 ! 'the size contained in the passed layout.'
683 ! end if
684 !
685 ! if( size(array,4) .lt. loc_size_x4 ) then
686 ! print *, 'sll_o_4d_parallel_array_initializer error: ', &
687 ! 'fourth dimension of passed array is inconsistent with ', &
688 ! 'the size contained in the passed layout.'
689 ! end if
690 !
691 ! if(present( transf_x1_x2 ) ) then
692 ! if( .not. associated(transf_x1_x2%get_cartesian_mesh(), mesh2d_eta1_eta2) ) then
693 ! print *, 'sll_o_4d_parallel_array_initializer warning: ', &
694 ! 'the mesh associated to the transf_x1_x2 transformation ', &
695 ! 'is not the same as the mesh2d_eta1_eta2 logical mesh. ', &
696 ! 'Unless the parameters of these meshes are the same, ', &
697 ! 'bad things will happen.'
698 ! end if
699 ! end if
700 !
701 ! if(present( transf_x3_x4 ) ) then
702 ! if( .not. associated(transf_x3_x4%get_cartesian_mesh(), mesh2d_eta3_eta4) ) then
703 ! print *, 'sll_o_4d_parallel_array_initializer warning: ', &
704 ! 'the mesh associated to the transf_x3_x4 transformation ', &
705 ! 'is not the same as the mesh2d_eta3_eta4 logical mesh. ', &
706 ! 'Unless the parameters of these meshes are the same, ', &
707 ! 'bad things will happen.'
708 ! end if
709 ! end if
710 !
711 ! case_selector = 0
712 !
713 ! if( present(transf_x1_x2) ) then
714 ! case_selector = case_selector + 1
715 ! end if
716 !
717 ! if( present(transf_x3_x4) ) then
718 ! case_selector = case_selector + 2
719 ! end if
720 ! eta1_min = mesh2d_eta1_eta2%eta1_min
721 ! eta2_min = mesh2d_eta1_eta2%eta2_min
722 ! eta3_min = mesh2d_eta3_eta4%eta1_min
723 ! eta4_min = mesh2d_eta3_eta4%eta2_min
724 ! delta1 = mesh2d_eta1_eta2%delta_eta1
725 ! delta2 = mesh2d_eta1_eta2%delta_eta2
726 ! delta3 = mesh2d_eta3_eta4%delta_eta1
727 ! delta4 = mesh2d_eta3_eta4%delta_eta2
728 !
729 ! ! This initializes a node-centered array. The loop should be repeated
730 ! ! below if cell-centered or if arbitrary positions are specified.
731 !
732 ! select case (case_selector)
733 !
734 ! case (0) ! none of the transformations was provided
735 ! do l=1,loc_size_x4
736 ! do k=1,loc_size_x3
737 ! do j=1,loc_size_x2
738 ! do i=1,loc_size_x1
739 ! gi(:) = sll_o_local_to_global( layout, (/i,j,k,l/) )
740 ! eta1 = eta1_min + real(gi(1)-1,f64)*delta1
741 ! eta2 = eta2_min + real(gi(2)-1,f64)*delta2
742 ! eta3 = eta3_min + real(gi(3)-1,f64)*delta3
743 ! eta4 = eta4_min + real(gi(4)-1,f64)*delta4
744 ! array(i,j,k,l) = func(eta1,eta2,eta3,eta4,func_params)
745 ! end do
746 ! end do
747 ! end do
748 ! end do
749 ! case(1) ! Only the x1,x2 transfomation was provided.
750 ! do l=1,loc_size_x4
751 ! do k=1,loc_size_x3
752 ! do j=1,loc_size_x2
753 ! do i=1,loc_size_x1
754 ! gi(:) = sll_o_local_to_global( layout, (/i,j,k,l/) )
755 ! eta1 = eta1_min + real(gi(1)-1,f64)*delta1
756 ! eta2 = eta2_min + real(gi(2)-1,f64)*delta2
757 ! eta3 = eta3_min + real(gi(3)-1,f64)*delta3
758 ! eta4 = eta4_min + real(gi(4)-1,f64)*delta4
759 ! x1 = transf_x1_x2%x1(eta1,eta2)
760 ! x2 = transf_x1_x2%x2(eta1,eta2)
761 ! array(i,j,k,l) = func(x1,x2,eta3,eta4,func_params)
762 ! end do
763 ! end do
764 ! end do
765 ! end do
766 ! case(2) ! Only the x3,x4 transformation was provided
767 ! do l=1,loc_size_x4
768 ! do k=1,loc_size_x3
769 ! do j=1,loc_size_x2
770 ! do i=1,loc_size_x1
771 ! gi(:) = sll_o_local_to_global( layout, (/i,j,k,l/) )
772 ! eta1 = eta1_min + real(gi(1)-1,f64)*delta1
773 ! eta2 = eta2_min + real(gi(2)-1,f64)*delta2
774 ! eta3 = eta3_min + real(gi(3)-1,f64)*delta3
775 ! eta4 = eta4_min + real(gi(4)-1,f64)*delta4
776 ! x3 = transf_x3_x4%x1(eta3,eta4)
777 ! x4 = transf_x3_x4%x2(eta3,eta4)
778 ! array(i,j,k,l) = func(eta1,eta2,x3,x4,func_params)
779 ! end do
780 ! end do
781 ! end do
782 ! end do
783 ! case(3) ! Both transformations were provided.
784 ! do l=1,loc_size_x4
785 ! do k=1,loc_size_x3
786 ! do j=1,loc_size_x2
787 ! do i=1,loc_size_x1
788 ! gi(:) = sll_o_local_to_global( layout, (/i,j,k,l/) )
789 ! eta1 = eta1_min + real(gi(1)-1,f64)*delta1
790 ! eta2 = eta2_min + real(gi(2)-1,f64)*delta2
791 ! eta3 = eta3_min + real(gi(3)-1,f64)*delta3
792 ! eta4 = eta4_min + real(gi(4)-1,f64)*delta4
793 ! x1 = transf_x1_x2%x1(eta1,eta2)
794 ! x2 = transf_x1_x2%x2(eta1,eta2)
795 ! x3 = transf_x3_x4%x1(eta3,eta4)
796 ! x4 = transf_x3_x4%x2(eta3,eta4)
797 ! array(i,j,k,l) = func(x1,x2,x3,x4,func_params)
798 ! end do
799 ! end do
800 ! end do
801 ! end do
802 ! end select
803 !
804 ! end subroutine sll_o_4d_parallel_array_initializer
805 
806 !in case of interpolation in the vx, vy, x, y direction
808  layout, &
809  mesh2d_eta1_eta2, &
810  mesh2d_eta3_eta4, &
811  array, &
812  func, &
813  func_params, &
814  transf_x1_x2, &
815  transf_x3_x4, &
816  ! in case of a mesh with cell subdivisions
817  subcells1, &
818  subcells2, &
819  subcells3, &
820  subcells4)
821 
822  type(sll_t_layout_4d), pointer :: layout
823  type(sll_t_cartesian_mesh_2d), pointer :: mesh2d_eta1_eta2
824  type(sll_t_cartesian_mesh_2d), pointer :: mesh2d_eta3_eta4
825  sll_real64, dimension(:, :, :, :), intent(out) :: array
826  procedure(sll_i_scalar_initializer_4d) :: func
827  sll_real64, dimension(:), intent(in) :: func_params
828  sll_int32, optional :: subcells1
829  sll_int32, optional :: subcells2
830  sll_int32, optional :: subcells3
831  sll_int32, optional :: subcells4
832 
833  class(sll_c_coordinate_transformation_2d_base), pointer, optional :: &
834  transf_x1_x2
835 
836  class(sll_c_coordinate_transformation_2d_base), pointer, optional :: &
837  transf_x3_x4
838 
839  sll_int32 :: i
840  sll_int32 :: j
841  sll_int32 :: k
842  sll_int32 :: l
843  sll_int32 :: loc_size_x1
844  sll_int32 :: loc_size_x2
845  sll_int32 :: loc_size_x3
846  sll_int32 :: loc_size_x4
847  sll_real64 :: delta1
848  sll_real64 :: delta2
849  sll_real64 :: delta3
850  sll_real64 :: delta4
851  sll_real64 :: eta1_min
852  sll_real64 :: eta2_min
853  sll_real64 :: eta3_min
854  sll_real64 :: eta4_min
855  sll_real64 :: eta1
856  sll_real64 :: eta2
857  sll_real64 :: eta3
858  sll_real64 :: eta4
859  sll_real64 :: x1
860  sll_real64 :: x2
861  sll_real64 :: x3
862  sll_real64 :: x4
863  sll_int32 :: case_selector
864  sll_int32 :: sub1, sub2, sub3, sub4
865  sll_int32, dimension(1:4) :: gi ! global indices in the distributed array
866 
867  if (.not. associated(layout)) then
868  print *, 'sll_o_4d_parallel_array_initializer error: ', &
869  'passed layout is uninitialized.'
870  end if
871 
872  if (.not. associated(mesh2d_eta1_eta2)) then
873  print *, 'sll_o_4d_parallel_array_initializer error: ', &
874  'passed mesh2d_eta1_eta2 argument is uninitialized.'
875  end if
876 
877  if (.not. associated(mesh2d_eta3_eta4)) then
878  print *, 'sll_o_4d_parallel_array_initializer error: ', &
879  'passed mesh2d_eta3_eta4 argument is uninitialized.'
880  end if
881 
882  call sll_o_compute_local_sizes(layout, loc_size_x1, loc_size_x2, loc_size_x3, &
883  loc_size_x4)
884 
885  if (size(array, 1) .lt. loc_size_x1) then
886  print *, 'sll_o_4d_parallel_array_initializer error: ', &
887  'first dimension of passed array is inconsistent with ', &
888  'the size contained in the passed layout.'
889  end if
890 
891  if (size(array, 2) .lt. loc_size_x2) then
892  print *, 'sll_o_4d_parallel_array_initializer error: ', &
893  'second dimension of passed array is inconsistent with ', &
894  'the size contained in the passed layout.'
895  end if
896 
897  if (size(array, 3) .lt. loc_size_x3) then
898  print *, 'sll_o_4d_parallel_array_initializer error: ', &
899  'third dimension of passed array is inconsistent with ', &
900  'the size contained in the passed layout.'
901  end if
902 
903  if (size(array, 4) .lt. loc_size_x4) then
904  print *, 'sll_o_4d_parallel_array_initializer error: ', &
905  'fourth dimension of passed array is inconsistent with ', &
906  'the size contained in the passed layout.'
907  end if
908 
909  ! if(present( transf_x1_x2 ) ) then
910  ! if(.not. associated(transf_x1_x2%get_cartesian_mesh(), mesh2d_eta1_eta2) )&
911  ! then
912 
913  ! print *, 'sll_o_4d_parallel_array_initializer warning: ', &
914  ! 'the mesh associated to the transf_x1_x2 transformation ', &
915  ! 'is not the same as the mesh2d_eta1_eta2 logical mesh. ', &
916  ! 'Unless the parameters of these meshes are the same, ', &
917  ! 'bad things will happen.'
918  ! end if
919  ! end if
920 
921  ! if(present( transf_x3_x4 ) ) then
922  ! if(.not. associated(transf_x3_x4%get_cartesian_mesh(),mesh2d_eta3_eta4) )&
923  ! then
924 
925  ! print *, 'sll_o_4d_parallel_array_initializer warning: ', &
926  ! 'the mesh associated to the transf_x3_x4 transformation ', &
927  ! 'is not the same as the mesh2d_eta3_eta4 logical mesh. ', &
928  ! 'Unless the parameters of these meshes are the same, ', &
929  ! 'bad things will happen.'
930  ! end if
931  ! end if
932 
933  if (.not. present(subcells1)) then
934  sub1 = 1
935  else
936  sub1 = subcells1
937  end if
938 
939  if (.not. present(subcells2)) then
940  sub2 = 1
941  else
942  sub2 = subcells2
943  end if
944 
945  if (.not. present(subcells3)) then
946  sub3 = 1
947  else
948  sub3 = subcells3
949  end if
950 
951  if (.not. present(subcells3)) then
952  sub4 = 1
953  else
954  sub4 = subcells4
955  end if
956 
957  case_selector = 0
958 
959  if (present(transf_x1_x2)) then
960  case_selector = case_selector + 1
961  end if
962 
963  if (present(transf_x3_x4)) then
964  case_selector = case_selector + 2
965  end if
966 
967  eta1_min = mesh2d_eta1_eta2%eta1_min
968  eta2_min = mesh2d_eta1_eta2%eta2_min
969  eta3_min = mesh2d_eta3_eta4%eta1_min
970  eta4_min = mesh2d_eta3_eta4%eta2_min
971  delta1 = mesh2d_eta1_eta2%delta_eta1/sub1
972  delta2 = mesh2d_eta1_eta2%delta_eta2/sub2
973  delta3 = mesh2d_eta3_eta4%delta_eta1/sub3
974  delta4 = mesh2d_eta3_eta4%delta_eta2/sub4
975 
976  ! This initializes a node-centered array. The loop should be repeated
977  ! below if cell-centered or if arbitrary positions are specified.
978 
979  select case (case_selector)
980 
981  case (0) ! none of the transformations was provided
982  do l = 1, loc_size_x4
983  do k = 1, loc_size_x3
984  do j = 1, loc_size_x2
985  do i = 1, loc_size_x1
986  gi(:) = sll_o_local_to_global(layout, (/i, j, k, l/))
987  eta1 = eta1_min + real(gi(1) - 1, f64)*delta1
988  eta2 = eta2_min + real(gi(2) - 1, f64)*delta2
989  eta3 = eta3_min + real(gi(3) - 1, f64)*delta3
990  eta4 = eta4_min + real(gi(4) - 1, f64)*delta4
991  array(i, j, k, l) = func(eta1, eta2, eta3, eta4, func_params)
992  end do
993  end do
994  end do
995  end do
996  case (1) ! Only the x1,x2 transfomation was provided.
997  do l = 1, loc_size_x4
998  do k = 1, loc_size_x3
999  do j = 1, loc_size_x2
1000  do i = 1, loc_size_x1
1001  gi(:) = sll_o_local_to_global(layout, (/i, j, k, l/))
1002  eta1 = eta1_min + real(gi(1) - 1, f64)*delta1
1003  eta2 = eta2_min + real(gi(2) - 1, f64)*delta2
1004  eta3 = eta3_min + real(gi(3) - 1, f64)*delta3
1005  eta4 = eta4_min + real(gi(4) - 1, f64)*delta4
1006  x1 = transf_x1_x2%x1(eta1, eta2)
1007  x2 = transf_x1_x2%x2(eta1, eta2)
1008  array(i, j, k, l) = func(x1, x2, eta3, eta4, func_params)
1009  end do
1010  end do
1011  end do
1012  end do
1013  case (2) ! Only the x3,x4 transformation was provided
1014  do l = 1, loc_size_x4
1015  do k = 1, loc_size_x3
1016  do j = 1, loc_size_x2
1017  do i = 1, loc_size_x1
1018  gi(:) = sll_o_local_to_global(layout, (/i, j, k, l/))
1019  eta1 = eta1_min + real(gi(1) - 1, f64)*delta1
1020  eta2 = eta2_min + real(gi(2) - 1, f64)*delta2
1021  eta3 = eta3_min + real(gi(3) - 1, f64)*delta3
1022  eta4 = eta4_min + real(gi(4) - 1, f64)*delta4
1023  x3 = transf_x3_x4%x1(eta3, eta4)
1024  x4 = transf_x3_x4%x2(eta3, eta4)
1025  array(i, j, k, l) = func(eta1, eta2, x3, x4, func_params)
1026  end do
1027  end do
1028  end do
1029  end do
1030  case (3) ! Both transformations were provided.
1031  do l = 1, loc_size_x4
1032  do k = 1, loc_size_x3
1033  do j = 1, loc_size_x2
1034  do i = 1, loc_size_x1
1035  gi(:) = sll_o_local_to_global(layout, (/i, j, k, l/))
1036  eta1 = eta1_min + real(gi(1) - 1, f64)*delta1
1037  eta2 = eta2_min + real(gi(2) - 1, f64)*delta2
1038  eta3 = eta3_min + real(gi(3) - 1, f64)*delta3
1039  eta4 = eta4_min + real(gi(4) - 1, f64)*delta4
1040  x1 = transf_x1_x2%x1(eta1, eta2)
1041  x2 = transf_x1_x2%x2(eta1, eta2)
1042  x3 = transf_x3_x4%x1(eta3, eta4)
1043  x4 = transf_x3_x4%x2(eta3, eta4)
1044  array(i, j, k, l) = func(x1, x2, x3, x4, func_params)
1045  end do
1046  end do
1047  end do
1048  end do
1049  end select
1050 
1052 
1053 ! subroutine sll_o_4d_parallel_array_initializer_cartesian( &
1054 ! layout, &
1055 ! mesh4d, &
1056 ! array, &
1057 ! func, &
1058 ! func_params)
1059 !
1060 ! type(sll_t_layout_4d), pointer :: layout
1061 ! type(sll_t_cartesian_mesh_4d), pointer :: mesh4d
1062 ! sll_real64, dimension(:,:,:,:), intent(out) :: array
1063 ! procedure(sll_i_scalar_initializer_4d) :: func
1064 ! sll_real64, dimension(:), optional :: func_params
1065 !
1066 !
1067 ! sll_int32 :: i
1068 ! sll_int32 :: j
1069 ! sll_int32 :: k
1070 ! sll_int32 :: l
1071 ! sll_int32 :: loc_size_x1
1072 ! sll_int32 :: loc_size_x2
1073 ! sll_int32 :: loc_size_x3
1074 ! sll_int32 :: loc_size_x4
1075 ! sll_real64 :: delta1
1076 ! sll_real64 :: delta2
1077 ! sll_real64 :: delta3
1078 ! sll_real64 :: delta4
1079 ! sll_real64 :: eta1_min
1080 ! sll_real64 :: eta2_min
1081 ! sll_real64 :: eta3_min
1082 ! sll_real64 :: eta4_min
1083 ! sll_real64 :: eta1
1084 ! sll_real64 :: eta2
1085 ! sll_real64 :: eta3
1086 ! sll_real64 :: eta4
1087 ! sll_int32, dimension(1:4) :: gi ! global indices in the distributed array
1088 !
1089 ! if( .not. associated(layout) ) then
1090 ! print *, 'sll_o_4d_parallel_array_initializer_cartesian error: ', &
1091 ! 'passed layout is uninitialized.'
1092 ! end if
1093 !
1094 ! if( .not. associated(mesh4d) ) then
1095 ! print *, 'sll_o_4d_parallel_array_initializer_cartesian error: ', &
1096 ! 'passed mesh2d_eta1_eta2 argument is uninitialized.'
1097 ! end if
1098 !
1099 !
1100 ! call sll_o_compute_local_sizes( layout, loc_size_x1, loc_size_x2, loc_size_x3, &
1101 ! loc_size_x4 )
1102 !
1103 ! if( size(array,1) .lt. loc_size_x1 ) then
1104 ! print *, 'sll_o_4d_parallel_array_initializer_cartesian error: ', &
1105 ! 'first dimension of passed array is inconsistent with ', &
1106 ! 'the size contained in the passed layout.'
1107 ! end if
1108 !
1109 ! if( size(array,2) .lt. loc_size_x2 ) then
1110 ! print *, 'sll_o_4d_parallel_array_initializer_cartesian error: ', &
1111 ! 'second dimension of passed array is inconsistent with ', &
1112 ! 'the size contained in the passed layout.'
1113 ! end if
1114 ! if( size(array,3) .lt. loc_size_x3 ) then
1115 ! print *, 'sll_o_4d_parallel_array_initializer_cartesian error: ', &
1116 ! 'third dimension of passed array is inconsistent with ', &
1117 ! 'the size contained in the passed layout.'
1118 ! end if
1119 !
1120 ! if( size(array,4) .lt. loc_size_x4 ) then
1121 ! print *, 'sll_o_4d_parallel_array_initializer_cartesian error: ', &
1122 ! 'fourth dimension of passed array is inconsistent with ', &
1123 ! 'the size contained in the passed layout.'
1124 ! end if
1125 !
1126 !
1127 !
1128 ! eta1_min = mesh4d%eta1_min
1129 ! eta2_min = mesh4d%eta2_min
1130 ! eta3_min = mesh4d%eta3_min
1131 ! eta4_min = mesh4d%eta4_min
1132 ! delta1 = mesh4d%delta_eta1
1133 ! delta2 = mesh4d%delta_eta2
1134 ! delta3 = mesh4d%delta_eta3
1135 ! delta4 = mesh4d%delta_eta4
1136 !
1137 ! ! This initializes a node-centered array. The loop should be repeated
1138 ! ! below if cell-centered or if arbitrary positions are specified.
1139 !
1140 ! do l=1,loc_size_x4
1141 ! do k=1,loc_size_x3
1142 ! do j=1,loc_size_x2
1143 ! do i=1,loc_size_x1
1144 ! gi(:) = sll_o_local_to_global( layout, (/i,j,k,l/) )
1145 ! eta1 = eta1_min + real(gi(1)-1,f64)*delta1
1146 ! eta2 = eta2_min + real(gi(2)-1,f64)*delta2
1147 ! eta3 = eta3_min + real(gi(3)-1,f64)*delta3
1148 ! eta4 = eta4_min + real(gi(4)-1,f64)*delta4
1149 ! array(i,j,k,l) = func(eta1,eta2,eta3,eta4,func_params)
1150 ! end do
1151 ! end do
1152 ! end do
1153 ! end do
1154 !
1155 ! end subroutine sll_o_4d_parallel_array_initializer_cartesian
1156 
Cartesian mesh basic types.
Abstract class for coordinate transformations.
subroutine sll_4d_parallel_array_initializer_cartesian_logical_1d_1d_1d_1d(layout, mesh1d_eta1, mesh1d_eta2, mesh1d_eta3, mesh1d_eta4, array, func, func_params)
subroutine sll_2d_times_2d_parallel_array_initializer(layout, mesh2d_eta1_eta2, mesh2d_eta3_eta4, array, func, func_params, transf_x1_x2, transf_x3_x4)
subroutine, public sll_s_4d_parallel_array_initializer_finite_volume(layout, mesh2d_eta1_eta2, mesh2d_eta3_eta4, array, func, func_params, transf_x1_x2, transf_x3_x4, subcells1, subcells2, subcells3, subcells4)
subroutine sll_2d_parallel_array_initializer_cartesian_logical_2d(layout, mesh2d, array, func, func_params)
subroutine sll_2d_parallel_array_initializer_cartesian_array_1d_1d(layout, x1_array, x2_array, array, func, func_params)
subroutine sll_4d_parallel_array_initializer_cartesian_logical_4d(layout, mesh4d, array, func, func_params)
subroutine sll_4d_parallel_array_initializer_cartesian_aux(layout, eta1_min, eta2_min, eta3_min, eta4_min, delta1, delta2, delta3, delta4, array, func, func_params)
    Report Typos and Errors