7 #include "sll_assert.h"
8 #include "sll_working_precision.h"
32 sll_int32 :: no_gridpts(2)
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(:,:)
47 logical :: rho_collected
72 sll_real64,
intent( in ) :: position(self%dim)
73 sll_real64,
intent( in ) :: marker_charge
75 call self%kernel%add_charge( position, marker_charge, self%rho_dofs_local )
83 sll_real64,
intent( in ) :: position(self%dim)
84 sll_real64,
intent(out) :: func_value
86 if (self%rho_collected .EQV. .false.)
then
87 self%rho_collected = .true.
88 self%rho_dofs = 0.0_f64
90 self%no_dofs, mpi_sum, self%rho_dofs)
93 call self%kernel%evaluate( position, self%rho_dofs, func_value)
100 sll_real64,
intent( in ) :: position(self%dim)
101 sll_real64,
intent(out) :: func_value
103 call self%kernel%evaluate( position, self%phi_dofs, func_value)
110 sll_int32,
intent(in ) :: components(:)
111 sll_real64,
intent( in ) :: position(self%dim)
112 sll_real64,
intent(out) :: func_value(:)
114 call self%kernel%evaluate_multiple( position, components, self%efield_dofs, &
124 call self%solve_phi()
125 call self%solve_fields()
133 if (self%rho_collected .EQV. .false.)
then
134 self%rho_collected = .true.
135 self%rho_dofs = 0.0_f64
137 self%no_dofs, mpi_sum, self%rho_dofs)
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])
149 if (self%rho_collected .EQV. .false.)
then
150 self%rho_collected = .true.
151 self%rho_dofs = 0.0_f64
153 self%no_dofs, mpi_sum, self%rho_dofs)
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])
166 self%rho_dofs_local = 0.0_f64
167 self%rho_dofs = 0.0_f64
168 self%rho_collected = .false.
175 sll_real64,
intent( in ) :: factor_present
176 sll_real64,
intent( in ) :: factor_analytic
178 self%rho_dofs = factor_present * self%rho_dofs + &
179 factor_analytic * self%rho_analyt_dofs
188 call self%poisson%compute_rhs_from_function(func, self%rho_analyt_dofs)
195 sll_int32,
intent( in ) :: component
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)
211 sll_int32,
intent(in) :: no_gridpts(2)
219 self%no_gridpts = no_gridpts
220 self%no_dofs = product(no_gridpts)
221 self%rho_collected = .false.
223 self%poisson => poisson
224 self%kernel => kernel
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 )
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()
269 sll_int32,
intent(in) :: no_gridpts(2)
277 select type( pic_poisson)
279 call pic_poisson%init( no_gridpts, poisson, kernel)
Combines values from all processes and distributes the result back to all processes.
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 a kernel smoother used for PIC simulations.
Basic type of Poisson solver for PIC simulations.