Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_lobatto_poisson.F90
Go to the documentation of this file.
2 
3 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
4 #include "sll_assert.h"
5 #include "sll_working_precision.h"
6 
9 
10  use sll_m_dg_fields, only: &
12 
13  use sll_m_lobalap, only: &
14  sll_s_assemb, &
19  sll_s_init, &
22 
23  use sll_m_map_function, only: &
25 
26  implicit none
27 
28  public :: &
30  sll_o_create, &
31  sll_o_delete, &
33 
34  private
35 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
36 
39  sll_int32 :: order
41 
42  interface sll_o_create
43  module procedure initialize_lobatto_poisson
44  end interface sll_o_create
45 
46  interface sll_o_solve
47  module procedure solve_lobatto_poisson
48  end interface sll_o_solve
49 
50  interface sll_o_delete
51  module procedure delete_lobatto_poisson
52  end interface sll_o_delete
53 
54 contains
55 
56  subroutine initialize_lobatto_poisson(this, tau, order)
57 
58  type(sll_t_lobatto_poisson_solver) :: this
59  class(sll_c_coordinate_transformation_2d_base), pointer :: tau
60  sll_int32, optional :: order
61 
62  sll_int32 :: nx0
63  sll_int32 :: ny0
64 
65  this%tau => tau
66  nx0 = tau%mesh%num_cells1
67  ny0 = tau%mesh%num_cells2
68 
69  call sll_s_set_map_function(tau)
70 
71  if (present(order)) then
72  call sll_s_init(nx0, ny0, order)
73  else
74  call sll_s_init(nx0, ny0, 2)
75  end if
76  call sll_s_assemb()
77  call sll_s_computelu()
78 
79  end subroutine initialize_lobatto_poisson
80 
81  subroutine solve_lobatto_poisson(this, rhs, ex, ey)
82 
83  type(sll_t_lobatto_poisson_solver) :: this
84  type(sll_t_dg_field_2d) :: rhs
85  type(sll_t_dg_field_2d) :: ex
86  type(sll_t_dg_field_2d) :: ey
87 
88  sll_assert(this%order > 0)
89  call sll_s_assemb_rhs(rhs%array)
90  call sll_s_compute_phi()
91  call sll_s_compute_electric_field(ex%array, ey%array)
92 
93  end subroutine solve_lobatto_poisson
94 
95  subroutine delete_lobatto_poisson(this)
96 
97  type(sll_t_lobatto_poisson_solver) :: this
98 
99  call sll_s_plotgmsh()
100  call sll_s_release()
101  sll_assert(this%order > 0)
102 
103  end subroutine delete_lobatto_poisson
104 
105 end module sll_m_lobatto_poisson
Abstract class for coordinate transformations.
Solve Maxwell equations on cartesian domain with Disconituous Galerkine method:
subroutine, public sll_s_assemb()
subroutine, public sll_s_assemb_rhs(dg_rho)
subroutine, public sll_s_compute_electric_field(dg_ex, dg_ey)
subroutine, public sll_s_computelu()
subroutine, public sll_s_release()
subroutine, public sll_s_compute_phi()
subroutine, public sll_s_init(nx0, ny0, order0)
subroutine, public sll_s_plotgmsh()
subroutine delete_lobatto_poisson(this)
subroutine solve_lobatto_poisson(this, rhs, ex, ey)
subroutine initialize_lobatto_poisson(this, tau, order)
subroutine, public sll_s_set_map_function(mytau)
Object to describe field data with DG numerical method.
    Report Typos and Errors