Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_linear_operator_curl_cl_3d.F90
Go to the documentation of this file.
2 #include "sll_working_precision.h"
3 
5 
11 
12  use sll_m_linear_operator_block, only : &
14 
15  implicit none
16 
18 
19  private
21  type(sll_t_linear_operator_block), pointer :: mass1
22  type(sll_t_linear_operator_block), pointer :: mass2
23  sll_int32 :: n_total0, n_total1
24  sll_int32 :: n_dofs(3)
25  sll_int32 :: s_deg_0(3)
26  sll_real64 :: delta_x(3)
27  sll_real64 :: epsilon = 1._f64
28 
29 
30  contains
31  procedure :: create => create_linear_operator_curl_cl_3d
32  procedure :: free => free_curl_cl_3d
33  procedure :: dot => dot_curl_cl_3d
34  procedure :: print_info => print_info_curl_cl_3d
35 
37 
38 
39 contains
40 
41  subroutine create_linear_operator_curl_cl_3d( self, mass1, mass2, n_dofs, delta_x, s_deg_0 )
42  class(sll_t_linear_operator_curl_cl_3d), intent( inout ) :: self
43  type(sll_t_linear_operator_block), target :: mass1
44  type(sll_t_linear_operator_block), target :: mass2
45  sll_int32 :: n_dofs(3)
46  sll_real64 :: delta_x(3)
47  sll_int32 :: s_deg_0(3)
48 
49  self%mass1 => mass1
50  self%mass2 => mass2
51 
52  self%n_dofs = n_dofs
53  self%delta_x = delta_x
54  self%s_deg_0 = s_deg_0
55  self%n_total0 = product(self%n_dofs)
56  self%n_total1 = (self%n_dofs(1)-1)*self%n_dofs(2)*self%n_dofs(3)
57 
58  self%n_rows = self%n_total1+2*self%n_total0
59  self%n_cols = self%n_total1+2*self%n_total0
60 
61  self%n_global_rows = self%n_rows
62  self%n_global_cols = self%n_cols
63 
65 
66  subroutine free_curl_cl_3d( self )
67  class(sll_t_linear_operator_curl_cl_3d), intent( inout ) :: self
68 
69  self%mass1=> null()
70  self%mass2=> null()
71 
72  end subroutine free_curl_cl_3d
73 
74 
75  subroutine dot_curl_cl_3d ( self, x, y )
76  class(sll_t_linear_operator_curl_cl_3d), intent( in ) :: self
77  sll_real64, intent( in ) :: x(:)
78  sll_real64, intent( out ) :: y(:)
79  !local variables
80  sll_real64 :: scratch0(self%n_total1+2*self%n_total0), scratch1(self%n_total1+2*self%n_total0), scratch2(self%n_total0+2*self%n_total1), scratch3(self%n_total0+2*self%n_total1)
81 
82  ! Compute C x
83  call sll_s_multiply_c_clamped( self%n_dofs, self%delta_x, self%s_deg_0, x, scratch2 )
84 
85  ! Compute M2 C x
86  call self%mass2%dot( scratch2, scratch3 )
87 
88  ! Compute C^T M2 C x
89  call sll_s_multiply_ct_clamped( self%n_dofs, self%delta_x, self%s_deg_0, scratch3, y )
90 !!$
91  ! Compute M1 x
92  call self%mass1%dot( x, scratch1 )
93 
94  ! Compute G^T M1 x
95  call sll_s_multiply_gt_clamped( self%n_dofs, self%delta_x, self%s_deg_0, scratch1, scratch2(1:self%n_total0) )
96 
97  ! Compute G G^T M1 x
98  call sll_s_multiply_g_clamped( self%n_dofs, self%delta_x, self%s_deg_0, scratch2(1:self%n_total0), scratch1 )
99 
100  ! Compute M1 G G^T M_1 x
101  call self%mass1%dot( scratch1, scratch0 )
102 
103  ! y = (C^T M2 C + epsilon M1 G G^T M_1) x
104  y = y + self%epsilon*scratch0
105 
106  end subroutine dot_curl_cl_3d
107 
108 
109  subroutine print_info_curl_cl_3d( self )
110  class(sll_t_linear_operator_curl_cl_3d), intent(in) :: self
111 
112  end subroutine print_info_curl_cl_3d
113 
114 
module for abstract linear operator
module for a block linear operator
subroutine create_linear_operator_curl_cl_3d(self, mass1, mass2, n_dofs, delta_x, s_deg_0)
Utilites for Maxwell solver's with spline finite elements using sparse matrices.
subroutine, public sll_s_multiply_ct_clamped(n_dofs, delta_x, s_deg_0, field_in, field_out)
Multiply by transpose of discrete curl matrix.
subroutine, public sll_s_multiply_g_clamped(n_dofs, delta_x, s_deg_0, field_in, field_out)
Multiply by dicrete gradient matrix.
subroutine, public sll_s_multiply_c_clamped(n_dofs, delta_x, s_deg_0, field_in, field_out)
Multiply by discrete curl matrix.
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