GNU Octave 10.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 
Loading...
Searching...
No Matches
datv.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 datv (NEQ, Y, TN, YPRIME, SAVR, V, WGHT, YPTEM, RES,
6 * IRES, PSOL, Z, VTEM, WP, IWP, CJ, EPLIN, IER, NRE, NPSL,
7 * RPAR,IPAR)
8C
9C***BEGIN PROLOGUE DATV
10C***DATE WRITTEN 890101 (YYMMDD)
11C***REVISION DATE 900926 (YYMMDD)
12C
13C
14C-----------------------------------------------------------------------
15C***DESCRIPTION
16C
17C This routine computes the product
18C
19C Z = (D-inverse)*(P-inverse)*(dF/dY)*(D*V),
20C
21C where F(Y) = G(T, Y, CJ*(Y-A)), CJ is a scalar proportional to 1/H,
22C and A involves the past history of Y. The quantity CJ*(Y-A) is
23C an approximation to the first derivative of Y and is stored
24C in the array YPRIME. Note that dF/dY = dG/dY + CJ*dG/dYPRIME.
25C
26C D is a diagonal scaling matrix, and P is the left preconditioning
27C matrix. V is assumed to have L2 norm equal to 1.
28C The product is stored in Z and is computed by means of a
29C difference quotient, a call to RES, and one call to PSOL.
30C
31C On entry
32C
33C NEQ = Problem size, passed to RES and PSOL.
34C
35C Y = Array containing current dependent variable vector.
36C
37C YPRIME = Array containing current first derivative of y.
38C
39C SAVR = Array containing current value of G(T,Y,YPRIME).
40C
41C V = Real array of length NEQ (can be the same array as Z).
42C
43C WGHT = Array of length NEQ containing scale factors.
44C 1/WGHT(I) are the diagonal elements of the matrix D.
45C
46C YPTEM = Work array of length NEQ.
47C
48C VTEM = Work array of length NEQ used to store the
49C unscaled version of V.
50C
51C WP = Real work array used by preconditioner PSOL.
52C
53C IWP = Integer work array used by preconditioner PSOL.
54C
55C CJ = Scalar proportional to current value of
56C 1/(step size H).
57C
58C
59C On return
60C
61C Z = Array of length NEQ containing desired scaled
62C matrix-vector product.
63C
64C IRES = Error flag from RES.
65C
66C IER = Error flag from PSOL.
67C
68C NRE = The number of calls to RES.
69C
70C NPSL = The number of calls to PSOL.
71C
72C-----------------------------------------------------------------------
73C***ROUTINES CALLED
74C RES, PSOL
75C
76C***END PROLOGUE DATV
77C
78 INTEGER NEQ, IRES, IWP, IER, NRE, NPSL, IPAR
79 DOUBLE PRECISION Y, TN, YPRIME, SAVR, V, WGHT, YPTEM, Z, VTEM,
80 1 wp, cj, rpar
81 dimension y(*), yprime(*), savr(*), v(*), wght(*), yptem(*),
82 1 z(*), vtem(*), wp(*), iwp(*), rpar(*), ipar(*)
83 INTEGER I
84 DOUBLE PRECISION EPLIN
85 EXTERNAL res, psol
86C
87 ires = 0
88C-----------------------------------------------------------------------
89C Set VTEM = D * V.
90C-----------------------------------------------------------------------
91 DO 10 i = 1,neq
92 10 vtem(i) = v(i)/wght(i)
93 ier = 0
94C-----------------------------------------------------------------------
95C Store Y in Z and increment Z by VTEM.
96C Store YPRIME in YPTEM and increment YPTEM by VTEM*CJ.
97C-----------------------------------------------------------------------
98 DO 20 i = 1,neq
99 yptem(i) = yprime(i) + vtem(i)*cj
100 20 z(i) = y(i) + vtem(i)
101C-----------------------------------------------------------------------
102C Call RES with incremented Y, YPRIME arguments
103C stored in Z, YPTEM. VTEM is overwritten with new residual.
104C-----------------------------------------------------------------------
105 CONTINUE
106 CALL res(tn,z,yptem,cj,vtem,ires,rpar,ipar)
107 nre = nre + 1
108 IF (ires .LT. 0) RETURN
109C-----------------------------------------------------------------------
110C Set Z = (dF/dY) * VBAR using difference quotient.
111C (VBAR is old value of VTEM before calling RES)
112C-----------------------------------------------------------------------
113 DO 70 i = 1,neq
114 70 z(i) = vtem(i) - savr(i)
115C-----------------------------------------------------------------------
116C Apply inverse of left preconditioner to Z.
117C-----------------------------------------------------------------------
118 CALL psol (neq, tn, y, yprime, savr, yptem, cj, wght, wp, iwp,
119 1 z, eplin, ier, rpar, ipar)
120 npsl = npsl + 1
121 IF (ier .NE. 0) RETURN
122C-----------------------------------------------------------------------
123C Apply D-inverse to Z and return.
124C-----------------------------------------------------------------------
125 DO 90 i = 1,neq
126 90 z(i) = z(i)*wght(i)
127 RETURN
128C
129C------END OF SUBROUTINE DATV-------------------------------------------
130 END
subroutine datv(neq, y, tn, yprime, savr, v, wght, yptem, res, ires, psol, z, vtem, wp, iwp, cj, eplin, ier, nre, npsl, rpar, ipar)
Definition datv.f:8