xref: /phasta/phSolver/compressible/itrfdi.f (revision 8c06e07d75283b78095e57c11329d999e11c1ba6)
159599516SKenneth E. Jansen        subroutine itrFDI (ypre,      y,    ac,     x,
259599516SKenneth E. Jansen     &                     rmes,      uBrg,      BDiag,
359599516SKenneth E. Jansen     &                     iBC,       BC,        iper,
459599516SKenneth E. Jansen     &                     ilwork,    shp,       shgl,
559599516SKenneth E. Jansen     &                     shpb,      shglb)
659599516SKenneth E. Jansenc
759599516SKenneth E. Jansenc----------------------------------------------------------------------
859599516SKenneth E. Jansenc
959599516SKenneth E. Jansenc  This subroutine computes the "optimum" finite difference
1059599516SKenneth E. Jansenc  interval eps for the forward difference scheme
1159599516SKenneth E. Jansenc
1259599516SKenneth E. Jansenc                 Rmod(y + eps u) - Rmod(y)
1359599516SKenneth E. Jansenc                ---------------------------
1459599516SKenneth E. Jansenc                            eps
1559599516SKenneth E. Jansenc
1659599516SKenneth E. Jansenc  where  u is the step and Rmod is the modified residual.
1759599516SKenneth E. Jansenc
1859599516SKenneth E. Jansenc Note: A good theoretical reference is 'Practical Optimization' by
1959599516SKenneth E. Jansenc        P.E. Gill, W. Murray and M.H. Wright  [1981].
2059599516SKenneth E. Jansenc
2159599516SKenneth E. Jansenc input:
2259599516SKenneth E. Jansenc  y      (nshg,ndof)           : Y-variables
2359599516SKenneth E. Jansenc  ypre   (nshg,ndof)           : preconditioned Y-variables
2459599516SKenneth E. Jansenc  x      (numnp,nsd)            : node coordinates
2559599516SKenneth E. Jansenc  rmes   (nshg,nflow)           : modified residual
2659599516SKenneth E. Jansenc  uBrg   (nshg,nflow)           : step
2759599516SKenneth E. Jansenc  BDiag  (nshg,nflow,nflow)         : block-diagonal preconditioner
2859599516SKenneth E. Jansenc  iBC    (nshg)                : BC codes
2959599516SKenneth E. Jansenc  BC     (nshg,ndofBC)         : BC constraint parameters
3059599516SKenneth E. Jansenc
3159599516SKenneth E. Jansenc
3259599516SKenneth E. Jansenc
3359599516SKenneth E. Jansenc Zdenek Johan, Winter 1989.
3459599516SKenneth E. Jansenc Zdenek Johan, Winter 1991.  (Fortran 90)
3559599516SKenneth E. Jansenc----------------------------------------------------------------------
3659599516SKenneth E. Jansenc
3759599516SKenneth E. Jansen        include "common.h"
3859599516SKenneth E. Jansenc
3959599516SKenneth E. Jansen        dimension y(nshg,ndof),               ypre(nshg,nflow),
4059599516SKenneth E. Jansen     &            x(numnp,nsd),               ac(nshg,ndof),
4159599516SKenneth E. Jansen     &            rmes(nshg,nflow),
4259599516SKenneth E. Jansen     &            uBrg(nshg,nflow),           BDiag(nshg,nflow,nflow),
4359599516SKenneth E. Jansen     &            iBC(nshg),                  BC(nshg,ndofBC),
4459599516SKenneth E. Jansen     &            ilwork(nlwork),
4559599516SKenneth E. Jansen     &            iper(nshg)
4659599516SKenneth E. Jansenc
4759599516SKenneth E. Jansen        dimension ytmp(nshg,nflow),            rtmp(nshg,nflow)
4859599516SKenneth E. Jansen        dimension tmpy(nshg,ndof)
4959599516SKenneth E. Jansenc
5059599516SKenneth E. Jansen        dimension shp(MAXTOP,maxsh,MAXQPT),
5159599516SKenneth E. Jansen     &            shgl(MAXTOP,nsd,maxsh,MAXQPT),
5259599516SKenneth E. Jansen     &            shpb(MAXTOP,maxsh,MAXQPT),
5359599516SKenneth E. Jansen     &            shglb(MAXTOP,nsd,maxsh,MAXQPT)
5459599516SKenneth E. Jansenc
5559599516SKenneth E. Jansenc.... compute the accuracy (cancellation error)  ->  epsA
5659599516SKenneth E. Jansenc
5759599516SKenneth E. Jansen        rtmp = zero
5859599516SKenneth E. Jansenc
5959599516SKenneth E. Jansen        ytmp = ypre
6059599516SKenneth E. Jansenc
6159599516SKenneth E. Jansenc        call yshuffle(ytmp, 'new2old ')
6259599516SKenneth E. Jansenc
6359599516SKenneth E. Jansen        call i3LU (BDiag, ytmp, 'backward')
6459599516SKenneth E. Jansenc
6559599516SKenneth E. Jansen        call yshuffle(ytmp, 'old2new ')
6659599516SKenneth E. Jansenc
6759599516SKenneth E. Jansen        iabres = 1
6859599516SKenneth E. Jansenc
6959599516SKenneth E. Jansen        call itrRes (ytmp,            y,
7059599516SKenneth E. Jansen     &               x,               shp,
7159599516SKenneth E. Jansen     &               shgl,            iBC,
7259599516SKenneth E. Jansen     &               BC,              shpb,
7359599516SKenneth E. Jansen     &               shglb,           rtmp,
7459599516SKenneth E. Jansen     &               iper,            ilwork,
7559599516SKenneth E. Jansen     &               ac)
7659599516SKenneth E. Jansenc
7759599516SKenneth E. Jansen        iabres = 0
7859599516SKenneth E. Jansenc
7959599516SKenneth E. Jansen        call i3LU (BDiag, rtmp, 'forward ')
8059599516SKenneth E. Jansenc
8159599516SKenneth E. Jansen        rtmp = rtmp**2
8259599516SKenneth E. Jansen        call sumgat (rtmp, nflow, summed, ilwork)
8359599516SKenneth E. Jansen        epsA = (epsM**2) * sqrt(summed)
8459599516SKenneth E. Jansenc
8559599516SKenneth E. Jansenc.... compute the norm of the second derivative (truncation error)
8659599516SKenneth E. Jansenc
8759599516SKenneth E. Jansenc.... set interval
8859599516SKenneth E. Jansenc
8959599516SKenneth E. Jansen        epsSD = sqrt(epsM)
9059599516SKenneth E. Jansenc
9159599516SKenneth E. Jansenc.... compute the first residual
9259599516SKenneth E. Jansenc
9359599516SKenneth E. Jansen        rtmp = zero
9459599516SKenneth E. Jansenc
9559599516SKenneth E. Jansenc        call yshuffle(ypre, 'new2old ')
9659599516SKenneth E. Jansenc
9759599516SKenneth E. Jansen        ytmp = ypre + epsSD * uBrg
9859599516SKenneth E. Jansenc
9959599516SKenneth E. Jansen        call i3LU (BDiag, ytmp, 'backward')
10059599516SKenneth E. Jansenc
10159599516SKenneth E. Jansen        call yshuffle(ytmp, 'old2new ')
10259599516SKenneth E. Jansenc
10359599516SKenneth E. Jansen        call itrRes (ytmp,            y,
10459599516SKenneth E. Jansen     &               x,               shp,
10559599516SKenneth E. Jansen     &               shgl,            iBC,
10659599516SKenneth E. Jansen     &               BC,              shpb,
10759599516SKenneth E. Jansen     &               shglb,           rtmp,
108513954efSKenneth E. Jansen     &               iper,            ilwork, ac)
109513954efSKenneth E. Jansen!Added ac to the end if itrRes, but not tested - Nicholas
11059599516SKenneth E. Jansenc
11159599516SKenneth E. Jansenc.... compute the second residual and add it to the first one
11259599516SKenneth E. Jansenc
11359599516SKenneth E. Jansen         ytmp = ypre - epsSD * uBrg
11459599516SKenneth E. Jansenc
11559599516SKenneth E. Jansenc         call yshuffle(ypre, 'old2new ')
11659599516SKenneth E. Jansenc
11759599516SKenneth E. Jansen        call i3LU (BDiag, ytmp, 'backward')
11859599516SKenneth E. Jansenc
11959599516SKenneth E. Jansen        call yshuffle(ytmp, 'old2new ')
12059599516SKenneth E. Jansen        call itrRes (ytmp,            y,
12159599516SKenneth E. Jansen     &               x,               shp,
12259599516SKenneth E. Jansen     &               shgl,            iBC,
12359599516SKenneth E. Jansen     &               BC,              shpb,
12459599516SKenneth E. Jansen     &               shglb,           rtmp,
125513954efSKenneth E. Jansen     &               iper,            ilwork, ac)
126513954efSKenneth E. Jansen!Added ac to the end if itrRes, but not tested - Nicholas
12759599516SKenneth E. Jansenc
12859599516SKenneth E. Jansen        call i3LU (BDiag, rtmp, 'forward ')
12959599516SKenneth E. Jansenc
13059599516SKenneth E. Jansenc.... compute the second derivative and its norm
13159599516SKenneth E. Jansenc
13259599516SKenneth E. Jansen        rtmp = (( rtmp - two * rmes ) / epsM)**2
13359599516SKenneth E. Jansenc
13459599516SKenneth E. Jansen        call sumgat (rtmp, nflow, summed, ilwork)
13559599516SKenneth E. Jansen        SDnrm = sqrt(summed)
13659599516SKenneth E. Jansenc
13759599516SKenneth E. Jansenc.... compute the 'optimum' interval
13859599516SKenneth E. Jansenc
13959599516SKenneth E. Jansen        eGMRES = two * sqrt( epsA / SDnrm )
14059599516SKenneth E. Jansenc
14159599516SKenneth E. Jansenc.... flop count
14259599516SKenneth E. Jansenc
143f4d0b58bSKenneth E. Jansen   !      flops = flops + 10*nflow*nshg+3*nshg
14459599516SKenneth E. Jansenc
14559599516SKenneth E. Jansenc.... end
14659599516SKenneth E. Jansenc
14759599516SKenneth E. Jansen        return
14859599516SKenneth E. Jansen        end
14959599516SKenneth E. Jansen
15059599516SKenneth E. Jansen
15159599516SKenneth E. Jansen        subroutine itrFDISclr (y,         ypre,      x,
15259599516SKenneth E. Jansen     &                         rmes,      uBrg,      BDiag,
15359599516SKenneth E. Jansen     &                         iBC,       BC,        engBC,     iper,
154*95eb1d92SBen Matthews     &                         ilwork,    shp,       shgl,      wght,
155*95eb1d92SBen Matthews     &                         shpb,      shglb,     wghtb)
15659599516SKenneth E. Jansenc
15759599516SKenneth E. Jansenc----------------------------------------------------------------------
15859599516SKenneth E. Jansenc
15959599516SKenneth E. Jansenc  This subroutine computes the "optimum" finite difference
16059599516SKenneth E. Jansenc  interval eps for the forward difference scheme
16159599516SKenneth E. Jansenc
16259599516SKenneth E. Jansenc                 Rmod(y + eps u) - Rmod(y)
16359599516SKenneth E. Jansenc                ---------------------------
16459599516SKenneth E. Jansenc                            eps
16559599516SKenneth E. Jansenc
16659599516SKenneth E. Jansenc  where  u is the step and Rmod is the modified residual.
16759599516SKenneth E. Jansenc
16859599516SKenneth E. Jansenc Note: A good theoretical reference is 'Practical Optimization' by
16959599516SKenneth E. Jansenc        P.E. Gill, W. Murray and M.H. Wright  [1981].
17059599516SKenneth E. Jansenc
17159599516SKenneth E. Jansenc input:
17259599516SKenneth E. Jansenc  y      (nshg,ndof)           : Y-variables
17359599516SKenneth E. Jansenc  ypre   (nshg,ndof)           : preconditioned Y-variables
17459599516SKenneth E. Jansenc  x      (numnp,nsd)            : node coordinates
17559599516SKenneth E. Jansenc  rmes   (nshg,nflow)           : modified residual
17659599516SKenneth E. Jansenc  uBrg   (nshg,nflow)           : step
17759599516SKenneth E. Jansenc  BDiag  (nshg,nflow,nflow)         : block-diagonal preconditioner
17859599516SKenneth E. Jansenc  iBC    (nshg)                : BC codes
17959599516SKenneth E. Jansenc  BC     (nshg,ndofBC)         : BC constraint parameters
18059599516SKenneth E. Jansenc  engBC  (nshg)                : energy for BC on density or pressure
18159599516SKenneth E. Jansenc
18259599516SKenneth E. Jansenc
18359599516SKenneth E. Jansenc Zdenek Johan, Winter 1989.
18459599516SKenneth E. Jansenc Zdenek Johan, Winter 1991.  (Fortran 90)
18559599516SKenneth E. Jansenc----------------------------------------------------------------------
18659599516SKenneth E. Jansenc
18759599516SKenneth E. Jansen        include "common.h"
18859599516SKenneth E. Jansenc
18959599516SKenneth E. Jansen        dimension y(nshg,ndof),               ypre(nshg,ndof),
19059599516SKenneth E. Jansen     &            x(numnp,nsd),
19159599516SKenneth E. Jansen     &            rmes(nshg,nflow),
19259599516SKenneth E. Jansen     &            uBrg(nshg,nflow),            BDiag(nshg,nflow,nflow),
19359599516SKenneth E. Jansen     &            iBC(nshg),                  BC(nshg,ndofBC),
19459599516SKenneth E. Jansen     &            engBC(nshg),                ilwork(nlwork),
19559599516SKenneth E. Jansen     &            iper(nshg)
19659599516SKenneth E. Jansenc
19759599516SKenneth E. Jansen        dimension ytmp(nshg,ndof),            rtmp(nshg,nflow)
19859599516SKenneth E. Jansenc
19959599516SKenneth E. Jansenc.... compute the accuracy (cancellation error)  ->  epsA
20059599516SKenneth E. Jansenc
20159599516SKenneth E. Jansen        rtmp = zero
20259599516SKenneth E. Jansenc
20359599516SKenneth E. Jansen        ytmp = ypre
20459599516SKenneth E. Jansenc
20559599516SKenneth E. Jansenc       call tnanq(ytmp,ndof,"ytmp       ")
20659599516SKenneth E. Jansen        iabres = 1
20759599516SKenneth E. Jansenc
20859599516SKenneth E. Jansen        call itrRes (ytmp,            y,
209*95eb1d92SBen Matthews     &               x,               shp,
210*95eb1d92SBen Matthews     &               shgl,            wght,      iBC,
211*95eb1d92SBen Matthews     &               BC,              engBC,         shpb,
212*95eb1d92SBen Matthews     &               shglb,           wghtb,     rtmp,
213513954efSKenneth E. Jansen     &               iper,            ilwork, ac)
214513954efSKenneth E. Jansen!Added ac to the end if itrRes, but not tested - Nicholas
21559599516SKenneth E. Jansenc
21659599516SKenneth E. Jansen        iabres = 0
21759599516SKenneth E. Jansenc
21859599516SKenneth E. Jansen        rtmp = rtmp**2
21959599516SKenneth E. Jansen        call sumgat (rtmp, nflow, summed, ilwork)
22059599516SKenneth E. Jansen        epsA = (epsM**2) * sqrt(summed)
22159599516SKenneth E. Jansenc
22259599516SKenneth E. Jansenc.... compute the norm of the second derivative (truncation error)
22359599516SKenneth E. Jansenc
22459599516SKenneth E. Jansenc.... set interval
22559599516SKenneth E. Jansenc
22659599516SKenneth E. Jansen        epsSD = sqrt(epsM)
22759599516SKenneth E. Jansenc
22859599516SKenneth E. Jansenc.... compute the first residual
22959599516SKenneth E. Jansenc
23059599516SKenneth E. Jansen        rtmp = zero
23159599516SKenneth E. Jansenc
23259599516SKenneth E. Jansen        ytmp = ypre + epsSD * uBrg
23359599516SKenneth E. Jansenc
23459599516SKenneth E. Jansen        call itrRes (ytmp,            y,
235*95eb1d92SBen Matthews     &               x,               shp,
236*95eb1d92SBen Matthews     &               shgl,            wght,      iBC,
237*95eb1d92SBen Matthews     &               BC,              engBC,         shpb,
238*95eb1d92SBen Matthews     &               shglb,           wghtb,     rtmp,
239513954efSKenneth E. Jansen     &               iper,            ilwork, ac)
240513954efSKenneth E. Jansen!Added ac to the end if itrRes, but not tested - Nicholas
24159599516SKenneth E. Jansenc
24259599516SKenneth E. Jansenc.... compute the second residual and add it to the first one
24359599516SKenneth E. Jansenc
24459599516SKenneth E. Jansen        ytmp = ypre - epsSD * uBrg
24559599516SKenneth E. Jansenc
24659599516SKenneth E. Jansen        call itrRes (ytmp,            y,
247*95eb1d92SBen Matthews     &               x,               shp,
248*95eb1d92SBen Matthews     &               shgl,            wght,      iBC,
249*95eb1d92SBen Matthews     &               BC,              engBC,         shpb,
250*95eb1d92SBen Matthews     &               shglb,           wghtb,     rtmp,
251513954efSKenneth E. Jansen     &               iper,            ilwork, ac)
252513954efSKenneth E. Jansen!Added ac to the end if itrRes, but not tested - Nicholas
25359599516SKenneth E. Jansenc
25459599516SKenneth E. Jansenc.... compute the second derivative and its norm
25559599516SKenneth E. Jansenc
25659599516SKenneth E. Jansen        rtmp = (( rtmp - two * rmes ) / epsM)**2
25759599516SKenneth E. Jansenc
25859599516SKenneth E. Jansen        call sumgat (rtmp, nflow, summed, ilwork)
25959599516SKenneth E. Jansen        SDnrm = sqrt(summed)
26059599516SKenneth E. Jansenc
26159599516SKenneth E. Jansenc.... compute the 'optimum' interval
26259599516SKenneth E. Jansenc
26359599516SKenneth E. Jansen        eGMRES = two * sqrt( epsA / SDnrm )
26459599516SKenneth E. Jansenc
26559599516SKenneth E. Jansenc.... flop count
26659599516SKenneth E. Jansenc
267f4d0b58bSKenneth E. Jansen   !      flops = flops + 10*nflow*nshg+3*nshg
26859599516SKenneth E. Jansenc
26959599516SKenneth E. Jansenc.... end
27059599516SKenneth E. Jansenc
27159599516SKenneth E. Jansen        return
27259599516SKenneth E. Jansen        end
273