Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_linear_operator_particle_mass_3d_od.F90
Go to the documentation of this file.
2 #include "sll_working_precision.h"
3 
4  use sll_m_particle_mass_3d_base, only: &
6 
7 
8  implicit none
9 
11 
12  private
13 
15  sll_int32 :: degree2(3)
16 
17  contains
19  procedure :: free => free_particle_mass_3d_od
20  procedure :: dot => dot_particle_mass_3d_od
21  procedure :: print_info => print_info_particle_mass_3d_od
22 
24 
25 
26 contains
27 
28  subroutine create_linear_operator_particle_mass_3d_od( self, n_cells, degree, degree2 )
29  class(sll_t_linear_operator_particle_mass_3d_od), intent( inout ) :: self
30  sll_int32, intent( in ) :: n_cells(3)
31  sll_int32, intent( in ) :: degree(3)
32  sll_int32, optional, intent( in ) :: degree2(3)
33 
34  self%degree = degree
35  self%degree2 = degree2
36  self%n_dofs = n_cells
37  self%n_total = product(n_cells)
38 
39  allocate( self%particle_mass( product(self%degree+self%degree2+1), self%n_total))
40  self%particle_mass = 0._f64
41 
42  self%n_rows = self%n_total
43  self%n_cols = self%n_total
44 
45  self%n_global_rows = self%n_rows
46  self%n_global_cols = self%n_cols
47 
49 
50  subroutine free_particle_mass_3d_od( self )
51  class(sll_t_linear_operator_particle_mass_3d_od), intent( inout ) :: self
52 
53  deallocate( self%particle_mass )
54 
55  end subroutine free_particle_mass_3d_od
56 
57 
58  subroutine dot_particle_mass_3d_od ( self, x, y )
59  class(sll_t_linear_operator_particle_mass_3d_od), intent( in ) :: self
60  sll_real64, intent( in ) :: x(:)
61  sll_real64, intent( out ) :: y(:)
62  !local variables
63  sll_int32 :: ind_i, ind_j, ind
64  sll_int32 :: i1, i2, i3, j1, j2, j3
65 
66  ind_i = 1
67  do i3 = 1, self%n_dofs(3)
68  do i2 = 1, self%n_dofs(2)
69  do i1 = 1, self%n_dofs(1)
70  y(ind_i) = 0.0_f64
71  ind = 1
72  do j3 = i3-self%degree2(3), i3+self%degree(3)
73  do j2 = i2-self%degree2(2), i2+self%degree(2)
74  do j1 = i1-self%degree2(1), i1+self%degree(1)
75  ind_j = index_3dto1d( self%n_dofs, j1, j2, j3)
76  y(ind_i) = y(ind_i) + &
77  self%sign * self%particle_mass( ind, ind_i) * x(ind_j)
78  ind = ind+1
79  end do
80  end do
81  end do
82  ind_i = ind_i + 1
83  end do
84  end do
85  end do
86 
87 
88  end subroutine dot_particle_mass_3d_od
89 
90 
91  function index_3dto1d( num_pts, ind1, ind2, ind3 ) result( ind1d)
92  sll_int32 :: num_pts(3)
93  sll_int32 :: ind1
94  sll_int32 :: ind2
95  sll_int32 :: ind3
96  sll_int32 :: ind1d
97 
98  ind1d = 1+ modulo(ind1-1,num_pts(1)) + &
99  modulo(ind2-1,num_pts(2)) * num_pts(1) + &
100  modulo(ind3-1,num_pts(3))* num_pts(1)*num_pts(2)
101 
102 
103  end function index_3dto1d
104 
105 
107  class(sll_t_linear_operator_particle_mass_3d_od), intent(in) :: self
108  end subroutine print_info_particle_mass_3d_od
109 
integer(kind=i32) function index_3dto1d(num_pts, ind1, ind2, ind3)
subroutine create_linear_operator_particle_mass_3d_od(self, n_cells, degree, degree2)
    Report Typos and Errors