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_3d.F90
Go to the documentation of this file.
1 
6 
7  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
8 #include "sll_working_precision.h"
9 
10  use sll_m_splines_pp, only : &
12 
13  implicit none
14 
15  public :: &
17 
18  private
19  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
20 
21  ! Define parameters to set if Galerkin or collocation scaling should be used in accumulation routines
22  sll_int32, parameter :: sll_p_galerkin = 0
23  sll_int32, parameter :: sll_p_collocation = 1
24 
27  sll_int32 :: dim = 3
28  sll_real64 :: delta_x(3)
29  sll_real64 :: domain(3,2)
30  sll_int32 :: n_total
31  sll_int32 :: spline_degree(3)
32  sll_int32 :: n_cells(3)
33  sll_int32 :: n_dofs(3)
34 
35  type(sll_t_spline_pp_3d) :: spline_pp_0
36  type(sll_t_spline_pp_3d) :: spline_pp_1
37 
38  contains
39 
40  procedure(add_single), deferred :: add_charge
41  procedure(add_array), deferred :: add_particle_mass
42  procedure(add_array_mixed), deferred :: add_particle_mass_od
43  procedure(eval_single), deferred :: evaluate
44  procedure(eval_multiple), deferred :: evaluate_multiple
45  procedure(add_current), deferred :: add_current
47  procedure(add_update), deferred :: add_current_update_v_component1
48  procedure(add_update), deferred :: add_current_update_v_component2
49  procedure(add_update), deferred :: add_current_update_v_component3
50  procedure(empty), deferred :: free
51 
52 
54 
55  !---------------------------------------------------------------------------!
56  abstract interface
57  subroutine add_single(self, position, marker_charge, degree, rho_dofs)
60  class(sll_c_particle_mesh_coupling_3d), intent( inout ) :: self
61  sll_real64, intent( in ) :: position(3)
62  sll_real64, intent( in ) :: marker_charge
63  sll_int32, intent( in ) :: degree(3)
64  sll_real64, intent( inout ) :: rho_dofs(:)
65 
66  end subroutine add_single
67  end interface
68 
69  !---------------------------------------------------------------------------!
70  abstract interface
71  subroutine add_array(self, position, marker_charge, degree, particle_mass )
74  class(sll_c_particle_mesh_coupling_3d), intent( inout ) :: self
75  sll_real64, intent( in ) :: position(3)
76  sll_real64, intent( in ) :: marker_charge
77  sll_int32, intent( in ) :: degree(3)
78  sll_real64, intent( inout ) :: particle_mass(:, :)
79 
80  end subroutine add_array
81  end interface
82 
83  !---------------------------------------------------------------------------!
84  abstract interface
85  subroutine add_array_mixed(self, position, marker_charge, degree1, degree2, particle_mass )
88  class(sll_c_particle_mesh_coupling_3d), intent( inout ) :: self
89  sll_real64, intent( in ) :: position(3)
90  sll_real64, intent( in ) :: marker_charge
91  sll_int32, intent( in ) :: degree1(3), degree2(3)
92  sll_real64, intent( inout ) :: particle_mass(:, :)
93 
94  end subroutine add_array_mixed
95  end interface
96 
97  !---------------------------------------------------------------------------!
98  abstract interface
99  subroutine eval_single(self, position, degree, field_dofs, field_value)
102  class(sll_c_particle_mesh_coupling_3d), intent( inout ) :: self
103  sll_real64, intent( in ) :: position(3)
104  sll_int32 , intent( in ) :: degree(3)
105  sll_real64, intent( in ) :: field_dofs(:)
106  sll_real64, intent( out ) :: field_value
107  end subroutine eval_single
108  end interface
109 
110  !---------------------------------------------------------------------------!
111  abstract interface
112  subroutine eval_multiple(self, position, components, field_dofs, field_value)
115  class(sll_c_particle_mesh_coupling_3d), intent( inout ) :: self
116  sll_real64, intent( in ) :: position(3)
117  sll_int32, intent(in) :: components(:)
118  sll_real64, intent( in ) :: field_dofs(:,:)
119  sll_real64, intent(out) :: field_value(:)
120  end subroutine eval_multiple
121  end interface
122 
123  !---------------------------------------------------------------------------!
124  abstract interface
125  subroutine add_current(self, position_old, position_new, xdot, j_dofs)
128  class(sll_c_particle_mesh_coupling_3d), intent( inout ) :: self
129  sll_real64, intent( in ) :: position_old(3)
130  sll_real64, intent( in ) :: position_new(3)
131  sll_real64, intent( in ) :: xdot(3)
132  sll_real64, intent( inout ) :: j_dofs(:)
133 
134  end subroutine add_current
135  end interface
136 
137  !---------------------------------------------------------------------------!
138  abstract interface
139  subroutine add_current_evaluate(self, position_old, position_new, xdot, efield_dofs, j_dofs, efield_val)
142  class(sll_c_particle_mesh_coupling_3d), intent( inout ) :: self
143  sll_real64, intent( in ) :: position_old(3)
144  sll_real64, intent( in ) :: position_new(3)
145  sll_real64, intent( in ) :: xdot(3)
146  sll_real64, intent( in ) :: efield_dofs(:)
147  sll_real64, intent( inout ) :: j_dofs(:)
148  sll_real64, intent( out ) :: efield_val(3)
149 
150  end subroutine add_current_evaluate
151  end interface
152 
153  !---------------------------------------------------------------------------!
154  abstract interface
155  subroutine add_update (self, position_old, position_new, marker_charge, &
156  qoverm, bfield_dofs, vi, j_dofs)
159  class(sll_c_particle_mesh_coupling_3d), intent(inout) :: self
160  sll_real64, intent(in) :: position_old(3)
161  sll_real64, intent(in) :: position_new
162  sll_real64, intent(in) :: marker_charge
163  sll_real64, intent(in) :: qoverm
164  sll_real64, intent(in) :: bfield_dofs(:)
165  sll_real64, intent(inout) :: vi(3)
166  sll_real64, intent(inout) :: j_dofs(:)
167 
168  end subroutine add_update
169  end interface
170 
171  !---------------------------------------------------------------------------!
172  abstract interface
173  subroutine empty(self)
175  class(sll_c_particle_mesh_coupling_3d), intent( inout ) :: self
176 
177  end subroutine empty
178  end interface
179 
Base class for kernel smoothers for accumulation and field evaluation in PIC.
Splines in pp form.
Module to select the kind parameter.
    Report Typos and Errors