8 #include "sll_assert.h"
9 #include "sll_memory.h"
10 #include "sll_working_precision.h"
11 #include "sll_errors.h"
17 sll_s_uniform_bsplines_eval_basis, &
18 sll_s_uniform_bsplines_eval_deriv
50 sll_int32 :: no_particles
53 sll_int32 :: n_quad_points
54 sll_int32 :: n_quadp1_points
55 sll_int32 :: n_quadp2_points
56 sll_int32 :: n_quads2_points
58 sll_real64,
allocatable :: spline_val(:)
59 sll_real64,
allocatable :: spline_pol_val(:)
60 sll_real64,
allocatable :: spline_val_more(:)
61 sll_real64,
allocatable :: quad_xw(:,:)
62 sll_real64,
allocatable :: quadp1_xw(:,:)
63 sll_real64,
allocatable :: quadp2_xw(:,:)
64 sll_real64,
allocatable :: quads1_xw(:,:)
65 sll_real64,
allocatable :: quads2_xw(:,:)
67 sll_int32 :: n_quad_points_line
68 sll_real64,
allocatable :: quad_xw_line(:,:)
105 sll_real64,
intent( in ) :: position(self%dim)
106 sll_real64,
intent( in ) :: marker_charge
107 sll_real64,
intent( inout ) :: rho_dofs(self%n_cells)
110 sll_int32 :: index1d, box
114 box = box - self%spline_degree
115 call sll_s_uniform_bsplines_eval_basis(self%spline_degree, xi, self%spline_val)
119 do i1 = 1, self%n_span
120 index1d = modulo(box+i1-2,self%n_cells)+1
121 rho_dofs(index1d) = rho_dofs(index1d) +&
122 (marker_charge * self%spline_val(i1)* self%scaling)
132 sll_real64,
intent( in ) :: position(self%dim)
133 sll_real64,
intent( in ) :: marker_charge
134 sll_real64,
intent( inout ) :: particle_mass(:,:)
137 sll_int32 :: i1, column
138 sll_int32 :: index1d, box
141 sll_assert(
size(particle_mass,1) == self%n_span )
142 sll_assert(
size(particle_mass,2) == self%n_cells )
145 box = box - self%spline_degree
146 call sll_s_uniform_bsplines_eval_basis(self%spline_degree, xi, self%spline_val)
150 do i1 = 1, self%n_span
151 index1d = modulo(box+i1-2,self%n_cells)+1
152 do column = i1, self%n_span
153 particle_mass(column-i1+1, index1d) = particle_mass( column-i1+1, index1d) + &
154 marker_charge * self%spline_val(i1) * self%spline_val(column)
165 sll_real64,
intent( in ) :: position(self%dim)
166 sll_real64,
intent( in ) :: marker_charge
167 sll_real64,
intent( inout ) :: particle_mass(:,:)
170 sll_int32 :: i1, column, ind
171 sll_int32 :: index1d, box
174 sll_assert(
size(particle_mass,1) == 2*self%spline_degree+1 )
175 sll_assert(
size(particle_mass,2) == self%n_cells )
178 box = box - self%spline_degree
179 call sll_s_uniform_bsplines_eval_basis(self%spline_degree, xi, self%spline_val)
182 do i1 = 1, self%n_span
183 index1d = modulo(box+i1-2,self%n_cells)+1
184 ind=1+(self%n_span-i1)
185 do column = 1, self%n_span
186 particle_mass(ind, index1d) = particle_mass( ind, index1d) + &
187 marker_charge * self%spline_val(i1) * self%spline_val(column)
199 sll_real64,
intent( in ) :: position(self%dim)
200 sll_real64,
intent( in ) :: field_dofs(self%n_cells)
201 sll_real64,
intent( out ) :: field_value
205 sll_int32 :: index1d, box
209 box = box - self%spline_degree
211 call sll_s_uniform_bsplines_eval_basis(self%spline_degree, xi, self%spline_val)
215 field_value = 0.0_f64
216 do i1 = 1, self%n_span
217 index1d = modulo(box+i1-2, self%n_cells)+1
218 field_value = field_value + &
219 field_dofs(index1d) * &
230 sll_real64,
intent( in ) :: position(self%dim)
231 sll_int32,
intent(in) :: components(:)
232 sll_real64,
intent( in ) :: field_dofs(:,:)
233 sll_real64,
intent(out) :: field_value(:)
237 sll_int32 :: index1d, box
240 sll_assert(
size(field_dofs,1) == self%n_cells )
241 sll_assert(
size(field_dofs,2) ==
size(field_value) )
245 box = box - self%spline_degree
246 call sll_s_uniform_bsplines_eval_basis(self%spline_degree, xi, self%spline_val)
250 field_value = 0.0_f64
251 do i1 = 1, self%n_span
252 index1d = modulo(box+i1-2, self%n_cells)+1
253 field_value = field_value + &
254 field_dofs(index1d,components) * &
264 sll_real64,
intent( in ) :: position_old(self%dim)
265 sll_real64,
intent( in ) :: position_new(self%dim)
266 sll_real64,
intent( in ) :: marker_charge
267 sll_real64,
intent( inout ) :: j_dofs(self%n_cells)
269 sll_real64 :: dx_new, dx_old
270 sll_int32 :: box_new, box_old, i, index, j, ind
285 if (box_old<box_new)
then
287 do ind = box_old-self%spline_degree, min(box_new-self%spline_degree-1, box_old )
288 index = modulo(ind-1,self%n_cells)+1
290 j_dofs(index) = j_dofs(index) + (marker_charge * &
291 (self%delta_x-self%spline_val(i)))
294 do ind = box_new-self%spline_degree, box_old
297 index = modulo(ind-1,self%n_cells)+1
298 j_dofs(index) = j_dofs(index) + marker_charge * &
299 ((self%spline_val_more(j) - self%spline_val(i)))
301 do ind = max(box_old+1, box_new-self%spline_degree), box_new
303 index = modulo(ind-1,self%n_cells)+1
304 j_dofs(index) = j_dofs(index) + (marker_charge * &
305 (self%spline_val_more(j) ))
308 do ind = box_old+1, box_new-self%spline_degree-1
309 index = modulo(ind-1, self%n_cells)+1
310 j_dofs(index) = j_dofs(index) + marker_charge * self%delta_x
315 do ind = box_new-self%spline_degree, min(box_old-self%spline_degree-1, box_new)
316 index = modulo(ind-1,self%n_cells)+1
318 j_dofs(index) = j_dofs(index) + (marker_charge * &
319 (-self%delta_x+self%spline_val_more(j)))
322 do ind = box_old-self%spline_degree, box_new
325 index = modulo(ind-1,self%n_cells)+1
326 j_dofs(index) = j_dofs(index) + (marker_charge * &
327 (self%spline_val_more(j) - self%spline_val(i)))
329 do ind = max(box_new+1, box_old-self%spline_degree), box_old
331 index = modulo(ind-1,self%n_cells)+1
332 j_dofs(index) = j_dofs(index) + (marker_charge * &
333 (-self%spline_val(i) ))
336 do ind = box_new+1, box_old-self%spline_degree-1
337 index = modulo(ind-1, self%n_cells)+1
338 j_dofs(index) = j_dofs(index) - marker_charge * self%delta_x
348 sll_real64,
intent( in ) :: position_old(self%dim)
349 sll_real64,
intent( in ) :: position_new(self%dim)
350 sll_real64,
intent( in ) :: marker_charge
351 sll_real64,
intent( in ) :: vbar
352 sll_real64,
intent( in ) :: field_dofs(self%n_dofs)
353 sll_real64,
intent( inout ) :: j_dofs(self%n_dofs)
354 sll_real64,
intent( out ) :: field
356 sll_real64 :: dx_new, dx_old
357 sll_int32 :: box_new, box_old, i, index, ind, j
374 if (box_old<box_new)
then
376 do ind = box_old-self%spline_degree, min(box_new-self%spline_degree-1, box_old)
377 index = modulo(ind-1,self%n_cells)+1
379 j_dofs(index) = j_dofs(index) + (marker_charge * &
380 (self%delta_x-self%spline_val(i)))
381 field = field + ((self%delta_x - self%spline_val(i))) * field_dofs(index)
384 do ind = box_new-self%spline_degree, box_old
387 index = modulo(ind-1,self%n_cells)+1
388 j_dofs(index) = j_dofs(index) + marker_charge * &
389 ((self%spline_val_more(j) - self%spline_val(i)))
390 field = field + ((self%spline_val_more(j) - self%spline_val(i))) * field_dofs(index)
393 do ind = max(box_old+1, box_new-self%spline_degree), box_new
395 index = modulo(ind-1,self%n_cells)+1
396 j_dofs(index) = j_dofs(index) + (marker_charge * &
397 (self%spline_val_more(j) ))
398 field = field + (self%spline_val_more(j) ) * field_dofs(index)
401 do ind = box_old+1, box_new-self%spline_degree-1
402 index = modulo(ind-1, self%n_cells)+1
403 j_dofs(index) = j_dofs(index) + marker_charge * self%delta_x
404 field = field + self%delta_x * field_dofs(index)
409 do ind = box_new-self%spline_degree, min(box_old-self%spline_degree-1, box_new)
410 index = modulo(ind-1,self%n_cells)+1
412 j_dofs(index) = j_dofs(index) + (marker_charge * &
413 (-self%delta_x+self%spline_val_more(j)))
414 field = field + (-self%delta_x+self%spline_val_more(j)) * field_dofs(index)
417 do ind = box_old-self%spline_degree, box_new
420 index = modulo(ind-1,self%n_cells)+1
421 j_dofs(index) = j_dofs(index) + (marker_charge * &
422 (self%spline_val_more(j) - self%spline_val(i)))
423 field = field + (self%spline_val_more(j) - self%spline_val(i)) * field_dofs(index)
425 do ind = max(box_new+1, box_old-self%spline_degree), box_old
427 index = modulo(ind-1,self%n_cells)+1
428 j_dofs(index) = j_dofs(index) + (marker_charge * &
429 (-self%spline_val(i) ))
430 field = field + (-self%spline_val(i) ) * field_dofs(index)
433 do ind = box_new+1, box_old-self%spline_degree-1
434 index = modulo(ind-1, self%n_cells)+1
435 j_dofs(index) = j_dofs(index) - marker_charge * self%delta_x
436 field = field - self%delta_x * field_dofs(index)
448 sll_real64,
intent( in ) :: position_old(self%dim)
449 sll_real64,
intent( in ) :: position_new(self%dim)
450 sll_real64,
intent( in ) :: vbar
451 sll_real64,
intent( in ) :: field_dofs(self%n_cells)
452 sll_real64,
intent( out ) :: field
454 sll_real64 :: dx_new, dx_old
455 sll_int32 :: box_new, box_old, i, index, ind, j
472 if (box_old<box_new)
then
474 do ind = box_old-self%spline_degree, min(box_new-self%spline_degree-1, box_old)
475 index = modulo(ind-1,self%n_cells)+1
477 field = field + ((self%delta_x - self%spline_val(i))) * field_dofs(index)
480 do ind = box_new-self%spline_degree, box_old
483 index = modulo(ind-1,self%n_cells)+1
484 field = field + ((self%spline_val_more(j) - self%spline_val(i))) * field_dofs(index)
486 do ind = max(box_old+1, box_new-self%spline_degree), box_new
488 index = modulo(ind-1,self%n_cells)+1
489 field = field + (self%spline_val_more(j) ) * field_dofs(index)
491 do ind = box_old+1, box_new-self%spline_degree-1
492 index = modulo(ind-1, self%n_cells)+1
493 field = field + self%delta_x * field_dofs(index)
498 do ind = box_new-self%spline_degree, min(box_old-self%spline_degree-1, box_new)
499 index = modulo(ind-1,self%n_cells)+1
501 field = field + (-self%delta_x+self%spline_val_more(j)) * field_dofs(index)
504 do ind = box_old-self%spline_degree, box_new
507 index = modulo(ind-1,self%n_cells)+1
508 field = field + (self%spline_val_more(j) - self%spline_val(i)) * field_dofs(index)
510 do ind = max(box_new+1,box_old-self%spline_degree), box_old
512 index = modulo(ind-1,self%n_cells)+1
513 field = field + (-self%spline_val(i) ) * field_dofs(index)
516 do ind = box_new+1, box_old-self%spline_degree-1
517 index = modulo(ind-1, self%n_cells)+1
518 field = field - self%delta_x * field_dofs(index)
530 sll_real64,
intent( in ) :: position_old(self%dim)
531 sll_real64,
intent( in ) :: position_new(self%dim)
532 sll_real64,
intent( in ) :: field_dofs(self%n_cells)
533 sll_real64,
intent( out ) :: field
535 sll_real64 :: dx_new, dx_old
536 sll_int32 :: box_new, box_old, ind
537 sll_real64 :: rinterval, lowert, uppert
542 rinterval = self%delta_x/(position_new(1)-position_old(1))
546 if ( box_old == box_new )
then
547 call update_v_partial(self, dx_new, dx_old, 1.0_f64, 0.0_f64, box_old, field_dofs, field )
548 elseif ( box_old < box_new )
then
550 uppert = (1.0_f64-dx_old)*rinterval
551 call update_v_partial( self, 1.0_f64, dx_old, uppert, lowert, box_old, field_dofs, field )
552 do ind = box_old+1, box_new-1
554 uppert = uppert + rinterval
555 call update_v_partial( self, 1.0_f64, 0.0_f64, uppert, lowert, ind, field_dofs, field )
557 call update_v_partial( self, dx_new, 0.0_f64, 1.0_f64, uppert, box_new, field_dofs, field )
560 uppert = 1.0_f64 + (rinterval-dx_new*rinterval)
561 call update_v_partial( self, 1.0_f64, dx_new, uppert, lowert, box_new, field_dofs, field )
562 do ind = box_new+1, box_old-1
564 uppert = lowert + rinterval
565 call update_v_partial( self, 1.0_f64, 0.0_f64, uppert, lowert, ind, field_dofs, field )
567 call update_v_partial( self, dx_old, 0.0_f64, 0.0_f64, uppert, box_old, field_dofs, field )
576 sll_real64,
intent( in ) :: position_old(self%dim)
577 sll_real64,
intent( in ) :: position_new(self%dim)
578 sll_real64,
intent( in ) :: field_dofs(self%n_cells)
579 sll_real64,
intent( out ) :: field
581 sll_real64 :: dx_new, dx_old
582 sll_int32 :: box_new, box_old, ind
583 sll_real64 :: rinterval, lowert, uppert
588 rinterval = self%delta_x/(position_new(1)-position_old(1))
592 if ( box_old == box_new )
then
594 elseif ( box_old < box_new )
then
596 uppert = (1.0_f64-dx_old)*rinterval
598 do ind = box_old+1, box_new-1
600 uppert = uppert + rinterval
606 uppert = 1.0_f64 + (rinterval-dx_new*rinterval)
608 do ind = box_new+1, box_old-1
610 uppert = lowert + rinterval
623 sll_real64,
intent(in) :: position_old(self%dim)
624 sll_real64,
intent(in) :: position_new(self%dim)
625 sll_real64,
intent(in) :: marker_charge
626 sll_real64,
intent(in) :: qoverm
627 sll_real64,
intent(in) :: bfield_dofs(:)
628 sll_real64,
intent(inout) :: vi(:)
629 sll_real64,
intent(inout) :: j_dofs(:)
631 sll_int32 :: index_old, index_new, ind
632 sll_real64 :: r_old, r_new
637 index_old = index_old - 1
638 index_new = index_new - 1
640 if (index_old == index_new)
then
641 call update_jv_pp( self, r_old, r_new, index_old, marker_charge, &
642 qoverm, vi(2), j_dofs, bfield_dofs)
643 elseif (index_old < index_new)
then
644 call update_jv_pp ( self, r_old, 1.0_f64, index_old, marker_charge, &
645 qoverm, vi(2), j_dofs, bfield_dofs)
646 call update_jv_pp ( self, 0.0_f64, r_new, index_new, marker_charge, &
647 qoverm, vi(2), j_dofs, bfield_dofs)
648 do ind = index_old+1, index_new-1
649 call update_jv_pp ( self, 0.0_f64, 1.0_f64, ind, marker_charge, &
650 qoverm, vi(2), j_dofs, bfield_dofs)
653 call update_jv_pp ( self, 1.0_f64,r_new, index_new, marker_charge, qoverm, &
654 vi(2), j_dofs, bfield_dofs)
655 call update_jv_pp ( self, r_old, 0.0_f64, index_old, marker_charge, qoverm, &
656 vi(2), j_dofs, bfield_dofs)
657 do ind = index_new+1, index_old-1
658 call update_jv_pp ( self, 1.0_f64,0.0_f64, ind, marker_charge, qoverm, &
659 vi(2), j_dofs, bfield_dofs)
666 subroutine update_jv_pp(self, lower, upper, index, marker_charge, qoverm, vi, j_dofs, bfield_dofs)
668 sll_real64,
intent(in) :: lower
669 sll_real64,
intent(in) :: upper
670 sll_int32,
intent(in) :: index
671 sll_real64,
intent(in) :: marker_charge
672 sll_real64,
intent(in) :: qoverm
673 sll_real64,
intent(inout) :: vi
674 sll_real64,
intent(in) :: bfield_dofs(self%n_cells)
675 sll_real64,
intent(inout) :: j_dofs(self%n_cells)
677 sll_int32 :: ind, i_grid, i_mod, n_cells
679 n_cells = self%n_cells
684 self%spline_val = (self%spline_val_more - self%spline_val) *(self%delta_x)
688 do i_grid = index - self%spline_degree , index
689 i_mod = modulo(i_grid, n_cells ) + 1
690 j_dofs(i_mod) = j_dofs(i_mod) + &
691 (marker_charge*self%spline_val(ind)* self%scaling)
692 vi = vi - qoverm* self%spline_val(ind)*bfield_dofs(i_mod)
702 sll_real64,
intent(in) :: position_old(self%dim)
703 sll_real64,
intent(in) :: position_new(self%dim)
704 sll_real64,
intent(in) :: marker_charge
705 sll_real64,
intent(in) :: qoverm
706 sll_real64,
intent(in) :: bfield_dofs(:)
707 sll_real64,
intent(inout) :: vi(:)
708 sll_real64,
intent(inout) :: j_dofs(:)
710 sll_int32 :: index_old, index_new, ind
711 sll_real64 :: r_old, r_new
716 index_old = index_old - 1
717 index_new = index_new - 1
719 if (index_old == index_new)
then
720 if (r_old < r_new)
then
721 call update_jv( self, r_old, r_new, index_old, marker_charge, &
722 qoverm, 1.0_f64, vi(2), j_dofs, bfield_dofs)
724 call update_jv( self, r_new, r_old, index_old, marker_charge, qoverm, &
725 -1.0_f64, vi(2), j_dofs, bfield_dofs)
727 elseif (index_old < index_new)
then
728 call update_jv ( self, r_old, 1.0_f64, index_old, marker_charge, &
729 qoverm, 1.0_f64, vi(2), j_dofs, bfield_dofs)
730 call update_jv ( self, 0.0_f64, r_new, index_new, marker_charge, &
731 qoverm, 1.0_f64, vi(2), j_dofs, bfield_dofs)
732 do ind = index_old+1, index_new-1
733 call update_jv ( self, 0.0_f64, 1.0_f64, ind, marker_charge, &
734 qoverm, 1.0_f64, vi(2), j_dofs, bfield_dofs)
737 call update_jv ( self, r_new, 1.0_f64, index_new, marker_charge, qoverm, &
738 -1.0_f64, vi(2), j_dofs, bfield_dofs)
739 call update_jv ( self, 0.0_f64, r_old, index_old, marker_charge, qoverm, &
740 -1.0_f64, vi(2), j_dofs, bfield_dofs)
741 do ind = index_new+1, index_old-1
742 call update_jv ( self, 0.0_f64, 1.0_f64, ind, marker_charge, qoverm, &
743 -1.0_f64, vi(2), j_dofs, bfield_dofs)
751 subroutine update_jv(self, lower, upper, index, marker_charge, qoverm, sign, vi, j_dofs, bfield_dofs)
753 sll_real64,
intent(in) :: lower
754 sll_real64,
intent(in) :: upper
755 sll_int32,
intent(in) :: index
756 sll_real64,
intent(in) :: marker_charge
757 sll_real64,
intent(in) :: qoverm
758 sll_real64,
intent(in) :: sign
759 sll_real64,
intent(inout) :: vi
760 sll_real64,
intent(in) :: bfield_dofs(self%n_cells)
761 sll_real64,
intent(inout) :: j_dofs(self%n_cells)
763 sll_int32 :: ind, i_grid, i_mod, n_cells, j
766 n_cells = self%n_cells
768 c1 = 0.5_f64*(upper-lower)
769 c2 = 0.5_f64*(upper+lower)
772 call sll_s_uniform_bsplines_eval_basis(self%spline_degree, c1*self%quad_xw(1,1)+c2, &
774 self%spline_val = self%spline_val * (self%quad_xw(2,1)*c1)
775 do j=2,self%n_quad_points
776 call sll_s_uniform_bsplines_eval_basis(self%spline_degree, c1*self%quad_xw(1,j)+c2, &
777 self%spline_val_more)
778 self%spline_val = self%spline_val + self%spline_val_more * (self%quad_xw(2,j)*c1)
780 self%spline_val = self%spline_val * (sign*self%delta_x)
783 do i_grid = index - self%spline_degree , index
784 i_mod = modulo(i_grid, n_cells ) + 1
785 j_dofs(i_mod) = j_dofs(i_mod) + &
786 (marker_charge*self%spline_val(ind)* self%scaling)
787 vi = vi - qoverm* self%spline_val(ind)*bfield_dofs(i_mod)
796 sll_real64,
intent( in ) :: upper
797 sll_real64,
intent( in ) :: lower
798 sll_real64,
intent( in ) :: uppert
799 sll_real64,
intent( in ) :: lowert
800 sll_int32,
intent( in ) :: index
801 sll_real64,
intent( in ) :: field_dofs(self%n_cells)
802 sll_real64,
intent( inout ) :: bint
804 sll_int32 :: j, ind, i_mod, i_grid, n_cells
805 sll_real64 :: c1, c2, c1t, c2t
808 n_cells = self%n_cells
810 c1 = 0.5_f64*(upper-lower)
811 c2 = 0.5_f64*(upper+lower)
813 c1t = 0.5_f64*(uppert-lowert)
814 c2t = 0.5_f64*(uppert+lowert)
816 call sll_s_uniform_bsplines_eval_basis ( self%spline_degree, c1*self%quadp1_xw(1,1) + c2, &
818 self%spline_val = self%spline_val * (self%quadp1_xw(2,1)*c1t) * ( c1t* self%quadp1_xw(1,1) + c2t )
820 do j=2, self%n_quadp1_points
821 call sll_s_uniform_bsplines_eval_basis( self%spline_degree, c1*self%quadp1_xw(1,j)+c2, &
822 self%spline_val_more )
823 self%spline_val = self%spline_val + self%spline_val_more * (self%quadp1_xw(2,j)*c1t) * ( c1t* self%quadp1_xw(1,j) + c2t )
825 self%spline_val = self%spline_val
828 do i_grid = index - self%spline_degree, index
829 i_mod = modulo(i_grid-1, n_cells ) + 1
830 bint = bint + self%spline_val(ind) * field_dofs( i_mod )
841 sll_real64,
intent( in ) :: upper
842 sll_real64,
intent( in ) :: lower
843 sll_real64,
intent( in ) :: uppert
844 sll_real64,
intent( in ) :: lowert
845 sll_int32,
intent( in ) :: index
846 sll_real64,
intent( in ) :: field_dofs(self%n_cells)
847 sll_real64,
intent( inout ) :: bint
849 sll_int32 :: j, ind, i_mod, i_grid, n_cells
850 sll_real64 :: c1, c2, c1t, c2t, tau
852 n_cells = self%n_cells
854 c1 = 0.5_f64*(upper-lower)
855 c2 = 0.5_f64*(upper+lower)
857 c1t = 0.5_f64*(uppert-lowert)
858 c2t = 0.5_f64*(uppert+lowert)
861 call sll_s_uniform_bsplines_eval_deriv ( self%spline_degree, c1*self%quadp1_xw(1,1) + c2, &
863 tau = ( c1t* self%quadp1_xw(1,1) + c2t )
864 self%spline_val = self%spline_val * (self%quadp1_xw(2,1)*c1t) * ( 1.0_f64 - tau ) * tau
866 do j=2, self%n_quadp1_points
867 call sll_s_uniform_bsplines_eval_deriv( self%spline_degree, c1*self%quadp1_xw(1,j)+c2, &
868 self%spline_val_more )
869 tau = ( c1t* self%quadp1_xw(1,j) + c2t )
870 self%spline_val = self%spline_val + self%spline_val_more * (self%quadp1_xw(2,j)*c1t) * ( 1.0_f64 - tau ) * tau
872 self%spline_val = self%spline_val
875 do i_grid = index - self%spline_degree, index
876 i_mod = modulo(i_grid-1, n_cells ) + 1
877 bint = bint + self%spline_val(ind) * field_dofs( i_mod )
887 sll_real64,
intent( in ) :: position(self%dim)
888 sll_real64,
intent( out ) :: xi
889 sll_int32,
intent( out ) :: box
891 xi = (position(1) - self%domain(1)) /self%delta_x
892 box = floor( xi ) + 1
893 xi = xi - real(box-1, f64)
899 subroutine init_spline_1d( self, domain, n_cells, no_particles, spline_degree, smoothing_type )
901 sll_int32,
intent(in) :: n_cells
902 sll_real64,
intent(in) :: domain(2)
903 sll_int32,
intent(in) :: no_particles
904 sll_int32,
intent(in) :: spline_degree
905 sll_int32,
intent(in) :: smoothing_type
914 self%n_cells = n_cells
915 self%n_dofs = n_cells
916 self%delta_x = (self%domain(2)-self%domain(1))/real(n_cells, f64)
919 self%no_particles = no_particles
922 self%spline_degree = spline_degree
923 self%n_span = spline_degree + 1
927 self%scaling = 1.0_f64/self%delta_x
929 self%scaling = 1.0_f64
931 print*,
'Smoothing Type ', smoothing_type,
' not implemented for kernel_smoother_spline_1d.'
934 self%n_quad_points = (self%spline_degree+2)/2
936 ALLOCATE( self%spline_val(self%n_span), stat = ierr)
937 sll_assert( ierr == 0 )
938 ALLOCATE( self%spline_val_more(self%n_span), stat = ierr )
939 sll_assert( ierr == 0 )
940 ALLOCATE( self%quad_xw(2,self%n_quad_points), stat = ierr )
941 sll_assert( ierr == 0 )
948 self%n_quadp1_points = (self%spline_degree+3)/2
949 allocate( self%quadp1_xw(2, self%n_quadp1_points) )
950 allocate( self%spline_pol_val(self%n_span) )
955 self%n_quadp2_points = (self%spline_degree+4)/2
956 allocate( self%quadp2_xw(2, self%n_quadp2_points) )
962 allocate( self%quads1_xw(2, self%n_quadp1_points) )
964 self%quads1_xw = self%quadp1_xw
965 self%quads1_xw(1,:) = (self%quads1_xw(1,:) + 1.0_f64)*0.5_f64
966 self%quads1_xw(2,:) = self%quads1_xw(2,:)*0.5_f64
971 self%n_quads2_points = (self%spline_degree+4)/2
972 allocate( self%quads2_xw(2, self%n_quads2_points) )
975 self%quads2_xw(1,:) = (self%quads2_xw(1,:) + 1.0_f64)*0.5_f64
976 self%quads2_xw(2,:) = self%quads2_xw(2,:)*0.5_f64
978 self%n_quad_points_line = (self%spline_degree+2)/2
980 allocate( self%quad_xw_line(2,self%n_quad_points_line) )
984 self%quad_xw_line(1,:) = 0.5_f64*(self%quad_xw_line(1,:)+1.0_f64)
985 self%quad_xw_line(2,:) = 0.5_f64*self%quad_xw_line(2,:)
993 deallocate(self%spline_val)
994 deallocate(self%spline_val_more)
995 deallocate(self%quad_xw)
1005 sll_int32,
intent(in) :: n_cells
1006 sll_real64,
intent(in) :: domain(2)
1007 sll_int32,
intent(in) :: no_particles
1008 sll_int32,
intent(in) :: spline_degree
1009 sll_int32,
intent(in) :: smoothing_type
1016 sll_assert( ierr == 0)
1018 select type( smoother )
1020 call smoother%init( domain, n_cells, no_particles, spline_degree, smoothing_type )
1029 sll_int32,
intent(in) :: n_cells
1030 sll_real64,
intent(in) :: domain(2)
1031 sll_int32,
intent(in) :: no_particles
1032 sll_int32,
intent(in) :: spline_degree
1033 sll_int32,
intent(in) :: smoothing_type
1040 sll_assert( ierr == 0)
1042 select type( smoother )
1044 call smoother%init( domain, n_cells, no_particles, spline_degree, smoothing_type )
1055 sll_real64,
intent( in ) :: position_old(self%dim)
1056 sll_real64,
intent( in ) :: position_new(self%dim)
1057 sll_int32,
intent( in ) :: iter
1058 sll_int32,
intent( in ) :: total_iter
1059 sll_real64,
intent( in ) :: marker_charge
1060 sll_real64,
intent( inout ) :: j_dofs(self%n_cells,2)
1062 sll_real64 :: dx_new, dx_old
1063 sll_int32 :: box_new, box_old, ind
1064 sll_real64 :: rinterval, lowert, uppert
1066 dx_new = (position_new(1)-self%domain(1))/self%delta_x
1067 box_new = floor(dx_new)+1
1068 dx_new = dx_new - real(box_new-1,f64)
1069 dx_old = (position_old(1)-self%domain(1))/self%delta_x
1070 box_old = floor(dx_old)+1
1071 dx_old = dx_old - real(box_old-1,f64)
1073 if ( box_old == box_new )
then
1075 iter, total_iter, marker_charge, j_dofs )
1076 elseif ( box_old < box_new )
then
1077 rinterval = self%delta_x/(position_new(1)-position_old(1))
1079 uppert = (1.0_f64-dx_old)*rinterval
1081 iter, total_iter, marker_charge, j_dofs )
1082 do ind = box_old+1, box_new-1
1084 uppert = uppert + rinterval
1086 iter, total_iter, marker_charge, j_dofs )
1089 iter, total_iter, marker_charge, j_dofs )
1091 rinterval = self%delta_x/(position_new(1)-position_old(1))
1093 uppert = 1.0_f64 + (rinterval - dx_new*rinterval)
1096 iter, total_iter, -marker_charge, j_dofs )
1097 do ind = box_new+1, box_old-1
1099 uppert = lowert + rinterval
1101 iter, total_iter, -marker_charge, j_dofs )
1104 iter, total_iter, -marker_charge, j_dofs )
1113 sll_real64,
intent( in ) :: upper
1114 sll_real64,
intent( in ) :: lower
1115 sll_real64,
intent( in ) :: uppert
1116 sll_real64,
intent( in ) :: lowert
1117 sll_int32,
intent( in ) :: index
1118 sll_int32,
intent( in ) :: iter
1119 sll_int32,
intent( in ) :: total_iter
1120 sll_real64,
intent( in ) :: marker_charge
1121 sll_real64,
intent( inout ) :: j_dofs(self%n_cells,2)
1124 sll_int32 :: j, ind, i_mod, i_grid, n_cells
1125 sll_real64 :: c1, c2, c1t, c2t
1128 n_cells = self%n_cells
1130 c1 = 0.5_f64*(upper-lower)
1131 c2 = 0.5_f64*(upper+lower)
1133 c1t = 0.5_f64*(uppert-lowert)
1134 c2t = 0.5_f64*(uppert+lowert)
1136 call sll_s_uniform_bsplines_eval_basis ( self%spline_degree, &
1137 c1*self%quadp1_xw(1,1) + c2, &
1139 time = (c1t* self%quadp1_xw(1,1) + c2t + real(iter,f64))/real(total_iter, f64)
1140 self%spline_pol_val = self%spline_val * (self%quadp1_xw(2,1)*c1t) * &
1142 self%spline_val = self%spline_val * (self%quadp1_xw(2,1)*c1t) * time
1144 do j=2, self%n_quadp1_points
1145 call sll_s_uniform_bsplines_eval_basis( self%spline_degree, &
1146 c1*self%quadp1_xw(1,j)+c2, &
1147 self%spline_val_more )
1148 time = (c1t* self%quadp1_xw(1,j) + c2t + real(iter,f64))/real(total_iter, f64)
1149 self%spline_val = self%spline_val + &
1150 self%spline_val_more * (self%quadp1_xw(2,j)*c1t) * time
1151 self%spline_pol_val = self%spline_pol_val + &
1152 self%spline_val_more * (self%quadp1_xw(2,j)*c1t) * ( 1.0_f64 - time )
1156 do i_grid = index - self%spline_degree, index
1157 i_mod = modulo(i_grid-1, n_cells ) + 1
1158 j_dofs(i_mod,1) = j_dofs(i_mod,1) + marker_charge * self%spline_val(ind)
1159 j_dofs(i_mod,2) = j_dofs(i_mod,2) + marker_charge * self%spline_pol_val(ind)
1168 sll_real64,
intent( in ) :: position_old(self%dim)
1169 sll_real64,
intent( in ) :: position_new(self%dim)
1170 sll_real64,
intent( in ) :: field_dofs_1(self%n_cells)
1171 sll_real64,
intent( in ) :: field_dofs_2(self%n_cells)
1172 sll_int32,
intent( in ) :: iter
1173 sll_int32,
intent( in ) :: total_iter
1174 sll_real64,
intent( out ) :: field(2)
1176 sll_real64 :: dx_new, dx_old
1177 sll_int32 :: box_new, box_old, ind
1178 sll_real64 :: rinterval, lowert, uppert
1180 dx_new = (position_new(1)-self%domain(1))/self%delta_x
1181 box_new = floor(dx_new)+1
1182 dx_new = dx_new - real(box_new-1,f64)
1183 dx_old = (position_old(1)-self%domain(1))/self%delta_x
1184 box_old = floor(dx_old)+1
1185 dx_old = dx_old - real(box_old-1,f64)
1187 if ( box_old == box_new )
then
1190 rinterval = self%delta_x/(position_new(1)-position_old(1))
1194 if ( box_old == box_new )
then
1196 field_dofs_1, field_dofs_2, iter, total_iter, field )
1197 elseif ( box_old < box_new )
then
1199 uppert = (1.0_f64-dx_old)*rinterval
1201 field_dofs_1, field_dofs_2, iter, total_iter, field )
1202 do ind = box_old+1, box_new-1
1204 uppert = uppert + rinterval
1206 field_dofs_1, field_dofs_2, iter, total_iter, field )
1209 field_dofs_1, field_dofs_2, iter, total_iter, field )
1212 uppert = 1.0_f64 + (rinterval-dx_new*rinterval)
1214 field_dofs_1, field_dofs_2, iter, total_iter, field )
1215 do ind = box_new+1, box_old-1
1217 uppert = lowert + rinterval
1219 field_dofs_1, field_dofs_2, iter, total_iter, field )
1222 field_dofs_1, field_dofs_2, iter, total_iter, field )
1230 subroutine update_field_linear_partial(self, upper, lower, uppert, lowert, index, field_dofs_1, field_dofs_2, iter, total_iter, bint )
1232 sll_real64,
intent( in ) :: upper
1233 sll_real64,
intent( in ) :: lower
1234 sll_real64,
intent( in ) :: uppert
1235 sll_real64,
intent( in ) :: lowert
1236 sll_int32,
intent( in ) :: index
1237 sll_real64,
intent( in ) :: field_dofs_1(self%n_cells)
1238 sll_real64,
intent( in ) :: field_dofs_2(self%n_cells)
1239 sll_int32,
intent( in ) :: iter
1240 sll_int32,
intent( in ) :: total_iter
1241 sll_real64,
intent( inout ) :: bint(2)
1243 sll_int32 :: j, ind, i_mod, i_grid, n_cells
1244 sll_real64 :: c1, c2, c1t, c2t
1245 sll_real64 :: tau, time, weight, field_dof_t
1247 n_cells = self%n_cells
1249 c1 = 0.5_f64*(upper-lower)
1250 c2 = 0.5_f64*(upper+lower)
1252 c1t = 0.5_f64*(uppert-lowert)
1253 c2t = 0.5_f64*(uppert+lowert)
1255 do j=1, self%n_quadp2_points
1256 tau = c1t* self%quadp2_xw(1,j) + c2t
1257 time = ( tau +real(iter,f64))/real(total_iter, f64)
1258 weight = (self%quadp2_xw(2,j)*c1t)
1259 call sll_s_uniform_bsplines_eval_basis( self%spline_degree, c1*self%quadp2_xw(1,j)+c2, &
1262 do i_grid = index - self%spline_degree, index
1263 i_mod = modulo(i_grid-1, n_cells ) + 1
1264 field_dof_t = field_dofs_1( i_mod ) + time * ( field_dofs_2( i_mod ) - field_dofs_1( i_mod ) )
1265 bint(1) = bint(1) + self%spline_val(ind) * field_dof_t * weight * tau
1266 bint(2) = bint(2) + self%spline_val(ind) * field_dof_t * weight * (1.0_f64 - tau )
1276 sll_real64,
intent( in ) :: position_old(self%dim)
1277 sll_real64,
intent( in ) :: position_new(self%dim)
1278 sll_real64,
intent( in ) :: field_dofs_1(self%n_cells)
1279 sll_real64,
intent( in ) :: field_dofs_2(self%n_cells)
1280 sll_int32,
intent( in ) :: iter
1281 sll_int32,
intent( in ) :: total_iter
1282 sll_real64,
intent( out ) :: field
1284 sll_real64 :: dx_new, dx_old
1285 sll_int32 :: box_new, box_old, ind
1286 sll_real64 :: rinterval, lowert, uppert
1288 dx_new = (position_new(1)-self%domain(1))/self%delta_x
1289 box_new = floor(dx_new)+1
1290 dx_new = dx_new - real(box_new-1,f64)
1291 dx_old = (position_old(1)-self%domain(1))/self%delta_x
1292 box_old = floor(dx_old)+1
1293 dx_old = dx_old - real(box_old-1,f64)
1295 rinterval = self%delta_x/(position_new(1)-position_old(1))
1299 if ( box_old == box_new )
then
1301 field_dofs_1, field_dofs_2, iter, total_iter, field )
1302 elseif ( box_old < box_new )
then
1304 uppert = (1.0_f64-dx_old)*rinterval
1306 field_dofs_1, field_dofs_2, iter, total_iter, field )
1307 do ind = box_old+1, box_new-1
1309 uppert = uppert + rinterval
1311 field_dofs_1, field_dofs_2, iter, total_iter, field )
1314 field_dofs_1, field_dofs_2, iter, total_iter, field )
1317 uppert = 1.0_f64 + (rinterval-dx_new*rinterval)
1319 field_dofs_1, field_dofs_2, iter, total_iter, field )
1320 do ind = box_new+1, box_old-1
1322 uppert = lowert + rinterval
1324 field_dofs_1, field_dofs_2, iter, total_iter, field )
1327 field_dofs_1, field_dofs_2, iter, total_iter, field )
1335 subroutine update_field_linear_partial_newton(self, upper, lower, uppert, lowert, index, field_dofs_1, field_dofs_2, iter, total_iter, bint )
1337 sll_real64,
intent( in ) :: upper
1338 sll_real64,
intent( in ) :: lower
1339 sll_real64,
intent( in ) :: uppert
1340 sll_real64,
intent( in ) :: lowert
1341 sll_int32,
intent( in ) :: index
1342 sll_real64,
intent( in ) :: field_dofs_1(self%n_cells)
1343 sll_real64,
intent( in ) :: field_dofs_2(self%n_cells)
1344 sll_int32,
intent( in ) :: iter
1345 sll_int32,
intent( in ) :: total_iter
1346 sll_real64,
intent( inout ) :: bint
1348 sll_int32 :: j, ind, i_mod, i_grid, n_cells
1349 sll_real64 :: c1, c2, c1t, c2t
1350 sll_real64 :: tau, time, weight, field_dof_t
1352 n_cells = self%n_cells
1354 c1 = 0.5_f64*(upper-lower)
1355 c2 = 0.5_f64*(upper+lower)
1357 c1t = 0.5_f64*(uppert-lowert)
1358 c2t = 0.5_f64*(uppert+lowert)
1360 do j=1, self%n_quadp2_points
1361 tau = c1t* self%quadp2_xw(1,j) + c2t
1362 time = ( tau +real(iter,f64))/real(total_iter, f64)
1363 weight = (self%quadp2_xw(2,j)*c1t)
1364 call sll_s_uniform_bsplines_eval_deriv ( self%spline_degree, c1*self%quadp2_xw(1,1) + c2, &
1367 do i_grid = index - self%spline_degree, index
1368 i_mod = modulo(i_grid-1, n_cells ) + 1
1369 field_dof_t = field_dofs_1( i_mod ) + time * ( field_dofs_2( i_mod ) - field_dofs_1( i_mod ) )
1370 bint = bint + self%spline_val(ind) * field_dof_t * weight * (1.0_f64 - tau ) * tau
1380 sll_real64,
intent( in ) :: position_old(self%dim)
1381 sll_real64,
intent( in ) :: position_new(self%dim)
1382 sll_real64,
intent( in ) :: marker_charge
1383 sll_real64,
intent( inout ) :: j_dofs(self%n_cells)
1385 sll_real64 :: dx_new, dx_old
1386 sll_int32 :: box_new, box_old, ind
1387 sll_real64 :: rinterval, lowert, uppert
1389 dx_new = (position_new(1)-self%domain(1))/self%delta_x
1390 box_new = floor(dx_new)+1
1391 dx_new = dx_new - real(box_new-1,f64)
1392 dx_old = (position_old(1)-self%domain(1))/self%delta_x
1393 box_old = floor(dx_old)+1
1394 dx_old = dx_old - real(box_old-1,f64)
1396 if ( box_old == box_new )
then
1398 marker_charge, j_dofs )
1399 elseif ( box_old < box_new )
then
1400 rinterval = self%delta_x/(position_new(1)-position_old(1))
1402 uppert = (1.0_f64-dx_old)*rinterval
1404 marker_charge, j_dofs )
1405 do ind = box_old+1, box_new-1
1407 uppert = uppert + rinterval
1409 marker_charge, j_dofs )
1412 marker_charge, j_dofs )
1414 rinterval = self%delta_x/(position_new(1)-position_old(1))
1416 uppert = 1.0_f64 + (rinterval - dx_new*rinterval)
1419 -marker_charge, j_dofs )
1420 do ind = box_new+1, box_old-1
1422 uppert = lowert + rinterval
1424 -marker_charge, j_dofs )
1427 -marker_charge, j_dofs )
1436 sll_real64,
intent( in ) :: upper
1437 sll_real64,
intent( in ) :: lower
1438 sll_real64,
intent( in ) :: uppert
1439 sll_real64,
intent( in ) :: lowert
1440 sll_int32,
intent( in ) :: index
1441 sll_real64,
intent( in ) :: marker_charge
1442 sll_real64,
intent( inout ) :: j_dofs(self%n_cells)
1444 sll_int32 :: j, ind, i_mod, i_grid, n_cells
1445 sll_real64 :: c1, c2, c1t, c2t
1447 n_cells = self%n_cells
1449 c1 = 0.5_f64*(upper-lower)
1450 c2 = 0.5_f64*(upper+lower)
1452 c1t = 0.5_f64*(uppert-lowert)
1453 c2t = 0.5_f64*(uppert+lowert)
1455 call sll_s_uniform_bsplines_eval_basis ( self%spline_degree, &
1456 c1*self%quadp1_xw(1,1) + c2, &
1458 self%spline_val = self%spline_val * (self%quadp1_xw(2,1)*c1t)
1460 do j=2, self%n_quadp1_points
1461 call sll_s_uniform_bsplines_eval_basis( self%spline_degree, &
1462 c1*self%quadp1_xw(1,j)+c2, &
1463 self%spline_val_more )
1464 self%spline_val = self%spline_val + &
1465 self%spline_val_more * (self%quadp1_xw(2,j)*c1t)
1469 do i_grid = index - self%spline_degree, index
1470 i_mod = modulo(i_grid-1, n_cells ) + 1
1471 j_dofs(i_mod) = j_dofs(i_mod) + marker_charge * self%spline_val(ind)
Gauss-Legendre integration.
real(kind=f64) function, dimension(2, 1:npoints), public sll_f_gauss_legendre_points_and_weights(npoints, a, b)
Returns a 2d array of size (2,npoints) containing gauss-legendre points and weights in the interval [...
Low level arbitrary degree splines.
Base class for kernel smoothers for accumulation and field evaluation in PIC.
integer(kind=i32), parameter, public sll_p_collocation
integer(kind=i32), parameter, public sll_p_galerkin
Kernel smoother for 1d with splines of arbitrary degree placed on a uniform mesh.
subroutine add_particle_mass_spline_1d(self, position, marker_charge, particle_mass)
Add the contribution of one particle to the approximate mass matrix.
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)
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)
subroutine add_particle_mass_spline_1d_full(self, position, marker_charge, particle_mass)
Add the contribution of one particle to the approximate mass matrix.
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 der...
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 ...
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 f...
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 u...
subroutine add_current_spline_1d(self, position_old, position_new, marker_charge, j_dofs)
Add current with integration over x.
subroutine evaluate_field_single_spline_1d(self, position, field_dofs, field_value)
Evaluate field at at position position.
subroutine update_jv(self, lower, upper, index, marker_charge, qoverm, sign, vi, j_dofs, bfield_dofs)
Helper function for add_current_update_v.
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)
subroutine, public sll_s_new_particle_mesh_coupling_spline_1d_ptr(smoother, domain, n_cells, no_particles, spline_degree, smoothing_type)
subroutine update_jv_pp(self, lower, upper, index, marker_charge, qoverm, vi, j_dofs, bfield_dofs)
Helper function for add_current_update_v.
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)
subroutine init_spline_1d(self, domain, n_cells, no_particles, spline_degree, smoothing_type)
Initializer.
subroutine evaluate_multiple_spline_1d(self, position, components, field_dofs, field_value)
Evaluate several fields at position position.
subroutine evaluate_int_spline_1d_subnewton(self, position_old, position_new, field_dofs, field)
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.
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.
subroutine free_spline_1d(self)
Destructor.
subroutine, public sll_s_new_particle_mesh_coupling_spline_1d(smoother, domain, n_cells, no_particles, spline_degree, smoothing_type)
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...
subroutine update_current_partial(self, upper, lower, uppert, lowert, index, iter, total_iter, marker_charge, j_dofs)
Helper function for add_current_split.
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.
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.
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)
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 evaluat...
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)
subroutine, public sll_s_spline_pp_horner_primitive_1d(val, degree, pp_coeffs, x)
Perform a 1d hornerschema on the pp_coeffs evaluate at x.
subroutine, public sll_s_spline_pp_init_1d(self, degree, n_cells, boundary)
Initialize sll_t_spline_pp_1d object (Set poly_coeffs depending on degree)
subroutine, public sll_s_spline_pp_free_1d(self)
Destructor 1d.
real(kind=f64) function, public sll_f_spline_pp_horner_1d(degree, pp_coeffs, x, index)
Perform a 1d hornerschema on the pp_coeffs at index.
Basic type of a kernel smoother used for PIC simulations.
Spline kernel smoother in1d.
arbitrary degree 1d spline