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_ev_1d.F90
Go to the documentation of this file.
2 #include "sll_working_precision.h"
3 
6 
7  use sll_m_maxwell_1d_base, only: &
9 
10  use sll_m_particle_mass_1d_base, only : &
12 
13  implicit none
14 
16 
17  private
19 
20  class(sll_c_maxwell_1d_base), pointer :: maxwell_solver
21  class(sll_c_particle_mass_1d_base), pointer :: particle_mass
22  sll_int32 :: n_dofs
23  sll_real64 :: sign = 1.0_f64
24  sll_int32 :: degree
25 
26  contains
27  procedure :: create => create_linear_operator_schur_ev_1d
28  procedure :: free => free_schur_ev_1d
29  procedure :: dot => dot_schur_ev_1d
30  procedure :: print_info => print_info_schur_ev_1d
31 
33 
34 
35 contains
36 
37  subroutine create_linear_operator_schur_ev_1d( self, maxwell_solver, particle_mass, n_dofs, degree)
38  class(sll_t_linear_operator_schur_ev_1d), intent( inout ) :: self
39  class(sll_c_maxwell_1d_base), target :: maxwell_solver
40  class( sll_c_particle_mass_1d_base), target :: particle_mass
41  sll_int32, intent(in) :: n_dofs
42  sll_int32, intent(in) :: degree
43 
44  self%degree=degree
45  self%particle_mass => particle_mass
46  self%maxwell_solver => maxwell_solver
47 
48  self%n_dofs = n_dofs
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
55 
57 
58  subroutine free_schur_ev_1d( self )
59  class(sll_t_linear_operator_schur_ev_1d), intent( inout ) :: self
60 
61  self%maxwell_solver => null()
62  self%particle_mass => null()
63 
64  end subroutine free_schur_ev_1d
65 
66 
67  subroutine dot_schur_ev_1d ( self, x, y )
68  class(sll_t_linear_operator_schur_ev_1d), intent( in ) :: self
69  sll_real64, intent( in ) :: x(:)
70  sll_real64, intent( out ) :: y(:)
71  !local variables
72  sll_real64 :: z(self%n_dofs)
73 
74  !M_1 E+ dt^2/4 q^2/m M^\star E
75  call self%particle_mass%dot( x, z )
76  call self%maxwell_solver%multiply_mass( x, y, self%degree )
77  y= y + self%sign*z
78 
79  end subroutine dot_schur_ev_1d
80 
81  subroutine print_info_schur_ev_1d( self )
82  class(sll_t_linear_operator_schur_ev_1d), intent(in) :: self
83 
84  end subroutine print_info_schur_ev_1d
85 
86 
module for abstract linear operator
subroutine create_linear_operator_schur_ev_1d(self, maxwell_solver, particle_mass, n_dofs, degree)
Module interface to solve Maxwell's equations in 1D.
    Report Typos and Errors