dfnrmk Subroutine

subroutine dfnrmk(neq, y, t, ydot, savr, r, cj, tscale, wght, sqrtn, rsqrtn, res, ires, psol, irin, ierr, rnorm, epslin, rwp, iwp, wk, rpar, ipar)

Uses

  • proc~~dfnrmk~~UsesGraph proc~dfnrmk dfnrmk module~blas_interfaces blas_interfaces proc~dfnrmk->module~blas_interfaces module~daskr daskr proc~dfnrmk->module~daskr module~daskr_kinds daskr_kinds proc~dfnrmk->module~daskr_kinds module~blas_interfaces->module~daskr_kinds module~daskr->module~daskr_kinds iso_fortran_env iso_fortran_env module~daskr_kinds->iso_fortran_env

This routine calculates the scaled preconditioned norm of the nonlinear function used in the nonlinear iteration for obtaining consistent initial conditions. Specifically, it calculates the weighted root-mean-square norm of the vector:

where is the preconditioner matrix and is the DAE equation vector.

Arguments

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

Problem size.

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

Solution vector.

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

Independent variable.

real(kind=rk), intent(in) :: ydot(neq)

Derivative of solution vector after successful step.

real(kind=rk), intent(inout) :: savr(neq)

Saved residual vector.

real(kind=rk), intent(out) :: r(neq)

Result vector.

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

Scalar used in forming the system Jacobian.

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

Scale factor in t; used for stopping tests if nonzero.

real(kind=rk), intent(inout) :: wght(neq)

Scaling factors.

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

Square root of neq.

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

Reciprocal of square root of neq.

procedure(res_t) :: res

User-defined residuals routine.

integer, intent(out) :: ires

Error flag from res.

procedure(psol_t) :: psol

User-defined preconditioner routine.

integer, intent(in) :: irin

Flag indicating whether the current residual vector is input in savr. 0: it is not. 1: it is.

integer, intent(out) :: ierr

Error flag from psol.

real(kind=rk), intent(out) :: rnorm

Weighted root-mean-square norm of r.

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

Tolerance for linear system.

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

Real work array used by preconditioner psol.

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

Integer work array used by preconditioner psol.

real(kind=rk), intent(inout) :: wk(neq)

Real work array used by psol.

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

User real workspace.

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

User integer workspace.


Calls

proc~~dfnrmk~~CallsGraph proc~dfnrmk dfnrmk ddwnrm ddwnrm proc~dfnrmk->ddwnrm interface~dcopy dcopy proc~dfnrmk->interface~dcopy interface~dscal dscal proc~dfnrmk->interface~dscal

Source Code

subroutine dfnrmk( &
   neq, y, t, ydot, savr, r, cj, tscale, wght, &
   sqrtn, rsqrtn, res, ires, psol, irin, ierr, &
   rnorm, epslin, rwp, iwp, wk, rpar, ipar)
!! This routine calculates the scaled preconditioned norm of the nonlinear function used
!! in the nonlinear iteration for obtaining consistent initial conditions. Specifically,
!! it calculates the weighted root-mean-square norm of the vector:
!!
!! $$  r = P^{-1} G(t,y,\dot{y}) $$
!!
!! where \(P\) is the preconditioner matrix and \(G\) is the DAE equation vector.

   use daskr_kinds, only: rk, zero
   use blas_interfaces, only: dcopy, dscal
   use daskr
   implicit none

   integer, intent(in) :: neq
      !! Problem size.
   real(rk), intent(in) :: y(neq)
      !! Solution vector.
   real(rk), intent(in) :: t
      !! Independent variable.
   real(rk), intent(in) :: ydot(neq)
      !! Derivative of solution vector after successful step.
   real(rk), intent(inout) :: savr(neq)
      !! Saved residual vector.
   real(rk), intent(out) :: r(neq)
      !! Result vector.
   real(rk), intent(in) :: cj
      !! Scalar used in forming the system Jacobian.
   real(rk), intent(in) :: tscale ! @todo: what is "t?
      !! Scale factor in `t`; used for stopping tests if nonzero.
   real(rk), intent(inout) :: wght(neq)
      !! Scaling factors.
   real(rk), intent(in) :: sqrtn
      !! Square root of `neq`.
   real(rk), intent(in) :: rsqrtn
      !! Reciprocal of square root of `neq`.
   procedure(res_t) :: res
      !! User-defined residuals routine.
   integer, intent(out) :: ires
      !! Error flag from `res`.
   procedure(psol_t) :: psol
      !! User-defined preconditioner routine.
   integer, intent(in) :: irin
      !! Flag indicating whether the current residual vector is input in `savr`.
      !! `0`: it is not.
      !! `1`: it is.
   integer, intent(out) :: ierr
      !! Error flag from `psol`.
   real(rk), intent(out) :: rnorm
      !! Weighted root-mean-square norm of `r`.
   real(rk), intent(in) :: epslin
      !! Tolerance for linear system.
   real(rk), intent(inout) :: rwp(*)
      !! Real work array used by preconditioner `psol`.
   integer, intent(inout) :: iwp(*)
      !! Integer work array used by preconditioner `psol`.
   real(rk), intent(inout) :: wk(neq)
      !! Real work array used by `psol`.
   real(rk), intent(inout) :: rpar(*)
      !! User real workspace.
   integer, intent(inout) :: ipar(*)
      !! User integer workspace.

   real(rk), external :: ddwnrm ! @todo: remove this once inside module

   ! Call RES routine if IRIN = 0.
   if (irin == 0) then
      ires = 0
      call res(t, y, ydot, cj, savr, ires, rpar, ipar)
      ! @todo: count missing here, or is it counted in the caller?
      if (ires < 0) return
   end if

   ! Apply inverse of left preconditioner to vector R.
   ! First scale WT array by 1/sqrt(N), and undo scaling afterward.
   call dcopy(neq, savr, 1, r, 1)
   call dscal(neq, rsqrtn, wght, 1)
   ierr = 0
   call psol(neq, t, y, ydot, savr, wk, cj, wght, rwp, iwp, r, epslin, ierr, rpar, ipar)
   ! @todo: count missing here, or is it counted in the caller?
   call dscal(neq, sqrtn, wght, 1)
   if (ierr /= 0) return

   ! Calculate norm of R.
   rnorm = ddwnrm(neq, r, wght, rpar, ipar)
   if (tscale > zero) rnorm = rnorm*tscale*abs(cj)

end subroutine dfnrmk