Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_linear_operator_ecsim_eb.F90
Go to the documentation of this file.
2 #include "sll_working_precision.h"
4 
5  use sll_m_spline_fem_utilities, only : &
7 
8  implicit none
9 
11 
12  private
13 
15 
16  sll_int32 :: degree
17  sll_real64 :: dt
18  sll_real64 :: dx
19  sll_real64, allocatable :: mass_line_0(:)
20  sll_real64, allocatable :: mass_line_1(:)
21 
22  contains
23  procedure :: create => create_linear_operator_ecsim_eb
24  procedure :: free => free_ecsim_eb
25  procedure :: dot => dot_mono_r2r_ecsim_eb
26  procedure :: print_info => print_info_ecsim_eb
27 
29 
30 contains
31 
32  subroutine create_linear_operator_ecsim_eb( self, n_dof, degree, &
33  mass_line_0, mass_line_1)
34  class(sll_t_linear_operator_ecsim_eb), intent( inout ) :: self
35  sll_int32 :: n_dof
36  sll_int32 :: degree
37  sll_real64 :: mass_line_0(:)
38  sll_real64 :: mass_line_1(:)
39 
40 
41  self%n_dof=n_dof
42  self%degree=degree
43  self%dt=0.0_f64
44 
45  allocate(self%mass_line_0(degree+1))
46  allocate(self%mass_line_1(degree))
47  self%mass_line_0=mass_line_0
48  self%mass_line_1=mass_line_1
49 
50 
51  self%n_rows = 2*n_dof
52  self%n_cols = 2*n_dof
53 
54  self%n_global_rows = self%n_rows
55  self%n_global_cols = self%n_cols
56  end subroutine create_linear_operator_ecsim_eb
57 
58 
59  subroutine free_ecsim_eb( self )
60  class(sll_t_linear_operator_ecsim_eb), intent( inout ) :: self
61  deallocate(self%mass_line_0)
62  deallocate(self%mass_line_1)
63 
64  end subroutine free_ecsim_eb
65 
66  subroutine dot_mono_r2r_ecsim_eb( self, x, y )
67  class(sll_t_linear_operator_ecsim_eb), intent( in ) :: self
68  sll_real64, intent( in ) :: x(:)
69  sll_real64, intent( out ) :: y(:)
70  ! local variables
71  sll_int32 :: row, column
72  sll_real64 :: scratch(self%n_dof+1)
73  scratch=0.0_f64
74  y=0.0_f64
75 
76  do row=1, self%n_dof
77  ! Multiplication of the inputvector with entries of the mass and particle-mass matrices
78  y(row)= self%mass_line_0(1)*x(row)
79  scratch(row)=self%mass_line_1(1)*x(self%n_dof+row)
80  do column = 2, self%degree
81  scratch(row)=scratch(row)+&
82  self%mass_line_1(column) * &
83  (x(self%n_dof+ modulo(row+column-2,self%n_dof)+1) +&
84  x(self%n_dof+ modulo(row-column,self%n_dof)+1))
85  end do
86  do column = 2, self%degree+1
87  y(row) = y(row) +&
88  self%mass_line_0(column) * &
89  (x(modulo(row+column-2,self%n_dof)+1) +&
90  x(modulo(row-column,self%n_dof)+1))
91  end do
92  end do
93  ! Update of the outputvector with derivatives of the inputvector
94  scratch(self%n_dof+1)=scratch(1)
95  do row=1, self%n_dof
96  y(row)= y(row)-&
97  0.5_f64*self%dt/self%dx*(scratch(row)-scratch(row+1))
98 
99  y(self%n_dof+row)=x(self%n_dof+row)+ &
100  0.5_f64*self%dt/self%dx*(x(row)-x(modulo(row-2,self%n_dof)+1))
101  end do
102 
103  end subroutine dot_mono_r2r_ecsim_eb
104 
105 
106  subroutine print_info_ecsim_eb( self )
107  class(sll_t_linear_operator_ecsim_eb), intent(in) :: self
108  end subroutine print_info_ecsim_eb
109 
module for abstract linear operator
subroutine create_linear_operator_ecsim_eb(self, n_dof, degree, mass_line_0, mass_line_1)
Utilites for Maxwell solver's with spline finite elements.
subroutine, public sll_s_spline_fem_mass_line(degree, mass_line)
Computes the mass line for a mass matrix with degree splines.
    Report Typos and Errors