xref: /phasta/phSolver/compressible/i3pre.f (revision 712d3df0b59ebebaaeaea358162c8d2c043c6e08)
159599516SKenneth E. Jansen      subroutine i3pre (BDtmp,  Binv,  EGmass,  ilwork )
259599516SKenneth E. Jansenc
359599516SKenneth E. Jansenc---------------------------------------------------------------------
459599516SKenneth E. Jansenc This is the initialization routine for the EBE-GMRES solver.
559599516SKenneth E. Jansenc It pre-preconditions the LHS mass matrix and sets up the
659599516SKenneth E. Jansenc EBE preconditioners. (pre-preconditioning is block diagonal
759599516SKenneth E. Jansenc scaling).
859599516SKenneth E. Jansenc
959599516SKenneth E. Jansenc input:
1059599516SKenneth E. Jansenc     BDtmp  (nshg,nflow,nflow)  : block diagonal scaling matrix
1159599516SKenneth E. Jansenc                                  which is already LU factored.
1259599516SKenneth E. Jansenc     EGmass (numel,nedof,nedof) : element mass matrix
1359599516SKenneth E. Jansenc
1459599516SKenneth E. Jansenc output:
1559599516SKenneth E. Jansenc     EGmass (numel,nedof,nedof) : pre-preconditioned (scaled) mass matrix
1659599516SKenneth E. Jansenc     Binv   (numel,nedof,nedof) : EBE preconditioner for each element
1759599516SKenneth E. Jansenc
1859599516SKenneth E. Jansenc---------------------------------------------------------------------
1959599516SKenneth E. Jansenc
2059599516SKenneth E. Jansen      use pointer_data
2159599516SKenneth E. Jansen
2259599516SKenneth E. Jansen      include "common.h"
2359599516SKenneth E. Jansen
2459599516SKenneth E. Jansen      dimension  BDtmp(nshg,nflow,nflow),
2559599516SKenneth E. Jansen     &           EGmass(numel,nedof,nedof),
2659599516SKenneth E. Jansenctoomuchmemory     &           Binv(numel,nedof,nedof),
2759599516SKenneth E. Jansen     &           ilwork(nlwork)
2859599516SKenneth E. Jansenc
2959599516SKenneth E. Jansen      dimension  BDiagl(numel,nshape,nflow,nflow),
3059599516SKenneth E. Jansen     &           BDiag(nshg,nflow,nflow)
3159599516SKenneth E. Jansen
3259599516SKenneth E. Jansen      BDiag = BDtmp
3359599516SKenneth E. Jansen
3459599516SKenneth E. Jansen      if (numpe > 1) then
3559599516SKenneth E. Jansen        call commu (BDiag  , ilwork, nflow*nflow  , 'out')
3659599516SKenneth E. Jansen      endif
3759599516SKenneth E. Jansenc
3859599516SKenneth E. Jansenc.... block-diagonal pre-precondition LHS
3959599516SKenneth E. Jansenc
4059599516SKenneth E. Jansenc
4159599516SKenneth E. Jansenc.... loop over element blocks
4259599516SKenneth E. Jansenc
4359599516SKenneth E. Jansen      do iblk = 1, nelblk
4459599516SKenneth E. Jansen         iel    = lcblk(1,iblk)
4559599516SKenneth E. Jansen         nenl   = lcblk(5,iblk)
4659599516SKenneth E. Jansen         npro   = lcblk(1,iblk+1) - iel
4759599516SKenneth E. Jansen         n      = iel - 1
4859599516SKenneth E. Jansen         inum   = iel + npro - 1
4959599516SKenneth E. Jansen         nshl   = lcblk(10,iblk)
5059599516SKenneth E. Jansenc
5159599516SKenneth E. Jansenc.... localize block diagonal scaling matrices
5259599516SKenneth E. Jansenc
5359599516SKenneth E. Jansen         call local (BDiag,  BDiagl(iel:inum,:,:,:), abs(mien(iblk)%p),
5459599516SKenneth E. Jansen     &               nflow*nflow,  'gather  ' )
5559599516SKenneth E. Jansenc
5659599516SKenneth E. Jansenc.... loop through local nodes and reduce all columns and rows
5759599516SKenneth E. Jansenc
5859599516SKenneth E. Jansen         do inode = 1, nshl
5959599516SKenneth E. Jansen            i = (inode - 1) * nflow ! EGmass dof offset
6059599516SKenneth E. Jansenc
6159599516SKenneth E. Jansenc.... reduce by columns, (left block diagonal preconditioning)
6259599516SKenneth E. Jansenc
6359599516SKenneth E. Jansenc     EGmass <-- inverse(L^tilde) EGmass
6459599516SKenneth E. Jansenc
6559599516SKenneth E. Jansen            do j = 1, nedof
6659599516SKenneth E. Jansen               do iv = 1, npro
6759599516SKenneth E. Jansen                  EGmass(n+iv,i+1,j) = EGmass(n+iv,i+1,j)
6859599516SKenneth E. Jansenc
6959599516SKenneth E. Jansen                  EGmass(n+iv,i+2,j) = (EGmass(n+iv,i+2,j)
7059599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,2,1) * EGmass(n+iv,i+1,j))
7159599516SKenneth E. Jansenc
7259599516SKenneth E. Jansen                  EGmass(n+iv,i+3,j) = (EGmass(n+iv,i+3,j)
7359599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,3,1) * EGmass(n+iv,i+1,j)
7459599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,3,2) * EGmass(n+iv,i+2,j))
7559599516SKenneth E. Jansenc
7659599516SKenneth E. Jansen                  EGmass(n+iv,i+4,j) = (EGmass(n+iv,i+4,j)
7759599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,4,1) * EGmass(n+iv,i+1,j)
7859599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,4,2) * EGmass(n+iv,i+2,j)
7959599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,4,3) * EGmass(n+iv,i+3,j))
8059599516SKenneth E. Jansenc
8159599516SKenneth E. Jansen                  EGmass(n+iv,i+5,j) = (EGmass(n+iv,i+5,j)
8259599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,5,1) * EGmass(n+iv,i+1,j)
8359599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,5,2) * EGmass(n+iv,i+2,j)
8459599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,5,3) * EGmass(n+iv,i+3,j)
8559599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,5,4) * EGmass(n+iv,i+4,j))
8659599516SKenneth E. Jansen               enddo
8759599516SKenneth E. Jansen            enddo
8859599516SKenneth E. Jansen         enddo
8959599516SKenneth E. Jansen
9059599516SKenneth E. Jansen         do inode = 1, nshl
9159599516SKenneth E. Jansen            i = (inode - 1) * nflow ! EGmass dof offset
9259599516SKenneth E. Jansen
9359599516SKenneth E. Jansenc
9459599516SKenneth E. Jansenc.... reduce by rows, (right block diagonal preconditioning)
9559599516SKenneth E. Jansenc
9659599516SKenneth E. Jansenc     EGmass <-- EGmass inverse(U^tilde)
9759599516SKenneth E. Jansenc
9859599516SKenneth E. Jansen            do j = 1, nedof
9959599516SKenneth E. Jansen               do iv = 1, npro
10059599516SKenneth E. Jansen                  EGmass(n+iv,j,i+1) = BDiagl(n+iv,inode,1,1) *
10159599516SKenneth E. Jansen     &                 (EGmass(n+iv,j,i+1))
10259599516SKenneth E. Jansenc
10359599516SKenneth E. Jansen                  EGmass(n+iv,j,i+2) = BDiagl(n+iv,inode,2,2) * (
10459599516SKenneth E. Jansen     &                 EGmass(n+iv,j,i+2)
10559599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,1,2) * EGmass(n+iv,j,i+1))
10659599516SKenneth E. Jansenc
10759599516SKenneth E. Jansen                  EGmass(n+iv,j,i+3) = BDiagl(n+iv,inode,3,3) * (
10859599516SKenneth E. Jansen     &                 EGmass(n+iv,j,i+3)
10959599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,1,3) * EGmass(n+iv,j,i+1)
11059599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,2,3) * EGmass(n+iv,j,i+2))
11159599516SKenneth E. Jansenc
11259599516SKenneth E. Jansen                  EGmass(n+iv,j,i+4) = BDiagl(n+iv,inode,4,4) * (
11359599516SKenneth E. Jansen     &                 EGmass(n+iv,j,i+4)
11459599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,1,4) * EGmass(n+iv,j,i+1)
11559599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,2,4) * EGmass(n+iv,j,i+2)
11659599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,3,4) * EGmass(n+iv,j,i+3))
11759599516SKenneth E. Jansenc
11859599516SKenneth E. Jansen                  EGmass(n+iv,j,i+5) = BDiagl(n+iv,inode,5,5) * (
11959599516SKenneth E. Jansen     &                 EGmass(n+iv,j,i+5)
12059599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,1,5) * EGmass(n+iv,j,i+1)
12159599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,2,5) * EGmass(n+iv,j,i+2)
12259599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,3,5) * EGmass(n+iv,j,i+3)
12359599516SKenneth E. Jansen     &                 - BDiagl(n+iv,inode,4,5) * EGmass(n+iv,j,i+4))
12459599516SKenneth E. Jansen               enddo
12559599516SKenneth E. Jansen            enddo
12659599516SKenneth E. Jansenc
12759599516SKenneth E. Jansenc.... end loops over row and column nodes
12859599516SKenneth E. Jansenc
12959599516SKenneth E. Jansen         enddo
13059599516SKenneth E. Jansenc
13159599516SKenneth E. Jansenc.... end of element blocks loop
13259599516SKenneth E. Jansenc
13359599516SKenneth E. Jansen      enddo
13459599516SKenneth E. Jansenc
13559599516SKenneth E. Jansenc.... calculate non-symmetric Cholesky EBE preconditioners
13659599516SKenneth E. Jansenc
13759599516SKenneth E. Jansenctoomuchmemory      Binv = EGmass
13859599516SKenneth E. Jansenc$$$      if (iPcond .eq. 2) then
13959599516SKenneth E. Jansenc$$$         call itrPr2 (ieneg, lcblk, Binv, ubBgl, ubBgl, 'LDU_Fact')
14059599516SKenneth E. Jansenc$$$      endif
14159599516SKenneth E. Jansenc
14259599516SKenneth E. Jansenc.... end of (pre process), return
14359599516SKenneth E. Jansenc
14459599516SKenneth E. Jansen      return
14559599516SKenneth E. Jansen
14659599516SKenneth E. Jansen      end
14759599516SKenneth E. Jansenc
14859599516SKenneth E. Jansenc
14959599516SKenneth E. Jansenc
15059599516SKenneth E. Jansen      subroutine i3preSclr (Diag,  Dinv,  EGmassT,  ilwork )
15159599516SKenneth E. Jansenc
15259599516SKenneth E. Jansenc---------------------------------------------------------------------
15359599516SKenneth E. Jansenc This is the initialization routine for the EBE-GMRES solver.
15459599516SKenneth E. Jansenc It pre-preconditions the LHS mass matrix and sets up the
15559599516SKenneth E. Jansenc EBE preconditioners. (pre-preconditioning is block diagonal
15659599516SKenneth E. Jansenc scaling).
15759599516SKenneth E. Jansenc
15859599516SKenneth E. Jansenc input:
15959599516SKenneth E. Jansenc     Diag  (numnp,nflow,nflow)    : diagonal scaling matrix
16059599516SKenneth E. Jansenc     EGmass (numel,nedof,nedof) : element mass matrix
16159599516SKenneth E. Jansenc
16259599516SKenneth E. Jansenc output:
16359599516SKenneth E. Jansenc     EGmass (numel,nedof,nedof) : pre-preconditioned (scaled) mass matrix
16459599516SKenneth E. Jansenc     Dinv   (numel,nedof,nedof) : EBE preconditioner for each element
16559599516SKenneth E. Jansenc
16659599516SKenneth E. Jansenc---------------------------------------------------------------------
16759599516SKenneth E. Jansenc
16859599516SKenneth E. Jansen      use pointer_data
16959599516SKenneth E. Jansen
17059599516SKenneth E. Jansen      include "common.h"
17159599516SKenneth E. Jansen
17259599516SKenneth E. Jansen      dimension  EGmassT(numel,nshape,nshape),
17359599516SKenneth E. Jansen     &           Dinv(nshg), ilwork(nlwork)
17459599516SKenneth E. Jansenc
175*513954efSKenneth E. Jansen      dimension  Dinvl(numel,nshape), Diag(nshg)
17659599516SKenneth E. Jansenc
17759599516SKenneth E. Jansen      Dinv = one/Diag
17859599516SKenneth E. Jansenc
17959599516SKenneth E. Jansen      if (numpe > 1) then
18059599516SKenneth E. Jansen        call commu (Dinv  , ilwork, 1  , 'out')
18159599516SKenneth E. Jansen      endif
18259599516SKenneth E. Jansenc
18359599516SKenneth E. Jansenc.... loop over element blocks
18459599516SKenneth E. Jansenc
18559599516SKenneth E. Jansen      do iblk = 1, nelblk
18659599516SKenneth E. Jansen         iel    = lcblk(1,iblk)
18759599516SKenneth E. Jansen         nenl   = lcblk(5,iblk)
18859599516SKenneth E. Jansen         npro   = lcblk(1,iblk+1) - iel
18959599516SKenneth E. Jansen         n      = iel - 1
19059599516SKenneth E. Jansen         inum   = iel + npro - 1
19159599516SKenneth E. Jansen         nshl   = lcblk(10,iblk)
19259599516SKenneth E. Jansenc
19359599516SKenneth E. Jansenc.... localize diagonal scaling matrices
19459599516SKenneth E. Jansenc
19559599516SKenneth E. Jansen         call local (Dinv,  Dinvl(iel:inum,:), mien(iblk)%p,
19659599516SKenneth E. Jansen     &               1,  'gather  ' )
19759599516SKenneth E. Jansenc
19859599516SKenneth E. Jansenc.... loop through and reduce all columns
19959599516SKenneth E. Jansenc
20059599516SKenneth E. Jansen         do icol = 1, nshl
20159599516SKenneth E. Jansen            do irow = 1, nshl
20259599516SKenneth E. Jansen               do iv = 1, npro
20359599516SKenneth E. Jansen                  EGmassT(n+iv,irow,icol) = EGmassT(n+iv,irow,icol)
20459599516SKenneth E. Jansen     &                                      *Dinvl(n+iv,icol)
20559599516SKenneth E. Jansen               enddo
20659599516SKenneth E. Jansen            enddo
20759599516SKenneth E. Jansen         enddo
20859599516SKenneth E. Jansenc
20959599516SKenneth E. Jansenc.... end of element blocks loop
21059599516SKenneth E. Jansenc
21159599516SKenneth E. Jansen      enddo
21259599516SKenneth E. Jansenc
21359599516SKenneth E. Jansenc.... end of (pre process), return
21459599516SKenneth E. Jansenc
21559599516SKenneth E. Jansen      return
21659599516SKenneth E. Jansen
21759599516SKenneth E. Jansen      end
21859599516SKenneth E. Jansen
21959599516SKenneth E. Jansen
22059599516SKenneth E. Jansen
22159599516SKenneth E. Jansen
222