xref: /phasta/phSolver/compressible/au1mfg.f (revision 712d3df0b59ebebaaeaea358162c8d2c043c6e08)
159599516SKenneth E. Jansen        subroutine Au1MFG (ypre,      y,    ac,     x,
259599516SKenneth E. Jansen     &                     rmes,      res,       uBrg,
359599516SKenneth E. Jansen     &                     BDiag,     iBC,       BC,
459599516SKenneth E. Jansen     &                     iper,      ilwork,    shp,
559599516SKenneth E. Jansen     &                     shgl,      shpb,
659599516SKenneth E. Jansen     &                     shglb)
759599516SKenneth E. Jansenc
859599516SKenneth E. Jansenc----------------------------------------------------------------------
959599516SKenneth E. Jansenc
1059599516SKenneth E. Jansenc This routine performs a matrix-vector product for the Matrix-Free
1159599516SKenneth E. Jansenc  Implicit/Iterative solver using a one-sided scheme.
1259599516SKenneth E. Jansenc
1359599516SKenneth E. Jansenc input:
1459599516SKenneth E. Jansenc  y      (nshg,ndof)           : Y-variables
1559599516SKenneth E. Jansenc  ypre   (nshg,nflow)           : preconditioned Y-variables
1659599516SKenneth E. Jansenc                                  (perturbed, no-scalars)
1759599516SKenneth E. Jansenc  x      (numnp,nsd)            : node coordinates
1859599516SKenneth E. Jansenc  rmes   (nshg,nflow)           : modified residual
1959599516SKenneth E. Jansenc  res    (nshg,nflow)           : residual
2059599516SKenneth E. Jansenc  uBrg   (nshg,nflow)           : Krylov space vector
2159599516SKenneth E. Jansenc  BDiag  (nshg,nflow,nflow)         : block-diagonal preconditioner
2259599516SKenneth E. Jansenc  iBC    (nshg)                : BC codes
2359599516SKenneth E. Jansenc  BC     (nshg,ndofBC)         : BC constraint parameters
2459599516SKenneth E. Jansenc  engBC  (nshg)                : energy for BC on density or pressure
2559599516SKenneth E. Jansenc  shp(b) (nshape,ngauss)        : element shape functions (boundary)
2659599516SKenneth E. Jansenc  shgl(b)(nsd,nshape,ngauss)    : local gradients of shape functions
2759599516SKenneth E. Jansenc
2859599516SKenneth E. Jansenc output:
2959599516SKenneth E. Jansenc  uBrg   (nshg,nflow)           : Krylov space vectors
3059599516SKenneth E. Jansenc
3159599516SKenneth E. Jansenc
3259599516SKenneth E. Jansenc Zdenek Johan, Winter 1991.  (Fortran 90)
3359599516SKenneth E. Jansenc----------------------------------------------------------------------
3459599516SKenneth E. Jansenc
3559599516SKenneth E. Jansen        include "common.h"
3659599516SKenneth E. Jansen        include "mpif.h"
3759599516SKenneth E. Jansen        include "auxmpi.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),          ytmp(nshg,nflow),
4259599516SKenneth E. Jansen     &            res(nshg,nflow),           uBrg(nshg,nflow),
4359599516SKenneth E. Jansen     &            BDiag(nshg,nflow,nflow),       iBC(nshg),
4459599516SKenneth E. Jansen     &            BC(nshg,ndofBC),          iper(nshg)
4559599516SKenneth E. Jansenc
4659599516SKenneth E. Jansen        dimension uBtmp(nshg,nflow),         tmpBC(nshg),
4759599516SKenneth E. Jansen     &            ilwork(nlwork)
4859599516SKenneth E. Jansenc
4959599516SKenneth E. Jansen
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. Jansen
5559599516SKenneth E. Jansenc$$$        dimension shp(nshl,ngauss),  shgl(nsd,nshl,ngauss),
5659599516SKenneth E. Jansenc$$$     &            shpb(nshl,ngaussb),
5759599516SKenneth E. Jansenc$$$     &            shglb(nsd,nshl,ngaussb)
5859599516SKenneth E. Jansenc
5959599516SKenneth E. Jansenc.... calculate Rmod(y + eps u_i)
6059599516SKenneth E. Jansenc
6159599516SKenneth E. Jansen        uBtmp = zero
6259599516SKenneth E. Jansenc
6359599516SKenneth E. Jansenc        call yshuffle(ypre, 'new2old ')
6459599516SKenneth E. Jansen        uBrg = ypre + eGMRES * uBrg
6559599516SKenneth E. Jansenc
6659599516SKenneth E. Jansen        call i3LU (BDiag,  uBrg,  'backward')
6759599516SKenneth E. Jansenc
6859599516SKenneth E. Jansen        call yshuffle(uBrg, 'old2new ')
6959599516SKenneth E. Jansenc
7059599516SKenneth E. Jansen        call itrBC (uBrg, uBrg, iBC,  BC, iper, ilwork)
7159599516SKenneth E. Jansenc
7259599516SKenneth E. Jansen        call itrRes (uBrg,            y,
7359599516SKenneth E. Jansen     &               x,               shp,
7459599516SKenneth E. Jansen     &               shgl,            iBC,
7559599516SKenneth E. Jansen     &               BC,              shpb,
7659599516SKenneth E. Jansen     &               shglb,           uBtmp,
7759599516SKenneth E. Jansen     &               iper,            ilwork, ac)
7859599516SKenneth E. Jansenc
7959599516SKenneth E. Jansen        call i3LU (BDiag,  uBtmp,  'forward ')
8059599516SKenneth E. Jansenc
8159599516SKenneth E. Jansenc.... calculate ( Rmod(y + eps u_i) - Rmod(y) ) / eps
8259599516SKenneth E. Jansenc
8359599516SKenneth E. Jansen        uBrg = ( uBtmp - rmes ) / eGMRES
8459599516SKenneth E. Jansenc ... before returning lets put ypre back in the new format
8559599516SKenneth E. Jansenc         call yshuffle(ypre, 'old2new ')
8659599516SKenneth E. Jansenc
8759599516SKenneth E. Jansenc.... flop count
8859599516SKenneth E. Jansenc
89*f4d0b58bSKenneth E. Jansen   !      flops = flops + 4*nflow*nshg
9059599516SKenneth E. Jansenc
9159599516SKenneth E. Jansenc.... end
9259599516SKenneth E. Jansenc
9359599516SKenneth E. Jansen        return
9459599516SKenneth E. Jansen        end
9559599516SKenneth E. Jansen
9659599516SKenneth E. Jansen
9759599516SKenneth E. Jansen
9859599516SKenneth E. Jansen
99