Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_cartesian_meshes.F90
Go to the documentation of this file.
1 !**************************************************************
2 ! Copyright INRIA
3 ! Authors :
4 ! CALVI project team
5 !
6 ! This code SeLaLib (for Semi-Lagrangian-Library)
7 ! is a parallel library for simulating the plasma turbulence
8 ! in a tokamak.
9 !
10 ! This software is governed by the CeCILL-B license
11 ! under French law and abiding by the rules of distribution
12 ! of free software. You can use, modify and redistribute
13 ! the software under the terms of the CeCILL-B license as
14 ! circulated by CEA, CNRS and INRIA at the following URL
15 ! "http://www.cecill.info".
16 !**************************************************************
17 
21 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
22 #include "sll_assert.h"
23 #include "sll_memory.h"
24 #include "sll_working_precision.h"
25 #include "sll_errors.h"
26 
27  use sll_m_meshes_base, only: &
31 
32  implicit none
33 
34  public :: &
40  operator(*), &
50  sll_o_cell, &
52  sll_o_delete, &
53  sll_o_display, &
55  sll_o_new, &
57 
58  private
59 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
60 
63  sll_int32 :: num_cells
64  sll_real64 :: eta_min
65  sll_real64 :: eta_max
66  sll_real64 :: delta_eta
67  contains
68  procedure, pass(mesh) :: eta1_node => eta1_node_1d
69  procedure, pass(mesh) :: eta1_cell => eta1_cell_1d
70  procedure, pass(mesh) :: display => display_cartesian_mesh_1d
71  procedure, pass(mesh) :: delete => sll_s_cartesian_mesh_1d_free
72  procedure, pass(mesh) :: eta1_nodes => nodes_cartesian_mesh_1d
73  procedure, pass(mesh) :: num_nodes => num_nodes_cartesian_mesh_1d
74  procedure, pass(mesh) :: cellnum => cell_cartesian_mesh_1d
75  procedure, pass(mesh) :: area => length_cartesian_mesh_1d
76  procedure, pass(mesh) :: cell_margin => cell_margin_cartesian_mesh_1d
77  procedure, pass(mesh) :: period => period_cartesian_mesh_1d
79 
82  sll_int32 :: num_cells1
83  sll_int32 :: num_cells2
84  sll_real64 :: eta1_min
85  sll_real64 :: eta1_max
86  sll_real64 :: eta2_min
87  sll_real64 :: eta2_max
88  sll_real64 :: delta_eta1
89  sll_real64 :: delta_eta2
90  contains
91  procedure, pass(mesh) :: eta1_node => eta1_node_2d
92  procedure, pass(mesh) :: eta2_node => eta2_node_2d
93  procedure, pass(mesh) :: eta1_cell_one_arg => eta1_cell_2d_one_arg
94  procedure, pass(mesh) :: eta1_cell_two_arg => eta1_cell_2d_two_arg
95  procedure, pass(mesh) :: eta2_cell_one_arg => eta2_cell_2d_one_arg
96  procedure, pass(mesh) :: eta2_cell_two_arg => eta2_cell_2d_two_arg
97  procedure, pass(mesh) :: display => display_cartesian_mesh_2d
98  procedure, pass(mesh) :: delete => sll_s_cartesian_mesh_2d_free
100 
103  type(sll_t_cartesian_mesh_2d), pointer :: lm
105 
108  sll_int32 :: num_cells1
109  sll_int32 :: num_cells2
110  sll_int32 :: num_cells3
111  sll_real64 :: eta1_min
112  sll_real64 :: eta1_max
113  sll_real64 :: eta2_min
114  sll_real64 :: eta2_max
115  sll_real64 :: eta3_min
116  sll_real64 :: eta3_max
117  sll_real64 :: delta_eta1
118  sll_real64 :: delta_eta2
119  sll_real64 :: delta_eta3
120  contains
121  procedure, pass(mesh) :: eta1_node => eta1_node_3d
122  procedure, pass(mesh) :: eta2_node => eta2_node_3d
123  procedure, pass(mesh) :: eta3_node => eta3_node_3d
124  procedure, pass(mesh) :: eta1_cell => eta1_cell_3d
125  procedure, pass(mesh) :: eta2_cell => eta2_cell_3d
126  procedure, pass(mesh) :: eta3_cell => eta3_cell_3d
127  procedure, pass(mesh) :: display => display_cartesian_mesh_3d
128  procedure, pass(mesh) :: delete => sll_s_cartesian_mesh_3d_free
129  end type sll_t_cartesian_mesh_3d
130 
133  sll_int32 :: num_cells1
134  sll_int32 :: num_cells2
135  sll_int32 :: num_cells3
136  sll_int32 :: num_cells4
137  sll_real64 :: eta1_min
138  sll_real64 :: eta1_max
139  sll_real64 :: eta2_min
140  sll_real64 :: eta2_max
141  sll_real64 :: eta3_min
142  sll_real64 :: eta3_max
143  sll_real64 :: eta4_min
144  sll_real64 :: eta4_max
145  sll_real64 :: delta_eta1
146  sll_real64 :: delta_eta2
147  sll_real64 :: delta_eta3
148  sll_real64 :: delta_eta4
149  contains
150  procedure, pass(mesh) :: eta1_node => eta1_node_4d
151  procedure, pass(mesh) :: eta2_node => eta2_node_4d
152  procedure, pass(mesh) :: eta3_node => eta3_node_4d
153  procedure, pass(mesh) :: eta1_cell => eta1_cell_4d
154  procedure, pass(mesh) :: eta2_cell => eta2_cell_4d
155  procedure, pass(mesh) :: eta3_cell => eta3_cell_4d
156  procedure, pass(mesh) :: display => display_cartesian_mesh_4d
157  procedure, pass(mesh) :: delete => sll_s_cartesian_mesh_4d_free
158  end type sll_t_cartesian_mesh_4d
159 
161  type, public :: sll_t_cartesian_mesh_6d
162  sll_int32 :: num_cells(6)
163  sll_real64 :: eta_min(6)
164  sll_real64 :: eta_max(6)
165  sll_real64 :: delta_eta(6)
166  sll_real64 :: volume
167  sll_real64 :: volume_eta123
168  sll_real64 :: volume_eta456
169 
170  contains
171  procedure :: init => sll_s_cartesian_mesh_6d_init
172  procedure :: display => display_cartesian_mesh_6d
173  procedure :: delete => sll_s_cartesian_mesh_6d_free
174  end type sll_t_cartesian_mesh_6d
175 
178  interface sll_o_delete
179  module procedure sll_s_cartesian_mesh_1d_free
180  module procedure sll_s_cartesian_mesh_2d_free
181  module procedure sll_s_cartesian_mesh_3d_free
182  module procedure sll_s_cartesian_mesh_4d_free
183  end interface sll_o_delete
184 
186  interface operator(*)
187  module procedure sll_f_tensor_product_1d_1d
188  module procedure tensor_product_2d_2d
189  end interface operator(*)
190 
192  interface sll_o_display
193  module procedure display_cartesian_mesh_1d
194  module procedure display_cartesian_mesh_2d
195  module procedure display_cartesian_mesh_3d
196  module procedure display_cartesian_mesh_4d
197  end interface sll_o_display
198 
201  module procedure get_node_positions_1d
202  module procedure get_node_positions_2d
203  end interface sll_o_get_node_positions
204 
214  interface sll_o_new
215  module procedure sll_f_new_cartesian_mesh_1d
216  module procedure sll_f_new_cartesian_mesh_2d
217  module procedure sll_f_new_cartesian_mesh_3d
218  module procedure sll_f_new_cartesian_mesh_4d
219  end interface sll_o_new
220 
221  interface sll_o_cell
222  module procedure cell_cartesian_mesh_1d
223  end interface sll_o_cell
224 
226  module procedure cell_margin_cartesian_mesh_1d
227  end interface sll_o_cell_margin
228 
229  interface sll_o_mesh_area
230  module procedure length_cartesian_mesh_1d
231  module procedure area_cartesian_mesh_2d
232  module procedure area_cartesian_mesh_3d
233  module procedure area_cartesian_mesh_4d
234  end interface
235 
236 contains
237 
238 #define TEST_PRESENCE_AND_ASSIGN_VAL( obj, arg, slot, default_val ) \
239  if (present(arg)) then; \
240  obj%slot = arg; \
241  else; \
242  obj%slot = default_val; \
243  end if
244 
255  num_cells, &
256  eta_min, &
257  eta_max) result(m)
258 
259  type(sll_t_cartesian_mesh_1d), pointer :: m
260  sll_int32, intent(in) :: num_cells
261  sll_real64, optional, intent(in) :: eta_min
262  sll_real64, optional, intent(in) :: eta_max
263  !sll_real64 :: delta
264  sll_int32 :: ierr
265  sll_allocate(m, ierr)
266  call sll_s_cartesian_mesh_1d_init(m, num_cells, eta_min, eta_max)
267  end function sll_f_new_cartesian_mesh_1d
268 
276  subroutine sll_s_cartesian_mesh_1d_init(m, num_cells, eta_min, eta_max)
277  class(sll_t_cartesian_mesh_1d), intent(inout) :: m
278  sll_int32, intent(in) :: num_cells
279  sll_real64, optional, intent(in) :: eta_min
280  sll_real64, optional, intent(in) :: eta_max
281 
282  m%num_cells = num_cells
283 
284  test_presence_and_assign_val(m, eta_min, eta_min, 0.0_f64)
285  test_presence_and_assign_val(m, eta_max, eta_max, 1.0_f64)
286 
287  m%delta_eta = (m%eta_max - m%eta_min)/real(num_cells, f64)
288 
289  if (m%eta_max <= m%eta_min) then
290  print *, 'ERROR, sll_s_cartesian_mesh_1d_init(): ', &
291  'Problem to construct the mesh 1d '
292  print *, 'because eta_max <= eta_min'
293  end if
294  end subroutine sll_s_cartesian_mesh_1d_init
295 
297  function sll_f_tensor_product_1d_1d(m_a, m_b) result(m_c)
298  type(sll_t_cartesian_mesh_1d), intent(in), pointer :: m_a
299  type(sll_t_cartesian_mesh_1d), intent(in), pointer :: m_b
300  type(sll_t_cartesian_mesh_2d), pointer :: m_c
301 
303  m_a%num_cells, &
304  m_b%num_cells, &
305  m_a%eta_min, &
306  m_a%eta_max, &
307  m_b%eta_min, &
308  m_b%eta_max)
309 
310  end function sll_f_tensor_product_1d_1d
311 
313  function tensor_product_2d_2d(m_a, m_b) result(m_c)
314  type(sll_t_cartesian_mesh_2d), intent(in), pointer :: m_a
315  type(sll_t_cartesian_mesh_2d), intent(in), pointer :: m_b
316  type(sll_t_cartesian_mesh_4d), pointer :: m_c
317 
319  m_a%num_cells1, &
320  m_a%num_cells2, &
321  m_b%num_cells1, &
322  m_b%num_cells2, &
323  m_a%eta1_min, &
324  m_a%eta1_max, &
325  m_a%eta2_min, &
326  m_a%eta2_max, &
327  m_b%eta1_min, &
328  m_b%eta1_max, &
329  m_b%eta2_min, &
330  m_b%eta2_max)
331 
332  end function tensor_product_2d_2d
333 
334  subroutine get_node_positions_1d(m, eta1_node)
335  type(sll_t_cartesian_mesh_1d) :: m
336  sll_real64, dimension(:), allocatable :: eta1_node
337  sll_int32 :: num_cells
338  sll_real64 :: eta_min
339  sll_real64 :: delta_eta
340  sll_int32 :: i
341  sll_int32 :: ierr
342 
343  num_cells = m%num_cells
344  eta_min = m%eta_min
345  delta_eta = m%delta_eta
346  sll_allocate(eta1_node(num_cells + 1), ierr)
347  do i = 1, num_cells + 1
348  eta1_node(i) = eta_min + real(i - 1, f64)*delta_eta
349  end do
350  end subroutine get_node_positions_1d
351 
352  function eta1_node_1d(mesh, i) result(res)
353  class(sll_t_cartesian_mesh_1d), intent(in) :: mesh
354  sll_int32, intent(in) :: i
355  sll_real64 :: res
356  sll_real64 :: eta_min
357  sll_real64 :: delta_eta
358 
359  eta_min = mesh%eta_min
360  delta_eta = mesh%delta_eta
361  res = eta_min + real(i - 1, f64)*delta_eta
362  end function eta1_node_1d
363 
364  function eta1_cell_1d(mesh, i) result(res)
365  class(sll_t_cartesian_mesh_1d), intent(in) :: mesh
366  sll_int32, intent(in) :: i
367  sll_real64 :: res
368  sll_real64 :: eta_min
369  sll_real64 :: delta_eta
370 
371  eta_min = mesh%eta_min
372  delta_eta = mesh%delta_eta
373  res = eta_min + (real(i - 1, f64) + 0.5_f64)*delta_eta
374  end function eta1_cell_1d
375 
376  subroutine get_node_positions_2d(m, eta1, eta2)
377  class(sll_t_cartesian_mesh_2d) :: m
378  sll_real64, dimension(:, :), allocatable :: eta1
379  sll_real64, dimension(:, :), allocatable :: eta2
380  sll_int32 :: num_cells1
381  sll_int32 :: num_cells2
382  sll_real64 :: eta1_min
383  sll_real64 :: eta2_min
384  sll_real64 :: delta_eta1
385  sll_real64 :: delta_eta2
386  sll_int32 :: i
387  sll_int32 :: j
388  sll_int32 :: ierr
389 
390  num_cells1 = m%num_cells1
391  num_cells2 = m%num_cells2
392  eta1_min = m%eta1_min
393  delta_eta1 = m%delta_eta1
394  eta2_min = m%eta2_min
395  delta_eta2 = m%delta_eta2
396  sll_allocate(eta1(num_cells1 + 1, num_cells2 + 1), ierr)
397  sll_allocate(eta2(num_cells1 + 1, num_cells2 + 1), ierr)
398  do j = 1, num_cells2 + 1
399  do i = 1, num_cells1 + 1
400  eta1(i, j) = eta1_min + real(i - 1, f64)*delta_eta1
401  eta2(i, j) = eta2_min + real(j - 1, f64)*delta_eta2
402  end do
403  end do
404 
405  end subroutine get_node_positions_2d
406 
422  num_cells1, &
423  num_cells2, &
424  eta1_min, &
425  eta1_max, &
426  eta2_min, &
427  eta2_max) result(m)
428 
429  type(sll_t_cartesian_mesh_2d), pointer :: m
430  sll_int32, intent(in) :: num_cells1
431  sll_int32, intent(in) :: num_cells2
432  sll_real64, optional, intent(in) :: eta1_min
433  sll_real64, optional, intent(in) :: eta1_max
434  sll_real64, optional, intent(in) :: eta2_min
435  sll_real64, optional, intent(in) :: eta2_max
436  sll_int32 :: ierr
437 
438  sll_allocate(m, ierr)
440  m, &
441  num_cells1, &
442  num_cells2, &
443  eta1_min, &
444  eta1_max, &
445  eta2_min, &
446  eta2_max)
447 
448  end function sll_f_new_cartesian_mesh_2d
449 
464  m, &
465  num_cells1, &
466  num_cells2, &
467  eta1_min, &
468  eta1_max, &
469  eta2_min, &
470  eta2_max)
471 
472  class(sll_t_cartesian_mesh_2d), intent(inout) :: m
473  sll_int32, intent(in) :: num_cells1
474  sll_int32, intent(in) :: num_cells2
475  sll_real64, optional, intent(in) :: eta1_min
476  sll_real64, optional, intent(in) :: eta1_max
477  sll_real64, optional, intent(in) :: eta2_min
478  sll_real64, optional, intent(in) :: eta2_max
479 
480  m%num_cells1 = num_cells1
481  m%num_cells2 = num_cells2
482 
483  test_presence_and_assign_val(m, eta1_min, eta1_min, 0.0_f64)
484  test_presence_and_assign_val(m, eta1_max, eta1_max, 1.0_f64)
485  test_presence_and_assign_val(m, eta2_min, eta2_min, 0.0_f64)
486  test_presence_and_assign_val(m, eta2_max, eta2_max, 1.0_f64)
487 
488  m%delta_eta1 = (m%eta1_max - m%eta1_min)/real(num_cells1, f64)
489  m%delta_eta2 = (m%eta2_max - m%eta2_min)/real(num_cells2, f64)
490 
491  if (m%eta1_max <= m%eta1_min) then
492  sll_error('sll_s_cartesian_mesh_2d_init():', 'Problem to construct the mesh 2d because eta1_max <= eta1_min')
493  end if
494  if (m%eta2_max <= m%eta2_min) then
495  sll_error('sll_s_cartesian_mesh_2d_init():', 'Problem to construct the mesh 2d because eta2_max <= eta2_min')
496  end if
497  end subroutine sll_s_cartesian_mesh_2d_init
498 
499  function eta1_node_2d(mesh, i, j) result(res)
500  class(sll_t_cartesian_mesh_2d), intent(in) :: mesh
501  sll_int32, intent(in) :: i
502  sll_int32, intent(in) :: j
503  sll_real64 :: res
504  sll_real64 :: eta1_min
505  sll_real64 :: delta_eta1
506 #ifdef DEBUG
507  sll_int32 :: dummy
508  dummy = i + j
509 #endif
510 
511  eta1_min = mesh%eta1_min
512  delta_eta1 = mesh%delta_eta1
513  res = eta1_min + real(i - 1, f64)*delta_eta1
514  end function eta1_node_2d
515 
516  function eta2_node_2d(mesh, i, j) result(res)
517  class(sll_t_cartesian_mesh_2d), intent(in) :: mesh
518  sll_int32, intent(in) :: i
519  sll_int32, intent(in) :: j
520  sll_real64 :: res
521  sll_real64 :: eta2_min
522  sll_real64 :: delta_eta2
523 #ifdef DEBUG
524  sll_int32 :: dummy
525  dummy = i + j
526 #endif
527 
528  eta2_min = mesh%eta2_min
529  delta_eta2 = mesh%delta_eta2
530  res = eta2_min + real(j - 1, f64)*delta_eta2
531  end function eta2_node_2d
532 
533  function eta1_cell_2d_two_arg(mesh, i, j) result(res)
534  class(sll_t_cartesian_mesh_2d), intent(in) :: mesh
535  sll_int32, intent(in) :: i
536  sll_int32, intent(in) :: j
537  sll_real64 :: res
538  sll_real64 :: eta1_min
539  sll_real64 :: delta_eta1
540 #ifdef DEBUG
541  sll_int32 :: dummy
542  dummy = i + j
543 #endif
544 
545  eta1_min = mesh%eta1_min
546  delta_eta1 = mesh%delta_eta1
547  res = eta1_min + (real(i - 1, f64) + 0.5_f64)*delta_eta1
548  end function eta1_cell_2d_two_arg
549 
550  function eta2_cell_2d_two_arg(mesh, i, j) result(res)
551  class(sll_t_cartesian_mesh_2d), intent(in) :: mesh
552  sll_int32, intent(in) :: i
553  sll_int32, intent(in) :: j
554  sll_real64 :: res
555  sll_real64 :: eta2_min
556  sll_real64 :: delta_eta2
557 #ifdef DEBUG
558  sll_int32 :: dummy
559  dummy = i + j
560 #endif
561 
562  eta2_min = mesh%eta2_min
563  delta_eta2 = mesh%delta_eta2
564  res = eta2_min + (real(j - 1, f64) + 0.5_f64)*delta_eta2
565  end function eta2_cell_2d_two_arg
566 
567  function eta1_cell_2d_one_arg(mesh, cell_num) result(res)
568  class(sll_t_cartesian_mesh_2d), intent(in) :: mesh
569  sll_int32, intent(in) :: cell_num
570  sll_real64 :: res
571  sll_real64 :: eta1_min
572  sll_real64 :: delta_eta1
573  sll_int32 :: i
574 
575  i = mod(cell_num, mesh%num_cells1)
576  eta1_min = mesh%eta1_min
577  delta_eta1 = mesh%delta_eta1
578  res = eta1_min + (real(i - 1, f64) + 0.5_f64)*delta_eta1
579  end function eta1_cell_2d_one_arg
580 
581  function eta2_cell_2d_one_arg(mesh, cell_num) result(res)
582  class(sll_t_cartesian_mesh_2d), intent(in) :: mesh
583  sll_int32, intent(in) :: cell_num
584  sll_real64 :: res
585  sll_real64 :: eta2_min
586  sll_real64 :: delta_eta2
587  sll_int32 :: j
588 
589  j = cell_num/(mesh%num_cells1 + 1) + 1
590  eta2_min = mesh%eta2_min
591  delta_eta2 = mesh%delta_eta2
592  res = eta2_min + (real(j - 1, f64) + 0.5_f64)*delta_eta2
593  end function eta2_cell_2d_one_arg
594 
615  num_cells1, &
616  num_cells2, &
617  num_cells3, &
618  eta1_min, &
619  eta1_max, &
620  eta2_min, &
621  eta2_max, &
622  eta3_min, &
623  eta3_max) result(m)
624 
625  type(sll_t_cartesian_mesh_3d), pointer :: m
626  sll_int32, intent(in) :: num_cells1
627  sll_int32, intent(in) :: num_cells2
628  sll_int32, intent(in) :: num_cells3
629  sll_real64, optional, intent(in) :: eta1_min
630  sll_real64, optional, intent(in) :: eta1_max
631  sll_real64, optional, intent(in) :: eta2_min
632  sll_real64, optional, intent(in) :: eta2_max
633  sll_real64, optional, intent(in) :: eta3_min
634  sll_real64, optional, intent(in) :: eta3_max
635  sll_int32 :: ierr
636 
637  sll_allocate(m, ierr)
639  num_cells1, &
640  num_cells2, &
641  num_cells3, &
642  eta1_min, &
643  eta1_max, &
644  eta2_min, &
645  eta2_max, &
646  eta3_min, &
647  eta3_max)
648 
649  end function sll_f_new_cartesian_mesh_3d
650 
670  num_cells1, &
671  num_cells2, &
672  num_cells3, &
673  eta1_min, &
674  eta1_max, &
675  eta2_min, &
676  eta2_max, &
677  eta3_min, &
678  eta3_max)
679 
680  type(sll_t_cartesian_mesh_3d), intent(out) :: m
681  sll_int32, intent(in) :: num_cells1
682  sll_int32, intent(in) :: num_cells2
683  sll_int32, intent(in) :: num_cells3
684  sll_real64, optional, intent(in) :: eta1_min
685  sll_real64, optional, intent(in) :: eta1_max
686  sll_real64, optional, intent(in) :: eta2_min
687  sll_real64, optional, intent(in) :: eta2_max
688  sll_real64, optional, intent(in) :: eta3_min
689  sll_real64, optional, intent(in) :: eta3_max
690 
691  m%num_cells1 = num_cells1
692  m%num_cells2 = num_cells2
693  m%num_cells3 = num_cells3
694  test_presence_and_assign_val(m, eta1_min, eta1_min, 0.0_f64)
695  test_presence_and_assign_val(m, eta1_max, eta1_max, 1.0_f64)
696  test_presence_and_assign_val(m, eta2_min, eta2_min, 0.0_f64)
697  test_presence_and_assign_val(m, eta2_max, eta2_max, 1.0_f64)
698  test_presence_and_assign_val(m, eta3_min, eta3_min, 0.0_f64)
699  test_presence_and_assign_val(m, eta3_max, eta3_max, 1.0_f64)
700  m%delta_eta1 = (m%eta1_max - m%eta1_min)/real(num_cells1, f64)
701  m%delta_eta2 = (m%eta2_max - m%eta2_min)/real(num_cells2, f64)
702  m%delta_eta3 = (m%eta3_max - m%eta3_min)/real(num_cells3, f64)
703 
704  if (m%eta1_max <= m%eta1_min) then
705  print *, 'Problem to construct the mesh 3d '
706  print *, 'because eta1_max <= eta1_min'
707  end if
708  if (m%eta2_max <= m%eta2_min) then
709  print *, 'Problem to construct the mesh 3d '
710  print *, 'because eta2_max <= eta2_min'
711  end if
712  if (m%eta3_max <= m%eta3_min) then
713  print *, 'Problem to construct the mesh 3d '
714  print *, 'because eta3_max <= eta3_min'
715  end if
716 
717  end subroutine sll_s_cartesian_mesh_3d_init
718 
719  function eta1_node_3d(mesh, i1, i2, i3) result(res)
720  class(sll_t_cartesian_mesh_3d), intent(in) :: mesh
721  sll_int32, intent(in) :: i1
722  sll_int32, intent(in) :: i2
723  sll_int32, intent(in) :: i3
724  sll_real64 :: res
725  sll_real64 :: eta1_min
726  sll_real64 :: delta_eta1
727 #ifdef DEBUG
728  sll_int32 :: dummy
729  dummy = i1 + i2 + i3
730 #endif
731 
732  eta1_min = mesh%eta1_min
733  delta_eta1 = mesh%delta_eta1
734  res = eta1_min + real(i1 - 1, f64)*delta_eta1
735  end function eta1_node_3d
736 
737  function eta2_node_3d(mesh, i1, i2, i3) result(res)
738  class(sll_t_cartesian_mesh_3d), intent(in) :: mesh
739  sll_int32, intent(in) :: i1
740  sll_int32, intent(in) :: i2
741  sll_int32, intent(in) :: i3
742  sll_real64 :: res
743  sll_real64 :: eta2_min
744  sll_real64 :: delta_eta2
745 #ifdef DEBUG
746  sll_int32 :: dummy
747  dummy = i1 + i2 + i3
748 #endif
749 
750  eta2_min = mesh%eta2_min
751  delta_eta2 = mesh%delta_eta2
752  res = eta2_min + real(i2 - 1, f64)*delta_eta2
753  end function eta2_node_3d
754 
755  function eta3_node_3d(mesh, i1, i2, i3) result(res)
756  class(sll_t_cartesian_mesh_3d), intent(in) :: mesh
757  sll_int32, intent(in) :: i1
758  sll_int32, intent(in) :: i2
759  sll_int32, intent(in) :: i3
760  sll_real64 :: res
761  sll_real64 :: eta3_min
762  sll_real64 :: delta_eta3
763 #ifdef DEBUG
764  sll_int32 :: dummy
765  dummy = i1 + i2 + i3
766 #endif
767 
768  eta3_min = mesh%eta3_min
769  delta_eta3 = mesh%delta_eta3
770  res = eta3_min + real(i3 - 1, f64)*delta_eta3
771  end function eta3_node_3d
772 
773  function eta1_cell_3d(mesh, i1, i2, i3) result(res)
774  class(sll_t_cartesian_mesh_3d), intent(in) :: mesh
775  sll_int32, intent(in) :: i1
776  sll_int32, intent(in) :: i2
777  sll_int32, intent(in) :: i3
778  sll_real64 :: res
779  sll_real64 :: eta1_min
780  sll_real64 :: delta_eta1
781 #ifdef DEBUG
782  sll_int32 :: dummy
783  dummy = i1 + i2 + i3
784 #endif
785 
786  eta1_min = mesh%eta1_min
787  delta_eta1 = mesh%delta_eta1
788  res = eta1_min + (real(i1 - 1, f64) + 0.5_f64)*delta_eta1
789  end function eta1_cell_3d
790 
791  function eta2_cell_3d(mesh, i1, i2, i3) result(res)
792  class(sll_t_cartesian_mesh_3d), intent(in) :: mesh
793  sll_int32, intent(in) :: i1
794  sll_int32, intent(in) :: i2
795  sll_int32, intent(in) :: i3
796  sll_real64 :: res
797  sll_real64 :: eta2_min
798  sll_real64 :: delta_eta2
799 #ifdef DEBUG
800  sll_int32 :: dummy
801  dummy = i1 + i2 + i3
802 #endif
803 
804  eta2_min = mesh%eta2_min
805  delta_eta2 = mesh%delta_eta2
806  res = eta2_min + (real(i2 - 1, f64) + 0.5_f64)*delta_eta2
807  end function eta2_cell_3d
808 
809  function eta3_cell_3d(mesh, i1, i2, i3) result(res)
810  class(sll_t_cartesian_mesh_3d), intent(in) :: mesh
811  sll_int32, intent(in) :: i1
812  sll_int32, intent(in) :: i2
813  sll_int32, intent(in) :: i3
814  sll_real64 :: res
815  sll_real64 :: eta3_min
816  sll_real64 :: delta_eta3
817 #ifdef DEBUG
818  sll_int32 :: dummy
819  dummy = i1 + i2 + i3
820 #endif
821 
822  eta3_min = mesh%eta3_min
823  delta_eta3 = mesh%delta_eta3
824  res = eta3_min + (real(i3 - 1, f64) + 0.5_f64)*delta_eta3
825  end function eta3_cell_3d
826 
854  num_cells1, &
855  num_cells2, &
856  num_cells3, &
857  num_cells4, &
858  eta1_min, &
859  eta1_max, &
860  eta2_min, &
861  eta2_max, &
862  eta3_min, &
863  eta3_max, &
864  eta4_min, &
865  eta4_max) result(m)
866 
867  type(sll_t_cartesian_mesh_4d), pointer :: m
868  sll_int32, intent(in) :: num_cells1
869  sll_int32, intent(in) :: num_cells2
870  sll_int32, intent(in) :: num_cells3
871  sll_int32, intent(in) :: num_cells4
872  sll_real64, optional, intent(in) :: eta1_min
873  sll_real64, optional, intent(in) :: eta1_max
874  sll_real64, optional, intent(in) :: eta2_min
875  sll_real64, optional, intent(in) :: eta2_max
876  sll_real64, optional, intent(in) :: eta3_min
877  sll_real64, optional, intent(in) :: eta3_max
878  sll_real64, optional, intent(in) :: eta4_min
879  sll_real64, optional, intent(in) :: eta4_max
880  sll_int32 :: ierr
881 
882  sll_allocate(m, ierr)
883  m%num_cells1 = num_cells1
884  m%num_cells2 = num_cells2
885  m%num_cells3 = num_cells3
886  m%num_cells4 = num_cells4
887  test_presence_and_assign_val(m, eta1_min, eta1_min, 0.0_f64)
888  test_presence_and_assign_val(m, eta1_max, eta1_max, 1.0_f64)
889  test_presence_and_assign_val(m, eta2_min, eta2_min, 0.0_f64)
890  test_presence_and_assign_val(m, eta2_max, eta2_max, 1.0_f64)
891  test_presence_and_assign_val(m, eta3_min, eta3_min, 0.0_f64)
892  test_presence_and_assign_val(m, eta3_max, eta3_max, 1.0_f64)
893  test_presence_and_assign_val(m, eta4_min, eta4_min, 0.0_f64)
894  test_presence_and_assign_val(m, eta4_max, eta4_max, 1.0_f64)
895  m%delta_eta1 = (m%eta1_max - m%eta1_min)/real(num_cells1, f64)
896  m%delta_eta2 = (m%eta2_max - m%eta2_min)/real(num_cells2, f64)
897  m%delta_eta3 = (m%eta3_max - m%eta3_min)/real(num_cells3, f64)
898  m%delta_eta4 = (m%eta4_max - m%eta4_min)/real(num_cells4, f64)
899 
900  if (m%eta1_max <= m%eta1_min) then
901  print *, 'Problem to construct the mesh 4d '
902  print *, 'because eta1_max <= eta1_min'
903  end if
904  if (m%eta2_max <= m%eta2_min) then
905  print *, 'Problem to construct the mesh 4d '
906  print *, 'because eta2_max <= eta2_min'
907  end if
908  if (m%eta3_max <= m%eta3_min) then
909  print *, 'Problem to construct the mesh 4d '
910  print *, 'because eta3_max <= eta3_min'
911  end if
912  if (m%eta4_max <= m%eta4_min) then
913  print *, 'Problem to construct the mesh 4d '
914  print *, 'because eta4_max <= eta4_min'
915  end if
916  end function sll_f_new_cartesian_mesh_4d
917 
918  function eta1_node_4d(mesh, i1, i2, i3, i4) result(res)
919  class(sll_t_cartesian_mesh_4d), intent(in) :: mesh
920  sll_int32, intent(in) :: i1
921  sll_int32, intent(in) :: i2
922  sll_int32, intent(in) :: i3
923  sll_int32, intent(in) :: i4
924  sll_real64 :: res
925  sll_real64 :: eta1_min
926  sll_real64 :: delta_eta1
927 #ifdef DEBUG
928  sll_int32 :: dummy
929  dummy = i1 + i2 + i3 + i4
930 #endif
931 
932  eta1_min = mesh%eta1_min
933  delta_eta1 = mesh%delta_eta1
934  res = eta1_min + real(i1 - 1, f64)*delta_eta1
935  end function eta1_node_4d
936 
937  function eta2_node_4d(mesh, i1, i2, i3, i4) result(res)
938  class(sll_t_cartesian_mesh_4d), intent(in) :: mesh
939  sll_int32, intent(in) :: i1
940  sll_int32, intent(in) :: i2
941  sll_int32, intent(in) :: i3
942  sll_int32, intent(in) :: i4
943  sll_real64 :: res
944  sll_real64 :: eta2_min
945  sll_real64 :: delta_eta2
946 #ifdef DEBUG
947  sll_int32 :: dummy
948  dummy = i1 + i2 + i3 + i4
949 #endif
950 
951  eta2_min = mesh%eta2_min
952  delta_eta2 = mesh%delta_eta2
953  res = eta2_min + real(i2 - 1, f64)*delta_eta2
954  end function eta2_node_4d
955 
956  function eta3_node_4d(mesh, i1, i2, i3, i4) result(res)
957  class(sll_t_cartesian_mesh_4d), intent(in) :: mesh
958  sll_int32, intent(in) :: i1
959  sll_int32, intent(in) :: i2
960  sll_int32, intent(in) :: i3
961  sll_int32, intent(in) :: i4
962  sll_real64 :: res
963  sll_real64 :: eta3_min
964  sll_real64 :: delta_eta3
965 #ifdef DEBUG
966  sll_int32 :: dummy
967  dummy = i1 + i2 + i3 + i4
968 #endif
969 
970  eta3_min = mesh%eta3_min
971  delta_eta3 = mesh%delta_eta3
972  res = eta3_min + real(i3 - 1, f64)*delta_eta3
973  end function eta3_node_4d
974 
975  function eta4_node_4d(mesh, i1, i2, i3, i4) result(res)
976  class(sll_t_cartesian_mesh_4d), intent(in) :: mesh
977  sll_int32, intent(in) :: i1
978  sll_int32, intent(in) :: i2
979  sll_int32, intent(in) :: i3
980  sll_int32, intent(in) :: i4
981  sll_real64 :: res
982  sll_real64 :: eta4_min
983  sll_real64 :: delta_eta4
984 #ifdef DEBUG
985  sll_int32 :: dummy
986  dummy = i1 + i2 + i3 + i4
987 #endif
988 
989  eta4_min = mesh%eta4_min
990  delta_eta4 = mesh%delta_eta4
991  res = eta4_min + real(i4 - 1, f64)*delta_eta4
992  end function eta4_node_4d
993 
994  function eta1_cell_4d(mesh, i1, i2, i3, i4) result(res)
995  class(sll_t_cartesian_mesh_4d), intent(in) :: mesh
996  sll_int32, intent(in) :: i1
997  sll_int32, intent(in) :: i2
998  sll_int32, intent(in) :: i3
999  sll_int32, intent(in) :: i4
1000  sll_real64 :: res
1001  sll_real64 :: eta1_min
1002  sll_real64 :: delta_eta1
1003 #ifdef DEBUG
1004  sll_int32 :: dummy
1005  dummy = i1 + i2 + i3 + i4
1006 #endif
1007 
1008  eta1_min = mesh%eta1_min
1009  delta_eta1 = mesh%delta_eta1
1010  res = eta1_min + (real(i1 - 1, f64) + 0.5_f64)*delta_eta1
1011  end function eta1_cell_4d
1012 
1013  function eta2_cell_4d(mesh, i1, i2, i3, i4) result(res)
1014  class(sll_t_cartesian_mesh_4d), intent(in) :: mesh
1015  sll_int32, intent(in) :: i1
1016  sll_int32, intent(in) :: i2
1017  sll_int32, intent(in) :: i3
1018  sll_int32, intent(in) :: i4
1019  sll_real64 :: res
1020  sll_real64 :: eta2_min
1021  sll_real64 :: delta_eta2
1022 #ifdef DEBUG
1023  sll_int32 :: dummy
1024  dummy = i1 + i2 + i3 + i4
1025 #endif
1026 
1027  eta2_min = mesh%eta2_min
1028  delta_eta2 = mesh%delta_eta2
1029  res = eta2_min + (real(i2 - 1, f64) + 0.5_f64)*delta_eta2
1030  end function eta2_cell_4d
1031 
1032  function eta3_cell_4d(mesh, i1, i2, i3, i4) result(res)
1033  class(sll_t_cartesian_mesh_4d), intent(in) :: mesh
1034  sll_int32, intent(in) :: i1
1035  sll_int32, intent(in) :: i2
1036  sll_int32, intent(in) :: i3
1037  sll_int32, intent(in) :: i4
1038  sll_real64 :: res
1039  sll_real64 :: eta3_min
1040  sll_real64 :: delta_eta3
1041 #ifdef DEBUG
1042  sll_int32 :: dummy
1043  dummy = i1 + i2 + i3 + i4
1044 #endif
1045 
1046  eta3_min = mesh%eta3_min
1047  delta_eta3 = mesh%delta_eta3
1048  res = eta3_min + (real(i3 - 1, f64) + 0.5_f64)*delta_eta3
1049  end function eta3_cell_4d
1050 
1051  function eta4_cell_4d(mesh, i1, i2, i3, i4) result(res)
1052  class(sll_t_cartesian_mesh_4d), intent(in) :: mesh
1053  sll_int32, intent(in) :: i1
1054  sll_int32, intent(in) :: i2
1055  sll_int32, intent(in) :: i3
1056  sll_int32, intent(in) :: i4
1057  sll_real64 :: res
1058  sll_real64 :: eta4_min
1059  sll_real64 :: delta_eta4
1060 #ifdef DEBUG
1061  sll_int32 :: dummy
1062  dummy = i1 + i2 + i3 + i4
1063 #endif
1064 
1065  eta4_min = mesh%eta4_min
1066  delta_eta4 = mesh%delta_eta4
1067  res = eta4_min + (real(i4 - 1, f64) + 0.5_f64)*delta_eta4
1068  end function eta4_cell_4d
1069 
1073  subroutine display_cartesian_mesh_1d(mesh)
1074  class(sll_t_cartesian_mesh_1d), intent(in) :: mesh
1075 
1076  write (*, "(/,(a))") '1D mesh : num_cell eta_min eta_max delta_eta'
1077  write (*, "(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells, &
1078  mesh%eta_min, &
1079  mesh%eta_max, &
1080  mesh%delta_eta
1081 
1082  end subroutine display_cartesian_mesh_1d
1083 
1087  subroutine display_cartesian_mesh_2d(mesh)
1088  class(sll_t_cartesian_mesh_2d), intent(in) :: mesh
1089 
1090  write (*, "(/,(a))") '2D mesh : num_cell eta_min eta_max delta_eta'
1091  write (*, "(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells1, &
1092  mesh%eta1_min, &
1093  mesh%eta1_max, &
1094  mesh%delta_eta1
1095  write (*, "(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells2, &
1096  mesh%eta2_min, &
1097  mesh%eta2_max, &
1098  mesh%delta_eta2
1099  end subroutine display_cartesian_mesh_2d
1100 
1104  subroutine display_cartesian_mesh_3d(mesh)
1105  class(sll_t_cartesian_mesh_3d), intent(in) :: mesh
1106 
1107  write (*, "(/,(a))") '3D mesh : num_cell eta_min eta_max delta_eta'
1108  write (*, "(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells1, &
1109  mesh%eta1_min, &
1110  mesh%eta1_max, &
1111  mesh%delta_eta1
1112  write (*, "(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells2, &
1113  mesh%eta2_min, &
1114  mesh%eta2_max, &
1115  mesh%delta_eta2
1116  write (*, "(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells3, &
1117  mesh%eta3_min, &
1118  mesh%eta3_max, &
1119  mesh%delta_eta3
1120  end subroutine display_cartesian_mesh_3d
1121 
1125  subroutine display_cartesian_mesh_4d(mesh)
1126  class(sll_t_cartesian_mesh_4d), intent(in) :: mesh
1127 
1128  write (*, "(/,(a))") '4D mesh : num_cell eta_min eta_max delta_eta'
1129  write (*, "(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells1, &
1130  mesh%eta1_min, &
1131  mesh%eta1_max, &
1132  mesh%delta_eta1
1133  write (*, "(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells2, &
1134  mesh%eta2_min, &
1135  mesh%eta2_max, &
1136  mesh%delta_eta2
1137  write (*, "(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells3, &
1138  mesh%eta3_min, &
1139  mesh%eta3_max, &
1140  mesh%delta_eta3
1141  write (*, "(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells4, &
1142  mesh%eta4_min, &
1143  mesh%eta4_max, &
1144  mesh%delta_eta4
1145  end subroutine display_cartesian_mesh_4d
1146 
1151  class(sll_t_cartesian_mesh_1d), intent(inout) :: mesh
1152  ! sll_int32 :: ierr
1153  ! if(.not. associated(mesh))then
1154  ! print *, 'sll_s_cartesian_mesh_1d_free, ERROR: passed argument is not ', &
1155  ! 'associated. Crash imminent...'
1156  ! end if
1157  ! SLL_DEALLOCATE(mesh, ierr)
1158  sll_assert(mesh%num_cells > 0)
1159  end subroutine sll_s_cartesian_mesh_1d_free
1160 
1165  class(sll_t_cartesian_mesh_2d), intent(inout) :: mesh
1166  ! sll_int32 :: ierr
1167  ! if(.not. associated(mesh))then
1168  ! print *, 'sll_s_cartesian_mesh_2d_free, ERROR: passed argument is not ', &
1169  ! 'associated. Crash imminent...'
1170  ! end if
1171  ! SLL_DEALLOCATE(mesh, ierr)
1172  sll_assert(mesh%num_cells1 > 0)
1173  end subroutine sll_s_cartesian_mesh_2d_free
1174 
1179  class(sll_t_cartesian_mesh_3d), intent(inout) :: mesh
1180  ! sll_int32 :: ierr
1181  ! if(.not. associated(mesh))then
1182  ! print *, 'sll_s_cartesian_mesh_3d_free, ERROR: passed argument is not ', &
1183  ! 'associated. Crash imminent...'
1184  ! end if
1185  ! SLL_DEALLOCATE(mesh, ierr)
1186  sll_assert(mesh%num_cells1 > 0)
1187  end subroutine sll_s_cartesian_mesh_3d_free
1188 
1193  class(sll_t_cartesian_mesh_4d), intent(inout) :: mesh
1194  ! sll_int32 :: ierr
1195  ! if(.not. associated(mesh))then
1196  ! print *, 'sll_s_cartesian_mesh_4d_free, ERROR: passed argument is not ', &
1197  ! 'associated. Crash imminent...'
1198  ! end if
1199  ! SLL_DEALLOCATE(mesh, ierr)
1200  sll_assert(mesh%num_cells1 > 0)
1201  end subroutine sll_s_cartesian_mesh_4d_free
1202 
1205  function nodes_cartesian_mesh_1d(mesh) result(nodes)
1206  class(sll_t_cartesian_mesh_1d), intent(in) :: mesh
1207  sll_real64, dimension(mesh%num_cells + 1) ::nodes
1208  sll_int32 :: idx
1209  sll_int32 :: nknots
1210  nknots = mesh%num_cells + 1
1211  do idx = 1, nknots
1212  nodes(idx) = mesh%eta_min + real(idx - 1, f64)*mesh%delta_eta
1213  end do
1214  nodes(1) = mesh%eta_min
1215  nodes(nknots) = mesh%eta_max
1216  end function
1217 
1222  function cell_cartesian_mesh_1d(mesh, point) result(cell)
1223  class(sll_t_cartesian_mesh_1d), intent(in) :: mesh
1224  sll_real64, dimension(:), intent(in) ::point
1225  sll_int32, dimension(size(point)) :: cell
1226  cell = floor((point - mesh%eta_min)/mesh%delta_eta) + 1
1227  !Last knot belongs to last cell
1228  where (cell == mesh%num_cells + 1) cell = mesh%num_cells
1229  end function
1230 
1236  function cell_margin_cartesian_mesh_1d(mesh, cell) result(margin)
1237  class(sll_t_cartesian_mesh_1d), intent(in) :: mesh
1238  sll_int32, intent(in) :: cell
1239  sll_real64, dimension(2) :: margin
1240 
1241  margin(1) = mesh%eta_min + real(cell - 1, f64)*mesh%delta_eta
1242  margin(2) = mesh%eta_min + real(cell, f64)*mesh%delta_eta
1243  !!SLL_ASSERT(margin(2)<=mesh%eta_max)
1244  end function
1245 
1248  function num_nodes_cartesian_mesh_1d(mesh) result(num_nodes)
1249  class(sll_t_cartesian_mesh_1d), intent(in) :: mesh
1250  sll_int32 :: num_nodes
1251  num_nodes = mesh%num_cells + 1
1252  end function num_nodes_cartesian_mesh_1d
1253 
1256  function length_cartesian_mesh_1d(mesh) result(length)
1257  class(sll_t_cartesian_mesh_1d), intent(in) :: mesh
1258  sll_real64 :: length
1259  length = mesh%eta_max - mesh%eta_min
1260  sll_assert(length >= 0)
1261  end function length_cartesian_mesh_1d
1262 
1265  function area_cartesian_mesh_2d(mesh) result(area)
1266  class(sll_t_cartesian_mesh_2d), intent(in) :: mesh
1267  sll_real64 :: area
1268  area = (mesh%eta1_max - mesh%eta1_min)* &
1269  (mesh%eta2_max - mesh%eta2_min)
1270 
1271  sll_assert(area >= 0)
1272  end function area_cartesian_mesh_2d
1273 
1276  function area_cartesian_mesh_3d(mesh) result(area)
1277  class(sll_t_cartesian_mesh_3d), intent(in) :: mesh
1278  sll_real64 :: area
1279  area = (mesh%eta1_max - mesh%eta1_min)* &
1280  (mesh%eta2_max - mesh%eta2_min)* &
1281  (mesh%eta3_max - mesh%eta3_min)
1282 
1283  sll_assert(area >= 0)
1284  end function area_cartesian_mesh_3d
1285 
1288  function area_cartesian_mesh_4d(mesh) result(area)
1289  class(sll_t_cartesian_mesh_4d), intent(in) :: mesh
1290  sll_real64 :: area
1291  area = (mesh%eta1_max - mesh%eta1_min)* &
1292  (mesh%eta2_max - mesh%eta2_min)* &
1293  (mesh%eta3_max - mesh%eta3_min)* &
1294  (mesh%eta4_max - mesh%eta4_min)
1295 
1296  sll_assert(area >= 0)
1297  end function area_cartesian_mesh_4d
1298 
1302  function period_cartesian_mesh_1d(mesh, point) result(posmod)
1303  class(sll_t_cartesian_mesh_1d), intent(in) :: mesh
1304  sll_real64, intent(in) :: point
1305  sll_real64 :: posmod
1306 
1307  posmod = modulo(point - mesh%eta_min, mesh%eta_max - mesh%eta_min) + mesh%eta_min
1308 
1309  sll_assert(posmod >= mesh%eta_min)
1310  sll_assert(posmod <= mesh%eta_max)
1311  end function period_cartesian_mesh_1d
1312 
1313  ! 6D procedures
1314 
1321  this, &
1322  num_cells, &
1323  eta_min, &
1324  eta_max)
1325  class(sll_t_cartesian_mesh_6d), intent(inout) :: this
1326  sll_int32, intent(in) :: num_cells(6)
1327  sll_real64, intent(in) :: eta_min(6)
1328  sll_real64, intent(in) :: eta_max(6)
1329 
1330  this%num_cells = num_cells
1331  this%eta_min = eta_min
1332  this%eta_max = eta_max
1333 
1334  this%delta_eta = (eta_max - eta_min)/real(num_cells, f64)
1335 
1336  this%volume_eta123 = product(this%delta_eta(1:3))
1337  this%volume_eta456 = product(this%delta_eta(4:6))
1338  this%volume = this%volume_eta123*this%volume_eta456
1339 
1340  end subroutine sll_s_cartesian_mesh_6d_init
1341 
1345  subroutine display_cartesian_mesh_6d(mesh)
1346  class(sll_t_cartesian_mesh_6d), intent(in) :: mesh
1347 
1348  write (*, "(/,(a))") '6D mesh : num_cell eta_min eta_max delta_eta'
1349  write (*, "(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells(1), &
1350  mesh%eta_min(1), &
1351  mesh%eta_max(1), &
1352  mesh%delta_eta(1)
1353  write (*, "(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells(2), &
1354  mesh%eta_min(2), &
1355  mesh%eta_max(2), &
1356  mesh%delta_eta(2)
1357  write (*, "(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells(3), &
1358  mesh%eta_min(3), &
1359  mesh%eta_max(3), &
1360  mesh%delta_eta(3)
1361  write (*, "(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells(4), &
1362  mesh%eta_min(4), &
1363  mesh%eta_max(4), &
1364  mesh%delta_eta(4)
1365  write (*, "(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells(5), &
1366  mesh%eta_min(5), &
1367  mesh%eta_max(5), &
1368  mesh%delta_eta(5)
1369  write (*, "(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells(6), &
1370  mesh%eta_min(6), &
1371  mesh%eta_max(6), &
1372  mesh%delta_eta(6)
1373  end subroutine display_cartesian_mesh_6d
1374 
1378  class(sll_t_cartesian_mesh_6d), intent(inout) :: this
1379 
1380  end subroutine sll_s_cartesian_mesh_6d_free
1381 
1382 #undef TEST_PRESENCE_AND_ASSIGN_VAL
1383 
1384 end module sll_m_cartesian_meshes
Deallocates memory for the cartesian mesh.
allocates the memory space for a new cartesian mesh on the heap, initializes it with the given argume...
Cartesian mesh basic types.
real(kind=f64) function, dimension(2) cell_margin_cartesian_mesh_1d(mesh, cell)
Returns the margin (a,b) for a given cell.
real(kind=f64) function eta3_node_3d(mesh, i1, i2, i3)
real(kind=f64) function eta1_node_1d(mesh, i)
subroutine display_cartesian_mesh_1d(mesh)
display contents of a 1D cartesian mesh. Recommended access through the generic interface sll_o_displ...
integer(kind=i32) function, dimension(size(point)) cell_cartesian_mesh_1d(mesh, point)
Returns cell number(s) for given point(s) in cartesian mesh.
subroutine sll_s_cartesian_mesh_4d_free(mesh)
deallocates memory for the 4D cartesian mesh. Recommended access through the generic interface delete...
real(kind=f64) function eta1_cell_2d_one_arg(mesh, cell_num)
real(kind=f64) function eta1_cell_4d(mesh, i1, i2, i3, i4)
subroutine sll_s_cartesian_mesh_2d_free(mesh)
deallocates memory for the 2D cartesian mesh. Recommended access through the generic interface delete...
subroutine display_cartesian_mesh_6d(mesh)
display contents of a 6d cartesian mesh. Recommended access through the generic interface sll_display...
subroutine display_cartesian_mesh_4d(mesh)
display contents of a 4d cartesian mesh. Recommended access through the generic interface sll_o_displ...
subroutine, public sll_s_cartesian_mesh_3d_free(mesh)
deallocates memory for the 3D cartesian mesh. Recommended access through the generic interface delete...
real(kind=f64) function eta2_cell_3d(mesh, i1, i2, i3)
real(kind=f64) function eta2_node_2d(mesh, i, j)
subroutine, public sll_s_cartesian_mesh_3d_init(m, num_cells1, num_cells2, num_cells3, eta1_min, eta1_max, eta2_min, eta2_max, eta3_min, eta3_max)
allocates the memory space for a new 3D cartesian mesh, initializes it with the given arguments and r...
subroutine sll_s_cartesian_mesh_6d_free(this)
destructor of 6D cartesian mesh
real(kind=f64) function, dimension(mesh%num_cells+1) nodes_cartesian_mesh_1d(mesh)
Returns all nodes for the 1D cartesian mesh.
type(sll_t_cartesian_mesh_2d) function, pointer, public sll_f_new_cartesian_mesh_2d(num_cells1, num_cells2, eta1_min, eta1_max, eta2_min, eta2_max)
allocates the memory space for a new 2D cartesian mesh on the heap, initializes it with the given arg...
real(kind=f64) function eta2_cell_2d_two_arg(mesh, i, j)
real(kind=f64) function eta3_cell_4d(mesh, i1, i2, i3, i4)
real(kind=f64) function eta4_node_4d(mesh, i1, i2, i3, i4)
subroutine get_node_positions_2d(m, eta1, eta2)
subroutine sll_s_cartesian_mesh_1d_free(mesh)
deallocates memory for the 1D cartesian mesh. Recommended access through the generic interface delete...
real(kind=f64) function period_cartesian_mesh_1d(mesh, point)
Returns coordinate in a periodic mesh.
type(sll_t_cartesian_mesh_4d) function, pointer, public sll_f_new_cartesian_mesh_4d(num_cells1, num_cells2, num_cells3, num_cells4, eta1_min, eta1_max, eta2_min, eta2_max, eta3_min, eta3_max, eta4_min, eta4_max)
allocates the memory space for a new 4D cartesian mesh on the heap,
type(sll_t_cartesian_mesh_3d) function, pointer, public sll_f_new_cartesian_mesh_3d(num_cells1, num_cells2, num_cells3, eta1_min, eta1_max, eta2_min, eta2_max, eta3_min, eta3_max)
allocates the memory space for a new 3D cartesian mesh on the heap, initializes it with the given arg...
subroutine sll_s_cartesian_mesh_6d_init(this, num_cells, eta_min, eta_max)
initialize a 6D cartesian mesh
real(kind=f64) function area_cartesian_mesh_4d(mesh)
Returns the area size for a 3D cartesian mesh.
real(kind=f64) function eta1_cell_2d_two_arg(mesh, i, j)
real(kind=f64) function length_cartesian_mesh_1d(mesh)
Returns the interval length for a 1D cartesian mesh.
integer(kind=i32) function num_nodes_cartesian_mesh_1d(mesh)
Returns the number of nodes for the 1D cartesian mesh.
type(sll_t_cartesian_mesh_1d) function, pointer, public sll_f_new_cartesian_mesh_1d(num_cells, eta_min, eta_max)
allocates the memory space for a new 1D cartesian mesh on the heap, initializes it with the given arg...
real(kind=f64) function eta2_node_3d(mesh, i1, i2, i3)
subroutine display_cartesian_mesh_3d(mesh)
display contents of a 3d cartesian mesh. Recommended access through the generic interface sll_o_displ...
subroutine, public sll_s_cartesian_mesh_1d_init(m, num_cells, eta_min, eta_max)
Initializes a previously allocated 1D cartesian mesh object.
real(kind=f64) function eta1_cell_1d(mesh, i)
real(kind=f64) function eta3_node_4d(mesh, i1, i2, i3, i4)
real(kind=f64) function eta2_cell_2d_one_arg(mesh, cell_num)
real(kind=f64) function eta1_cell_3d(mesh, i1, i2, i3)
real(kind=f64) function eta1_node_4d(mesh, i1, i2, i3, i4)
subroutine, public sll_s_cartesian_mesh_2d_init(m, num_cells1, num_cells2, eta1_min, eta1_max, eta2_min, eta2_max)
initializes a cartesian mesh 2D object that has been already allocated.
real(kind=f64) function area_cartesian_mesh_2d(mesh)
Returns the area size for a 2D cartesian mesh.
type(sll_t_cartesian_mesh_2d) function, pointer, public sll_f_tensor_product_1d_1d(m_a, m_b)
Create a 2d mesh from two 1d meshes.
real(kind=f64) function eta2_node_4d(mesh, i1, i2, i3, i4)
real(kind=f64) function area_cartesian_mesh_3d(mesh)
Returns the area size for a 3D cartesian mesh.
subroutine display_cartesian_mesh_2d(mesh)
display contents of a 2d cartesian mesh. Recommended access through the generic interface sll_o_displ...
real(kind=f64) function eta1_node_2d(mesh, i, j)
real(kind=f64) function eta4_cell_4d(mesh, i1, i2, i3, i4)
real(kind=f64) function eta1_node_3d(mesh, i1, i2, i3)
type(sll_t_cartesian_mesh_4d) function, pointer tensor_product_2d_2d(m_a, m_b)
Create a 4d mesh from two 2d meshes.
real(kind=f64) function eta3_cell_3d(mesh, i1, i2, i3)
real(kind=f64) function eta2_cell_4d(mesh, i1, i2, i3, i4)
subroutine get_node_positions_1d(m, eta1_node)
    Report Typos and Errors