Find the trapped pitch angle weights using an old (finite-difference) integration scheme
Here we find the trapped pitch angle weights corresponding to a trapezium rule integration in v||/v. With this method we can write: Int(f) ~ Sum_i((f(l_{i+1}) - f(l_{i}))*(l_{i+1}-l_{i}))/2 with some handling required for the upper and lower boundaries. Here we find the v||/v spacing. Note there is no factor 1/2 due to the fact the lambda weight is the sum of the (identical) positive and negative v||/v weights, which cancels this factor of 1/2. Also note that there is no weight for the bouncing lambda.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real, | intent(in), | dimension(:) | :: | lambda_grid |
The lambda grid points. |
|
real, | intent(inout), | dimension(-ntgrid:, :) | :: | weights |
The integration weights for the lambda grid. Note this is
intent |
subroutine setup_trapped_lambda_grids_old_finite_difference(lambda_grid, weights)
use theta_grid, only: ntgrid, bmag
implicit none
!> The lambda grid points.
real, dimension(:), intent(in) :: lambda_grid
!> The integration weights for the lambda grid. Note this is
!> intent `in out` as we only set a portion of the grid so may
!> want to keep the rest of the input values.
real, dimension(-ntgrid:, :), intent(in out) :: weights
real :: wwo
integer :: ig, il
do il = ng2+1, nlambda-1
do ig = -ntgrid, ntgrid
if (1.0-lambda_grid(il)*bmag(ig) > -bouncefuzz &
.and. 1.0-lambda_grid(il+1)*bmag(ig) > -bouncefuzz) &
then
wwo = sqrt(max(1.0 - lambda_grid(il)*bmag(ig),0.0)) - &
sqrt(max(1.0 - lambda_grid(il+1)*bmag(ig),0.0))
weights(ig,il) = weights(ig,il) + wwo
weights(ig,il+1) = weights(ig,il+1) + wwo
end if
end do
end do
end subroutine setup_trapped_lambda_grids_old_finite_difference