psol_rbdpre Subroutine

public pure subroutine psol_rbdpre(b, bd, ipbd)

Uses

  • proc~~psol_rbdpre~~UsesGraph proc~psol_rbdpre psol_rbdpre module~linpack linpack proc~psol_rbdpre->module~linpack module~daskr_kinds daskr_kinds module~linpack->module~daskr_kinds iso_fortran_env iso_fortran_env module~daskr_kinds->iso_fortran_env

This routine solves a linear system , using the LU factors of the diagonal blocks computed in jac_rbdpre and mesh parameters passed as module variables. The solution is carried out by dgesl from LINPACK.

Arguments

Type IntentOptional Attributes Name
real(kind=rk), intent(inout) :: b(*)

Right-hand side vector on entry and solution vector on return.

real(kind=rk), intent(in) :: bd(*)

LU factors of the diagonal blocks. bd corresponds to the segment rwp of rwork.

integer, intent(in) :: ipbd(*)

Pivots for the LU factorizations. ipbd corresponds to the segment iwp of iwork.


Calls

proc~~psol_rbdpre~~CallsGraph proc~psol_rbdpre psol_rbdpre proc~dgesl dgesl proc~psol_rbdpre->proc~dgesl

Called by

proc~~psol_rbdpre~~CalledByGraph proc~psol_rbdpre psol_rbdpre proc~psol psol proc~psol->proc~psol_rbdpre

Source Code

   pure subroutine psol_rbdpre(b, bd, ipbd)
   !! This routine solves a linear system \(A_R x = b\), using the LU factors of the diagonal 
   !! blocks computed in [[jac_rbdpre]] and mesh parameters passed as module variables. The
   !! solution is carried out by [[dgesl]] from LINPACK.

      use linpack, only: dgesl

      real(rk), intent(inout) :: b(*)
         !! Right-hand side vector on entry and solution vector on return.
      real(rk), intent(in) :: bd(*)
         !! LU factors of the diagonal blocks. `bd` corresponds to the segment `rwp` of `rwork`.
      integer, intent(in) :: ipbd(*)
         !! Pivots for the LU factorizations. `ipbd` corresponds to the segment `iwp` of `iwork`.

      integer :: ib, ibd, jx, jy

      ib = 1
      ibd = 1
      do jy = 1, meshy
         do jx = 1, meshx
            call dgesl(bd(ibd), mp, mp, ipbd(ib), b(ib), 0)
            ib = ib + mp
            ibd = ibd + mpsq
         end do
      end do

   end subroutine psol_rbdpre