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_3d.F90
Go to the documentation of this file.
2 #include "sll_working_precision.h"
4 
5  use sll_m_maxwell_3d_base, only: &
7 
8  use sll_m_linear_operator_block, only : &
10 
11  implicit none
12 
14 
15  private
17 
18  class(sll_c_maxwell_3d_base), pointer :: maxwell_solver
19  type(sll_t_linear_operator_block), pointer :: particle_mass
20  sll_int32 :: n_total0, n_total1
21  sll_real64 :: sign = 1.0_f64
22 
23  contains
24  procedure :: create => create_linear_operator_schur_phiv_3d
25  procedure :: free => free_schur_phiv_3d
26  procedure :: dot => dot_schur_phiv_3d
27  procedure :: print_info => print_info_schur_phiv_3d
28 
30 
31 
32 contains
33 
34  subroutine create_linear_operator_schur_phiv_3d( self, maxwell_solver, particle_mass)
35  class(sll_t_linear_operator_schur_phiv_3d), intent( inout ) :: self
36  class(sll_c_maxwell_3d_base), target :: maxwell_solver
37  type(sll_t_linear_operator_block), target :: particle_mass
38 
39 
40  self%particle_mass => particle_mass
41  self%maxwell_solver => maxwell_solver
42 
43  self%n_total0 = maxwell_solver%n_total0
44  self%n_total1 = maxwell_solver%n_total1
45 
46  self%n_rows = self%n_total0
47  self%n_cols = self%n_total0
48 
49  self%n_global_rows = self%n_rows
50  self%n_global_cols = self%n_cols
51 
53 
54  subroutine free_schur_phiv_3d( self )
55  class(sll_t_linear_operator_schur_phiv_3d), intent( inout ) :: self
56 
57  self%maxwell_solver => null()
58  self%particle_mass => null()
59 
60  end subroutine free_schur_phiv_3d
61 
62 
63  subroutine dot_schur_phiv_3d ( self, x, y )
64  class(sll_t_linear_operator_schur_phiv_3d), intent( in ) :: self
65  sll_real64, intent( in ) :: x(:)
66  sll_real64, intent( out ) :: y(:)
67 
68  !local variables
69  sll_real64 :: f(self%n_total1+2*self%n_total0), g(self%n_total1+2*self%n_total0), z(self%n_total0)
70  !M_0 phi+ dt^2/4 q^2/m G^T M^\star G phi
71  call self%maxwell_solver%multiply_g( x, f )
72  call self%particle_mass%dot( f, g )
73  call self%maxwell_solver%multiply_gt( g, z )
74  call self%maxwell_solver%multiply_mass([0], x, y)
75  y= y + self%sign*z
76 
77  end subroutine dot_schur_phiv_3d
78 
79  subroutine print_info_schur_phiv_3d( self )
80  class(sll_t_linear_operator_schur_phiv_3d), intent(in) :: self
81 
82  end subroutine print_info_schur_phiv_3d
83 
84 
module for abstract linear operator
module for a block linear operator
subroutine create_linear_operator_schur_phiv_3d(self, maxwell_solver, particle_mass)
Module interface to solve Maxwell's equations in 3D.
    Report Typos and Errors