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_3d.F90
Go to the documentation of this file.
2 #include "sll_working_precision.h"
3 
5 
6  use sll_m_spline_fem_utilities, only : &
9 
10  use sll_m_linear_operator_block, only : &
12 
13  implicit none
14 
16 
17  private
19  type(sll_t_linear_operator_block), pointer :: mass1
20  type(sll_t_linear_operator_block), pointer :: mass2
21  sll_int32 :: n_dofs(3)
22  sll_real64 :: delta_x(3)
23  sll_int32 :: n_total
24  sll_real64 :: sign = 1.0_f64
25 
26  contains
27  procedure :: create => create_linear_operator_schur_eb_3d
28  procedure :: free => free_schur_eb_3d
29  procedure :: dot => dot_schur_eb_3d
30  procedure :: print_info => print_info_schur_eb_3d
31 
33 
34 
35 contains
36 
37  subroutine create_linear_operator_schur_eb_3d( self, mass1, mass2, n_total, n_dofs, delta_x )
38  class(sll_t_linear_operator_schur_eb_3d), intent( inout ) :: self
39  type(sll_t_linear_operator_block), target :: mass1
40  type(sll_t_linear_operator_block), target :: mass2
41  sll_int32 :: n_total
42  sll_int32 :: n_dofs(3)
43  sll_real64 :: delta_x(3)
44 
45  self%mass1 => mass1
46  self%mass2 => mass2
47 
48  self%n_total = n_total
49  self%n_dofs = n_dofs
50  self%delta_x = delta_x
51 
52  self%n_rows = 3*self%n_total
53  self%n_cols = 3*self%n_total
54 
55  self%n_global_rows = self%n_rows
56  self%n_global_cols = self%n_cols
57 
59 
60  subroutine free_schur_eb_3d( self )
61  class(sll_t_linear_operator_schur_eb_3d), intent( inout ) :: self
62 
63  self%mass1=> null()
64  self%mass2=> null()
65 
66  end subroutine free_schur_eb_3d
67 
68 
69  subroutine dot_schur_eb_3d ( self, x, y )
70  class(sll_t_linear_operator_schur_eb_3d), intent( in ) :: self
71  sll_real64, intent( in ) :: x(:)
72  sll_real64, intent( out ) :: y(:)
73  !local variables
74  sll_real64 :: scratch1(3*self%n_total), scratch2(3*self%n_total)
75 
76  ! Compute C x
77  call sll_s_multiply_c( self%n_dofs, self%delta_x, x, scratch1 )
78 
79  ! Compute M2 C x
80  call self%mass2%dot( scratch1, scratch2 )
81 
82  ! Compute C^T M2 C x
83  call sll_s_multiply_ct( self%n_dofs, self%delta_x, scratch2, scratch1 )
84 
85  ! Compute M1 x
86  call self%mass1%dot( x, y )
87 
88  ! Sum up the two parts
89  y = y + self%sign * scratch1
90 
91  end subroutine dot_schur_eb_3d
92 
93 
94  subroutine print_info_schur_eb_3d( self )
95  class(sll_t_linear_operator_schur_eb_3d), intent(in) :: self
96 
97  end subroutine print_info_schur_eb_3d
98 
99 
module for abstract linear operator
module for a block linear operator
subroutine create_linear_operator_schur_eb_3d(self, mass1, mass2, n_total, n_dofs, delta_x)
Utilites for Maxwell solver's with spline finite elements.
subroutine, public sll_s_multiply_c(n_dofs, delta_x, field_in, field_out)
Multiply by discrete curl matrix.
subroutine, public sll_s_multiply_ct(n_dofs, delta_x, field_in, field_out)
Multiply by transpose of discrete curl matrix.
    Report Typos and Errors