xref: /phasta/phSolver/compressible/solgmr.f (revision d06028c1ed09954bf598b0145fe1ac8aecabad11)
159599516SKenneth E. Jansen        subroutine SolGMRe (y,         ac,        yold,      acold,
259599516SKenneth E. Jansen     &			   x,         iBC,       BC,        EGmass,
359599516SKenneth E. Jansen     &                     res,       BDiag,     HBrg,      eBrg,
459599516SKenneth E. Jansen     &                     yBrg,      Rcos,      Rsin,      iper,
559599516SKenneth E. Jansen     &                     ilwork,    shp,       shgl,      shpb,
659599516SKenneth E. Jansen     &                     shglb,     Dy, rerr)
759599516SKenneth E. Jansenc
859599516SKenneth E. Jansenc----------------------------------------------------------------------
959599516SKenneth E. Jansenc
1059599516SKenneth E. Jansenc  This is the preconditioned GMRES driver routine.
1159599516SKenneth E. Jansenc
1259599516SKenneth E. Jansenc input:
1359599516SKenneth E. Jansenc  y      (nshg,ndof)           : Y-variables at n+alpha_v
1459599516SKenneth E. Jansenc  ac     (nshg,ndof)           : Primvar. accel. variable n+alpha_m
1559599516SKenneth E. Jansenc  yold   (nshg,ndof)           : Y-variables at beginning of step
1659599516SKenneth E. Jansenc  acold  (nshg,ndof)           : Primvar. accel. variable at begng step
1759599516SKenneth E. Jansenc  x      (numnp,nsd)            : node coordinates
1859599516SKenneth E. Jansenc  iBC    (nshg)                : BC codes
1959599516SKenneth E. Jansenc  BC     (nshg,ndofBC)         : BC constraint parameters
2059599516SKenneth E. Jansenc  HBrg   (Kspace+1,Kspace)      : Hessenberg matrix (LHS matrix)
2159599516SKenneth E. Jansenc  eBrg   (Kspace+1)             : RHS      of Hessenberg minim. problem
2259599516SKenneth E. Jansenc  yBrg   (Kspace)               : solution of Hessenberg minim. problem
2359599516SKenneth E. Jansenc  Rcos   (Kspace)               : Rotational cosine of QR algorithm
2459599516SKenneth E. Jansenc  Rsin   (Kspace)               : Rotational sine   of QR algorithm
2559599516SKenneth E. Jansenc  shp(b) (nen,maxsh,melCat)     : element shape functions (boundary)
2659599516SKenneth E. Jansenc  shgl(b)(nsd,nen,maxsh,melCat) : local gradients of shape functions
2759599516SKenneth E. Jansenc
2859599516SKenneth E. Jansenc output:
2959599516SKenneth E. Jansenc  res    (nshg,nflow)           : preconditioned residual
3059599516SKenneth E. Jansenc  BDiag  (nshg,nflow,nflow)      : block-diagonal preconditioner
3159599516SKenneth E. Jansenc
3259599516SKenneth E. Jansenc
3359599516SKenneth E. Jansenc Zdenek Johan,  Winter 1991.  (Fortran 90)
3459599516SKenneth E. Jansenc----------------------------------------------------------------------
3559599516SKenneth E. Jansenc
3659599516SKenneth E. Jansen        use pointer_data
3759599516SKenneth E. Jansen
3859599516SKenneth E. Jansen        include "common.h"
3959599516SKenneth E. Jansen        include "mpif.h"
4059599516SKenneth E. Jansen        include "auxmpi.h"
4159599516SKenneth E. Jansenc
4259599516SKenneth E. Jansen        dimension y(nshg,ndof),             ac(nshg,ndof),
4359599516SKenneth E. Jansen     &            yold(nshg,ndof),          acold(nshg,ndof),
4459599516SKenneth E. Jansen     &            x(numnp,nsd),
4559599516SKenneth E. Jansen     &            iBC(nshg),                BC(nshg,ndofBC),
4659599516SKenneth E. Jansen     &            res(nshg,nflow),
4759599516SKenneth E. Jansen     &            BDiag(nshg,nflow,nflow),
4859599516SKenneth E. Jansen     &            HBrg(Kspace+1,*),         eBrg(*),
4959599516SKenneth E. Jansen     &            yBrg(*),                  Rcos(*),
5059599516SKenneth E. Jansen     &            Rsin(*),                  ilwork(nlwork),
5159599516SKenneth E. Jansen     &            iper(nshg),               EGmass(numel,nedof,nedof)!,
5259599516SKenneth E. Jansenctoomuchmem     &            Binv(numel,nedof,nedof)
5359599516SKenneth E. Jansenc
5459599516SKenneth E. Jansen        dimension Dy(nshg,nflow),            rmes(nshg,nflow),
5559599516SKenneth E. Jansen     &            temp(nshg,nflow),
5659599516SKenneth E. Jansen     &            uBrg(nshg,nflow,Kspace+1)
5759599516SKenneth E. Jansenc
5859599516SKenneth E. Jansen        dimension shp(MAXTOP,maxsh,MAXQPT),
5959599516SKenneth E. Jansen     &            shgl(MAXTOP,nsd,maxsh,MAXQPT),
6059599516SKenneth E. Jansen     &            shpb(MAXTOP,maxsh,MAXQPT),
6159599516SKenneth E. Jansen     &            shglb(MAXTOP,nsd,maxsh,MAXQPT)
6259599516SKenneth E. Jansen      real*8    rerr(nshg,10)
6359599516SKenneth E. Jansenc
6459599516SKenneth E. Jansenc.... *******************>> Element Data Formation <<******************
6559599516SKenneth E. Jansenc
6659599516SKenneth E. Jansenc
6759599516SKenneth E. Jansenc.... set the parameters for flux and surface tension calculations
6859599516SKenneth E. Jansenc
6959599516SKenneth E. Jansenc
7059599516SKenneth E. Jansen        idflx = 0
7159599516SKenneth E. Jansen        if(idiff >= 1)  idflx= idflx + (nflow-1) * nsd
7259599516SKenneth E. Jansen        if (isurf == 1) idflx=idflx + nsd
7359599516SKenneth E. Jansenc
7459599516SKenneth E. Jansenc.... form the LHS matrices, the residual vector, and the block
7559599516SKenneth E. Jansenc     diagonal preconditioner
7659599516SKenneth E. Jansenc
7759599516SKenneth E. Jansen        call ElmGMRe(y,             ac,            x,
7859599516SKenneth E. Jansen     &               shp,           shgl,          iBC,
7959599516SKenneth E. Jansen     &               BC,            shpb,
8059599516SKenneth E. Jansen     &               shglb,         res,
8159599516SKenneth E. Jansen     &               rmes,          BDiag,         iper,
8259599516SKenneth E. Jansen     &               ilwork,        EGmass,        rerr )
83*6d194905SKenneth E. Jansen      rmes=res  ! saving the b vector (residual)
8459599516SKenneth E. Jansenc
8559599516SKenneth E. Jansenc.... **********************>>    EBE - GMRES    <<********************
8659599516SKenneth E. Jansenc
8759599516SKenneth E. Jansen        call timer ('Solver  ')
8859599516SKenneth E. Jansenc
8959599516SKenneth E. Jansenc.... ------------------------> Initialization <-----------------------
9059599516SKenneth E. Jansenc
9159599516SKenneth E. Jansenc
9259599516SKenneth E. Jansenc.... LU decompose the block diagonals
9359599516SKenneth E. Jansenc
9459599516SKenneth E. Jansen        if (iprec .ne. 0)
9559599516SKenneth E. Jansen     &  call i3LU (BDiag, res,  'LU_Fact ')
9659599516SKenneth E. Jansen
9759599516SKenneth E. Jansenc
9859599516SKenneth E. Jansenc.... block diagonal precondition residual
9959599516SKenneth E. Jansenc
10059599516SKenneth E. Jansen        call i3LU (BDiag, res,  'forward ')
10159599516SKenneth E. Jansenc
10259599516SKenneth E. Jansenc.... initialize Dy
10359599516SKenneth E. Jansenc
10459599516SKenneth E. Jansen        Dy = zero
10559599516SKenneth E. Jansenc
10659599516SKenneth E. Jansenc.... Pre-precondition the LHS mass matrix and set up the element
10759599516SKenneth E. Jansenc     by element preconditioners
10859599516SKenneth E. Jansenc
10959599516SKenneth E. Jansenctoomuchmemory note that Binv is demoted from huge array to just one
11059599516SKenneth E. Jansenc        real*8 in i3pre because it takes too much memory
11159599516SKenneth E. Jansen
11259599516SKenneth E. Jansen        call i3pre (BDiag,    Binv,   EGmass,  ilwork)
11359599516SKenneth E. Jansenc
11459599516SKenneth E. Jansenc.... left EBE precondition the residual
11559599516SKenneth E. Jansenc
11659599516SKenneth E. Jansenctoomuchmem        call i3Pcond (Binv,  res,  ilwork,   'L_Pcond ')
11759599516SKenneth E. Jansenc
11859599516SKenneth E. Jansenc.... copy res in uBrg(1)
11959599516SKenneth E. Jansenc
12059599516SKenneth E. Jansen        uBrg(:,:,1) = res
12159599516SKenneth E. Jansenc
12259599516SKenneth E. Jansenc.... calculate norm of residual
12359599516SKenneth E. Jansenc
12459599516SKenneth E. Jansen        temp  = res**2
12559599516SKenneth E. Jansen
12659599516SKenneth E. Jansen        call sumgat (temp, nflow, summed, ilwork)
12759599516SKenneth E. Jansen        unorm = sqrt(summed)
12859599516SKenneth E. Jansenc
12959599516SKenneth E. Jansenc.... check if GMRES iterations are required
13059599516SKenneth E. Jansenc
13159599516SKenneth E. Jansen        iKs    = 0
13259599516SKenneth E. Jansen        lGMRES = 0
13359599516SKenneth E. Jansenc
13459599516SKenneth E. Jansenc.... if we are down to machine precision, don't bother solving
13559599516SKenneth E. Jansenc
13659599516SKenneth E. Jansen        if (unorm .lt. 100.*epsM**2) goto 3000
13759599516SKenneth E. Jansenc
13859599516SKenneth E. Jansenc.... set up tolerance of the Hessenberg's problem
13959599516SKenneth E. Jansenc
14059599516SKenneth E. Jansen        epsnrm = etol * unorm
14159599516SKenneth E. Jansenc
14259599516SKenneth E. Jansenc.... ------------------------>  GMRES Loop  <-------------------------
14359599516SKenneth E. Jansenc
14459599516SKenneth E. Jansenc.... loop through GMRES cycles
14559599516SKenneth E. Jansenc
14659599516SKenneth E. Jansen        do 2000 mGMRES = 1, nGMRES
14759599516SKenneth E. Jansen        lGMRES = mGMRES - 1
14859599516SKenneth E. Jansenc
14959599516SKenneth E. Jansen        if (lGMRES .gt. 0) then
15059599516SKenneth E. Jansenc
15159599516SKenneth E. Jansenc.... if GMRES restarts are necessary, calculate  R - A x
15259599516SKenneth E. Jansenc
15359599516SKenneth E. Jansenc
15459599516SKenneth E. Jansenc.... right precondition Dy
15559599516SKenneth E. Jansenc
15659599516SKenneth E. Jansen           temp = Dy
15759599516SKenneth E. Jansen
15859599516SKenneth E. Jansenctoomuchmem           call i3Pcond (Binv,  temp,  ilwork,  'R_Pcond ')
15959599516SKenneth E. Jansenc
16059599516SKenneth E. Jansenc.... perform the A x product
16159599516SKenneth E. Jansenc
16259599516SKenneth E. Jansen           call Au1GMR (EGmass,  temp,  ilwork, iBC,iper)
16359599516SKenneth E. Jansenc
16459599516SKenneth E. Jansenc.... periodic nodes have to assemble results to their partners
16559599516SKenneth E. Jansenc
16659599516SKenneth E. Jansen           call bc3per (iBC,  temp,  iper, ilwork, nflow)
16759599516SKenneth E. Jansenc
16859599516SKenneth E. Jansenc.... left preconditioning
16959599516SKenneth E. Jansenc
17059599516SKenneth E. Jansenctoomuchmem           call i3Pcond (Binv,  temp,  ilwork,  'L_Pcond ')
17159599516SKenneth E. Jansenc
17259599516SKenneth E. Jansenc.... subtract A x from residual and calculate the norm
17359599516SKenneth E. Jansenc
17459599516SKenneth E. Jansen           temp = res - temp
17559599516SKenneth E. Jansen           uBrg(:,:,1) = temp
17659599516SKenneth E. Jansenc
17759599516SKenneth E. Jansenc.... calculate the norm
17859599516SKenneth E. Jansenc
17959599516SKenneth E. Jansen           temp  = temp**2
18059599516SKenneth E. Jansen           call sumgat (temp, nflow, summed, ilwork)
18159599516SKenneth E. Jansen           unorm = sqrt(summed)
18259599516SKenneth E. Jansenc
18359599516SKenneth E. Jansenc.... flop count
18459599516SKenneth E. Jansenc
185f4d0b58bSKenneth E. Jansen      !      flops = flops + nflow*nshg+nshg
18659599516SKenneth E. Jansenc
18759599516SKenneth E. Jansen        endif
18859599516SKenneth E. Jansenc
18959599516SKenneth E. Jansenc.... set up RHS of the Hessenberg's problem
19059599516SKenneth E. Jansenc
19159599516SKenneth E. Jansen        call clear (eBrg, Kspace+1)
19259599516SKenneth E. Jansen        eBrg(1) = unorm
19359599516SKenneth E. Jansenc
19459599516SKenneth E. Jansenc.... normalize the first Krylov vector
19559599516SKenneth E. Jansenc
19659599516SKenneth E. Jansen        uBrg(:,:,1) = uBrg(:,:,1) / unorm
19759599516SKenneth E. Jansenc
19859599516SKenneth E. Jansenc.... loop through GMRES iterations
19959599516SKenneth E. Jansenc
20059599516SKenneth E. Jansen        do 1000 iK = 1, Kspace
20159599516SKenneth E. Jansen           iKs = iK
20259599516SKenneth E. Jansen
20359599516SKenneth E. Jansen           uBrg(:,:,iKs+1) = uBrg(:,:,iKs)
20459599516SKenneth E. Jansenc
20559599516SKenneth E. Jansenc.... right EBE precondition the LHS ( u_{i+1} <-- inverse(U) u_i )
20659599516SKenneth E. Jansenc
20759599516SKenneth E. Jansenctoomuchmem           call i3Pcond (Binv,  uBrg(:,:,iKs+1), ilwork,  'R_Pcond ')
20859599516SKenneth E. Jansenc
20959599516SKenneth E. Jansenc.... Au product  ( u_{i+1} <-- EGmass u_{i+1} )
21059599516SKenneth E. Jansenc
21159599516SKenneth E. Jansen           call Au1GMR ( EGmass, uBrg(:,:,iKs+1),  ilwork, iBC,iper)
21259599516SKenneth E. Jansenc
21359599516SKenneth E. Jansenc.... periodic nodes have to assemble results to their partners
21459599516SKenneth E. Jansenc
21559599516SKenneth E. Jansen           call bc3per (iBC,  uBrg(:,:,iKs+1),  iper, ilwork, nflow)
21659599516SKenneth E. Jansen
21759599516SKenneth E. Jansenc
21859599516SKenneth E. Jansenc.... left EBE precondition the LHS ( u_{i+1} <-- inverse(L) u_{i+1} )
21959599516SKenneth E. Jansenc
22059599516SKenneth E. Jansenctoomuchmem           call i3Pcond (Binv,  uBrg(:,:,iKs+1), ilwork, 'L_Pcond ')
22159599516SKenneth E. Jansenc
22259599516SKenneth E. Jansenc.... orthogonalize and get the norm
22359599516SKenneth E. Jansenc
22459599516SKenneth E. Jansen          do jK = 1, iKs+1
22559599516SKenneth E. Jansenc
22659599516SKenneth E. Jansen            if (jK .eq. 1) then
22759599516SKenneth E. Jansenc
22859599516SKenneth E. Jansen              temp = uBrg(:,:,iKs+1) * uBrg(:,:,1)  ! {u_{i+1}*u_1} vector
22959599516SKenneth E. Jansen              call sumgat (temp, nflow, beta, ilwork) ! sum vector=(u_{i+1},u_1)
23059599516SKenneth E. Jansenc
23159599516SKenneth E. Jansen            else
23259599516SKenneth E. Jansenc
23359599516SKenneth E. Jansenc project off jK-1 vector
23459599516SKenneth E. Jansenc
23559599516SKenneth E. Jansen              uBrg(:,:,iKs+1) = uBrg(:,:,iKs+1) - beta * uBrg(:,:,jK-1)
23659599516SKenneth E. Jansenc
23759599516SKenneth E. Jansen              temp = uBrg(:,:,iKs+1) * uBrg(:,:,jK) !{u_{i+1}*u_j} vector
23859599516SKenneth E. Jansen              call sumgat (temp, nflow, beta, ilwork) ! sum vector=(u_{i+1},u_j)
23959599516SKenneth E. Jansenc
24059599516SKenneth E. Jansen            endif
24159599516SKenneth E. Jansenc
24259599516SKenneth E. Jansen            HBrg(jK,iKs) = beta   ! put this in the Hessenberg Matrix
24359599516SKenneth E. Jansenc
24459599516SKenneth E. Jansen        enddo
24559599516SKenneth E. Jansenc
246f4d0b58bSKenneth E. Jansen   !      flops = flops + (3*iKs+1)*nflow*numnp+(iKs+1)*numnp
24759599516SKenneth E. Jansenc
24859599516SKenneth E. Jansenc  the last inner product was with what was left of the vector (after
24959599516SKenneth E. Jansenc  projecting off all of the previous vectors
25059599516SKenneth E. Jansenc
25159599516SKenneth E. Jansen        unorm           = sqrt(beta)
25259599516SKenneth E. Jansen        HBrg(iKs+1,iKs) = unorm   ! this fills the 1 sub diagonal band
25359599516SKenneth E. Jansenc
25459599516SKenneth E. Jansenc.... normalize the Krylov vector
25559599516SKenneth E. Jansenc
25659599516SKenneth E. Jansen        uBrg(:,:,iKs+1) = uBrg(:,:,iKs+1) / unorm  ! normalize the next Krylov
25759599516SKenneth E. Jansenc vector
25859599516SKenneth E. Jansenc
25959599516SKenneth E. Jansenc.... construct and reduce the Hessenberg Matrix
26059599516SKenneth E. Jansenc  since there is only one subdiagonal we can use a Givens rotation to
26159599516SKenneth E. Jansenc  rotate off each subdiagonal AS IT IS FORMED.   We do this because it
26259599516SKenneth E. Jansenc  allows us to check progress of solution and quit when satisfied.  Note
26359599516SKenneth E. Jansenc  that all future K vects will put a subdiagonal in the next column so
26459599516SKenneth E. Jansenc  there is no penalty to work ahead as  the rotation for the next vector
26559599516SKenneth E. Jansenc  will be unaffected by this rotation.
26659599516SKenneth E. Jansen
26759599516SKenneth E. Jansenc
26859599516SKenneth E. Jansenc     H Y = E ========>   R_i H Y = R_i E
26959599516SKenneth E. Jansenc
27059599516SKenneth E. Jansen           do jK = 1, iKs-1
27159599516SKenneth E. Jansen              tmp            =  Rcos(jK) * HBrg(jK,  iKs) +
27259599516SKenneth E. Jansen     &                          Rsin(jK) * HBrg(jK+1,iKs)
27359599516SKenneth E. Jansen              HBrg(jK+1,iKs) = -Rsin(jK) * HBrg(jK,  iKs) +
27459599516SKenneth E. Jansen     &                          Rcos(jK) * HBrg(jK+1,iKs)
27559599516SKenneth E. Jansen              HBrg(jK,  iKs) =  tmp
27659599516SKenneth E. Jansen           enddo
27759599516SKenneth E. Jansenc
27859599516SKenneth E. Jansen           tmp             = sqrt(HBrg(iKs,iKs)**2 + HBrg(iKs+1,iKs)**2)
27959599516SKenneth E. Jansen           Rcos(iKs)       = HBrg(iKs,  iKs) / tmp
28059599516SKenneth E. Jansen           Rsin(iKs)       = HBrg(iKs+1,iKs) / tmp
28159599516SKenneth E. Jansen           HBrg(iKs,  iKs) = tmp
28259599516SKenneth E. Jansen           HBrg(iKs+1,iKs) = zero
28359599516SKenneth E. Jansenc
28459599516SKenneth E. Jansenc.... rotate eBrg    R_i E
28559599516SKenneth E. Jansenc
28659599516SKenneth E. Jansen           tmp         = Rcos(iKs) * eBrg(iKs) + Rsin(iKs) * eBrg(iKs+1)
28759599516SKenneth E. Jansen           eBrg(iKs+1) =-Rsin(iKs) * eBrg(iKs) + Rcos(iKs) * eBrg(iKs+1)
28859599516SKenneth E. Jansen           eBrg(iKs)   = tmp
28959599516SKenneth E. Jansenc
29059599516SKenneth E. Jansenc.... check for convergence
29159599516SKenneth E. Jansenc
29259599516SKenneth E. Jansen           ntotGM = ntotGM + 1
29359599516SKenneth E. Jansen           echeck=abs(eBrg(iKs+1))
29459599516SKenneth E. Jansen           if (echeck .le. epsnrm) exit
29559599516SKenneth E. Jansenc
29659599516SKenneth E. Jansenc.... end of GMRES iteration loop
29759599516SKenneth E. Jansenc
29859599516SKenneth E. Jansen 1000   continue
29959599516SKenneth E. Jansenc
30059599516SKenneth E. Jansenc.... ------------------------->   Solution   <------------------------
30159599516SKenneth E. Jansenc
30259599516SKenneth E. Jansenc.... if converged or end of Krylov space
30359599516SKenneth E. Jansenc
30459599516SKenneth E. Jansenc.... solve for yBrg
30559599516SKenneth E. Jansenc
30659599516SKenneth E. Jansen        do jK = iKs, 1, -1
30759599516SKenneth E. Jansen           yBrg(jK) = eBrg(jK) / HBrg(jK,jK)
30859599516SKenneth E. Jansen           do lK = 1, jK-1
30959599516SKenneth E. Jansen              eBrg(lK) = eBrg(lK) - yBrg(jK) * HBrg(lK,jK)
31059599516SKenneth E. Jansen           enddo
31159599516SKenneth E. Jansen        enddo
31259599516SKenneth E. Jansenc
31359599516SKenneth E. Jansenc.... update Dy
31459599516SKenneth E. Jansenc
31559599516SKenneth E. Jansen        do jK = 1, iKs
31659599516SKenneth E. Jansen           Dy = Dy + yBrg(jK) * uBrg(:,:,jK)
31759599516SKenneth E. Jansen        enddo
31859599516SKenneth E. Jansenc
31959599516SKenneth E. Jansenc.... flop count
32059599516SKenneth E. Jansenc
321f4d0b58bSKenneth E. Jansen   !      flops = flops + (3*iKs+1)*nflow*nshg
32259599516SKenneth E. Jansenc
32359599516SKenneth E. Jansenc.... check for convergence
32459599516SKenneth E. Jansenc
32559599516SKenneth E. Jansen
32659599516SKenneth E. Jansen        echeck=abs(eBrg(iKs+1))
32759599516SKenneth E. Jansen        if (echeck .le. epsnrm) exit
32859599516SKenneth E. Jansen        if(myrank.eq.master) write(*,*)'solver tolerance %satisfaction',
32959599516SKenneth E. Jansen     &  (one-echeck/unorm)/(one-etol)*100
33059599516SKenneth E. Jansenc
33159599516SKenneth E. Jansenc.... end of mGMRES loop
33259599516SKenneth E. Jansenc
33359599516SKenneth E. Jansen 2000 continue
33459599516SKenneth E. Jansenc
33559599516SKenneth E. Jansenc.... ------------------------>   Converged   <------------------------
33659599516SKenneth E. Jansenc
33759599516SKenneth E. Jansenc.... if converged
33859599516SKenneth E. Jansenc
33959599516SKenneth E. Jansen 3000 continue
34059599516SKenneth E. Jansenc
34159599516SKenneth E. Jansenc.... back EBE precondition the results
34259599516SKenneth E. Jansenc
34359599516SKenneth E. Jansenctoomuchmem      call i3Pcond (Binv,  Dy, ilwork, 'R_Pcond ')
34459599516SKenneth E. Jansenc
34559599516SKenneth E. Jansenc.... back block-diagonal precondition the results
34659599516SKenneth E. Jansenc
34759599516SKenneth E. Jansen      call i3LU (BDiag, Dy, 'backward')
34859599516SKenneth E. Jansenc
34959599516SKenneth E. Jansenc
35059599516SKenneth E. Jansenc.... output the statistics
35159599516SKenneth E. Jansenc
352*6d194905SKenneth E. Jansen              call rstat (res, ilwork,rmes)
35359599516SKenneth E. Jansenc
35459599516SKenneth E. Jansenc.... stop the timer
35559599516SKenneth E. Jansenc
35659599516SKenneth E. Jansen 3002 continue                  ! no solve just res.
35759599516SKenneth E. Jansen      call timer ('Back    ')
35859599516SKenneth E. Jansenc
35959599516SKenneth E. Jansenc.... end
36059599516SKenneth E. Jansenc
36159599516SKenneth E. Jansen      return
36259599516SKenneth E. Jansen      end
36359599516SKenneth E. Jansen
36459599516SKenneth E. Jansen
36559599516SKenneth E. Jansen
36659599516SKenneth E. Jansen
36759599516SKenneth E. Jansen
36859599516SKenneth E. Jansen      subroutine SolGMRs(y,         ac,        yold,      acold,
36959599516SKenneth E. Jansen     &			 x,         iBC,       BC,
37059599516SKenneth E. Jansen     &                   col,       row,       lhsk,
37159599516SKenneth E. Jansen     &                   res,       BDiag,     HBrg,      eBrg,
37259599516SKenneth E. Jansen     &                   yBrg,      Rcos,      Rsin,      iper,
37359599516SKenneth E. Jansen     &                   ilwork,    shp,       shgl,      shpb,
37459599516SKenneth E. Jansen     &                   shglb,     Dy,        rerr)
37559599516SKenneth E. Jansenc
37659599516SKenneth E. Jansenc----------------------------------------------------------------------
37759599516SKenneth E. Jansenc
37859599516SKenneth E. Jansenc  This is the preconditioned GMRES driver routine.
37959599516SKenneth E. Jansenc
38059599516SKenneth E. Jansenc input:
38159599516SKenneth E. Jansenc  y      (nshg,ndof)           : Y-variables at n+alpha_v
38259599516SKenneth E. Jansenc  ac     (nshg,ndof)           : Primvar. accel. variable n+alpha_m
38359599516SKenneth E. Jansenc  yold   (nshg,ndof)           : Y-variables at beginning of step
38459599516SKenneth E. Jansenc  acold  (nshg,ndof)           : Primvar. accel. variable at begng step
38559599516SKenneth E. Jansenc  x      (numnp,nsd)            : node coordinates
38659599516SKenneth E. Jansenc  iBC    (nshg)                : BC codes
38759599516SKenneth E. Jansenc  BC     (nshg,ndofBC)         : BC constraint parameters
38859599516SKenneth E. Jansenc  HBrg   (Kspace+1,Kspace)      : Hessenberg matrix (LHS matrix)
38959599516SKenneth E. Jansenc  eBrg   (Kspace+1)             : RHS      of Hessenberg minim. problem
39059599516SKenneth E. Jansenc  yBrg   (Kspace)               : solution of Hessenberg minim. problem
39159599516SKenneth E. Jansenc  Rcos   (Kspace)               : Rotational cosine of QR algorithm
39259599516SKenneth E. Jansenc  Rsin   (Kspace)               : Rotational sine   of QR algorithm
39359599516SKenneth E. Jansenc  shp(b) (nen,maxsh,melCat)     : element shape functions (boundary)
39459599516SKenneth E. Jansenc  shgl(b)(nsd,nen,maxsh,melCat) : local gradients of shape functions
39559599516SKenneth E. Jansenc
39659599516SKenneth E. Jansenc output:
39759599516SKenneth E. Jansenc  res    (nshg,nflow)           : preconditioned residual
39859599516SKenneth E. Jansenc  BDiag  (nshg,nflow,nflow)      : block-diagonal preconditioner
39959599516SKenneth E. Jansenc
40059599516SKenneth E. Jansenc
40159599516SKenneth E. Jansenc Zdenek Johan,  Winter 1991.  (Fortran 90)
40259599516SKenneth E. Jansenc----------------------------------------------------------------------
40359599516SKenneth E. Jansenc
40459599516SKenneth E. Jansen      use pointer_data
40559599516SKenneth E. Jansen
40659599516SKenneth E. Jansen      include "common.h"
40759599516SKenneth E. Jansen      include "mpif.h"
40859599516SKenneth E. Jansen      include "auxmpi.h"
40959599516SKenneth E. Jansenc
41059599516SKenneth E. Jansen      integer col(nshg+1), row(nnz*nshg)
41159599516SKenneth E. Jansen      real*8 lhsK(nflow*nflow,nnz_tot)
41259599516SKenneth E. Jansen
41359599516SKenneth E. Jansen
41459599516SKenneth E. Jansen      dimension y(nshg,ndof),             ac(nshg,ndof),
41559599516SKenneth E. Jansen     &          yold(nshg,ndof),          acold(nshg,ndof),
41659599516SKenneth E. Jansen     &          x(numnp,nsd),
41759599516SKenneth E. Jansen     &          iBC(nshg),                BC(nshg,ndofBC),
41859599516SKenneth E. Jansen     &          res(nshg,nflow),
41959599516SKenneth E. Jansen     &          BDiag(nshg,nflow,nflow),
42059599516SKenneth E. Jansen     &          HBrg(Kspace+1,Kspace),    eBrg(Kspace+1),
42159599516SKenneth E. Jansen     &          yBrg(Kspace),             Rcos(Kspace),
42259599516SKenneth E. Jansen     &          Rsin(Kspace),             ilwork(nlwork),
42359599516SKenneth E. Jansen     &          iper(nshg)
42459599516SKenneth E. Jansenc
42559599516SKenneth E. Jansen      dimension Dy(nshg,nflow),            rmes(nshg,nflow),
42659599516SKenneth E. Jansen     &          temp(nshg,nflow),
42759599516SKenneth E. Jansen     &          uBrg(nshg,nflow,Kspace+1)
42859599516SKenneth E. Jansenc
42959599516SKenneth E. Jansen      dimension shp(MAXTOP,maxsh,MAXQPT),
43059599516SKenneth E. Jansen     &          shgl(MAXTOP,nsd,maxsh,MAXQPT),
43159599516SKenneth E. Jansen     &          shpb(MAXTOP,maxsh,MAXQPT),
43259599516SKenneth E. Jansen     &          shglb(MAXTOP,nsd,maxsh,MAXQPT)
43359599516SKenneth E. Jansen      real*8    rerr(nshg,10)
43459599516SKenneth E. Jansenc
43559599516SKenneth E. Jansenc
43659599516SKenneth E. Jansenc.... *******************>> Element Data Formation <<******************
43759599516SKenneth E. Jansenc
43859599516SKenneth E. Jansenc
43959599516SKenneth E. Jansenc.... set the parameters for flux and surface tension calculations
44059599516SKenneth E. Jansenc
44159599516SKenneth E. Jansenc
44259599516SKenneth E. Jansen      idflx = 0
44359599516SKenneth E. Jansen      if(idiff >= 1)  idflx= idflx + (nflow-1) * nsd
44459599516SKenneth E. Jansen      if (isurf == 1) idflx=idflx + nsd
44559599516SKenneth E. Jansenc
44659599516SKenneth E. Jansenc.... form the LHS matrices, the residual vector, and the block
44759599516SKenneth E. Jansenc     diagonal preconditioner
44859599516SKenneth E. Jansenc
44959599516SKenneth E. Jansen      call ElmGMRs(y,             ac,            x,
45059599516SKenneth E. Jansen     &             shp,           shgl,          iBC,
45159599516SKenneth E. Jansen     &             BC,            shpb,
45259599516SKenneth E. Jansen     &             shglb,         res,
45359599516SKenneth E. Jansen     &             rmes,          BDiag,         iper,
45459599516SKenneth E. Jansen     &             ilwork,        lhsK,          col,
45559599516SKenneth E. Jansen     &             row,           rerr )
456*6d194905SKenneth E. Jansen      rmes=res  ! saving the b vector (residual)
45713c4f35aSKenneth E. Jansenc
45859599516SKenneth E. Jansen
45959599516SKenneth E. Jansen	call tnanq(res,5, 'res_egmr')
46059599516SKenneth E. Jansen	call tnanq(BDiag,25, 'bdg_egmr')
46159599516SKenneth E. Jansenc
46259599516SKenneth E. Jansenc.... **********************>>    EBE - GMRES    <<********************
46359599516SKenneth E. Jansenc
46459599516SKenneth E. Jansen      call timer ('Solver  ')
46559599516SKenneth E. Jansenc
46659599516SKenneth E. Jansenc.... ------------------------> Initialization <-----------------------
46759599516SKenneth E. Jansenc
46859599516SKenneth E. Jansenc
46959599516SKenneth E. Jansenc.... LU decompose the block diagonals
47059599516SKenneth E. Jansenc
47159599516SKenneth E. Jansen      if (iprec .ne. 0) then
47259599516SKenneth E. Jansen         call i3LU (BDiag, res,  'LU_Fact ')
47359599516SKenneth E. Jansen         if (numpe > 1) then
47459599516SKenneth E. Jansen            call commu (BDiag  , ilwork, nflow*nflow  , 'out')
47559599516SKenneth E. Jansen         endif
47659599516SKenneth E. Jansen      endif
47759599516SKenneth E. Jansenc
47859599516SKenneth E. Jansenc.... block diagonal precondition residual
47959599516SKenneth E. Jansenc
48059599516SKenneth E. Jansen      call i3LU (BDiag, res,  'forward ')
481*6d194905SKenneth E. Jansen!  from this point forward b is btilde (Preconditioned residual)
48259599516SKenneth E. Jansenc
48359599516SKenneth E. Jansenc Check the residual for divering trend
48459599516SKenneth E. Jansenc
48559599516SKenneth E. Jansen	call rstatCheck(res,ilwork,y,ac)
48659599516SKenneth E. Jansenc
48759599516SKenneth E. Jansenc.... initialize Dy
48859599516SKenneth E. Jansenc
48959599516SKenneth E. Jansen      Dy = zero
49059599516SKenneth E. Jansenc
49159599516SKenneth E. Jansenc.... Pre-precondition the LHS mass matrix and set up the sparse
49259599516SKenneth E. Jansenc     preconditioners
49359599516SKenneth E. Jansenc
49459599516SKenneth E. Jansen
49559599516SKenneth E. Jansen      if(lhs.eq.1) call Spsi3pre (BDiag,    lhsK,  col, row)
49659599516SKenneth E. Jansenc
49759599516SKenneth E. Jansenc.... copy res in uBrg(1)
49859599516SKenneth E. Jansenc
49959599516SKenneth E. Jansen      uBrg(:,:,1) = res
50059599516SKenneth E. Jansenc
50159599516SKenneth E. Jansenc.... calculate norm of residual
50259599516SKenneth E. Jansenc
50359599516SKenneth E. Jansen      temp  = res**2
50459599516SKenneth E. Jansen
50559599516SKenneth E. Jansen      call sumgat (temp, nflow, summed, ilwork)
50659599516SKenneth E. Jansen      unorm = sqrt(summed)
50759599516SKenneth E. Jansenc
50859599516SKenneth E. Jansenc.... check if GMRES iterations are required
50959599516SKenneth E. Jansenc
51059599516SKenneth E. Jansen      iKs    = 0
511513954efSKenneth E. Jansen      lGMRESs = 0
51259599516SKenneth E. Jansenc
51359599516SKenneth E. Jansenc.... if we are down to machine precision, don't bother solving
51459599516SKenneth E. Jansenc
51559599516SKenneth E. Jansen      if (unorm .lt. 100.*epsM**2) goto 3000
51659599516SKenneth E. Jansenc
51759599516SKenneth E. Jansenc.... set up tolerance of the Hessenberg's problem
51859599516SKenneth E. Jansenc
51959599516SKenneth E. Jansen      epsnrm = etol * unorm
52059599516SKenneth E. Jansenc
52159599516SKenneth E. Jansenc.... ------------------------>  GMRES Loop  <-------------------------
52259599516SKenneth E. Jansenc
52359599516SKenneth E. Jansenc.... loop through GMRES cycles
52459599516SKenneth E. Jansenc
52559599516SKenneth E. Jansen      do 2000 mGMRES = 1, nGMRES
526513954efSKenneth E. Jansen         lGMRESs = mGMRES - 1
52759599516SKenneth E. Jansenc
52859599516SKenneth E. Jansen         if (lGMRES .gt. 0) then
52959599516SKenneth E. Jansenc
53059599516SKenneth E. Jansenc.... if GMRES restarts are necessary, calculate  R - A x
53159599516SKenneth E. Jansenc
53259599516SKenneth E. Jansenc
53359599516SKenneth E. Jansenc.... right precondition Dy
53459599516SKenneth E. Jansenc
53559599516SKenneth E. Jansen            temp = Dy
53659599516SKenneth E. Jansen
53759599516SKenneth E. Jansenc
53859599516SKenneth E. Jansenc.... perform the A x product
53959599516SKenneth E. Jansenc
54059599516SKenneth E. Jansen            call SparseAp (iper,ilwork,iBC, col, row, lhsK,  temp)
541513954efSKenneth E. Jansenc           call tnanq(temp,5, 'q_spAPrs')
54259599516SKenneth E. Jansen
54359599516SKenneth E. Jansenc
54459599516SKenneth E. Jansenc.... periodic nodes have to assemble results to their partners
54559599516SKenneth E. Jansenc
54659599516SKenneth E. Jansen            call bc3per (iBC,  temp,  iper, ilwork, nflow)
547513954efSKenneth E. Jansenc           call tnanq(temp,5, 'q_BCprs')
54859599516SKenneth E. Jansenc
54959599516SKenneth E. Jansenc.... subtract A x from residual and calculate the norm
55059599516SKenneth E. Jansenc
55159599516SKenneth E. Jansen            temp = res - temp
55259599516SKenneth E. Jansen            uBrg(:,:,1) = temp
55359599516SKenneth E. Jansenc
55459599516SKenneth E. Jansenc.... calculate the norm
55559599516SKenneth E. Jansenc
55659599516SKenneth E. Jansen            temp  = temp**2
55759599516SKenneth E. Jansen            call sumgat (temp, nflow, summed, ilwork)
55859599516SKenneth E. Jansen            unorm = sqrt(summed)
55959599516SKenneth E. Jansenc
56059599516SKenneth E. Jansenc.... flop count
56159599516SKenneth E. Jansenc
562f4d0b58bSKenneth E. Jansen       !      flops = flops + nflow*nshg+nshg
56359599516SKenneth E. Jansenc
56459599516SKenneth E. Jansen         endif
56559599516SKenneth E. Jansenc
56659599516SKenneth E. Jansenc.... set up RHS of the Hessenberg's problem
56759599516SKenneth E. Jansenc
56859599516SKenneth E. Jansen         call clear (eBrg, Kspace+1)
56959599516SKenneth E. Jansen         eBrg(1) = unorm
57059599516SKenneth E. Jansenc
57159599516SKenneth E. Jansenc.... normalize the first Krylov vector
57259599516SKenneth E. Jansenc
57359599516SKenneth E. Jansen         uBrg(:,:,1) = uBrg(:,:,1) / unorm
57459599516SKenneth E. Jansenc
57559599516SKenneth E. Jansenc.... loop through GMRES iterations
57659599516SKenneth E. Jansenc
57759599516SKenneth E. Jansen         do 1000 iK = 1, Kspace
57859599516SKenneth E. Jansen            iKs = iK
57959599516SKenneth E. Jansen
58059599516SKenneth E. Jansen            uBrg(:,:,iKs+1) = uBrg(:,:,iKs)
58159599516SKenneth E. Jansenc
58259599516SKenneth E. Jansenc.... Au product  ( u_{i+1} <-- EGmass u_{i+1} )
58359599516SKenneth E. Jansenc
58459599516SKenneth E. Jansen            call SparseAp (iper, ilwork, iBC,
58559599516SKenneth E. Jansen     &                     col,  row,    lhsK,
58659599516SKenneth E. Jansen     &                     uBrg(:,:,iKs+1) )
587513954efSKenneth E. Jansenc           call tnanq(uBrg(:,:,iKS+1),5, 'q_spAP')
58859599516SKenneth E. Jansen
58959599516SKenneth E. Jansenc
59059599516SKenneth E. Jansenc.... periodic nodes have to assemble results to their partners
59159599516SKenneth E. Jansenc
59259599516SKenneth E. Jansen            call bc3per (iBC,  uBrg(:,:,iKs+1),  iper, ilwork, nflow)
593513954efSKenneth E. Jansenc           call tnanq(uBrg(:,:,iKS+1),5, 'q_bc')
59459599516SKenneth E. Jansen
59559599516SKenneth E. Jansenc
59659599516SKenneth E. Jansenc.... orthogonalize and get the norm
59759599516SKenneth E. Jansenc
59859599516SKenneth E. Jansen            do jK = 1, iKs+1
59959599516SKenneth E. Jansenc
60059599516SKenneth E. Jansen               if (jK .eq. 1) then
60159599516SKenneth E. Jansenc
60259599516SKenneth E. Jansen                  temp = uBrg(:,:,iKs+1) * uBrg(:,:,1) ! {u_{i+1}*u_1} vector
60359599516SKenneth E. Jansen                  call sumgat (temp, nflow, beta, ilwork) ! sum vector=(u_{i+1},u_1)
60459599516SKenneth E. Jansenc
60559599516SKenneth E. Jansen               else
60659599516SKenneth E. Jansenc
60759599516SKenneth E. Jansenc project off jK-1 vector
60859599516SKenneth E. Jansenc
60959599516SKenneth E. Jansen                  uBrg(:,:,iKs+1)=uBrg(:,:,iKs+1)-beta * uBrg(:,:,jK-1)
61059599516SKenneth E. Jansenc
61159599516SKenneth E. Jansen                  temp = uBrg(:,:,iKs+1) * uBrg(:,:,jK) !{u_{i+1}*u_j} vector
61259599516SKenneth E. Jansen                  call sumgat (temp, nflow, beta, ilwork) ! sum vector=(u_{i+1},u_j)
61359599516SKenneth E. Jansenc
61459599516SKenneth E. Jansen               endif
61559599516SKenneth E. Jansenc
61659599516SKenneth E. Jansen               HBrg(jK,iKs) = beta ! put this in the Hessenberg Matrix
61759599516SKenneth E. Jansenc
61859599516SKenneth E. Jansen            enddo
61959599516SKenneth E. Jansenc
620f4d0b58bSKenneth E. Jansen       !      flops = flops + (3*iKs+1)*nflow*numnp+(iKs+1)*numnp
62159599516SKenneth E. Jansenc
62259599516SKenneth E. Jansenc  the last inner product was with what was left of the vector (after
62359599516SKenneth E. Jansenc  projecting off all of the previous vectors
62459599516SKenneth E. Jansenc
625513954efSKenneth E. Jansen        if(beta.le.0) write(*,*) 'beta in solgmr non-positive'
62659599516SKenneth E. Jansen            unorm           = sqrt(beta)
62759599516SKenneth E. Jansen            HBrg(iKs+1,iKs) = unorm ! this fills the 1 sub diagonal band
62859599516SKenneth E. Jansenc
62959599516SKenneth E. Jansenc.... normalize the Krylov vector
63059599516SKenneth E. Jansenc
63159599516SKenneth E. Jansen            uBrg(:,:,iKs+1) = uBrg(:,:,iKs+1) / unorm ! normalize the next Krylov
63259599516SKenneth E. Jansenc vector
63359599516SKenneth E. Jansenc
63459599516SKenneth E. Jansenc.... construct and reduce the Hessenberg Matrix
63559599516SKenneth E. Jansenc  since there is only one subdiagonal we can use a Givens rotation to
63659599516SKenneth E. Jansenc  rotate off each subdiagonal AS IT IS FORMED.   We do this because it
63759599516SKenneth E. Jansenc  allows us to check progress of solution and quit when satisfied.  Note
63859599516SKenneth E. Jansenc  that all future K vects will put a subdiagonal in the next column so
63959599516SKenneth E. Jansenc  there is no penalty to work ahead as  the rotation for the next vector
64059599516SKenneth E. Jansenc  will be unaffected by this rotation.
64159599516SKenneth E. Jansen
64259599516SKenneth E. Jansenc
64359599516SKenneth E. Jansenc     H Y = E ========>   R_i H Y = R_i E
64459599516SKenneth E. Jansenc
64559599516SKenneth E. Jansen            do jK = 1, iKs-1
64659599516SKenneth E. Jansen               tmp            =  Rcos(jK) * HBrg(jK,  iKs) +
64759599516SKenneth E. Jansen     &                           Rsin(jK) * HBrg(jK+1,iKs)
64859599516SKenneth E. Jansen               HBrg(jK+1,iKs) = -Rsin(jK) * HBrg(jK,  iKs) +
64959599516SKenneth E. Jansen     &                           Rcos(jK) * HBrg(jK+1,iKs)
65059599516SKenneth E. Jansen               HBrg(jK,  iKs) =  tmp
65159599516SKenneth E. Jansen            enddo
65259599516SKenneth E. Jansenc
65359599516SKenneth E. Jansen            tmp            = sqrt(HBrg(iKs,iKs)**2 + HBrg(iKs+1,iKs)**2)
65459599516SKenneth E. Jansen            Rcos(iKs)      = HBrg(iKs,  iKs) / tmp
65559599516SKenneth E. Jansen            Rsin(iKs)      = HBrg(iKs+1,iKs) / tmp
65659599516SKenneth E. Jansen            HBrg(iKs,  iKs)= tmp
65759599516SKenneth E. Jansen            HBrg(iKs+1,iKs)= zero
65859599516SKenneth E. Jansenc
65959599516SKenneth E. Jansenc.... rotate eBrg    R_i E
66059599516SKenneth E. Jansenc
66159599516SKenneth E. Jansen            tmp        = Rcos(iKs) * eBrg(iKs) + Rsin(iKs) * eBrg(iKs+1)
66259599516SKenneth E. Jansen            eBrg(iKs+1)=-Rsin(iKs) * eBrg(iKs) + Rcos(iKs) * eBrg(iKs+1)
66359599516SKenneth E. Jansen            eBrg(iKs)  = tmp
66459599516SKenneth E. Jansenc
66559599516SKenneth E. Jansenc.... check for convergence
66659599516SKenneth E. Jansenc
66759599516SKenneth E. Jansen            ntotGM = ntotGM + 1
66859599516SKenneth E. Jansen            echeck=abs(eBrg(iKs+1))
66959599516SKenneth E. Jansen            if (echeck .le. epsnrm.and. iKs .ge. minIters) exit
67059599516SKenneth E. Jansenc
67159599516SKenneth E. Jansenc.... end of GMRES iteration loop
67259599516SKenneth E. Jansenc
67359599516SKenneth E. Jansen 1000    continue
67459599516SKenneth E. Jansenc
67559599516SKenneth E. Jansenc.... ------------------------->   Solution   <------------------------
67659599516SKenneth E. Jansenc
67759599516SKenneth E. Jansenc.... if converged or end of Krylov space
67859599516SKenneth E. Jansenc
67959599516SKenneth E. Jansenc.... solve for yBrg
68059599516SKenneth E. Jansenc
68159599516SKenneth E. Jansen         do jK = iKs, 1, -1
68259599516SKenneth E. Jansen            yBrg(jK) = eBrg(jK) / HBrg(jK,jK)
68359599516SKenneth E. Jansen            do lK = 1, jK-1
68459599516SKenneth E. Jansen               eBrg(lK) = eBrg(lK) - yBrg(jK) * HBrg(lK,jK)
68559599516SKenneth E. Jansen            enddo
68659599516SKenneth E. Jansen         enddo
68759599516SKenneth E. Jansenc
68859599516SKenneth E. Jansenc.... update Dy
68959599516SKenneth E. Jansenc
69059599516SKenneth E. Jansen         do jK = 1, iKs
69159599516SKenneth E. Jansen            Dy = Dy + yBrg(jK) * uBrg(:,:,jK)
69259599516SKenneth E. Jansen         enddo
69359599516SKenneth E. Jansenc
69459599516SKenneth E. Jansenc.... flop count
69559599516SKenneth E. Jansenc
696f4d0b58bSKenneth E. Jansen    !      flops = flops + (3*iKs+1)*nflow*nshg
69759599516SKenneth E. Jansenc
69859599516SKenneth E. Jansenc.... check for convergence
69959599516SKenneth E. Jansenc
70059599516SKenneth E. Jansen        echeck=abs(eBrg(iKs+1))
70159599516SKenneth E. Jansen        if (echeck .le. epsnrm) exit
702513954efSKenneth E. Jansen!        if(myrank.eq.master) write(*,*)'solver tolerance %satisfaction',
703513954efSKenneth E. Jansen!     &  (one-echeck*etol/epsnrm)/(one-etol)*100
70459599516SKenneth E. Jansen
70559599516SKenneth E. Jansenc
70659599516SKenneth E. Jansenc.... end of mGMRES loop
70759599516SKenneth E. Jansenc
70859599516SKenneth E. Jansen 2000 continue
70959599516SKenneth E. Jansenc
71059599516SKenneth E. Jansenc.... ------------------------>   Converged   <------------------------
71159599516SKenneth E. Jansenc
71259599516SKenneth E. Jansenc.... if converged
71359599516SKenneth E. Jansenc
71459599516SKenneth E. Jansen 3000 continue
71559599516SKenneth E. Jansen
71659599516SKenneth E. Jansenc
71759599516SKenneth E. Jansenc
71859599516SKenneth E. Jansenc.... back block-diagonal precondition the results
71959599516SKenneth E. Jansenc
72059599516SKenneth E. Jansen      call i3LU (BDiag, Dy, 'backward')
72159599516SKenneth E. Jansenc
72259599516SKenneth E. Jansenc
72359599516SKenneth E. Jansenc.... output the statistics
72459599516SKenneth E. Jansenc
725*6d194905SKenneth E. Jansen      call rstat (res, ilwork,rmes)
726513954efSKenneth E. Jansen
727513954efSKenneth E. Jansen      if(myrank.eq.master) then
728513954efSKenneth E. Jansen        if (echeck .le. epsnrm) then
729513954efSKenneth E. Jansen            write(*,*)
730513954efSKenneth E. Jansen        else
731513954efSKenneth E. Jansen            write(*,*)'solver tolerance %satisfaction',
732513954efSKenneth E. Jansen     &  (one-echeck*etol/epsnrm)/(one-etol)*100
733513954efSKenneth E. Jansen        endif
734513954efSKenneth E. Jansen      endif
73559599516SKenneth E. Jansenc
73659599516SKenneth E. Jansenc.... stop the timer
73759599516SKenneth E. Jansenc
73859599516SKenneth E. Jansen 3002 continue                  ! no solve just res.
73959599516SKenneth E. Jansen      call timer ('Back    ')
74059599516SKenneth E. Jansenc
74159599516SKenneth E. Jansenc.... end
74259599516SKenneth E. Jansenc
74359599516SKenneth E. Jansen      return
74459599516SKenneth E. Jansen      end
74559599516SKenneth E. Jansen
74659599516SKenneth E. Jansen        subroutine SolGMRSclr(y,       ac,      yold,
74759599516SKenneth E. Jansen     &                        acold,   EGmasst,
74859599516SKenneth E. Jansen     &                        x,       elDw,
74959599516SKenneth E. Jansen     &                        iBC,      BC,
75059599516SKenneth E. Jansen     &                        rest,     HBrg,     eBrg,
75159599516SKenneth E. Jansen     &                        yBrg,     Rcos,     Rsin,    iper,
75259599516SKenneth E. Jansen     &                        ilwork,
75359599516SKenneth E. Jansen     &                        shp,      shgl,
75459599516SKenneth E. Jansen     &                        shpb,     shglb,  Dyt)
75559599516SKenneth E. Jansenc
75659599516SKenneth E. Jansenc----------------------------------------------------------------------
75759599516SKenneth E. Jansenc
75859599516SKenneth E. Jansenc  This is the preconditioned GMRES driver routine.
75959599516SKenneth E. Jansenc
76059599516SKenneth E. Jansenc input:
76159599516SKenneth E. Jansenc  y      (nshg,ndof)           : Y-variables at n+alpha_f
76259599516SKenneth E. Jansenc  ac     (nshg,ndof)           : Primvar. accel. variable n+alpha_m
76359599516SKenneth E. Jansenc  yold   (nshg,ndof)           : Y-variables at beginning of step
76459599516SKenneth E. Jansenc  acold  (nshg,ndof)           : Primvar. accel. variable at begng step
76559599516SKenneth E. Jansenc  x      (numnp,nsd)            : node coordinates
76659599516SKenneth E. Jansenc  iBC    (nshg)                : BC codes
76759599516SKenneth E. Jansenc  BC     (nshg,ndofBC)         : BC constraint parameters
76859599516SKenneth E. Jansenc  iper   (nshg)                : periodic nodal information
76959599516SKenneth E. Jansenc
77059599516SKenneth E. Jansenc output:
77159599516SKenneth E. Jansenc  y      (nshg,ndof)           : Y-variables at n+alpha_f
77259599516SKenneth E. Jansenc  ac     (nshg,ndof)           : Primvar. accel. variable n+alpha_m
77359599516SKenneth E. Jansenc  HBrg   (Kspace+1,Kspace)      : Hessenberg matrix (LHS matrix)
77459599516SKenneth E. Jansenc  eBrg   (Kspace+1)             : RHS      of Hessenberg minim. problem
77559599516SKenneth E. Jansenc  yBrg   (Kspace)               : solution of Hessenberg minim. problem
77659599516SKenneth E. Jansenc  Rcos   (Kspace)               : Rotational cosine of QR algorithm
77759599516SKenneth E. Jansenc  Rsin   (Kspace)               : Rotational sine   of QR algorithm
77859599516SKenneth E. Jansenc output:
77959599516SKenneth E. Jansenc  rest   (numnp)           : preconditioned residual
78059599516SKenneth E. Jansenc
78159599516SKenneth E. Jansenc
78259599516SKenneth E. Jansenc Zdenek Johan,  Winter 1991.  (Fortran 90)
78359599516SKenneth E. Jansenc----------------------------------------------------------------------
78459599516SKenneth E. Jansenc
78559599516SKenneth E. Jansen        use pointer_data
78659599516SKenneth E. Jansen
78759599516SKenneth E. Jansen        include "common.h"
78859599516SKenneth E. Jansen        include "mpif.h"
78959599516SKenneth E. Jansen        include "auxmpi.h"
79059599516SKenneth E. Jansenc
79159599516SKenneth E. Jansen        dimension y(nshg,ndof),      ac(nshg,ndof),
79259599516SKenneth E. Jansen     &            yold(nshg,ndof),   acold(nshg,ndof),
79359599516SKenneth E. Jansen     &            x(numnp,nsd),
79459599516SKenneth E. Jansen     &            iBC(nshg),         BC(nshg,ndofBC),
79559599516SKenneth E. Jansen     &            rest(nshg),
79659599516SKenneth E. Jansen     &            Diag(nshg),
79759599516SKenneth E. Jansen     &            HBrg(Kspace+1,*),  eBrg(*),
79859599516SKenneth E. Jansen     &            yBrg(*),           Rcos(*),
79959599516SKenneth E. Jansen     &            Rsin(*),           ilwork(nlwork),
80059599516SKenneth E. Jansen     &            iper(nshg),        EGmasst(numel,nshape,nshape)
80159599516SKenneth E. Jansenc
80259599516SKenneth E. Jansen        dimension Dyt(nshg),         rmest(nshg),
80359599516SKenneth E. Jansen     &            tempt(nshg),       Dinv(nshg),
80459599516SKenneth E. Jansen     &            uBrgt(nshg,Kspace+1)
80559599516SKenneth E. Jansenc
80659599516SKenneth E. Jansen        dimension shp(MAXTOP,maxsh,MAXQPT),
80759599516SKenneth E. Jansen     &            shgl(MAXTOP,nsd,maxsh,MAXQPT),
80859599516SKenneth E. Jansen     &            shpb(MAXTOP,maxsh,MAXQPT),
80959599516SKenneth E. Jansen     &            shglb(MAXTOP,nsd,maxsh,MAXQPT)
81059599516SKenneth E. Jansen        real*8    elDw(numel)
81159599516SKenneth E. Jansenc
81259599516SKenneth E. Jansenc.... *******************>> Element Data Formation <<******************
81359599516SKenneth E. Jansenc
81459599516SKenneth E. Jansenc
81559599516SKenneth E. Jansenc.... form the LHS matrices, the residual vector, and the block
81659599516SKenneth E. Jansenc     diagonal preconditioner
81759599516SKenneth E. Jansenc
81859599516SKenneth E. Jansen        call ElmGMRSclr(y,             ac,
81959599516SKenneth E. Jansen     &                  x,             elDw,      shp,        shgl,
82059599516SKenneth E. Jansen     &                  iBC,           BC,
82159599516SKenneth E. Jansen     &                  shpb,          shglb,
82259599516SKenneth E. Jansen     &                  rest,
82359599516SKenneth E. Jansen     &                  rmest,         Diag,       iper,
82459599516SKenneth E. Jansen     &                  ilwork,        EGmasst)
82559599516SKenneth E. Jansenc
82659599516SKenneth E. Jansenc.... **********************>>    EBE - GMRES    <<********************
82759599516SKenneth E. Jansenc
82859599516SKenneth E. Jansen        call timer ('Solver  ')
82959599516SKenneth E. Jansenc
83059599516SKenneth E. Jansenc.... ------------------------> Initialization <-----------------------
83159599516SKenneth E. Jansenc
83259599516SKenneth E. Jansenc
83359599516SKenneth E. Jansen      id = isclr+5
83459599516SKenneth E. Jansenc.... initialize Dy
83559599516SKenneth E. Jansenc
83659599516SKenneth E. Jansen        Dyt = zero
83759599516SKenneth E. Jansenc
83859599516SKenneth E. Jansenc.... Right preconditioner
83959599516SKenneth E. Jansenc
84059599516SKenneth E. Jansen        call i3preSclr(Diag, Dinv, EGmassT, ilwork)
84159599516SKenneth E. Jansenc
84259599516SKenneth E. Jansenc Check the residual for divering trend
84359599516SKenneth E. Jansenc
844513954efSKenneth E. Jansen
84559599516SKenneth E. Jansen        call rstatCheckSclr(rest,ilwork,y,ac)
84659599516SKenneth E. Jansen
84759599516SKenneth E. Jansenc
84859599516SKenneth E. Jansenc.... copy rest in uBrgt(1)
84959599516SKenneth E. Jansenc
85059599516SKenneth E. Jansen        uBrgt(:,1) = rest
85159599516SKenneth E. Jansenc
85259599516SKenneth E. Jansenc.... calculate norm of residual
85359599516SKenneth E. Jansenc
85459599516SKenneth E. Jansen        tempt  = rest**2
85559599516SKenneth E. Jansen
85659599516SKenneth E. Jansen        call sumgat (tempt, 1, summed, ilwork)
85759599516SKenneth E. Jansen        unorm = sqrt(summed)
85859599516SKenneth E. Jansenc
85959599516SKenneth E. Jansenc.... check if GMRES iterations are required
86059599516SKenneth E. Jansenc
861513954efSKenneth E. Jansen        iKss    = 0
86259599516SKenneth E. Jansen        lGMRESt = 0
86359599516SKenneth E. Jansenc
86459599516SKenneth E. Jansenc.... if we are down to machine precision, don't bother solving
86559599516SKenneth E. Jansenc
86659599516SKenneth E. Jansen        if (unorm .lt. 100.*epsM**2) goto 3000
86759599516SKenneth E. Jansenc
86859599516SKenneth E. Jansenc.... set up tolerance of the Hessenberg's problem
86959599516SKenneth E. Jansenc
87059599516SKenneth E. Jansen        epsnrm = etol * unorm
87159599516SKenneth E. Jansenc
87259599516SKenneth E. Jansenc.... ------------------------>  GMRES Loop  <-------------------------
87359599516SKenneth E. Jansenc
87459599516SKenneth E. Jansenc.... loop through GMRES cycles
87559599516SKenneth E. Jansenc
87659599516SKenneth E. Jansen        do 2000 mGMRES = 1, nGMRES
87759599516SKenneth E. Jansen        lGMRESt = mGMRES - 1
87859599516SKenneth E. Jansenc
87959599516SKenneth E. Jansen        if (lGMRESt .gt. 0) then
88059599516SKenneth E. Jansenc
88159599516SKenneth E. Jansenc.... if GMRES restarts are necessary, calculate  R - A x
88259599516SKenneth E. Jansenc
88359599516SKenneth E. Jansenc
88459599516SKenneth E. Jansen
88559599516SKenneth E. Jansenc.... perform the A x product
88659599516SKenneth E. Jansenc
88759599516SKenneth E. Jansen           call Au1GMRSclr (EGmasst,  tempt,  ilwork, iper)
88859599516SKenneth E. Jansenc
88959599516SKenneth E. Jansenc.... periodic nodes have to assemble results to their partners
89059599516SKenneth E. Jansenc
89159599516SKenneth E. Jansenc          call bc3perSclr (iBC,  tempt,  iper)
89259599516SKenneth E. Jansenc
89359599516SKenneth E. Jansenc.... subtract A x from residual and calculate the norm
89459599516SKenneth E. Jansenc
89559599516SKenneth E. Jansen           tempt = rest - tempt
89659599516SKenneth E. Jansen           uBrgt(:,1) = tempt
89759599516SKenneth E. Jansenc
89859599516SKenneth E. Jansenc.... calculate the norm
89959599516SKenneth E. Jansenc
90059599516SKenneth E. Jansen           tempt  = tempt**2
90159599516SKenneth E. Jansen           call sumgat (tempt, 1, summed, ilwork)
90259599516SKenneth E. Jansen           unorm = sqrt(summed)
90359599516SKenneth E. Jansenc
90459599516SKenneth E. Jansenc.... flop count
90559599516SKenneth E. Jansenc
906f4d0b58bSKenneth E. Jansen      !      flops = flops + ndof*numnp+numnp
90759599516SKenneth E. Jansenc
90859599516SKenneth E. Jansen        endif
90959599516SKenneth E. Jansenc
91059599516SKenneth E. Jansenc.... set up RHS of the Hessenberg's problem
91159599516SKenneth E. Jansenc
91259599516SKenneth E. Jansen        call clear (eBrg, Kspace+1)
91359599516SKenneth E. Jansen        call clear (HBrg, Kspace+1)
91459599516SKenneth E. Jansen        eBrg(1) = unorm
91559599516SKenneth E. Jansenc
91659599516SKenneth E. Jansenc.... normalize the first Krylov vector
91759599516SKenneth E. Jansenc
91859599516SKenneth E. Jansen        uBrgt(:,1) = uBrgt(:,1) / unorm
91959599516SKenneth E. Jansenc
92059599516SKenneth E. Jansenc.... loop through GMRES iterations
92159599516SKenneth E. Jansenc
92259599516SKenneth E. Jansen        do 1000 iK = 1, Kspace
923513954efSKenneth E. Jansen           iKss = iK
92459599516SKenneth E. Jansen
925513954efSKenneth E. Jansen           uBrgt(:,iKss+1) = uBrgt(:,iKss)
92659599516SKenneth E. Jansen
92759599516SKenneth E. Jansenc.... Au product  ( u_{i+1} <-- EGmass u_{i+1} )
92859599516SKenneth E. Jansenc
929513954efSKenneth E. Jansen           call Au1GMRSclr ( EGmasst, uBrgt(:,iKss+1),  ilwork, iper )
93059599516SKenneth E. Jansen
93159599516SKenneth E. Jansenc
93259599516SKenneth E. Jansenc.... periodic nodes have to assemble results to their partners
93359599516SKenneth E. Jansenc
934513954efSKenneth E. Jansen           call bc3perSclr (iBC,  uBrgt(:,iKss+1),  iper)
93559599516SKenneth E. Jansenc
93659599516SKenneth E. Jansenc.... orthogonalize and get the norm
93759599516SKenneth E. Jansenc
938513954efSKenneth E. Jansen          do jK = 1, iKss+1
93959599516SKenneth E. Jansenc
94059599516SKenneth E. Jansen            if (jK .eq. 1) then
94159599516SKenneth E. Jansenc
942513954efSKenneth E. Jansen              tempt = uBrgt(:,iKss+1) * uBrgt(:,1)  ! {u_{i+1}*u_1} vector
94359599516SKenneth E. Jansen              call sumgat (tempt, 1, beta, ilwork) ! sum vector=(u_{i+1},u_1)
94459599516SKenneth E. Jansenc
94559599516SKenneth E. Jansen            else
94659599516SKenneth E. Jansenc
94759599516SKenneth E. Jansenc project off jK-1 vector
94859599516SKenneth E. Jansenc
949513954efSKenneth E. Jansen          uBrgt(:,iKss+1) = uBrgt(:,iKss+1) - beta * uBrgt(:,jK-1)
95059599516SKenneth E. Jansenc
951513954efSKenneth E. Jansen              tempt = uBrgt(:,iKss+1) * uBrgt(:,jK) !{u_{i+1}*u_j} vector
95259599516SKenneth E. Jansen              call sumgat (tempt, 1, beta, ilwork) ! sum vector=(u_{i+1},u_j)
95359599516SKenneth E. Jansenc
95459599516SKenneth E. Jansen            endif
95559599516SKenneth E. Jansenc
956513954efSKenneth E. Jansen            HBrg(jK,iKss) = beta   ! put this in the Hessenberg Matrix
95759599516SKenneth E. Jansenc
95859599516SKenneth E. Jansen        enddo
95959599516SKenneth E. Jansenc
960f4d0b58bSKenneth E. Jansen   !      flops = flops + (3*iKss+1)*ndof*numnp+(iKss+1)*numnp
96159599516SKenneth E. Jansenc
96259599516SKenneth E. Jansenc  the last inner product was with what was left of the vector (after
96359599516SKenneth E. Jansenc  projecting off all of the previous vectors
96459599516SKenneth E. Jansenc
96559599516SKenneth E. Jansen        unorm           = sqrt(beta)
966513954efSKenneth E. Jansen        HBrg(iKss+1,iKss) = unorm   ! this fills the 1 sub diagonal band
96759599516SKenneth E. Jansenc
96859599516SKenneth E. Jansenc.... normalize the Krylov vector
96959599516SKenneth E. Jansenc
970513954efSKenneth E. Jansen        uBrgt(:,iKss+1) = uBrgt(:,iKss+1) / unorm  ! normalize the next Krylov
97159599516SKenneth E. Jansenc vector
97259599516SKenneth E. Jansenc
97359599516SKenneth E. Jansenc.... construct and reduce the Hessenberg Matrix
97459599516SKenneth E. Jansenc  since there is only one subdiagonal we can use a Givens rotation to
97559599516SKenneth E. Jansenc  rotate off each subdiagonal AS IT IS FORMED.   We do this because it
97659599516SKenneth E. Jansenc  allows us to check progress of solution and quit when satisfied.  Note
97759599516SKenneth E. Jansenc  that all future K vects will put a subdiagonal in the next column so
97859599516SKenneth E. Jansenc  there is no penalty to work ahead as  the rotation for the next vector
97959599516SKenneth E. Jansenc  will be unaffected by this rotation.
98059599516SKenneth E. Jansen
98159599516SKenneth E. Jansenc
98259599516SKenneth E. Jansenc     H Y = E ========>   R_i H Y = R_i E
98359599516SKenneth E. Jansenc
984513954efSKenneth E. Jansen           do jK = 1, iKss-1
985513954efSKenneth E. Jansen              tmp            =  Rcos(jK) * HBrg(jK,  iKss) +
986513954efSKenneth E. Jansen     &                          Rsin(jK) * HBrg(jK+1,iKss)
987513954efSKenneth E. Jansen              HBrg(jK+1,iKss) = -Rsin(jK) * HBrg(jK,  iKss) +
988513954efSKenneth E. Jansen     &                          Rcos(jK) * HBrg(jK+1,iKss)
989513954efSKenneth E. Jansen              HBrg(jK,  iKss) =  tmp
99059599516SKenneth E. Jansen           enddo
99159599516SKenneth E. Jansenc
992513954efSKenneth E. Jansen           tmp        = sqrt(HBrg(iKss,iKss)**2 + HBrg(iKss+1,iKss)**2)
993513954efSKenneth E. Jansen           Rcos(iKss) = HBrg(iKss,  iKss) / tmp
994513954efSKenneth E. Jansen           Rsin(iKss) = HBrg(iKss+1,iKss) / tmp
995513954efSKenneth E. Jansen           HBrg(iKss,  iKss) = tmp
996513954efSKenneth E. Jansen           HBrg(iKss+1,iKss) = zero
99759599516SKenneth E. Jansenc
99859599516SKenneth E. Jansenc.... rotate eBrg    R_i E
99959599516SKenneth E. Jansenc
1000513954efSKenneth E. Jansen           tmp         = Rcos(iKss)*eBrg(iKss) + Rsin(iKss)*eBrg(iKss+1)
1001513954efSKenneth E. Jansen           eBrg(iKss+1)=-Rsin(iKss)*eBrg(iKss) + Rcos(iKss)*eBrg(iKss+1)
1002513954efSKenneth E. Jansen           eBrg(iKss)  = tmp
100359599516SKenneth E. Jansenc
100459599516SKenneth E. Jansenc.... check for convergence
100559599516SKenneth E. Jansenc
1006513954efSKenneth E. Jansen           ercheck=eBrg(iKss+1)
1007513954efSKenneth E. Jansen           ntotGMs = ntotGMs + 1
1008513954efSKenneth E. Jansen           if (abs(eBrg(iKss+1)) .le. epsnrm) exit
100959599516SKenneth E. Jansenc
101059599516SKenneth E. Jansenc.... end of GMRES iteration loop
101159599516SKenneth E. Jansenc
101259599516SKenneth E. Jansen 1000   continue
101359599516SKenneth E. Jansenc
101459599516SKenneth E. Jansenc.... ------------------------->   Solution   <------------------------
101559599516SKenneth E. Jansenc
101659599516SKenneth E. Jansenc.... if converged or end of Krylov space
101759599516SKenneth E. Jansenc
101859599516SKenneth E. Jansenc.... solve for yBrg
101959599516SKenneth E. Jansenc
1020513954efSKenneth E. Jansen        do jK = iKss, 1, -1
102159599516SKenneth E. Jansen           yBrg(jK) = eBrg(jK) / HBrg(jK,jK)
102259599516SKenneth E. Jansen           do lK = 1, jK-1
102359599516SKenneth E. Jansen              eBrg(lK) = eBrg(lK) - yBrg(jK) * HBrg(lK,jK)
102459599516SKenneth E. Jansen           enddo
102559599516SKenneth E. Jansen        enddo
102659599516SKenneth E. Jansenc
102759599516SKenneth E. Jansenc.... update Dy
102859599516SKenneth E. Jansenc
1029513954efSKenneth E. Jansen        do jK = 1, iKss
103059599516SKenneth E. Jansen           Dyt = Dyt + yBrg(jK) * uBrgt(:,jK)
103159599516SKenneth E. Jansen        enddo
103259599516SKenneth E. Jansenc
103359599516SKenneth E. Jansenc.... flop count
103459599516SKenneth E. Jansenc
1035f4d0b58bSKenneth E. Jansen   !      flops = flops + (3*iKss+1)*ndof*numnp
103659599516SKenneth E. Jansenc
103759599516SKenneth E. Jansenc.... check for convergence
103859599516SKenneth E. Jansenc
1039513954efSKenneth E. Jansen        if (abs(eBrg(iKss+1)) .le. epsnrm) exit
104059599516SKenneth E. Jansenc
104159599516SKenneth E. Jansenc.... end of mGMRES loop
104259599516SKenneth E. Jansenc
104359599516SKenneth E. Jansen 2000 continue
104459599516SKenneth E. Jansenc
104559599516SKenneth E. Jansenc.... ------------------------>   Converged   <------------------------
104659599516SKenneth E. Jansenc
104759599516SKenneth E. Jansenc.... if converged
104859599516SKenneth E. Jansenc
104959599516SKenneth E. Jansen 3000 continue
105059599516SKenneth E. Jansenc
105159599516SKenneth E. Jansenc.... back precondition the result
105259599516SKenneth E. Jansenc
105359599516SKenneth E. Jansen      Dyt(:) = Dyt(:) * Dinv(:)
105459599516SKenneth E. Jansenc
105559599516SKenneth E. Jansenc.... output the statistics
105659599516SKenneth E. Jansenc
1057513954efSKenneth E. Jansen      call rstatSclr(rest, ilwork)
105859599516SKenneth E. Jansenc.... stop the timer
105959599516SKenneth E. Jansenc
106059599516SKenneth E. Jansen 3002 continue                  ! no solve just res.
106159599516SKenneth E. Jansen      call timer ('Back    ')
106259599516SKenneth E. Jansenc
106359599516SKenneth E. Jansenc.... end
106459599516SKenneth E. Jansenc
106559599516SKenneth E. Jansen      return
106659599516SKenneth E. Jansen      end
106759599516SKenneth E. Jansen
106859599516SKenneth E. Jansen
1069