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_phiv_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 
14  implicit none
15 
17 
18  private
20 
21  class(sll_c_maxwell_1d_base), pointer :: maxwell_solver
22  class(sll_c_particle_mass_1d_base), pointer :: particle_mass
23  sll_int32 :: n_dofs0, n_dofs1
24  sll_real64 :: sign = 1.0_f64
25  sll_int32 :: degree
26 
27  contains
28  procedure :: create => create_linear_operator_schur_phiv_1d
29  procedure :: free => free_schur_phiv_1d
30  procedure :: dot => dot_schur_phiv_1d
31  procedure :: print_info => print_info_schur_phiv_1d
32 
34 
35 
36 contains
37 
38  subroutine create_linear_operator_schur_phiv_1d( self, maxwell_solver, particle_mass )
39  class(sll_t_linear_operator_schur_phiv_1d), intent( inout ) :: self
40  class(sll_c_maxwell_1d_base), target :: maxwell_solver
41  class( sll_c_particle_mass_1d_base), target :: particle_mass
42 
43  self%particle_mass => particle_mass
44  self%maxwell_solver => maxwell_solver
45 
46  self%degree = self%maxwell_solver%s_deg_1
47  self%n_dofs0 = self%maxwell_solver%n_dofs0
48  self%n_dofs1 = self%maxwell_solver%n_dofs1
49 
50  self%n_rows = self%n_dofs0
51  self%n_cols = self%n_dofs0
52 
53  self%n_global_rows = self%n_rows
54  self%n_global_cols = self%n_cols
55 
57 
58  subroutine free_schur_phiv_1d( self )
59  class(sll_t_linear_operator_schur_phiv_1d), intent( inout ) :: self
60 
61  self%maxwell_solver => null()
62  self%particle_mass => null()
63 
64  end subroutine free_schur_phiv_1d
65 
66 
67  subroutine dot_schur_phiv_1d ( self, x, y )
68  class(sll_t_linear_operator_schur_phiv_1d), intent( in ) :: self
69  sll_real64, intent( in ) :: x(:)
70  sll_real64, intent( out ) :: y(:)
71  !local variables
72  sll_real64 :: g(self%n_dofs1), h(self%n_dofs1), z(self%n_dofs0)
73 
74  !M_0 phi+ dt^2/4 q^2/m G^T M^\star G phi
75  call self%maxwell_solver%multiply_g( x, g )
76  call self%particle_mass%dot( g, h )
77  call self%maxwell_solver%multiply_gt( h, z )
78  call self%maxwell_solver%multiply_mass( x, y, self%maxwell_solver%s_deg_0 )
79  y = y + self%sign*z
80 
81  end subroutine dot_schur_phiv_1d
82 
83  subroutine print_info_schur_phiv_1d( self )
84  class(sll_t_linear_operator_schur_phiv_1d), intent(in) :: self
85 
86  end subroutine print_info_schur_phiv_1d
87 
module for abstract linear operator
subroutine create_linear_operator_schur_phiv_1d(self, maxwell_solver, particle_mass)
Module interface to solve Maxwell's equations in 1D.
    Report Typos and Errors