Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_linear_operator_MG.F90
Go to the documentation of this file.
2 #include "sll_working_precision.h"
3 
6 
7  use sll_m_linear_operator_block, only : &
9 
10  use sll_m_spline_fem_utilities, only : &
12 
13  implicit none
14 
15  public :: sll_t_linear_operator_mg
16 
17  private
18 
20  type(sll_t_linear_operator_block), pointer :: mass
21  sll_int32 :: n_dofs(3)
22  sll_real64 :: delta_x(3)
23  sll_int32 :: n_total
24 
25 
26  contains
27  procedure :: create => create_linear_operator_mg
28  procedure :: free => free_mg
29  procedure :: dot => dot_mg
30  procedure :: print_info => print_info_mg
31 
33 
34 contains
35 
36 
37  subroutine create_linear_operator_mg( self, mass, n_dofs, delta_x)
38  class(sll_t_linear_operator_mg), intent( inout ) :: self
39  type(sll_t_linear_operator_block), target :: mass
40  sll_int32 :: n_dofs(3)
41  sll_real64 :: delta_x(3)
42 
43  self%mass => mass
44  self%n_dofs = n_dofs
45  self%n_total = product(n_dofs)
46  self%delta_x = delta_x
47 
48  self%n_rows = 3*self%n_total
49  self%n_cols = self%n_total
50 
51  self%n_global_rows = self%n_rows
52  self%n_global_cols = self%n_cols
53  end subroutine create_linear_operator_mg
54 
55 
56  subroutine free_mg( self )
57  class(sll_t_linear_operator_mg), intent( inout ) :: self
58  self%mass => null()
59 
60  end subroutine free_mg
61 
62 
63  subroutine dot_mg( self, x, y )
64  class(sll_t_linear_operator_mg), intent( in ) :: self
65  sll_real64, intent( in ) :: x(:)
66  sll_real64, intent( out ) :: y(:)
67  !local variable
68  sll_real64 :: scratch(3*self%n_total)
69 
70  call sll_s_multiply_g(self%n_dofs, self%delta_x, x, scratch)
71  call self%mass%dot(scratch, y)
72 
73  end subroutine dot_mg
74 
75  subroutine print_info_mg( self )
76  class(sll_t_linear_operator_mg), intent(in) :: self
77  end subroutine print_info_mg
78 
79 
80 end module sll_m_linear_operator_mg
module for abstract linear operator
module for a block linear operator
subroutine create_linear_operator_mg(self, mass, n_dofs, delta_x)
Utilites for Maxwell solver's with spline finite elements.
subroutine, public sll_s_multiply_g(n_dofs, delta_x, field_in, field_out)
Multiply by dicrete gradient matrix (3d version)
    Report Typos and Errors