Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_particle_mesh_coupling_base_1d.F90
Go to the documentation of this file.
1 
6 
7  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
8 #include "sll_working_precision.h"
9 
10  implicit none
11 
12  public :: &
16 
17  private
18  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
19 
20  ! Define parameters to set if Galerkin or collocation scaling should be used in accumulation routines
21  sll_int32, parameter :: sll_p_galerkin = 0
22  sll_int32, parameter :: sll_p_collocation = 1
23 
26  sll_int32 :: dim
27  sll_real64 :: delta_x
28  sll_int32 :: n_cells
29  sll_int32 :: n_dofs
30  sll_real64 :: domain(2)
31  sll_int32 :: spline_degree
32 
33  contains
34 
35  procedure(add_single), deferred :: add_charge
36  procedure(add_array), deferred :: add_particle_mass
37  procedure(add_array), deferred :: add_particle_mass_full
38  procedure(eval_single), deferred :: evaluate
39  procedure(eval_multiple), deferred :: evaluate_multiple
40  procedure(add_update) , deferred :: add_current_update_v
41  procedure(add_current), deferred :: add_current
43  procedure(empty), deferred :: free
44 
45 
47 
48  !---------------------------------------------------------------------------!
49  abstract interface
50  subroutine add_single(self, position, marker_charge, rho_dofs)
53  class(sll_c_particle_mesh_coupling_1d), intent( inout ) :: self
54  sll_real64, intent( in ) :: position(self%dim)
55  sll_real64, intent( in ) :: marker_charge
56  sll_real64, intent( inout ) :: rho_dofs(self%n_dofs)
57 
58  end subroutine add_single
59  end interface
60 
61  !---------------------------------------------------------------------------!
62  abstract interface
63  subroutine add_array(self, position, marker_charge, particle_mass)
66  class(sll_c_particle_mesh_coupling_1d), intent( inout ) :: self
67  sll_real64, intent( in ) :: position(self%dim)
68  sll_real64, intent( in ) :: marker_charge
69  sll_real64, intent( inout ) :: particle_mass(:, :)
70 
71  end subroutine add_array
72  end interface
73 
74  !---------------------------------------------------------------------------!
75  abstract interface
76  subroutine eval_single(self, position, field_dofs, field_value)
79  class(sll_c_particle_mesh_coupling_1d), intent( inout ) :: self
80  sll_real64, intent( in ) :: position(self%dim)
81  sll_real64, intent( in ) :: field_dofs(self%n_dofs)
82  sll_real64, intent(out) :: field_value
83  end subroutine eval_single
84  end interface
85 
86  !---------------------------------------------------------------------------!
87  abstract interface
88  subroutine eval_multiple(self, position, components, field_dofs, field_value)
91  class(sll_c_particle_mesh_coupling_1d), intent( inout ) :: self
92  sll_real64, intent( in ) :: position(self%dim)
93  sll_int32, intent(in) :: components(:)
94  sll_real64, intent( in ) :: field_dofs(:,:)
95  sll_real64, intent(out) :: field_value(:)
96  end subroutine eval_multiple
97  end interface
98 
99  !---------------------------------------------------------------------------!
100  abstract interface
101  subroutine add_current(self, position_old, position_new, marker_charge, j_dofs)
104  class(sll_c_particle_mesh_coupling_1d), intent( inout ) :: self
105  sll_real64, intent( in ) :: position_old(self%dim)
106  sll_real64, intent( in ) :: position_new(self%dim)
107  sll_real64, intent( in ) :: marker_charge
108  sll_real64, intent( inout ) :: j_dofs(self%n_dofs)
109 
110  end subroutine add_current
111  end interface
112 
113  !---------------------------------------------------------------------------!
114  abstract interface
115  subroutine add_current_evaluate(self, position_old, position_new, marker_charge, vbar, field_dofs, j_dofs, field)
118  class(sll_c_particle_mesh_coupling_1d), intent( inout ) :: self
119  sll_real64, intent( in ) :: position_old(self%dim)
120  sll_real64, intent( in ) :: position_new(self%dim)
121  sll_real64, intent( in ) :: marker_charge
122  sll_real64, intent( in ) :: vbar
123  sll_real64, intent( in ) :: field_dofs(self%n_dofs)
124  sll_real64, intent( inout ) :: j_dofs(self%n_dofs)
125  sll_real64, intent( out ) :: field
126 
127  end subroutine add_current_evaluate
128  end interface
129 
130  !---------------------------------------------------------------------------!
131  abstract interface
132  subroutine add_update (self, &
133  position_old, &
134  position_new, &
135  marker_charge, &
136  qoverm, &
137  bfield_dofs, &
138  vi, &
139  j_dofs)
142  class(sll_c_particle_mesh_coupling_1d), intent(inout) :: self
143  sll_real64, intent(in) :: position_old(self%dim)
144  sll_real64, intent(in) :: position_new(self%dim)
145  sll_real64, intent(in) :: marker_charge
146  sll_real64, intent(in) :: qoverm
147  sll_real64, intent(in) :: bfield_dofs(:)
148  sll_real64, intent(inout) :: vi(:)
149  sll_real64, intent(inout) :: j_dofs(:)
150 
151  end subroutine add_update
152  end interface
153 
154  !---------------------------------------------------------------------------!
155  abstract interface
156  subroutine empty(self)
158  class(sll_c_particle_mesh_coupling_1d), intent( inout ) :: self
159 
160  end subroutine empty
161  end interface
162 
Base class for kernel smoothers for accumulation and field evaluation in PIC.
integer(kind=i32), parameter, public sll_p_collocation
Module to select the kind parameter.
    Report Typos and Errors