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_clamped_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 
13 
14  implicit none
15 
17 
18  private
19 
21  type(sll_t_matrix_csr), pointer :: mass
22  sll_int32 :: s_deg_0
23  sll_int32 :: n_dofs
24  sll_real64 :: delta_x
25 
26  contains
28  procedure :: free => free_poisson_clamped_1d
29  procedure :: dot => dot_poisson_clamped_1d
30  procedure :: print_info => print_info_poisson_clamped_1d
31 
33 
34 
35 contains
36 
37 
38  subroutine create_linear_operator_poisson_clamped_1d( self, mass, s_deg_0, n_dofs, delta_x)
39  class(sll_t_linear_operator_poisson_clamped_1d), intent( inout ) :: self
40  type(sll_t_matrix_csr), target :: mass
41  sll_int32, intent(in) :: s_deg_0
42  sll_int32, intent(in) :: n_dofs
43  sll_real64, intent(in) :: delta_x
44 
45  self%mass => mass
46  self%s_deg_0 = s_deg_0
47  self%n_dofs = n_dofs
48  self%delta_x = delta_x
49 
50  self%n_rows = self%n_dofs
51  self%n_cols = self%n_dofs
52 
53  self%n_global_rows = self%n_rows
54  self%n_global_cols = self%n_cols
56 
57 
58  subroutine free_poisson_clamped_1d( self )
59  class(sll_t_linear_operator_poisson_clamped_1d), intent( inout ) :: self
60  self%mass => null()
61 
62  end subroutine free_poisson_clamped_1d
63 
64 
65  subroutine dot_poisson_clamped_1d( self, x, y )
66  class(sll_t_linear_operator_poisson_clamped_1d), intent( in ) :: self
67  sll_real64, intent( in ) :: x(:)
68  sll_real64, intent( out ) :: y(:)
69  !local variable
70  sll_real64 :: scratch1(self%n_dofs-1), scratch2(self%n_dofs-1)
71 
72  call sll_s_multiply_g_clamped_1d(self%n_dofs, self%delta_x, self%s_deg_0, x, scratch1)
73  call self%mass%dot(scratch1, scratch2)
74  call sll_s_multiply_gt_clamped_1d(self%n_dofs, self%delta_x, self%s_deg_0, scratch2, y)
75  end subroutine dot_poisson_clamped_1d
76 
77 
78  subroutine print_info_poisson_clamped_1d( self )
79  class(sll_t_linear_operator_poisson_clamped_1d), intent(in) :: self
80  end subroutine print_info_poisson_clamped_1d
81 
module for abstract linear operator
subroutine create_linear_operator_poisson_clamped_1d(self, mass, s_deg_0, n_dofs, delta_x)
module for Compressed Sparse Row Matrix (CSR)
Utilites for Maxwell solver's with spline finite elements using sparse matrices.
subroutine, public sll_s_multiply_gt_clamped_1d(n_dofs, delta_x, s_deg_0, in, out)
Multiplication of the input vector in by the transposed clamped derivative matrix G^T.
subroutine, public sll_s_multiply_g_clamped_1d(n_dofs, delta_x, s_deg_0, in, out)
Multiplication of the input vector in by the clamped derivative matrix G.
    Report Typos and Errors