Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Kernel smoother for 2d with splines of arbitrary degree placed on a uniform mesh.
Spline with index i starts at point i
Derived types and interfaces | |
type | sll_t_particle_mesh_coupling_spline_smooth_1d |
Spline kernel smoother in1d. More... | |
Functions/Subroutines | |
subroutine | add_charge_single_smoothened_spline_1d (self, position, marker_charge, rho_dofs) |
Evaluates the integral int_{poisition_old}^{position_new} field(x) d x. More... | |
subroutine | add_current_update_v_smoothened_spline_1d (self, position_old, position_new, marker_charge, qoverm, bfield_dofs, vi, j_dofs) |
Add current for one particle and update v (according to H_p1 part in Hamiltonian splitting) More... | |
subroutine | update_smoothened_single_spline_1d (self, position_normalized, box, marker_charge, qoverm, rho_dofs, bfield_dofs, vi) |
Integration over the prime function of the smoothing kernel S (for v update and j accumulation) More... | |
subroutine | evaluate_int_prime_smoothened_single_spline_1d (self, position_normalized, box, field_dofs, field_int) |
Integration over the prime function of the smoothing kernel S (for evaluation only) More... | |
subroutine | evaluate_int_simple_smoothened_spline_1d (self, lower, upper, index, sign, field_dofs, field_int) |
Helper function to evaluate the integrated efield in one interval (without smoothing kernel) (similar to update_jv). More... | |
subroutine | evaluate_field_single_smoothened_spline_1d (self, position, field_dofs, field_value) |
Evaluate field with smoothing at at position position. More... | |
subroutine | add_current_evaluate_smoothened_spline_1d (self, position_old, position_new, marker_charge, vbar, field_dofs, j_dofs, field) |
Add current for one particle and evaluate the field (according to iterative part in discrete gradient method) More... | |
subroutine | evaluate_int_smoothened_spline_1d (self, position_old, position_new, vbar, field_dofs, field) |
Add current for one particle and evaluate the field (according to iterative part in discrete gradient method) More... | |
subroutine, public | sll_s_new_particle_mesh_coupling_spline_smooth_1d_ptr (smoother, domain, n_cells, no_particles, spline_degree, smoothing_type) |
subroutine, public | sll_s_new_particle_mesh_coupling_spline_smooth_1d (smoother, domain, n_cells, no_particles, spline_degree, smoothing_type) |
subroutine | add_particle_mass_spline_smooth_1d (self, position, marker_charge, particle_mass) |
Add particle mass for one particle (upper triangular version) More... | |
subroutine | add_particle_mass_full_spline_smooth_1d (self, position, marker_charge, particle_mass) |
Add particle mass for one particle (full matrix) WARNING: THIS FUNCTION SEEMS TO HAVE BUGS. More... | |
subroutine | add_charge_box (self, position, rho_dofs) |
Add charge of one particle with smoothing. More... | |
|
private |
Add charge of one particle with smoothing.
[in,out] | self | kernel smoother object |
[in] | position | Position of the particle |
[in,out] | rho_dofs | Coefficient vector of the charge distribution |
Definition at line 734 of file sll_m_particle_mesh_coupling_spline_smooth_1d.F90.
|
private |
Evaluates the integral int_{poisition_old}^{position_new} field(x) d x.
Add charge of one particle with smoothing
[in,out] | self | kernel smoother object |
[in] | position | Position of the particle |
[in] | marker_charge | Particle weights time charge |
[in,out] | rho_dofs | Coefficient vector of the charge distribution |
Definition at line 65 of file sll_m_particle_mesh_coupling_spline_smooth_1d.F90.
|
private |
Add current for one particle and evaluate the field (according to iterative part in discrete gradient method)
[in,out] | self | kernel smoother object |
[in] | position_old | Position of the particle |
[in] | position_new | Position of the particle |
[in] | marker_charge | Particle weights time charge |
[in] | vbar | Particle weights time charge |
[in] | field_dofs | Coefficient vector of the current density |
[in,out] | j_dofs | Coefficient vector of the current density |
Definition at line 502 of file sll_m_particle_mesh_coupling_spline_smooth_1d.F90.
|
private |
Add current for one particle and update v (according to H_p1 part in Hamiltonian splitting)
[in,out] | self | kernel smoother object |
[in] | position_old | Position at time t |
[in] | position_new | Position at time t + \Delta t |
[in] | marker_charge | Particle weight time charge |
[in] | qoverm | charge to mass ration |
[in] | bfield_dofs | Coefficient of B-field expansion |
[in,out] | vi | Velocity of the particles |
[in,out] | j_dofs | Coefficients of current expansion |
Definition at line 140 of file sll_m_particle_mesh_coupling_spline_smooth_1d.F90.
|
private |
Add particle mass for one particle (full matrix) WARNING: THIS FUNCTION SEEMS TO HAVE BUGS.
[in,out] | self | kernel smoother object |
[in] | position | Position of the particle |
[in] | marker_charge | Particle weights time charge |
[in,out] | particle_mass | Coefficient vector of the charge distribution |
Definition at line 693 of file sll_m_particle_mesh_coupling_spline_smooth_1d.F90.
|
private |
Add particle mass for one particle (upper triangular version)
[in,out] | self | kernel smoother object |
[in] | position | Position of the particle |
[in] | marker_charge | Particle weights time charge |
[in,out] | particle_mass | Coefficient vector of the charge distribution |
Definition at line 655 of file sll_m_particle_mesh_coupling_spline_smooth_1d.F90.
|
private |
Evaluate field with smoothing at at position position.
[in,out] | self | Kernel smoother object |
[in] | position | Position of the particle |
[in] | field_dofs | Coefficient vector for the field DoFs |
[out] | field_value | Value(s) of the electric fields at given position |
Definition at line 441 of file sll_m_particle_mesh_coupling_spline_smooth_1d.F90.
|
private |
Integration over the prime function of the smoothing kernel S (for evaluation only)
[in,out] | self | kernel smoother object |
[in] | position_normalized | Position of the particle |
Definition at line 321 of file sll_m_particle_mesh_coupling_spline_smooth_1d.F90.
|
private |
Helper function to evaluate the integrated efield in one interval (without smoothing kernel) (similar to update_jv).
[in,out] | self | time splitting object |
Definition at line 396 of file sll_m_particle_mesh_coupling_spline_smooth_1d.F90.
|
private |
Add current for one particle and evaluate the field (according to iterative part in discrete gradient method)
[in,out] | self | kernel smoother object |
[in] | position_old | Position of the particle |
[in] | position_new | Position of the particle |
[in] | vbar | Particle weights time charge |
[in] | field_dofs | Coefficient vector of the current density |
Definition at line 526 of file sll_m_particle_mesh_coupling_spline_smooth_1d.F90.
subroutine, public sll_m_particle_mesh_coupling_spline_smooth_1d::sll_s_new_particle_mesh_coupling_spline_smooth_1d | ( | class( sll_c_particle_mesh_coupling_1d), intent(out), allocatable | smoother, |
real(kind=f64), dimension(2), intent(in) | domain, | ||
integer(kind=i32), intent(in) | n_cells, | ||
integer(kind=i32), intent(in) | no_particles, | ||
integer(kind=i32), intent(in) | spline_degree, | ||
integer(kind=i32), intent(in) | smoothing_type | ||
) |
[out] | smoother | kernel smoother object |
[in] | n_cells | number of DoFs (spline coefficients) |
[in] | domain | x_min and x_max of the domain |
[in] | no_particles | number of particles |
[in] | spline_degree | Degree of smoothing kernel spline |
[in] | smoothing_type | Define if Galerkin or collocation smoothing for right scaling in accumulation routines |
Definition at line 632 of file sll_m_particle_mesh_coupling_spline_smooth_1d.F90.
subroutine, public sll_m_particle_mesh_coupling_spline_smooth_1d::sll_s_new_particle_mesh_coupling_spline_smooth_1d_ptr | ( | class( sll_c_particle_mesh_coupling_1d), intent(out), pointer | smoother, |
real(kind=f64), dimension(2), intent(in) | domain, | ||
integer(kind=i32), intent(in) | n_cells, | ||
integer(kind=i32), intent(in) | no_particles, | ||
integer(kind=i32), intent(in) | spline_degree, | ||
integer(kind=i32), intent(in) | smoothing_type | ||
) |
[out] | smoother | kernel smoother object |
[in] | n_cells | number of DoFs (spline coefficients) |
[in] | domain | x_min and x_max of the domain |
[in] | no_particles | number of particles |
[in] | spline_degree | Degree of smoothing kernel spline |
[in] | smoothing_type | Define if Galerkin or collocation smoothing for right scaling in accumulation routines |
Definition at line 608 of file sll_m_particle_mesh_coupling_spline_smooth_1d.F90.
|
private |
Integration over the prime function of the smoothing kernel S (for v update and j accumulation)
[in,out] | self | kernel smoother object |
[in] | position_normalized | Position of the particle |
[in] | marker_charge | Particle weights time charge |
[in,out] | rho_dofs | Coefficient vector of the charge distribution |
Definition at line 237 of file sll_m_particle_mesh_coupling_spline_smooth_1d.F90.