banps Subroutine

public subroutine banps(neq, t, y, yprime, savr, wk, cj, wght, wp, iwp, b, eplin, ier, rpar, ipar)

This subroutine uses the factors produced by banja to solve linear systems for the banded preconditioner , given a vector . It calls the LINPACK routine dgbsl for this.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: neq

Problem size.

real(kind=rk), intent(in) :: t

Independent variable (not used).

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

Current dependent variables (not used).

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

Current derivatives of dependent variables (not used).

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

Current residual evaluated at (t, y, yprime) (not used).

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

Real work space of length neq (not used).

real(kind=rk), intent(in) :: cj

Scalar proportional to 1/h (not used).

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

Error weights for computing norms (not used).

real(kind=rk), intent(inout) :: wp(*)

Real work array containing the LU decomposition of P.

integer, intent(inout) :: iwp(*)

Integer array containing matrix pivot information.

real(kind=rk), intent(inout) :: b(*)

Right-hand side vector on input; solution on output.

real(kind=rk), intent(in) :: eplin

Tolerance for linear system (not used).

integer, intent(inout) :: ier

Output error flag (not used; assumed 0 on input).

real(kind=rk), intent(inout) :: rpar(*)

Real array used for communication between the calling program and external user routines (not used).

integer, intent(inout) :: ipar(*)

Integer array used for communication between the calling program and external user routines (not used). ipar(1) and ipar(2) must contain ml and mu, respectively.


Calls

proc~~banps~~CallsGraph proc~banps banps dgbsl dgbsl proc~banps->dgbsl

Variables

Type Visibility Attributes Name Initial
integer, public :: meband
integer, public :: ml
integer, public :: mu

Source Code

   subroutine banps(neq, t, y, yprime, savr, wk, cj, wght, wp, iwp, b, eplin, ier, rpar, ipar)
   !! This subroutine uses the factors produced by [[banja]] to solve linear systems \(P x = b\)
   !! for the banded preconditioner \(P\), given a vector \(b\). It calls the LINPACK routine
   !! [[DGBSL]] for this.
      integer, intent(in) :: neq
        !! Problem size.
      real(rk), intent(in) :: t
        !! Independent variable (not used).
      real(rk), intent(in) :: y(*)
        !! Current dependent variables (not used).
      real(rk), intent(in) :: yprime(*)
        !! Current derivatives of dependent variables (not used).
      real(rk), intent(in) :: savr(*)
        !! Current residual evaluated at `(t, y, yprime)` (not used).
      real(rk), intent(in) :: wk(*)
        !! Real work space of length `neq` (not used).
      real(rk), intent(in) :: cj
        !! Scalar proportional to `1/h` (not used).
      real(rk), intent(in) :: wght(*)
        !! Error weights for computing norms (not used).
      real(rk), intent(inout) :: wp(*)
        !! Real work array containing the LU decomposition of P.
      integer, intent(inout) :: iwp(*)
        !! Integer array containing matrix pivot information.
      real(rk), intent(inout) :: b(*)
        !! Right-hand side vector on input; solution on output.
      real(rk), intent(in) :: eplin
        !! Tolerance for linear system (not used).
      integer, intent(inout) :: ier
        !! Output error flag (not used; assumed 0 on input).
      real(rk), intent(inout) :: rpar(*)
        !! Real array used for communication between the calling program and external user
        !! routines (not used).
      integer, intent(inout) :: ipar(*)
        !! Integer array used for communication between the calling program and external user
        !! routines (not used). `ipar(1)` and `ipar(2)` must contain `ml` and `mu`, respectively.

      external :: dgbsl
      integer :: meband, ml, mu

      ml = ipar(1)
      mu = ipar(2)
      meband = 2*ml + mu + 1
      call dgbsl(wp, meband, neq, ml, mu, iwp, b, 0)

   end subroutine banps