Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_linear_operator_poisson_1d.F90
Go to the documentation of this file.
2 #include "sll_working_precision.h"
3 
6 
7  use sll_m_matrix_csr, only: &
9 
10  use sll_m_spline_fem_utilities, only : &
13 
14  implicit none
15 
17 
18  private
19 
21  type(sll_t_matrix_csr), pointer :: mass
22  sll_int32 :: n_dofs
23  sll_real64 :: delta_x
24 
25  contains
26  procedure :: create => create_linear_operator_poisson_1d
27  procedure :: dot => dot_poisson_1d
28  procedure :: print_info => print_info_poisson_1d
29  procedure :: free => free_poisson_1d
30 
32 
33 
34 contains
35 
36 
37  subroutine create_linear_operator_poisson_1d( self, mass, n_dofs, delta_x)
38  class(sll_t_linear_operator_poisson_1d), intent( inout ) :: self
39  type(sll_t_matrix_csr), target :: mass
40  sll_int32, intent(in) :: n_dofs
41  sll_real64, intent(in) :: delta_x
42 
43  self%mass => mass
44  self%n_dofs = n_dofs
45  self%delta_x = delta_x
46 
47  self%n_rows = self%n_dofs
48  self%n_cols = self%n_dofs
49 
50  self%n_global_rows = self%n_rows
51  self%n_global_cols = self%n_cols
52 
54 
55  subroutine dot_poisson_1d( self, x, y )
56  class(sll_t_linear_operator_poisson_1d), intent( in ) :: self
57  sll_real64, intent( in ) :: x(:)
58  sll_real64, intent( out ) :: y(:)
59  !local variable
60  sll_real64 :: scratch(self%n_dofs), scratch1(self%n_dofs)
61 
62  call sll_s_multiply_g_1d( self%n_dofs, self%delta_x, x, scratch )
63  call self%mass%dot( scratch, scratch1 )
64  call sll_s_multiply_gt_1d( self%n_dofs, self%delta_x, scratch1, y )
65 
66  end subroutine dot_poisson_1d
67 
68  subroutine print_info_poisson_1d( self )
69  class(sll_t_linear_operator_poisson_1d), intent(in) :: self
70  end subroutine print_info_poisson_1d
71 
72 
73  subroutine free_poisson_1d( self )
74  class(sll_t_linear_operator_poisson_1d), intent( inout ) :: self
75  self%mass => null()
76 
77  end subroutine free_poisson_1d
78 
module for abstract linear operator
subroutine create_linear_operator_poisson_1d(self, mass, n_dofs, delta_x)
module for Compressed Sparse Row Matrix (CSR)
Utilites for Maxwell solver's with spline finite elements.
subroutine, public sll_s_multiply_g_1d(n_dofs, delta_x, in, out)
Multiplication of the input vector in by the derivative matrix G.
subroutine, public sll_s_multiply_gt_1d(n_dofs, delta_x, in, out)
Multiplication of the input vector in by the transposed derivative matrix G^T
    Report Typos and Errors