Report Typos and Errors
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
particle_methods
pic_basic
pic_time_integration
sll_m_linear_operator_particle_mass_1d.F90
Go to the documentation of this file.
1
module
sll_m_linear_operator_particle_mass_1d
2
#include "sll_working_precision.h"
3
4
use
sll_m_particle_mass_1d_base
,
only
: &
5
sll_c_particle_mass_1d_base
6
7
8
implicit none
9
10
public
::
sll_t_linear_operator_particle_mass_1d
11
12
private
13
14
type
,
extends
(
sll_c_particle_mass_1d_base
) ::
sll_t_linear_operator_particle_mass_1d
15
16
17
contains
18
procedure
:: create =>
create_linear_operator_particle_mass_1d
19
procedure
:: free =>
free_particle_mass_1d
20
procedure
:: dot =>
dot_particle_mass_1d
21
procedure
:: print_info =>
print_info_particle_mass_1d
22
23
end type
sll_t_linear_operator_particle_mass_1d
24
25
26
contains
27
28
subroutine
create_linear_operator_particle_mass_1d
( self, degree, n_dofs )
29
class
(
sll_t_linear_operator_particle_mass_1d
),
intent( inout )
:: self
30
sll_int32,
intent( in )
:: degree
31
sll_int32,
intent( in )
:: n_dofs
32
33
self%degree = degree
34
self%n_dofs = n_dofs
35
36
allocate
( self%particle_mass( 2*self%degree+1, self%n_dofs) )
37
self%particle_mass = 0._f64
38
self%n_rows = self%n_dofs
39
self%n_cols = self%n_dofs
40
41
self%n_global_rows = self%n_rows
42
self%n_global_cols = self%n_cols
43
44
end subroutine
create_linear_operator_particle_mass_1d
45
46
subroutine
free_particle_mass_1d
( self )
47
class
(
sll_t_linear_operator_particle_mass_1d
),
intent( inout )
:: self
48
49
deallocate
( self%particle_mass )
50
51
end subroutine
free_particle_mass_1d
52
53
54
subroutine
dot_particle_mass_1d
( self, x, y )
55
class
(
sll_t_linear_operator_particle_mass_1d
),
intent( in )
:: self
56
sll_real64,
intent( in )
:: x(:)
57
sll_real64,
intent( out )
:: y(:)
58
!local variables
59
sll_int32 :: ind
60
sll_int32 :: i, j
61
62
do
i = 1, self%n_dofs
63
y(i) = 0.0_f64
64
ind = 1
65
do
j = i-self%degree, i+self%degree
66
y(i) = y(i) + self%sign * self%particle_mass( ind, i) * x(1 + modulo(j-1,self%n_dofs))
67
ind = ind+1
68
end do
69
end do
70
71
72
73
end subroutine
dot_particle_mass_1d
74
75
76
subroutine
print_info_particle_mass_1d
( self )
77
class
(
sll_t_linear_operator_particle_mass_1d
),
intent(in)
:: self
78
end subroutine
print_info_particle_mass_1d
79
80
end module
sll_m_linear_operator_particle_mass_1d
sll_m_linear_operator_particle_mass_1d
Definition:
sll_m_linear_operator_particle_mass_1d.F90:1
sll_m_linear_operator_particle_mass_1d::free_particle_mass_1d
subroutine free_particle_mass_1d(self)
Definition:
sll_m_linear_operator_particle_mass_1d.F90:47
sll_m_linear_operator_particle_mass_1d::create_linear_operator_particle_mass_1d
subroutine create_linear_operator_particle_mass_1d(self, degree, n_dofs)
Definition:
sll_m_linear_operator_particle_mass_1d.F90:29
sll_m_linear_operator_particle_mass_1d::print_info_particle_mass_1d
subroutine print_info_particle_mass_1d(self)
Definition:
sll_m_linear_operator_particle_mass_1d.F90:77
sll_m_linear_operator_particle_mass_1d::dot_particle_mass_1d
subroutine dot_particle_mass_1d(self, x, y)
Definition:
sll_m_linear_operator_particle_mass_1d.F90:55
sll_m_particle_mass_1d_base
Definition:
sll_m_particle_mass_1d_base.F90:1
sll_m_linear_operator_particle_mass_1d::sll_t_linear_operator_particle_mass_1d
Definition:
sll_m_linear_operator_particle_mass_1d.F90:14
sll_m_particle_mass_1d_base::sll_c_particle_mass_1d_base
Definition:
sll_m_particle_mass_1d_base.F90:13
Report Typos and Errors
Generated on Mon Oct 23 2023 19:15:40 for Semi-Lagrangian Library by
1.9.1