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_diag.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  implicit none
8 
10 
11  private
12 
14 
15  contains
17  procedure :: free => free_particle_mass_3d_diag
18  procedure :: dot => dot_particle_mass_3d_diag
19  procedure :: print_info => print_info_particle_mass_3d_diag
20 
22 
23 
24 contains
25 
26  subroutine create_linear_operator_particle_mass_3d_diag( self, n_cells, degree, degree2 )
27  class(sll_t_linear_operator_particle_mass_3d_diag), intent( inout ) :: self
28  sll_int32, intent( in ) :: n_cells(3)
29  sll_int32, intent( in ) :: degree(3)
30  sll_int32, optional, intent( in ) :: degree2(3)
31 
32  self%degree = degree
33  self%n_dofs = n_cells
34  self%n_total = product(n_cells)
35 
36  allocate( self%particle_mass( product(2*self%degree+1), self%n_total) )
37  self%particle_mass = 0._f64
38  self%n_rows = self%n_total
39  self%n_cols = self%n_total
40 
41  self%n_global_rows = self%n_rows
42  self%n_global_cols = self%n_cols
43 
45 
46  subroutine free_particle_mass_3d_diag( self )
47  class(sll_t_linear_operator_particle_mass_3d_diag), intent( inout ) :: self
48 
49  deallocate( self%particle_mass )
50 
51  end subroutine free_particle_mass_3d_diag
52 
53 
54  subroutine dot_particle_mass_3d_diag ( self, x, y )
55  class(sll_t_linear_operator_particle_mass_3d_diag), intent( in ) :: self
56  sll_real64, intent( in ) :: x(:)
57  sll_real64, intent( out ) :: y(:)
58  !local variables
59  sll_int32 :: ind_i, ind_j, ind
60  sll_int32 :: i1, i2, i3, j1, j2, j3
61 
62  ind_i = 1
63  do i3 = 1, self%n_dofs(3)
64  do i2 = 1, self%n_dofs(2)
65  do i1 = 1, self%n_dofs(1)
66  y(ind_i) = 0.0_f64
67  ind = 1
68  do j3 = i3-self%degree(3), i3+self%degree(3)
69  do j2 = i2-self%degree(2), i2+self%degree(2)
70  do j1 = i1-self%degree(1), i1+self%degree(1)
71  ind_j = index_3dto1d( self%n_dofs, j1, j2, j3)
72  y(ind_i) = y(ind_i) + &
73  self%sign * self%particle_mass( ind, ind_i) * x(ind_j)
74  ind = ind+1
75  end do
76  end do
77  end do
78  ind_i = ind_i + 1
79  end do
80  end do
81  end do
82 
83 
84  end subroutine dot_particle_mass_3d_diag
85 
86 
87  function index_3dto1d( num_pts, ind1, ind2, ind3 ) result( ind1d)
88  sll_int32 :: num_pts(3)
89  sll_int32 :: ind1
90  sll_int32 :: ind2
91  sll_int32 :: ind3
92  sll_int32 :: ind1d
93 
94  ind1d = 1+ modulo(ind1-1,num_pts(1)) + &
95  modulo(ind2-1,num_pts(2)) * num_pts(1) + &
96  modulo(ind3-1,num_pts(3))* num_pts(1)*num_pts(2)
97 
98 
99  end function index_3dto1d
100 
101 
103  class(sll_t_linear_operator_particle_mass_3d_diag), intent(in) :: self
104  end subroutine print_info_particle_mass_3d_diag
105 
subroutine create_linear_operator_particle_mass_3d_diag(self, n_cells, degree, degree2)
integer(kind=i32) function index_3dto1d(num_pts, ind1, ind2, ind3)
    Report Typos and Errors