GNU Octave 10.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 
Loading...
Searching...
No Matches
dfnrmk.f
Go to the documentation of this file.
1C Work performed under the auspices of the U.S. Department of Energy
2C by Lawrence Livermore National Laboratory under contract number
3C W-7405-Eng-48.
4C
5 SUBROUTINE dfnrmk (NEQ, Y, T, YPRIME, SAVR, R, CJ, WT,
6 * SQRTN, RSQRTN, RES, IRES, PSOL, IRIN, IER,
7 * FNORM, EPLIN, WP, IWP, PWK, RPAR, IPAR)
8C
9C***BEGIN PROLOGUE DFNRMK
10C***REFER TO DLINSK
11C***DATE WRITTEN 940830 (YYMMDD)
12C***REVISION DATE 951006 (SQRTN, RSQRTN, and scaling of WT added.)
13C
14C
15C-----------------------------------------------------------------------
16C***DESCRIPTION
17C
18C DFNRMK calculates the scaled preconditioned norm of the nonlinear
19C function used in the nonlinear iteration for obtaining consistent
20C initial conditions. Specifically, DFNRMK calculates the weighted
21C root-mean-square norm of the vector (P-inverse)*G(T,Y,YPRIME),
22C where P is the preconditioner matrix.
23C
24C In addition to the parameters described in the calling program
25C DLINSK, the parameters represent
26C
27C IRIN -- Flag showing whether the current residual vector is
28C input in SAVR. 1 means it is, 0 means it is not.
29C R -- Array of length NEQ that contains
30C (P-inverse)*G(T,Y,YPRIME) on return.
31C FNORM -- Scalar containing the weighted norm of R on return.
32C-----------------------------------------------------------------------
33C
34C***ROUTINES CALLED
35C RES, DCOPY, DSCAL, PSOL, DDWNRM
36C
37C***END PROLOGUE DFNRMK
38C
39C
40 IMPLICIT DOUBLE PRECISION (a-h,o-z)
41 EXTERNAL res, psol
42 dimension y(*), yprime(*), wt(*), savr(*), r(*), pwk(*)
43 dimension wp(*), iwp(*), rpar(*), ipar(*)
44C-----------------------------------------------------------------------
45C Call RES routine if IRIN = 0.
46C-----------------------------------------------------------------------
47 IF (irin .EQ. 0) THEN
48 ires = 0
49 CALL res (t, y, yprime, cj, savr, ires, rpar, ipar)
50 IF (ires .LT. 0) RETURN
51 ENDIF
52C-----------------------------------------------------------------------
53C Apply inverse of left preconditioner to vector R.
54C First scale WT array by 1/sqrt(N), and undo scaling afterward.
55C-----------------------------------------------------------------------
56 CALL dcopy(neq, savr, 1, r, 1)
57 CALL dscal (neq, rsqrtn, wt, 1)
58 ier = 0
59 CALL psol (neq, t, y, yprime, savr, pwk, cj, wt, wp, iwp,
60 * r, eplin, ier, rpar, ipar)
61 CALL dscal (neq, sqrtn, wt, 1)
62 IF (ier .NE. 0) RETURN
63C-----------------------------------------------------------------------
64C Calculate norm of R.
65C-----------------------------------------------------------------------
66 fnorm = ddwnrm(neq, r, wt, rpar, ipar)
67C
68 RETURN
69C----------------------- END OF SUBROUTINE DFNRMK ----------------------
70 END
double precision function ddwnrm(neq, v, rwt, rpar, ipar)
Definition ddwnrm.f:6
subroutine dfnrmk(neq, y, t, yprime, savr, r, cj, wt, sqrtn, rsqrtn, res, ires, psol, irin, ier, fnorm, eplin, wp, iwp, pwk, rpar, ipar)
Definition dfnrmk.f:8