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_clamped_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 
13 
14  implicit none
15 
17 
18  private
19 
21  type(sll_t_linear_operator_block) , pointer :: mass
22  sll_int32 :: s_deg_0(3)
23  sll_int32 :: n_dofs(3)
24  sll_real64 :: delta_x(3)
25  sll_int32 :: n_total0, n_total1
26 
27 
28 
29  contains
31  procedure :: free => free_poisson_clamped_3d
32  procedure :: dot => dot_poisson_clamped_3d
33  procedure :: print_info => print_info_poisson_clamped_3d
34 
36 
37 contains
38 
39 
40  subroutine create_linear_operator_poisson_clamped_3d( self, mass, s_deg_0, n_dofs, delta_x)
41  class(sll_t_linear_operator_poisson_clamped_3d), intent( inout ) :: self
42  type(sll_t_linear_operator_block) , target :: mass
43  sll_int32, intent(in) :: s_deg_0(3)
44  sll_int32, intent(in) :: n_dofs(3)
45  sll_real64, intent(in) :: delta_x(3)
46 
47  self%mass => mass
48  self%s_deg_0 = s_deg_0
49  self%n_dofs = n_dofs
50  self%n_total0 = product(n_dofs)
51  self%n_total1 = (n_dofs(1)-1)*n_dofs(2)*n_dofs(3)
52  self%delta_x = delta_x
53 
54  self%n_rows = self%n_total0
55  self%n_cols = self%n_total0
56 
57  self%n_global_rows = self%n_rows
58  self%n_global_cols = self%n_cols
60 
61 
62  subroutine free_poisson_clamped_3d( self )
63  class(sll_t_linear_operator_poisson_clamped_3d), intent( inout ) :: self
64  self%mass => null()
65 
66  end subroutine free_poisson_clamped_3d
67 
68 
69  subroutine dot_poisson_clamped_3d( self, x, y )
70  class(sll_t_linear_operator_poisson_clamped_3d), intent( in ) :: self
71  sll_real64, intent( in ) :: x(:)
72  sll_real64, intent( out ) :: y(:)
73  !local variable
74  sll_real64 :: scratch(2*self%n_total0+self%n_total1), scratch1(2*self%n_total0+self%n_total1)
75 
76 
77  call sll_s_multiply_g_clamped(self%n_dofs, self%delta_x, self%s_deg_0, x, scratch)
78  call self%mass%dot(scratch, scratch1)
79  call sll_s_multiply_gt_clamped(self%n_dofs, self%delta_x, self%s_deg_0, scratch1, y)
80 
81  end subroutine dot_poisson_clamped_3d
82 
83  subroutine print_info_poisson_clamped_3d( self )
84  class(sll_t_linear_operator_poisson_clamped_3d), intent(in) :: self
85  end subroutine print_info_poisson_clamped_3d
86 
87 
module for abstract linear operator
module for a block linear operator
subroutine create_linear_operator_poisson_clamped_3d(self, mass, s_deg_0, n_dofs, delta_x)
Utilites for Maxwell solver's with spline finite elements using sparse matrices.
subroutine, public sll_s_multiply_g_clamped(n_dofs, delta_x, s_deg_0, field_in, field_out)
Multiply by dicrete gradient matrix.
subroutine, public sll_s_multiply_gt_clamped(n_dofs, delta_x, s_deg_0, field_in, field_out)
Multiply by transpose of dicrete gradient matrix.
    Report Typos and Errors