Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_pic_poisson_base.F90
Go to the documentation of this file.
1 
9 
10 #include "sll_working_precision.h"
11 
12  use sll_m_poisson_2d_base, only : &
14 
15  implicit none
16 
17  public :: sll_c_pic_poisson
18 
19  private
20 
22  type, abstract :: sll_c_pic_poisson
23  sll_int32 :: dim
24  sll_int32 :: no_weights = 1
25 
26  contains
27  procedure(add_single), deferred :: add_charge_single
28  procedure :: add_charge_vector
29  procedure(eval_component_single), deferred :: evaluate_field_single
30  procedure :: evaluate_field_vector
31  procedure(eval_single), deferred :: evaluate_rho_single
32  procedure :: evaluate_rho_vector
33  procedure(eval_single), deferred :: evaluate_phi_single
34  procedure :: evaluate_phi_vector
35  procedure(empty), deferred :: reset
36  procedure(empty), deferred :: solve
37  procedure :: solve_phi
38  procedure :: solve_fields
39  procedure(compute_energy), deferred :: compute_field_energy
40  procedure(linear_combination), deferred :: add_analytic_charge
41  procedure(update_dofs_function), deferred :: set_analytic_charge
42 
43  procedure(empty), deferred :: free
44 
45 
46  generic :: add_charge => add_charge_single, add_charge_vector
47  generic :: evaluate_field => evaluate_field_single, evaluate_field_vector
48  generic :: evaluate_rho => evaluate_rho_single, evaluate_rho_vector
49  generic :: evaluate_phi => evaluate_phi_single, evaluate_phi_vector
50 
51  end type sll_c_pic_poisson
52 
53 
54  !---------------------------------------------------------------------------!
55  abstract interface
56  subroutine add_single(self, position, marker_charge)
58  import sll_c_pic_poisson
59  class(sll_c_pic_poisson), intent( inout ) :: self
60  sll_real64, intent( in ) :: position(self%dim)
61  sll_real64, intent( in ) :: marker_charge
62  end subroutine add_single
63  end interface
64 
65  !---------------------------------------------------------------------------!
66  abstract interface
67  subroutine eval_single(self, position, func_value)
69  import sll_c_pic_poisson
70  class(sll_c_pic_poisson), intent( inout ) :: self
71  sll_real64, intent( in ) :: position(self%dim)
72  sll_real64, intent( out) :: func_value
73  end subroutine eval_single
74  end interface
75 
76  !---------------------------------------------------------------------------!
77  abstract interface
78  subroutine eval_component_single(self, position, components, func_value)
80  import sll_c_pic_poisson
81  class(sll_c_pic_poisson), intent( inout ) :: self
82  sll_real64, intent( in ) :: position(self%dim)
83  sll_int32, intent( in ) :: components(:)
84  sll_real64, intent( out) :: func_value(:)
85  end subroutine eval_component_single
86  end interface
87 
88 
89  !---------------------------------------------------------------------------!
90  abstract interface
91  subroutine empty(self)
92  import sll_c_pic_poisson
93  class(sll_c_pic_poisson), intent( inout ) :: self
94  end subroutine empty
95  end interface
96 
97 
98  !---------------------------------------------------------------------------!
99  abstract interface
100  function compute_energy(self, component) result(energy)
102  import sll_c_pic_poisson
103  class(sll_c_pic_poisson), intent( in ) :: self
104  sll_int32, intent( in ) :: component
105  sll_real64 :: energy
106  end function compute_energy
107  end interface
108 
109 
110  !---------------------------------------------------------------------------!
111  abstract interface
112  subroutine update_dofs_function(self, func)
114  import sll_c_pic_poisson
116  class( sll_c_pic_poisson ), intent( inout ) :: self
117  procedure(sll_i_function_of_position) :: func
118  end subroutine update_dofs_function
119  end interface
120 
121  !---------------------------------------------------------------------------!
122  abstract interface
123  subroutine linear_combination(self, factor_present, factor_analytic)
125  import sll_c_pic_poisson
126  class( sll_c_pic_poisson ), intent( inout ) :: self
127  sll_real64, intent( in ) :: factor_present
128  sll_real64, intent( in ) :: factor_analytic
129 
130  end subroutine linear_combination
131  end interface
132 
133 contains
134 
135  !---------------------------------------------------------------------------!
136 
137  subroutine add_charge_vector(self,n_part, position, marker_charge)
138  class(sll_c_pic_poisson), intent( inout ) :: self
139  sll_int32, intent( in ) :: n_part
140  sll_real64, intent( in ) :: position(self%dim, n_part)
141  sll_real64, intent( in ) :: marker_charge(n_part)
142 
143  !local variables
144  sll_int32 :: i_part
145 
146  do i_part = 1, n_part
147  call self%add_charge_single( position(:,i_part), marker_charge(i_part))
148  end do
149 
150  end subroutine add_charge_vector
151 
152 
153  !---------------------------------------------------------------------------!
154 
155  subroutine evaluate_rho_vector(self, position, n_part, func_value)
156  class(sll_c_pic_poisson), intent( inout ) :: self
157  sll_real64, intent( in ) :: position(:, :)
158  sll_int32, intent( in ) :: n_part
159  sll_real64, intent( out) :: func_value(n_part)
160 
161  !local variables
162  sll_int32 :: i_part
163 
164  do i_part = 1, n_part
165  call self%evaluate_rho_single&
166  (position(:, i_part), func_value(i_part))
167  end do
168 
169 
170  end subroutine evaluate_rho_vector
171 
172  !---------------------------------------------------------------------------!
173 
174  subroutine evaluate_phi_vector(self, position, n_part, func_value)
175  class(sll_c_pic_poisson), intent( inout ) :: self
176  sll_real64, intent( in ) :: position(:,:)
177  sll_int32, intent( in ) :: n_part
178  sll_real64, intent( out) :: func_value(n_part)
179 
180  !local variables
181  sll_int32 :: i_part
182 
183  do i_part = 1, n_part
184  call self%evaluate_phi_single&
185  (position(:, i_part), func_value(i_part))
186  end do
187 
188 
189  end subroutine evaluate_phi_vector
190 
191 
192  !---------------------------------------------------------------------------!
193  subroutine evaluate_field_vector(self, position, n_part, components,func_value)
194  class(sll_c_pic_poisson), intent( inout ) :: self
195  sll_real64, intent( in ) :: position(:,:)
196  sll_int32, intent( in ) :: n_part
197  sll_int32, intent( in ) :: components(:)
198  sll_real64, intent( out) :: func_value(:,:)
199 
200  !local variables
201  sll_int32 :: i_part
202 
203  do i_part = 1, n_part
204  call self%evaluate_field_single&
205  (position(:, i_part), components, func_value(:,i_part))
206  end do
207 
208 
209  end subroutine evaluate_field_vector
210 
212  subroutine solve_phi( self )
213  class(sll_c_pic_poisson), intent( inout ) :: self
214 
215  call self%solve()
216 
217  end subroutine solve_phi
218 
220  subroutine solve_fields( self )
221  class(sll_c_pic_poisson), intent( inout ) :: self
222 
223  call self%solve()
224 
225  end subroutine solve_fields
226 
227 end module sll_m_pic_poisson_base
Base class for Poisson solver for particle methods.
subroutine add_charge_vector(self, n_part, position, marker_charge)
subroutine evaluate_rho_vector(self, position, n_part, func_value)
subroutine solve_fields(self)
Just compute the efield, by default it calls the combined procedure solve.
subroutine evaluate_field_vector(self, position, n_part, components, func_value)
subroutine evaluate_phi_vector(self, position, n_part, func_value)
subroutine solve_phi(self)
Just compute phi, by default it calls the combined procedure solve.
Module interface to solve Poisson equation in 2D.
Module to select the kind parameter.
Basic type of Poisson solver for PIC simulations.
    Report Typos and Errors