Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_linear_operator_GTM_cl.F90
Go to the documentation of this file.
2 #include "sll_working_precision.h"
3 
6 
7  use sll_m_linear_operator_block, only : &
9 
12 
13  implicit none
14 
16 
17  private
18 
20  type(sll_t_linear_operator_block), pointer :: mass
21  sll_int32 :: n_dofs(3)
22  sll_real64 :: delta_x(3)
23  sll_int32 :: s_deg_0(3)
24  sll_int32 :: n_total0, n_total1
25 
26 
27  contains
28  procedure :: create => create_linear_operator_gtm_cl
29  procedure :: free => free_gtm_cl
30  procedure :: dot => dot_gtm_cl
31  procedure :: print_info => print_info_gtm_cl
32 
34 
35 contains
36 
37 
38  subroutine create_linear_operator_gtm_cl( self, mass, n_dofs, delta_x, s_deg_0 )
39  class(sll_t_linear_operator_gtm_cl), intent( inout ) :: self
40  type(sll_t_linear_operator_block), target :: mass
41  sll_int32 :: n_dofs(3)
42  sll_real64 :: delta_x(3)
43  sll_int32 :: s_deg_0(3)
44 
45  self%mass => mass
46  self%n_dofs = n_dofs
47  self%s_deg_0 = s_deg_0
48  self%n_total0 = product(self%n_dofs)
49  self%n_total1 = (self%n_dofs(1)-1)*self%n_dofs(2)*self%n_dofs(3)
50  self%delta_x = delta_x
51 
52  self%n_rows = self%n_total0
53  self%n_cols = self%n_total1+2*self%n_total0
54 
55  self%n_global_rows = self%n_rows
56  self%n_global_cols = self%n_cols
57  end subroutine create_linear_operator_gtm_cl
58 
59 
60  subroutine free_gtm_cl( self )
61  class(sll_t_linear_operator_gtm_cl), intent( inout ) :: self
62  self%mass => null()
63 
64  end subroutine free_gtm_cl
65 
66 
67  subroutine dot_gtm_cl( self, x, y )
68  class(sll_t_linear_operator_gtm_cl), intent( in ) :: self
69  sll_real64, intent( in ) :: x(:)
70  sll_real64, intent( out ) :: y(:)
71  !local variable
72  sll_real64 :: scratch(self%n_total1+2*self%n_total0)
73 
74  call self%mass%dot(x, scratch)
75  call sll_s_multiply_gt_clamped(self%n_dofs, self%delta_x, self%s_deg_0, scratch, y)
76 
77  end subroutine dot_gtm_cl
78 
79  subroutine print_info_gtm_cl( self )
80  class(sll_t_linear_operator_gtm_cl), intent(in) :: self
81  end subroutine print_info_gtm_cl
82 
83 
module for abstract linear operator
module for a block linear operator
subroutine create_linear_operator_gtm_cl(self, mass, n_dofs, delta_x, s_deg_0)
Utilites for Maxwell solver's with spline finite elements using sparse matrices.
subroutine, public sll_s_multiply_gt_clamped(n_dofs, delta_x, s_deg_0, field_in, field_out)
Multiply by transpose of dicrete gradient matrix.
    Report Typos and Errors