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_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_cl_3d_diag
18  procedure :: dot => dot_particle_mass_cl_3d_diag
19  procedure :: print_info => print_info_particle_mass_cl_3d_diag
20 
22 
23 
24 contains
25 
26  subroutine create_linear_operator_particle_mass_cl_3d_diag( self, n_cells, degree, degree2 )
27  class(sll_t_linear_operator_particle_mass_cl_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_dofs(1) = n_cells(1)+degree(1)
35  self%n_total = product(self%n_dofs)
36 
37  allocate( self%particle_mass( product(2*self%degree+1), self%n_total) )
38  self%particle_mass = 0._f64
39  self%n_rows = self%n_total
40  self%n_cols = self%n_total
41 
42  self%n_global_rows = self%n_rows
43  self%n_global_cols = self%n_cols
44 
46 
47  subroutine free_particle_mass_cl_3d_diag( self )
48  class(sll_t_linear_operator_particle_mass_cl_3d_diag), intent( inout ) :: self
49 
50  deallocate( self%particle_mass )
51 
52  end subroutine free_particle_mass_cl_3d_diag
53 
54 
55  subroutine dot_particle_mass_cl_3d_diag ( self, x, y )
56  class(sll_t_linear_operator_particle_mass_cl_3d_diag), intent( in ) :: self
57  sll_real64, intent( in ) :: x(:)
58  sll_real64, intent( out ) :: y(:)
59  !local variables
60  sll_int32 :: ind_i, ind_j, ind
61  sll_int32 :: i1, i2, i3, j1, j2, j3
62 
63  ind_i = 1
64  do i3 = 1, self%n_dofs(3)
65  do i2 = 1, self%n_dofs(2)
66  do i1 = 1, self%degree(1)
67  ind = 1
68  y(ind_i) = 0.0_f64
69  do j3 = i3-self%degree(3), i3+self%degree(3)
70  do j2 = i2-self%degree(2), i2+self%degree(2)
71  ind = ind + self%degree(1)+1-i1
72  do j1 = 1, i1+self%degree(1)
73  ind_j = j1 + modulo(j2-1,self%n_dofs(2))*self%n_dofs(1) + modulo(j3-1,self%n_dofs(3))*self%n_dofs(1)*self%n_dofs(2)
74  y(ind_i) = y(ind_i) + &
75  self%sign * self%particle_mass( ind, ind_i) * x(ind_j)
76  ind = ind+1
77  end do
78  end do
79  end do
80  ind_i = ind_i + 1
81  end do
82  do i1 = self%degree(1)+1, self%n_dofs(1)-self%degree(1)
83  y(ind_i) = 0.0_f64
84  ind = 1
85  do j3 = i3-self%degree(3), i3+self%degree(3)
86  do j2 = i2-self%degree(2), i2+self%degree(2)
87  do j1 = i1-self%degree(1), i1+self%degree(1)
88  ind_j = j1 + modulo(j2-1,self%n_dofs(2))*self%n_dofs(1) + modulo(j3-1,self%n_dofs(3))*self%n_dofs(1)*self%n_dofs(2)
89  y(ind_i) = y(ind_i) + &
90  self%sign * self%particle_mass( ind, ind_i) * x(ind_j)
91  ind = ind+1
92  end do
93  end do
94  end do
95  ind_i = ind_i + 1
96  end do
97  do i1 = self%n_dofs(1)-self%degree(1)+1, self%n_dofs(1)
98  y(ind_i) = 0.0_f64
99  ind = 1
100  do j3 = i3-self%degree(3), i3+self%degree(3)
101  do j2 = i2-self%degree(2), i2+self%degree(2)
102  do j1 = i1-self%degree(1), self%n_dofs(1)
103  ind_j = j1 + modulo(j2-1,self%n_dofs(2))*self%n_dofs(1) + modulo(j3-1,self%n_dofs(3))*self%n_dofs(1)*self%n_dofs(2)
104  y(ind_i) = y(ind_i) + &
105  self%sign * self%particle_mass( ind, ind_i) * x(ind_j)
106  ind = ind +1
107  end do
108  ind = ind + self%degree(1)-self%n_dofs(1)+i1
109  end do
110  end do
111  ind_i = ind_i + 1
112  end do
113  end do
114  end do
115 
116 
117  end subroutine dot_particle_mass_cl_3d_diag
118 
119 
121  class(sll_t_linear_operator_particle_mass_cl_3d_diag), intent(in) :: self
123 
subroutine create_linear_operator_particle_mass_cl_3d_diag(self, n_cells, degree, degree2)
    Report Typos and Errors