sdc_qmats Subroutine

public subroutine sdc_qmats(qmat, smat, dst, src, flags, ndst, nsrc)

Subroutine to compute the quadrature matrices

Arguments

Type IntentOptional AttributesName
real(kind=pfdp), intent(out) :: qmat(ndst-1,nsrc)

O to dst quadrature weights

real(kind=pfdp), intent(out) :: smat(ndst-1,nsrc)

dst(m) to dst(m+1) quadrature weights

real(kind=pfqp), intent(in) :: dst(ndst)

Destination points

real(kind=pfqp), intent(in) :: src(nsrc)

Source points

integer, intent(in) :: flags(nsrc)
integer, intent(in), value:: ndst

Number of destination points

integer, intent(in), value:: nsrc

Number of source points


Calls

proc~~sdc_qmats~~CallsGraph proc~sdc_qmats sdc_qmats proc~not_proper not_proper proc~sdc_qmats->proc~not_proper interface~poly_eval poly_eval proc~sdc_qmats->interface~poly_eval proc~poly_int poly_int proc~sdc_qmats->proc~poly_int interface~poly_eval->interface~poly_eval proc~poly_eval_complex poly_eval_complex interface~poly_eval->proc~poly_eval_complex

Called by

proc~~sdc_qmats~~CalledByGraph proc~sdc_qmats sdc_qmats proc~pf_quadrature pf_quadrature proc~pf_quadrature->proc~sdc_qmats proc~pf_init_sdcmats pf_init_sdcmats proc~pf_init_sdcmats->proc~pf_quadrature proc~pf_level_setup pf_level_setup proc~pf_level_setup->proc~pf_init_sdcmats proc~pf_pfasst_setup pf_pfasst_setup proc~pf_pfasst_setup->proc~pf_level_setup

Contents

Source Code


Source Code

  subroutine sdc_qmats(qmat, smat, dst, src, flags, ndst, nsrc)
    integer ,  intent(in), value  :: ndst   !!  Number of destination points
    integer ,   intent(in), value  :: nsrc   !!  Number of source points
    real(pfqp), intent(in)  :: dst(ndst)     !!  Destination points
    real(pfqp), intent(in)  :: src(nsrc)     !!  Source points
    real(pfdp),      intent(out) :: qmat(ndst-1, nsrc)  !!  O to dst quadrature weights
    real(pfdp),      intent(out) :: smat(ndst-1, nsrc)  !! dst(m) to dst(m+1) quadrature weights
    integer ,      intent(in)  :: flags(nsrc)     

    integer  :: i, j, m
    real(pfqp) :: q, s, den, p(0:nsrc)

    qmat = 0.0_pfqp
    smat = 0.0_pfqp

    ! construct qmat and smat
    do i = 1, nsrc

       if (not_proper(flags, i)) cycle

       ! construct interpolating polynomial coefficients
       p    = 0.0_pfdp
       p(0) = 1.0_pfdp
       do m = 1, nsrc
          if (not_proper(flags, m) .or. m == i) cycle
          p = eoshift(p, -1) - src(m) * p
       end do

       den = poly_eval(p, nsrc, src(i))

       call poly_int(p, nsrc)

       ! evaluate integrals
       do j = 2, ndst
          q = poly_eval(p, nsrc, dst(j)) - poly_eval(p, nsrc,   0.0_pfqp)
          s = poly_eval(p, nsrc, dst(j)) - poly_eval(p, nsrc, dst(j-1))

          qmat(j-1,i) = real(q / den, pfqp)
          smat(j-1,i) = real(s / den, pfqp)
       end do
    end do

  end subroutine sdc_qmats