Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_linear_operator_schur_eb_cl_1d.F90
Go to the documentation of this file.
2 #include "sll_working_precision.h"
3 
5 
6  use sll_m_matrix_csr, only: &
8 
12 
13 
14  implicit none
15 
17 
18  private
20 
21  type(sll_t_matrix_csr), pointer :: mass0
22  type(sll_t_matrix_csr), pointer :: mass1
23  sll_int32 :: n_dofs
24  sll_real64 :: delta_x
25  sll_real64 :: sign = 1.0_f64
26  sll_int32 :: s_deg_0
27 
28  contains
30  procedure :: free => free_schur_eb_cl_1d
31  procedure :: dot => dot_schur_eb_cl_1d
32  procedure :: print_info => print_info_schur_eb_cl_1d
33 
35 
36 
37 contains
38 
39  subroutine create_linear_operator_schur_eb_cl_1d( self, mass0, mass1, n_dofs, delta_x, s_deg_0 )
40  class(sll_t_linear_operator_schur_eb_cl_1d), intent( inout ) :: self
41  type(sll_t_matrix_csr), target :: mass0
42  type(sll_t_matrix_csr), target :: mass1
43  sll_int32 :: n_dofs
44  sll_real64 :: delta_x
45  sll_int32 :: s_deg_0
46 
47  self%mass0 => mass0
48  self%mass1 => mass1
49  self%n_dofs = n_dofs
50  self%delta_x = delta_x
51  self%s_deg_0 = s_deg_0
52 
53  self%n_rows = self%n_dofs
54  self%n_cols = self%n_dofs
55 
56  self%n_global_rows = self%n_rows
57  self%n_global_cols = self%n_cols
58 
60 
61  subroutine free_schur_eb_cl_1d( self )
62  class(sll_t_linear_operator_schur_eb_cl_1d), intent( inout ) :: self
63 
64  self%mass0=> null()
65  self%mass1=> null()
66 
67  end subroutine free_schur_eb_cl_1d
68 
69 
70  subroutine dot_schur_eb_cl_1d ( self, x, y )
71  class(sll_t_linear_operator_schur_eb_cl_1d), intent( in ) :: self
72  sll_real64, intent( in ) :: x(:)
73  sll_real64, intent( out ) :: y(:)
74 
75  !local variables
76  sll_real64 :: scratch0(self%n_dofs), scratch1(self%n_dofs-1), scratch2(self%n_dofs-1)
77 
78  ! Compute D x
79  call sll_s_multiply_g_clamped_1d( self%n_dofs, self%delta_x, self%s_deg_0, x, scratch1)
80 
81  ! Compute M1 D x
82  call self%mass1%dot( scratch1, scratch2 )
83 
84  ! Compute D^T M1 D x + M_b D x
85  call sll_s_multiply_gt_clamped_1d( self%n_dofs, self%delta_x, self%s_deg_0, scratch2, scratch0 )
86 
87  !scratch0(1)= scratch0(1)+scratch1(1)
88  !scratch0(self%n_dofs)= scratch0(self%n_dofs)-scratch1(self%n_dofs-1)
89 
90  ! Compute M0 x
91  call self%mass0%dot( x, y )
92 
93  ! Sum up the two parts
94  y = y + self%sign * scratch0
95 
96  end subroutine dot_schur_eb_cl_1d
97 
98  subroutine print_info_schur_eb_cl_1d( self )
99  class(sll_t_linear_operator_schur_eb_cl_1d), intent(in) :: self
100 
101  end subroutine print_info_schur_eb_cl_1d
102 
103 
module for abstract linear operator
subroutine create_linear_operator_schur_eb_cl_1d(self, mass0, mass1, n_dofs, delta_x, s_deg_0)
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