Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_linear_operator_poisson_3d.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 : &
13 
14  implicit none
15 
17 
18  private
19 
21  type(sll_t_linear_operator_block), pointer :: mass
22  sll_int32 :: n_dofs(3)
23  sll_real64 :: delta_x(3)
24  sll_int32 :: n_total
25 
26 
27  contains
28  procedure :: create => create_linear_operator_poisson_3d
29  procedure :: free => free_poisson_3d
30  procedure :: dot => dot_poisson_3d
31  procedure :: print_info => print_info_poisson_3d
32 
34 
35 contains
36 
37 
38  subroutine create_linear_operator_poisson_3d( self, mass, n_dofs, delta_x)
39  class(sll_t_linear_operator_poisson_3d), intent( inout ) :: self
40  type(sll_t_linear_operator_block), target :: mass
41  sll_int32 :: n_dofs(3)
42  sll_real64 :: delta_x(3)
43 
44  self%mass => mass
45  self%n_dofs = n_dofs
46  self%n_total = product(n_dofs)
47  self%delta_x = delta_x
48 
49  self%n_rows = self%n_total
50  self%n_cols = self%n_total
51 
52  self%n_global_rows = self%n_rows
53  self%n_global_cols = self%n_cols
55 
56 
57  subroutine free_poisson_3d( self )
58  class(sll_t_linear_operator_poisson_3d), intent( inout ) :: self
59  self%mass => null()
60 
61  end subroutine free_poisson_3d
62 
63 
64  subroutine dot_poisson_3d( self, x, y )
65  class(sll_t_linear_operator_poisson_3d), intent( in ) :: self
66  sll_real64, intent( in ) :: x(:)
67  sll_real64, intent( out ) :: y(:)
68  !local variable
69  sll_real64 :: scratch(3*self%n_total), scratch1(3*self%n_total)
70 
71  call sll_s_multiply_g(self%n_dofs, self%delta_x, x, scratch)
72  call self%mass%dot(scratch, scratch1)
73  call sll_s_multiply_gt(self%n_dofs, self%delta_x, scratch1, y)
74 
75  end subroutine dot_poisson_3d
76 
77  subroutine print_info_poisson_3d( self )
78  class(sll_t_linear_operator_poisson_3d), intent(in) :: self
79  end subroutine print_info_poisson_3d
80 
81 
module for abstract linear operator
module for a block linear operator
subroutine create_linear_operator_poisson_3d(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)
subroutine, public sll_s_multiply_gt(n_dofs, delta_x, field_in, field_out)
Multiply by transpose of dicrete gradient matrix (3d version)
    Report Typos and Errors