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_cl_1d.F90
Go to the documentation of this file.
2 #include "sll_working_precision.h"
3 
4  use sll_m_particle_mass_1d_base, only: &
6 
7  implicit none
8 
10 
11  private
12 
14 
15  contains
17  procedure :: free => free_particle_mass_cl_1d
18  procedure :: dot => dot_particle_mass_cl_1d
19  procedure :: print_info => print_info_particle_mass_cl_1d
20 
22 
23 
24 contains
25 
26  subroutine create_linear_operator_particle_mass_cl_1d( self, degree, n_dofs )
27  class(sll_t_linear_operator_particle_mass_cl_1d), intent( inout ) :: self
28  sll_int32, intent( in ) :: degree
29  sll_int32, intent( in ) :: n_dofs
30 
31  self%degree = degree
32  self%n_dofs = n_dofs
33 
34  allocate( self%particle_mass( 2*self%degree+1, self%n_dofs) )
35  self%particle_mass = 0._f64
36  self%n_rows = self%n_dofs
37  self%n_cols = self%n_dofs
38 
39  self%n_global_rows = self%n_rows
40  self%n_global_cols = self%n_cols
41 
43 
44  subroutine free_particle_mass_cl_1d( self )
45  class(sll_t_linear_operator_particle_mass_cl_1d), intent( inout ) :: self
46 
47  deallocate( self%particle_mass )
48 
49  end subroutine free_particle_mass_cl_1d
50 
51 
52  subroutine dot_particle_mass_cl_1d ( self, x, y )
53  class(sll_t_linear_operator_particle_mass_cl_1d), intent( in ) :: self
54  sll_real64, intent( in ) :: x(:)
55  sll_real64, intent( out ) :: y(:)
56  !local variables
57  sll_int32 :: ind
58  sll_int32 :: i, j
59 
60  do i = 1, self%degree
61  y(i) = 0.0_f64
62  ind = 1
63  do j = 1, i+self%degree
64  y(i) = y(i) + self%sign * self%particle_mass( ind, i) * x(j)
65  ind = ind+1
66  end do
67  end do
68 
69  do i= self%degree+1, self%n_dofs-self%degree
70  y(i) = 0.0_f64
71  ind = 1
72  do j = i-self%degree, i+self%degree
73  y(i) = y(i) + self%sign * self%particle_mass( ind, i) * x(j)
74  ind = ind+1
75  end do
76  end do
77 
78  do i= self%n_dofs-self%degree+1, self%n_dofs
79  y(i) = 0.0_f64
80  ind = 1
81  do j = i-self%degree, self%n_dofs
82  y(i) = y(i) + self%sign * self%particle_mass( ind, i) * x(j)
83  ind = ind+1
84  end do
85  end do
86 
87  end subroutine dot_particle_mass_cl_1d
88 
89  subroutine print_info_particle_mass_cl_1d( self )
90  class(sll_t_linear_operator_particle_mass_cl_1d), intent(in) :: self
91  end subroutine print_info_particle_mass_cl_1d
92 
    Report Typos and Errors