Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
Derived types and interfaces | Functions/Subroutines
sll_m_particle_mesh_coupling_spline_1d Module Reference

Description

Kernel smoother for 1d with splines of arbitrary degree placed on a uniform mesh.

Author
Katharina Kormann, IPP

Spline with index i starts at point i

Derived types and interfaces

type  sll_t_particle_mesh_coupling_spline_1d
 Spline kernel smoother in1d. More...
 

Functions/Subroutines

subroutine add_charge_single_spline_1d (self, position, marker_charge, rho_dofs)
 helper function to compute the integral of j using Gauss quadrature (needs to be provided for the derived type with smoothing) More...
 
subroutine add_particle_mass_spline_1d (self, position, marker_charge, particle_mass)
 Add the contribution of one particle to the approximate mass matrix. More...
 
subroutine add_particle_mass_spline_1d_full (self, position, marker_charge, particle_mass)
 Add the contribution of one particle to the approximate mass matrix. More...
 
subroutine evaluate_field_single_spline_1d (self, position, field_dofs, field_value)
 Evaluate field at at position position. More...
 
subroutine evaluate_multiple_spline_1d (self, position, components, field_dofs, field_value)
 Evaluate several fields at position position. More...
 
subroutine add_current_spline_1d (self, position_old, position_new, marker_charge, j_dofs)
 Add current with integration over x. More...
 
subroutine add_current_evaluate_spline_1d (self, position_old, position_new, marker_charge, vbar, field_dofs, j_dofs, field)
 Combines add_current and evaluate_int. More...
 
subroutine evaluate_int_spline_1d (self, position_old, position_new, vbar, field_dofs, field)
 Evaluates the integral int_{poisition_old}^{position_new} field(x) d x. More...
 
subroutine evaluate_int_quad_spline_1d (self, position_old, position_new, field_dofs, field)
 Evaluates an integrated field (between position_old and position_new) but with a quadrature formula for the integration. More...
 
subroutine evaluate_int_spline_1d_subnewton (self, position_old, position_new, field_dofs, field)
 
subroutine add_current_update_v_spline_pp_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_jv_pp (self, lower, upper, index, marker_charge, qoverm, vi, j_dofs, bfield_dofs)
 Helper function for add_current_update_v. More...
 
subroutine add_current_update_v_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_jv (self, lower, upper, index, marker_charge, qoverm, sign, vi, j_dofs, bfield_dofs)
 Helper function for add_current_update_v. More...
 
subroutine update_v_partial (self, upper, lower, uppert, lowert, index, field_dofs, bint)
 Helper function for evaluate_int_quad (despite the name it evaluates a field and does not update v) More...
 
subroutine update_v_partial_newton (self, upper, lower, uppert, lowert, index, field_dofs, bint)
 Helper function for evaluate_int_quad_subnewton (despite the name it evaluates a field and does not update v) More...
 
subroutine convert_x_to_xbox (self, position, xi, box)
 Helper function that identifies in which box the particle is found and its normalized poistion in the box. More...
 
subroutine init_spline_1d (self, domain, n_cells, no_particles, spline_degree, smoothing_type)
 Initializer. More...
 
subroutine free_spline_1d (self)
 Destructor. More...
 
subroutine, public sll_s_new_particle_mesh_coupling_spline_1d_ptr (smoother, domain, n_cells, no_particles, spline_degree, smoothing_type)
 
subroutine, public sll_s_new_particle_mesh_coupling_spline_1d (smoother, domain, n_cells, no_particles, spline_degree, smoothing_type)
 
subroutine add_current_split_spline_1d (self, position_old, position_new, iter, total_iter, marker_charge, j_dofs)
 Add current version with tau and (1-tau) which is needed for the implicit subcycling scheme. More...
 
subroutine update_current_partial (self, upper, lower, uppert, lowert, index, iter, total_iter, marker_charge, j_dofs)
 Helper function for add_current_split. More...
 
subroutine evaluate_int_linear_quad_spline_1d (self, position_old, position_new, iter, total_iter, field_dofs_1, field_dofs_2, field)
 Evaluates the integrals with tau function for implicit subcycling. More...
 
subroutine update_field_linear_partial (self, upper, lower, uppert, lowert, index, field_dofs_1, field_dofs_2, iter, total_iter, bint)
 Helper function for evaluate_int_linear_quad (implements part within one cell) More...
 
subroutine evaluate_int_linear_quad_subnewton_spline_1d (self, position_old, position_new, iter, total_iter, field_dofs_1, field_dofs_2, field)
 For implicit subcycling propagator, evaluates the integral needed for Newtons's method in the evaluate_int_linear_quad part. More...
 
subroutine update_field_linear_partial_newton (self, upper, lower, uppert, lowert, index, field_dofs_1, field_dofs_2, iter, total_iter, bint)
 Helper function for evaluate_int_linear_quad_subnewton (takes care of the per cell computations) More...
 
subroutine add_charge_int_spline_1d (self, position_old, position_new, marker_charge, j_dofs)
 For explicit version of the subcycling propagator (to evaluate the charge as it comes out integrated from the formulation) More...
 
subroutine update_charge_int_partial (self, upper, lower, uppert, lowert, index, marker_charge, j_dofs)
 Helper function for add_charge_int (takes care of the per cell computations) More...
 

Function/Subroutine Documentation

◆ add_charge_int_spline_1d()

subroutine sll_m_particle_mesh_coupling_spline_1d::add_charge_int_spline_1d ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
dimension(self%dim), intent(in)  position_old,
dimension(self%dim), intent(in)  position_new,
intent(in)  marker_charge,
dimension(self%n_cells), intent(inout)  j_dofs 
)
private

For explicit version of the subcycling propagator (to evaluate the charge as it comes out integrated from the formulation)

Parameters
[in,out]selfkernel smoother object

Definition at line 1378 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ add_charge_single_spline_1d()

subroutine sll_m_particle_mesh_coupling_spline_1d::add_charge_single_spline_1d ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
dimension(self%dim), intent(in)  position,
intent(in)  marker_charge,
dimension(self%n_cells), intent(inout)  rho_dofs 
)
private

helper function to compute the integral of j using Gauss quadrature (needs to be provided for the derived type with smoothing)

Add charge of one particle

Parameters
[in,out]selfkernel smoother object

Definition at line 103 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ add_current_evaluate_spline_1d()

subroutine sll_m_particle_mesh_coupling_spline_1d::add_current_evaluate_spline_1d ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
dimension(self%dim), intent(in)  position_old,
dimension(self%dim), intent(in)  position_new,
intent(in)  marker_charge,
intent(in)  vbar,
dimension(self%n_dofs), intent(in)  field_dofs,
dimension(self%n_dofs), intent(inout)  j_dofs,
intent(out)  field 
)
private

Combines add_current and evaluate_int.

Parameters
[in,out]selfkernel smoother object

Definition at line 346 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ add_current_spline_1d()

subroutine sll_m_particle_mesh_coupling_spline_1d::add_current_spline_1d ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
dimension(self%dim), intent(in)  position_old,
dimension(self%dim), intent(in)  position_new,
intent(in)  marker_charge,
dimension(self%n_cells), intent(inout)  j_dofs 
)
private

Add current with integration over x.

Parameters
[in,out]selfkernel smoother object

Definition at line 262 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ add_current_split_spline_1d()

subroutine sll_m_particle_mesh_coupling_spline_1d::add_current_split_spline_1d ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
real(kind=f64), dimension(self%dim), intent(in)  position_old,
real(kind=f64), dimension(self%dim), intent(in)  position_new,
integer(kind=i32), intent(in)  iter,
integer(kind=i32), intent(in)  total_iter,
real(kind=f64), intent(in)  marker_charge,
real(kind=f64), dimension(self%n_cells,2), intent(inout)  j_dofs 
)
private

Add current version with tau and (1-tau) which is needed for the implicit subcycling scheme.

Parameters
[in,out]selfkernel smoother object
[in]position_oldPosition of the particle
[in]position_newPosition of the particle
[in]marker_chargeParticle weights time charge
[in,out]j_dofsCoefficient vector of the current density

Definition at line 1053 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ add_current_update_v_spline_1d()

subroutine sll_m_particle_mesh_coupling_spline_1d::add_current_update_v_spline_1d ( class(sll_t_particle_mesh_coupling_spline_1d), intent(inout)  self,
real(kind=f64), dimension(self%dim), intent(in)  position_old,
real(kind=f64), dimension(self%dim), intent(in)  position_new,
real(kind=f64), intent(in)  marker_charge,
real(kind=f64), intent(in)  qoverm,
real(kind=f64), dimension(:), intent(in)  bfield_dofs,
real(kind=f64), dimension(:), intent(inout)  vi,
real(kind=f64), dimension(:), intent(inout)  j_dofs 
)
private

Add current for one particle and update v (according to H_p1 part in Hamiltonian splitting)

Parameters
[in,out]selfkernel smoother object
[in]position_oldPosition at time t
[in]position_newPosition at time t + \Delta t
[in]marker_chargeParticle weight time charge
[in]qovermcharge to mass ration
[in]bfield_dofsCoefficient of B-field expansion
[in,out]viVelocity of the particles
[in,out]j_dofsCoefficients of current expansion

Definition at line 700 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ add_current_update_v_spline_pp_1d()

subroutine sll_m_particle_mesh_coupling_spline_1d::add_current_update_v_spline_pp_1d ( class(sll_t_particle_mesh_coupling_spline_1d), intent(inout)  self,
real(kind=f64), dimension(self%dim), intent(in)  position_old,
real(kind=f64), dimension(self%dim), intent(in)  position_new,
real(kind=f64), intent(in)  marker_charge,
real(kind=f64), intent(in)  qoverm,
real(kind=f64), dimension(:), intent(in)  bfield_dofs,
real(kind=f64), dimension(:), intent(inout)  vi,
real(kind=f64), dimension(:), intent(inout)  j_dofs 
)
private

Add current for one particle and update v (according to H_p1 part in Hamiltonian splitting)

Parameters
[in,out]selfkernel smoother object
[in]position_oldPosition at time t
[in]position_newPosition at time t + \Delta t
[in]marker_chargeParticle weight time charge
[in]qovermcharge to mass ration
[in]bfield_dofsCoefficient of B-field expansion
[in,out]viVelocity of the particles
[in,out]j_dofsCoefficients of current expansion

Definition at line 621 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ add_particle_mass_spline_1d()

subroutine sll_m_particle_mesh_coupling_spline_1d::add_particle_mass_spline_1d ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
dimension(self%dim), intent(in)  position,
intent(in)  marker_charge,
dimension(:,:), intent(inout)  particle_mass 
)
private

Add the contribution of one particle to the approximate mass matrix.

Parameters
[in,out]selfkernel smoother object

Definition at line 130 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ add_particle_mass_spline_1d_full()

subroutine sll_m_particle_mesh_coupling_spline_1d::add_particle_mass_spline_1d_full ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
dimension(self%dim), intent(in)  position,
intent(in)  marker_charge,
dimension(:,:), intent(inout)  particle_mass 
)
private

Add the contribution of one particle to the approximate mass matrix.

Parameters
[in,out]selfkernel smoother object

Definition at line 163 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ convert_x_to_xbox()

subroutine sll_m_particle_mesh_coupling_spline_1d::convert_x_to_xbox ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
real(kind=f64), dimension(self%dim), intent(in)  position,
real(kind=f64), intent(out)  xi,
integer(kind=i32), intent(out)  box 
)
private

Helper function that identifies in which box the particle is found and its normalized poistion in the box.

Parameters
[in,out]selfkernel smoother object
[in]positionPosition of the particle
[out]xiPosition of the particle
[out]boxPosition of the particle

Definition at line 885 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the caller graph for this function:

◆ evaluate_field_single_spline_1d()

subroutine sll_m_particle_mesh_coupling_spline_1d::evaluate_field_single_spline_1d ( class (sll_t_particle_mesh_coupling_spline_1d), intent(inout)  self,
dimension(self%dim), intent(in)  position,
dimension(self%n_cells), intent(in)  field_dofs,
intent(out)  field_value 
)
private

Evaluate field at at position position.

Parameters
[in,out]selfKernel smoother object

Definition at line 197 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ evaluate_int_linear_quad_spline_1d()

subroutine sll_m_particle_mesh_coupling_spline_1d::evaluate_int_linear_quad_spline_1d ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
real(kind=f64), dimension(self%dim), intent(in)  position_old,
real(kind=f64), dimension(self%dim), intent(in)  position_new,
integer(kind=i32), intent(in)  iter,
integer(kind=i32), intent(in)  total_iter,
real(kind=f64), dimension(self%n_cells), intent(in)  field_dofs_1,
real(kind=f64), dimension(self%n_cells), intent(in)  field_dofs_2,
real(kind=f64), dimension(2), intent(out)  field 
)
private

Evaluates the integrals with tau function for implicit subcycling.

Parameters
[in,out]selfkernel smoother object
[in]position_oldPosition of the particle
[in]position_newPosition of the particle
[in]field_dofs_1Coefficient vector of the current density
[in]field_dofs_2Coefficient vector of the current density

Definition at line 1166 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ evaluate_int_linear_quad_subnewton_spline_1d()

subroutine sll_m_particle_mesh_coupling_spline_1d::evaluate_int_linear_quad_subnewton_spline_1d ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
dimension(self%dim), intent(in)  position_old,
dimension(self%dim), intent(in)  position_new,
intent(in)  iter,
intent(in)  total_iter,
dimension(self%n_cells), intent(in)  field_dofs_1,
dimension(self%n_cells), intent(in)  field_dofs_2,
intent(out)  field 
)
private

For implicit subcycling propagator, evaluates the integral needed for Newtons's method in the evaluate_int_linear_quad part.

Parameters
[in,out]selfkernel smoother object

Definition at line 1274 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ evaluate_int_quad_spline_1d()

subroutine sll_m_particle_mesh_coupling_spline_1d::evaluate_int_quad_spline_1d ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
dimension(self%dim), intent(in)  position_old,
dimension(self%dim), intent(in)  position_new,
dimension(self%n_cells), intent(in)  field_dofs,
intent(out)  field 
)
private

Evaluates an integrated field (between position_old and position_new) but with a quadrature formula for the integration.

Parameters
[in,out]selfkernel smoother object

Definition at line 528 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ evaluate_int_spline_1d()

subroutine sll_m_particle_mesh_coupling_spline_1d::evaluate_int_spline_1d ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
dimension(self%dim), intent(in)  position_old,
dimension(self%dim), intent(in)  position_new,
intent(in)  vbar,
dimension(self%n_cells), intent(in)  field_dofs,
intent(out)  field 
)
private

Evaluates the integral int_{poisition_old}^{position_new} field(x) d x.

Parameters
[in,out]selfkernel smoother object

Definition at line 446 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ evaluate_int_spline_1d_subnewton()

subroutine sll_m_particle_mesh_coupling_spline_1d::evaluate_int_spline_1d_subnewton ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
real(kind=f64), dimension(self%dim), intent(in)  position_old,
real(kind=f64), dimension(self%dim), intent(in)  position_new,
real(kind=f64), dimension(self%n_cells), intent(in)  field_dofs,
real(kind=f64), intent(out)  field 
)
private
Parameters
[in,out]selfkernel smoother object
[in]position_oldPosition of the particle
[in]position_newPosition of the particle
[in]field_dofsCoefficient vector of the current density

Definition at line 574 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ evaluate_multiple_spline_1d()

subroutine sll_m_particle_mesh_coupling_spline_1d::evaluate_multiple_spline_1d ( class (sll_t_particle_mesh_coupling_spline_1d), intent(inout)  self,
dimension(self%dim), intent(in)  position,
dimension(:), intent(in)  components,
dimension(:,:), intent(in)  field_dofs,
dimension(:), intent(out)  field_value 
)
private

Evaluate several fields at position position.

Parameters
[in,out]selfKernel smoother object

Definition at line 228 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ free_spline_1d()

subroutine sll_m_particle_mesh_coupling_spline_1d::free_spline_1d ( class (sll_t_particle_mesh_coupling_spline_1d), intent(inout)  self)
private

Destructor.

Parameters
[in,out]selfKernel smoother object

Definition at line 990 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ init_spline_1d()

subroutine sll_m_particle_mesh_coupling_spline_1d::init_spline_1d ( class( sll_t_particle_mesh_coupling_spline_1d), intent(out)  self,
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 
)
private

Initializer.

Parameters
[out]selfkernel smoother object
[in]n_cellsnumber of DoFs (spline coefficients)
[in]domainx_min and x_max of the domain
[in]no_particlesnumber of particles
[in]spline_degreeDegree of smoothing kernel spline
[in]smoothing_typeDefine if Galerkin or collocation smoothing for right scaling in accumulation routines

Definition at line 899 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:

◆ sll_s_new_particle_mesh_coupling_spline_1d()

subroutine, public sll_m_particle_mesh_coupling_spline_1d::sll_s_new_particle_mesh_coupling_spline_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 
)
Parameters
[out]smootherkernel smoother object
[in]n_cellsnumber of DoFs (spline coefficients)
[in]domainx_min and x_max of the domain
[in]no_particlesnumber of particles
[in]spline_degreeDegree of smoothing kernel spline
[in]smoothing_typeDefine if Galerkin or collocation smoothing for right scaling in accumulation routines

Definition at line 1027 of file sll_m_particle_mesh_coupling_spline_1d.F90.

◆ sll_s_new_particle_mesh_coupling_spline_1d_ptr()

subroutine, public sll_m_particle_mesh_coupling_spline_1d::sll_s_new_particle_mesh_coupling_spline_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 
)
Parameters
[out]smootherkernel smoother object
[in]n_cellsnumber of DoFs (spline coefficients)
[in]domainx_min and x_max of the domain
[in]no_particlesnumber of particles
[in]spline_degreeDegree of smoothing kernel spline
[in]smoothing_typeDefine if Galerkin or collocation smoothing for right scaling in accumulation routines

Definition at line 1003 of file sll_m_particle_mesh_coupling_spline_1d.F90.

◆ update_charge_int_partial()

subroutine sll_m_particle_mesh_coupling_spline_1d::update_charge_int_partial ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
intent(in)  upper,
intent(in)  lower,
intent(in)  uppert,
intent(in)  lowert,
intent(in)  index,
intent(in)  marker_charge,
dimension(self%n_cells), intent(inout)  j_dofs 
)
private

Helper function for add_charge_int (takes care of the per cell computations)

Parameters
[in,out]selfkernel smoother object

Definition at line 1434 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the caller graph for this function:

◆ update_current_partial()

subroutine sll_m_particle_mesh_coupling_spline_1d::update_current_partial ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
real(kind=f64), intent(in)  upper,
real(kind=f64), intent(in)  lower,
real(kind=f64), intent(in)  uppert,
real(kind=f64), intent(in)  lowert,
integer(kind=i32), intent(in)  index,
integer(kind=i32), intent(in)  iter,
integer(kind=i32), intent(in)  total_iter,
real(kind=f64), intent(in)  marker_charge,
real(kind=f64), dimension(self%n_cells,2), intent(inout)  j_dofs 
)
private

Helper function for add_current_split.

Parameters
[in,out]selfkernel smoother object
[in]marker_chargeParticle weights time charge
[in,out]j_dofsCoefficient vector of the current density

Definition at line 1111 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the caller graph for this function:

◆ update_field_linear_partial()

subroutine sll_m_particle_mesh_coupling_spline_1d::update_field_linear_partial ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
real(kind=f64), intent(in)  upper,
real(kind=f64), intent(in)  lower,
real(kind=f64), intent(in)  uppert,
real(kind=f64), intent(in)  lowert,
integer(kind=i32), intent(in)  index,
real(kind=f64), dimension(self%n_cells), intent(in)  field_dofs_1,
real(kind=f64), dimension(self%n_cells), intent(in)  field_dofs_2,
integer(kind=i32), intent(in)  iter,
integer(kind=i32), intent(in)  total_iter,
real(kind=f64), dimension(2), intent(inout)  bint 
)
private

Helper function for evaluate_int_linear_quad (implements part within one cell)

Parameters
[in,out]selfkernel smoother object
[in]field_dofs_1Coefficient vector of the current density
[in]field_dofs_2Coefficient vector of the current density

Definition at line 1230 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the caller graph for this function:

◆ update_field_linear_partial_newton()

subroutine sll_m_particle_mesh_coupling_spline_1d::update_field_linear_partial_newton ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
intent(in)  upper,
intent(in)  lower,
intent(in)  uppert,
intent(in)  lowert,
intent(in)  index,
dimension(self%n_cells), intent(in)  field_dofs_1,
dimension(self%n_cells), intent(in)  field_dofs_2,
intent(in)  iter,
intent(in)  total_iter,
intent(inout)  bint 
)
private

Helper function for evaluate_int_linear_quad_subnewton (takes care of the per cell computations)

Parameters
[in,out]selfkernel smoother object

Definition at line 1335 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the caller graph for this function:

◆ update_jv()

subroutine sll_m_particle_mesh_coupling_spline_1d::update_jv ( class(sll_t_particle_mesh_coupling_spline_1d), intent(inout)  self,
real(kind=f64), intent(in)  lower,
real(kind=f64), intent(in)  upper,
integer(kind=i32), intent(in)  index,
real(kind=f64), intent(in)  marker_charge,
real(kind=f64), intent(in)  qoverm,
real(kind=f64), intent(in)  sign,
real(kind=f64), intent(inout)  vi,
real(kind=f64), dimension(self%n_cells), intent(inout)  j_dofs,
real(kind=f64), dimension(self%n_cells), intent(in)  bfield_dofs 
)
private

Helper function for add_current_update_v.

Parameters
[in,out]selftime splitting object

Definition at line 751 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the caller graph for this function:

◆ update_jv_pp()

subroutine sll_m_particle_mesh_coupling_spline_1d::update_jv_pp ( class(sll_t_particle_mesh_coupling_spline_1d), intent(inout)  self,
real(kind=f64), intent(in)  lower,
real(kind=f64), intent(in)  upper,
integer(kind=i32), intent(in)  index,
real(kind=f64), intent(in)  marker_charge,
real(kind=f64), intent(in)  qoverm,
real(kind=f64), intent(inout)  vi,
real(kind=f64), dimension(self%n_cells), intent(inout)  j_dofs,
real(kind=f64), dimension(self%n_cells), intent(in)  bfield_dofs 
)
private

Helper function for add_current_update_v.

Parameters
[in,out]selftime splitting object

Definition at line 666 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ update_v_partial()

subroutine sll_m_particle_mesh_coupling_spline_1d::update_v_partial ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
real(kind=f64), intent(in)  upper,
real(kind=f64), intent(in)  lower,
real(kind=f64), intent(in)  uppert,
real(kind=f64), intent(in)  lowert,
integer(kind=i32), intent(in)  index,
real(kind=f64), dimension(self%n_cells), intent(in)  field_dofs,
real(kind=f64), intent(inout)  bint 
)
private

Helper function for evaluate_int_quad (despite the name it evaluates a field and does not update v)

Parameters
[in,out]selfkernel smoother object
[in]field_dofsCoefficient vector of the current density

Definition at line 794 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the caller graph for this function:

◆ update_v_partial_newton()

subroutine sll_m_particle_mesh_coupling_spline_1d::update_v_partial_newton ( class( sll_t_particle_mesh_coupling_spline_1d ), intent(inout)  self,
real(kind=f64), intent(in)  upper,
real(kind=f64), intent(in)  lower,
real(kind=f64), intent(in)  uppert,
real(kind=f64), intent(in)  lowert,
integer(kind=i32), intent(in)  index,
real(kind=f64), dimension(self%n_cells), intent(in)  field_dofs,
real(kind=f64), intent(inout)  bint 
)
private

Helper function for evaluate_int_quad_subnewton (despite the name it evaluates a field and does not update v)

Parameters
[in,out]selfkernel smoother object
[in]field_dofsCoefficient vector of the current density

Definition at line 839 of file sll_m_particle_mesh_coupling_spline_1d.F90.

Here is the caller graph for this function:
    Report Typos and Errors