Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_pic_poisson_2d.F90
Go to the documentation of this file.
1 
5 
7 #include "sll_assert.h"
8 #include "sll_working_precision.h"
9 
10  use sll_m_pic_poisson_base, only : &
14  use sll_m_poisson_2d_base, only : &
16  use sll_m_collective, only : &
18  use mpi, only: &
19  mpi_sum
20 
21 
22  implicit none
23 
24  public :: sll_t_pic_poisson_2d, &
26 
27  private
28 
31 
32  sll_int32 :: no_gridpts(2)
33  sll_int32 :: no_dofs
34 
35  class(sll_c_particle_mesh_coupling_1d), pointer :: kernel
36  class(sll_c_poisson_2d_base), pointer :: poisson
37  sll_real64, allocatable :: rho_dofs(:)
38  sll_real64, allocatable :: rho_dofs_local(:)
39  sll_real64, allocatable :: rho_analyt_dofs(:)
40  sll_real64, allocatable :: efield_dofs(:,:)
41  sll_real64, allocatable :: phi_dofs(:)
42  sll_real64, allocatable :: rho2d(:,:)
43  sll_real64, allocatable :: efield1(:,:)
44  sll_real64, allocatable :: efield2(:,:)
45  sll_real64, allocatable :: phi2d(:,:)
46 
47  logical :: rho_collected
48 
49  contains
50  procedure :: add_charge_single => add_charge_single_2d
51  procedure :: reset => reset_2d
52  procedure :: evaluate_rho_single => evaluate_rho_single_2d
53  procedure :: evaluate_phi_single => evaluate_phi_single_2d
54  procedure :: evaluate_field_single => evaluate_field_single_2d
55  procedure :: solve => solve_2d
56  procedure :: solve_phi => solve_phi_2d
57  procedure :: solve_fields => solve_fields_2d
58  procedure :: add_analytic_charge => add_analytic_charge_2d
59  procedure :: set_analytic_charge => set_analytic_charge_2d
60  procedure :: compute_field_energy => compute_field_energy_2d
61  procedure :: init => init_pic_poisson_2d
62  procedure :: free => free_pic_poisson_2d
63 
64 
65  end type sll_t_pic_poisson_2d
66 
67 contains
68 
70  subroutine add_charge_single_2d(self, position, marker_charge)
71  class(sll_t_pic_poisson_2d), intent( inout ) :: self
72  sll_real64, intent( in ) :: position(self%dim)
73  sll_real64, intent( in ) :: marker_charge
74 
75  call self%kernel%add_charge( position, marker_charge, self%rho_dofs_local )
76 
77 
78  end subroutine add_charge_single_2d
79 
81  subroutine evaluate_rho_single_2d(self, position, func_value)
82  class(sll_t_pic_poisson_2d), intent( inout ) :: self
83  sll_real64, intent( in ) :: position(self%dim)
84  sll_real64, intent(out) :: func_value
85 
86  if (self%rho_collected .EQV. .false.) then
87  self%rho_collected = .true.
88  self%rho_dofs = 0.0_f64
89  call sll_o_collective_allreduce( sll_v_world_collective, self%rho_dofs_local, &
90  self%no_dofs, mpi_sum, self%rho_dofs)
91  end if
92 
93  call self%kernel%evaluate( position, self%rho_dofs, func_value)
94 
95  end subroutine evaluate_rho_single_2d
96 
98  subroutine evaluate_phi_single_2d(self, position, func_value)
99  class(sll_t_pic_poisson_2d), intent( inout ) :: self
100  sll_real64, intent( in ) :: position(self%dim)
101  sll_real64, intent(out) :: func_value
102 
103  call self%kernel%evaluate( position, self%phi_dofs, func_value)
104 
105  end subroutine evaluate_phi_single_2d
106 
108  subroutine evaluate_field_single_2d(self, position, components, func_value)
109  class(sll_t_pic_poisson_2d), intent( inout ) :: self
110  sll_int32, intent(in ) :: components(:)
111  sll_real64, intent( in ) :: position(self%dim)
112  sll_real64, intent(out) :: func_value(:)
113 
114  call self%kernel%evaluate_multiple( position, components, self%efield_dofs, &
115  func_value)
116 
117  end subroutine evaluate_field_single_2d
118 
119 
121  subroutine solve_2d(self)
122  class(sll_t_pic_poisson_2d), intent( inout ) :: self
123 
124  call self%solve_phi()
125  call self%solve_fields()
126 
127  end subroutine solve_2d
128 
130  subroutine solve_phi_2d(self)
131  class(sll_t_pic_poisson_2d), intent( inout ) :: self
132 
133  if (self%rho_collected .EQV. .false.) then
134  self%rho_collected = .true.
135  self%rho_dofs = 0.0_f64
136  call sll_o_collective_allreduce( sll_v_world_collective, self%rho_dofs_local, &
137  self%no_dofs, mpi_sum, self%rho_dofs)
138  end if
139  self%rho2d = reshape(self%rho_dofs, self%no_gridpts)
140  call self%poisson%compute_phi_from_rho(self%phi2d, self%rho2d)
141  self%phi_dofs = reshape(self%phi2d, [self%no_dofs])
142 
143  end subroutine solve_phi_2d
144 
146  subroutine solve_fields_2d(self)
147  class(sll_t_pic_poisson_2d), intent( inout ) :: self
148 
149  if (self%rho_collected .EQV. .false.) then
150  self%rho_collected = .true.
151  self%rho_dofs = 0.0_f64
152  call sll_o_collective_allreduce( sll_v_world_collective, self%rho_dofs_local, &
153  self%no_dofs, mpi_sum, self%rho_dofs)
154  end if
155  self%rho2d = reshape(self%rho_dofs, self%no_gridpts)
156  call self%poisson%compute_E_from_rho(self%efield1, self%efield2, self%rho2d)
157  self%efield_dofs(:,1) = reshape(self%efield1, [self%no_dofs])
158  self%efield_dofs(:,2) = reshape(self%efield2, [self%no_dofs])
159 
160  end subroutine solve_fields_2d
161 
163  subroutine reset_2d(self)
164  class(sll_t_pic_poisson_2d), intent( inout ) :: self
165 
166  self%rho_dofs_local = 0.0_f64
167  self%rho_dofs = 0.0_f64
168  self%rho_collected = .false.
169 
170  end subroutine reset_2d
171 
173  subroutine add_analytic_charge_2d(self, factor_present, factor_analytic)
174  class(sll_t_pic_poisson_2d), intent( inout ) :: self
175  sll_real64, intent( in ) :: factor_present
176  sll_real64, intent( in ) :: factor_analytic
177 
178  self%rho_dofs = factor_present * self%rho_dofs + &
179  factor_analytic * self%rho_analyt_dofs
180 
181  end subroutine add_analytic_charge_2d
182 
184  subroutine set_analytic_charge_2d(self, func)
185  class( sll_t_pic_poisson_2d ), intent( inout ) :: self
186  procedure(sll_i_function_of_position) :: func
187 
188  call self%poisson%compute_rhs_from_function(func, self%rho_analyt_dofs)
189 
190  end subroutine set_analytic_charge_2d
191 
193  function compute_field_energy_2d(self, component) result(energy)
194  class(sll_t_pic_poisson_2d), intent( in ) :: self
195  sll_int32, intent( in ) :: component
196  sll_real64 :: energy
197 
198 
199  if (component == 1) then
200  energy = self%poisson%l2norm_squared(self%efield1)
201  elseif (component == 2) then
202  energy = self%poisson%l2norm_squared(self%efield2)
203  end if
204 
205  end function compute_field_energy_2d
206 
207  !-------------------------------------------------------------------------------------------
208 
209  subroutine init_pic_poisson_2d(self, no_gridpts, poisson, kernel)
210  class( sll_t_pic_poisson_2d), intent(out) :: self
211  sll_int32, intent(in) :: no_gridpts(2)
212  class( sll_c_poisson_2d_base), pointer, intent(in) :: poisson
213  class( sll_c_particle_mesh_coupling_1d), pointer, intent(in) :: kernel
214 
215  !local variables
216  sll_int32 :: ierr
217 
218  self%dim = 1
219  self%no_gridpts = no_gridpts
220  self%no_dofs = product(no_gridpts)
221  self%rho_collected = .false.
222 
223  self%poisson => poisson
224  self%kernel => kernel
225 
226  allocate(self%rho_dofs(self%no_dofs), stat=ierr)
227  sll_assert( ierr == 0 )
228  allocate(self%rho_dofs_local(self%no_dofs), stat=ierr)
229  sll_assert( ierr == 0 )
230  allocate(self%rho_analyt_dofs(self%no_dofs), stat=ierr)
231  sll_assert( ierr == 0 )
232  allocate(self%efield_dofs(self%no_dofs, 2), stat=ierr)
233  sll_assert( ierr == 0 )
234  allocate(self%phi_dofs(self%no_dofs), stat=ierr)
235  sll_assert( ierr == 0 )
236  allocate(self%rho2d(self%no_gridpts(1), self%no_gridpts(2)), stat=ierr)
237  sll_assert( ierr == 0 )
238  allocate(self%efield1(self%no_gridpts(1), self%no_gridpts(2)), stat=ierr)
239  sll_assert( ierr == 0 )
240  allocate(self%efield2(self%no_gridpts(1), self%no_gridpts(2)), stat=ierr)
241  sll_assert( ierr == 0 )
242  allocate(self%phi2d(self%no_gridpts(1), self%no_gridpts(2)), stat=ierr)
243  sll_assert( ierr == 0 )
244 
245  end subroutine init_pic_poisson_2d
246 
248  subroutine free_pic_poisson_2d(self)
249  class( sll_t_pic_poisson_2d), intent(inout) :: self
250 
251  deallocate(self%rho_dofs)
252  deallocate(self%rho_dofs_local)
253  deallocate(self%rho_analyt_dofs)
254  deallocate(self%efield_dofs)
255  deallocate(self%phi_dofs)
256  deallocate(self%rho2d)
257  deallocate(self%efield1)
258  deallocate(self%efield2)
259  deallocate(self%phi2d)
260  self%poisson => null()
261  self%kernel => null()
262 
263  end subroutine free_pic_poisson_2d
264 
265 
267  subroutine sll_s_new_pic_poisson_2d(pic_poisson, no_gridpts, poisson, kernel)
268  class( sll_c_pic_poisson), pointer, intent(out) :: pic_poisson
269  sll_int32, intent(in) :: no_gridpts(2)
270  class( sll_c_poisson_2d_base), pointer, intent(in) :: poisson
271  class( sll_c_particle_mesh_coupling_1d), pointer, intent(in) :: kernel
272 
273  !local variables
274  sll_int32 :: ierr
275 
276  allocate( sll_t_pic_poisson_2d :: pic_poisson, stat=ierr)
277  select type( pic_poisson)
278  type is ( sll_t_pic_poisson_2d )
279  call pic_poisson%init( no_gridpts, poisson, kernel)
280  end select
281 
282  end subroutine sll_s_new_pic_poisson_2d
283 
284 end module sll_m_pic_poisson_2d
Combines values from all processes and distributes the result back to all processes.
Parallelizing facility.
type(sll_t_collective_t), pointer, public sll_v_world_collective
Control of subset assignment. Processes with the same color are in the same new communicator.
Base class for kernel smoothers for accumulation and field evaluation in PIC.
Factory method for Poisson solver for particle methods in 2d build from 2d Poisson solver and a kerne...
subroutine set_analytic_charge_2d(self, func)
Set analytic charge defined by a function func obeying the interface sll_i_function_of_position.
subroutine solve_2d(self)
Solve for phi and fields.
subroutine evaluate_field_single_2d(self, position, components, func_value)
Evaluate components components of the electric field as one position.
subroutine evaluate_phi_single_2d(self, position, func_value)
Evaluate potential phi at one position.
subroutine free_pic_poisson_2d(self)
Destructor.
subroutine add_charge_single_2d(self, position, marker_charge)
Add charge from one particle.
subroutine solve_fields_2d(self)
Solve efields from rho.
subroutine add_analytic_charge_2d(self, factor_present, factor_analytic)
Add analytic charge (set by set_analytic_charge ) to the accumulated charge.
real(kind=f64) function compute_field_energy_2d(self, component)
Compute the squared l2 norm of component component of the field.
subroutine, public sll_s_new_pic_poisson_2d(pic_poisson, no_gridpts, poisson, kernel)
Constructor for abstract type.
subroutine reset_2d(self)
Reset charge to zero.
subroutine solve_phi_2d(self)
Solve for potential.
subroutine init_pic_poisson_2d(self, no_gridpts, poisson, kernel)
subroutine evaluate_rho_single_2d(self, position, func_value)
Evaluate charge density rho at one position.
Base class for Poisson solver for particle methods.
Module interface to solve Poisson equation in 2D.
Basic type of Poisson solver for PIC simulations.
    Report Typos and Errors