xref: /phasta/phSolver/incompressible/solfar.f (revision d06028c1ed09954bf598b0145fe1ac8aecabad11)
159599516SKenneth E. Jansen      subroutine SolFlow(y,          ac,         u,
259599516SKenneth E. Jansen     &                   yold,       acold,      uold,
359599516SKenneth E. Jansen     &                   x,          iBC,
4*da6af917SKenneth E. Jansen     &                   BC,         res,        iper,
559599516SKenneth E. Jansen     &                   ilwork,     shp,        shgl,
659599516SKenneth E. Jansen     &                   shpb,       shglb,      rowp,
7*da6af917SKenneth E. Jansen     &                   colm,
859599516SKenneth E. Jansen     &                   solinc,     rerr,       tcorecp,
9bd36043dSBen Matthews     &                   GradV,       sumtime
10bd36043dSBen Matthews#ifdef HAVE_SVLS
11bd36043dSBen Matthews     &                   ,svLS_lhs,  svLS_ls,   svLS_nFaces)
12bd36043dSBen Matthews#else
13bd36043dSBen Matthews     &                   )
14bd36043dSBen Matthews#endif
1559599516SKenneth E. Jansenc
1659599516SKenneth E. Jansenc----------------------------------------------------------------------
1759599516SKenneth E. Jansenc
18ae8d68e4SKenneth E. Jansenc This is the 2nd interface routine to the  linear equation
1959599516SKenneth E. Jansenc solver library that uses the CGP and GMRES methods.
2059599516SKenneth E. Jansenc
2159599516SKenneth E. Jansenc input:
2259599516SKenneth E. Jansenc  y      (nshg,ndof)           : Y-variables at n+alpha_f
2359599516SKenneth E. Jansenc  ac     (nshg,ndof)           : Primvar. accel. variable n+alpha_m
2459599516SKenneth E. Jansenc  yold   (nshg,ndof)           : Y-variables at beginning of step
2559599516SKenneth E. Jansenc  acold   (nshg,ndof)          : Primvar. accel. at beginning of step
2659599516SKenneth E. Jansenc  x      (numnp,nsd)            : node coordinates
2759599516SKenneth E. Jansenc  iBC    (nshg)                : BC codes
2859599516SKenneth E. Jansenc  BC     (nshg,ndofBC)         : BC constraint parameters
2959599516SKenneth E. Jansenc  iper   (nshg)                : periodic nodal information
3059599516SKenneth E. Jansenc
3159599516SKenneth E. Jansenc output:
3259599516SKenneth E. Jansenc  res    (nshg,nflow)           : preconditioned residual
3359599516SKenneth E. Jansenc  y      (nshg,ndof)           : Y-variables at n+alpha_f
3459599516SKenneth E. Jansenc  ac     (nshg,ndof)           : Primvar. accel. variable n+alpha_m
3559599516SKenneth E. Jansenc
3659599516SKenneth E. Jansenc
3759599516SKenneth E. Jansenc The followings are preliminary steps required to use Farzin's
3859599516SKenneth E. Jansenc solver library.  New way of writing has to be used such as
3959599516SKenneth E. Jansenc
4059599516SKenneth E. Jansenc          |  K     G | | du |    | Rmom  |
4159599516SKenneth E. Jansenc          |          | |    | =  |       |
4259599516SKenneth E. Jansenc          | G^t    C | | dp |    | Rcon  |
4359599516SKenneth E. Jansenc
4459599516SKenneth E. Jansenc          |     E    | | dT | =  | Rtemp |
4559599516SKenneth E. Jansenc
4659599516SKenneth E. Jansenc     where
4759599516SKenneth E. Jansenc
4859599516SKenneth E. Jansenc      xKebe : K_ab = dRmom_a/du_b    xTe : E_ab = dRtemp_a/dT_b
4959599516SKenneth E. Jansenc
5059599516SKenneth E. Jansenc              G_ab = dRmom_a/dp_b
5159599516SKenneth E. Jansenc      xGoC  :
5259599516SKenneth E. Jansenc              C_ab = dRcon_a/dp_b
5359599516SKenneth E. Jansenc
5459599516SKenneth E. Jansenc              resf = Rmon Rcon       rest = Rtemp
5559599516SKenneth E. Jansenc
5659599516SKenneth E. Jansenc
5759599516SKenneth E. Jansenc Zdenek Johan,  Winter 1991.  (Fortran 90)
5859599516SKenneth E. Jansenc Juin Kim, Summer 1998. (Incompressible flow solver interface)
5959599516SKenneth E. Jansenc Alberto Figueroa.  CMM-FSI
6059599516SKenneth E. Jansenc----------------------------------------------------------------------
6159599516SKenneth E. Jansenc
6259599516SKenneth E. Jansen      use pointer_data
63*da6af917SKenneth E. Jansen      use solvedata
6459599516SKenneth E. Jansen#ifdef AMG
6559599516SKenneth E. Jansen      use ramg_data
6659599516SKenneth E. Jansen#endif
6759599516SKenneth E. Jansen
6859599516SKenneth E. Jansen      include "common.h"
6959599516SKenneth E. Jansen      include "mpif.h"
7059599516SKenneth E. Jansen      include "auxmpi.h"
71bd36043dSBen Matthews#ifdef HAVE_SVLS
72ae8d68e4SKenneth E. Jansen        include "svLS.h"
73bd36043dSBen Matthews#endif
74ae8d68e4SKenneth E. JansenC
75ae8d68e4SKenneth E. Jansen      REAL*8                rdtmp
76ae8d68e4SKenneth E. JansenC
77bd36043dSBen Matthews#ifdef HAVE_SVLS
78efb88323SKenneth E. Jansen      TYPE(svLS_lhsType), INTENT(INOUT) :: svLS_lhs
79efb88323SKenneth E. Jansen      TYPE(svLS_lsType), INTENT(INOUT) ::  svLS_ls
80bd36043dSBen Matthews#endif
81ae8d68e4SKenneth E. Jansen
8259599516SKenneth E. Jansen      real*8    y(nshg,ndof),             ac(nshg,ndof),
8359599516SKenneth E. Jansen     &          yold(nshg,ndof),          acold(nshg,ndof),
8459599516SKenneth E. Jansen     &          u(nshg,nsd),              uold(nshg,nsd),
8559599516SKenneth E. Jansen     &          x(numnp,nsd),             BC(nshg,ndofBC),
8659599516SKenneth E. Jansen     &          res(nshg,nflow),          tmpres(nshg,nflow),
8759599516SKenneth E. Jansen     &          flowDiag(nshg,4),
8859599516SKenneth E. Jansen     &          sclrDiag(nshg,1),
8959599516SKenneth E. Jansen     &          GradV(nshg,nsdsq)
9059599516SKenneth E. Jansenc
9159599516SKenneth E. Jansen      real*8    shp(MAXTOP,maxsh,MAXQPT),
9259599516SKenneth E. Jansen     &          shgl(MAXTOP,nsd,maxsh,MAXQPT),
9359599516SKenneth E. Jansen     &          shpb(MAXTOP,maxsh,MAXQPT),
9459599516SKenneth E. Jansen     &          shglb(MAXTOP,nsd,maxsh,MAXQPT)
9559599516SKenneth E. Jansenc
9659599516SKenneth E. Jansen      integer   usr(100),                 eqnType,temp,
9759599516SKenneth E. Jansen     &          rowp(nshg*nnz),           colm(nshg+1),
9859599516SKenneth E. Jansen     &          iBC(nshg),                ilwork(nlwork),
9959599516SKenneth E. Jansen     &          iper(nshg)
10059599516SKenneth E. Jansenc
10159599516SKenneth E. Jansen      real*8    yAlpha(nshg,ndof),        acAlpha(nshg,ndof),
10259599516SKenneth E. Jansen     &          uAlpha(nshg,nsd),
10359599516SKenneth E. Jansen     &          lesP(nshg,4),             lesQ(nshg,4),
10459599516SKenneth E. Jansen     &          solinc(nshg,ndof),        CGsol(nshg)
10559599516SKenneth E. Jansen
10659599516SKenneth E. Jansen      real*8     tcorecp(2)
10759599516SKenneth E. Jansen
10859599516SKenneth E. Jansen      real*8    rerr(nshg,10),            rtmp(nshg,4),rtemp
10959599516SKenneth E. Jansen
11059599516SKenneth E. Jansen      real*8    msum(4),mval(4),cpusec(10)
111ae8d68e4SKenneth E. Jansen      REAL*8 sumtime
112bd36043dSBen Matthews#ifdef HAVE_SVLS
113bd36043dSBen Matthews      INTEGER svLS_nFaces
114bd36043dSBen Matthews#endif
115bd36043dSBen Matthews      INTEGER dof, i, j, k, l
116ae8d68e4SKenneth E. Jansen      INTEGER, ALLOCATABLE :: incL(:)
117ae8d68e4SKenneth E. Jansen      REAL*8, ALLOCATABLE :: faceRes(:), Res4(:,:), Val4(:,:)
11859599516SKenneth E. Jansen
11959599516SKenneth E. Jansenc
12059599516SKenneth E. Jansenc.... *******************>> Element Data Formation <<******************
12159599516SKenneth E. Jansenc
12259599516SKenneth E. Jansenc
12359599516SKenneth E. Jansenc.... set the parameters for flux and surface tension calculations
12459599516SKenneth E. Jansenc
12559599516SKenneth E. Jansenc
12659599516SKenneth E. Jansen        temp = npro
12759599516SKenneth E. Jansen
12859599516SKenneth E. Jansen
12959599516SKenneth E. Jansen        idflx = 0
13059599516SKenneth E. Jansen        if(idiff >= 1 )  idflx= (nflow-1) * nsd
13159599516SKenneth E. Jansen        if (isurf == 1) idflx=nflow*nsd
13259599516SKenneth E. Jansenc.... compute solution at n+alpha
13359599516SKenneth E. Jansenc
13459599516SKenneth E. Jansen      call itrYAlpha( uold,    yold,    acold,
13559599516SKenneth E. Jansen     &                u,       y,       ac,
13659599516SKenneth E. Jansen     &                uAlpha,  yAlpha,  acAlpha)
13759599516SKenneth E. Jansen
13859599516SKenneth E. Jansenc
13959599516SKenneth E. Jansenc.... form the LHS matrices, the residual vector (at alpha)
14059599516SKenneth E. Jansenc
14159599516SKenneth E. Jansenc      call summary_start()
14259599516SKenneth E. Jansen      impistat=1
14359599516SKenneth E. Jansen      impistat2=1
14459599516SKenneth E. Jansen      telmcp1 = TMRC()
14559599516SKenneth E. Jansen      call ElmGMR (uAlpha,    yAlpha,     acAlpha,    x,
14659599516SKenneth E. Jansen     &             shp,       shgl,       iBC,
14759599516SKenneth E. Jansen     &             BC,        shpb,       shglb,
14859599516SKenneth E. Jansen     &             res,       iper,       ilwork,
14959599516SKenneth E. Jansen     &             rowp,      colm,       lhsK,
15059599516SKenneth E. Jansen     &             lhsP,      rerr,       GradV   )
15159599516SKenneth E. Jansen      telmcp2 = TMRC()
15259599516SKenneth E. Jansen      impistat=0
15359599516SKenneth E. Jansen      impistat2=0
15459599516SKenneth E. Jansenc      call summary_stop()
15559599516SKenneth E. Jansen
15659599516SKenneth E. Jansen
15759599516SKenneth E. Jansen            tmpres(:,:) = res(:,:)
15859599516SKenneth E. Jansen            iblk = 1
159bd36043dSBen Matthews#ifdef HAVE_SVLS
160ae8d68e4SKenneth E. Jansen      IF (svLSFlag .EQ. 1) THEN
16159599516SKenneth E. Jansen
162ae8d68e4SKenneth E. Jansenc####################################################################
163ae8d68e4SKenneth E. Jansen!     Here calling svLS
164ae8d68e4SKenneth E. Jansen
165ae8d68e4SKenneth E. Jansen      ALLOCATE(faceRes(svLS_nFaces), incL(svLS_nFaces))
166f4e2c78fSKenneth E. Jansen      faceRes=zero  ! function to compute this left out at this time but would be needed to enable adnvanced p vs. Q BC's
167ae8d68e4SKenneth E. Jansen      incL = 1
168ae8d68e4SKenneth E. Jansen      dof = 4
169ae8d68e4SKenneth E. Jansen      IF (.NOT.ALLOCATED(Res4)) THEN
170ae8d68e4SKenneth E. Jansen         ALLOCATE (Res4(dof,nshg), Val4(dof*dof,nnz_tot))
171ae8d68e4SKenneth E. Jansen      END IF
172ae8d68e4SKenneth E. Jansen
173ae8d68e4SKenneth E. Jansen      DO i=1, nshg
174ae8d68e4SKenneth E. Jansen         Res4(1:dof,i) = res(i,1:dof)
175ae8d68e4SKenneth E. Jansen      END DO
176ae8d68e4SKenneth E. Jansen
177ae8d68e4SKenneth E. Jansen      DO i=1, nnz_tot
178ae8d68e4SKenneth E. Jansen         Val4(1:3,i)   = lhsK(1:3,i)
179ae8d68e4SKenneth E. Jansen         Val4(5:7,i)   = lhsK(4:6,i)
180ae8d68e4SKenneth E. Jansen         Val4(9:11,i)  = lhsK(7:9,i)
181ae8d68e4SKenneth E. Jansen         Val4(13:15,i) = lhsP(1:3,i)
182ae8d68e4SKenneth E. Jansen         Val4(16,i)    = lhsP(4,i)
183ae8d68e4SKenneth E. Jansen      END DO
184ae8d68e4SKenneth E. Jansen
185ae8d68e4SKenneth E. Jansen      !Val4(4:12:4,:) = -lhsP(1:3,:)^t
186ae8d68e4SKenneth E. Jansen      DO i=1, nshg
187ae8d68e4SKenneth E. Jansen         Do j=colm(i), colm(i+1) - 1
188ae8d68e4SKenneth E. Jansen            k = rowp(j)
189ae8d68e4SKenneth E. Jansen            DO l=colm(k), colm(k+1) - 1
190ae8d68e4SKenneth E. Jansen               IF (rowp(l) .EQ. i) THEN
191ae8d68e4SKenneth E. Jansen                  Val4(4:12:4,l) = -lhsP(1:3,j)
192ae8d68e4SKenneth E. Jansen                  EXIT
193ae8d68e4SKenneth E. Jansen               END IF
194ae8d68e4SKenneth E. Jansen            END DO
195ae8d68e4SKenneth E. Jansen         END DO
196ae8d68e4SKenneth E. Jansen      END DO
197ae8d68e4SKenneth E. Jansen      CALL svLS_SOLVE(svLS_lhs, svLS_ls, dof, Res4, Val4, incL,
198ae8d68e4SKenneth E. Jansen     2   faceRes)
199ae8d68e4SKenneth E. Jansen
200ec121c45SKenneth E. Jansen      if(myrank.eq.master) write(*,*) 'svLS outer iterations', svLS_ls%RI%itr
201efb88323SKenneth E. Jansen      statsflow(1)=1.0*svLS_ls%GM%itr
202efb88323SKenneth E. Jansen      statsflow(4)=1.0*svLS_ls%CG%itr
203ae8d68e4SKenneth E. Jansen      DO i=1, nshg
204ae8d68e4SKenneth E. Jansen         solinc(i,1:dof) = Res4(1:dof,i)
205ae8d68e4SKenneth E. Jansen      END DO
206bd36043dSBen Matthews      ENDIF
207bd36043dSBen Matthews#endif
20879f1763eSKenneth E. Jansen
209956f980fSKenneth E. Jansen#ifdef HAVE_LESLIB
21079f1763eSKenneth E. Jansen      if(leslib.eq.1) then
211ae8d68e4SKenneth E. Jansenc
21259599516SKenneth E. Jansenc.... lesSolve : main matrix solver
21359599516SKenneth E. Jansenc
21459599516SKenneth E. Jansen      lesId   = numeqns(1)
21559599516SKenneth E. Jansen      eqnType = 1
21659599516SKenneth E. Jansen
21759599516SKenneth E. Jansenc      call summary_start()
21859599516SKenneth E. Jansen      impistat=1
21959599516SKenneth E. Jansen      impistat2=1
22059599516SKenneth E. Jansen      tlescp1 = TMRC()
22159599516SKenneth E. Jansen#ifdef AMG
22259599516SKenneth E. Jansen      ! Initial Time Profiling
22359599516SKenneth E. Jansen      call cpu_time(cpusec(1))
22459599516SKenneth E. Jansen      if (irun_amg_prec.gt.0) then
22559599516SKenneth E. Jansen          call ramg_control(colm,rowp,lhsK,lhsP,
22659599516SKenneth E. Jansen     &         ilwork,BC,iBC,iper)
22759599516SKenneth E. Jansen      end if
22859599516SKenneth E. Jansen
22959599516SKenneth E. Jansen      call cpu_time(cpusec(6))
23059599516SKenneth E. Jansen      if (irun_amg_prec.gt.0) then
23159599516SKenneth E. Jansen      ramg_flag = 1
23259599516SKenneth E. Jansen      if (irun_amg_prec.eq.2) then ! Some setup work (mode a)
23359599516SKenneth E. Jansen        ramg_window = 1.0
23459599516SKenneth E. Jansen        ramg_redo = 0
23559599516SKenneth E. Jansen      endif
23659599516SKenneth E. Jansen      do while (ramg_flag.le.irun_amg_prec)
23759599516SKenneth E. Jansen      ! if smart solve, possible run solve twice
23859599516SKenneth E. Jansen      ! restart only if meets plateau
23959599516SKenneth E. Jansen#endif
24059599516SKenneth E. Jansen
24159599516SKenneth E. Jansenc
24259599516SKenneth E. Jansenc.... setup the linear algebra solver
24359599516SKenneth E. Jansenc
24459599516SKenneth E. Jansen      rtmp = res(:,1:4)
24559599516SKenneth E. Jansen      call usrNew ( usr,        eqnType,          aperm,
24659599516SKenneth E. Jansen     &              atemp,      rtmp,             solinc,
24759599516SKenneth E. Jansen     &              flowDiag,   sclrDiag,         lesP,
24859599516SKenneth E. Jansen     &              lesQ,       iBC,              BC,
24959599516SKenneth E. Jansen     &              iper,       ilwork,           numpe,
25059599516SKenneth E. Jansen     &              nshg,       nshl,             nPermDims,
25159599516SKenneth E. Jansen     &              nTmpDims,   rowp,             colm,
25259599516SKenneth E. Jansen     &              lhsK,       lhsP,             rdtmp,
25359599516SKenneth E. Jansen     &              nnz_tot,    CGsol )
25459599516SKenneth E. Jansenc
25559599516SKenneth E. Jansenc.... solve linear system
25659599516SKenneth E. Jansenc
25759599516SKenneth E. Jansen      call myfLesSolve ( lesId, usr )
25859599516SKenneth E. Jansen#ifdef AMG
25959599516SKenneth E. Jansen      ramg_flag = ramg_flag + 2 ! Default no second run in mode a
26059599516SKenneth E. Jansen      if (irun_amg_prec.eq.3) then
26159599516SKenneth E. Jansen          if (maxIters.gt.int(statsflow(4))) then
26259599516SKenneth E. Jansen          ramg_flag = ramg_flag + 1 ! Default no second run in mode b
26359599516SKenneth E. Jansen          endif
26459599516SKenneth E. Jansen      endif
26559599516SKenneth E. Jansen      enddo
26659599516SKenneth E. Jansen      else
26759599516SKenneth E. Jansenc
26859599516SKenneth E. Jansenc.... setup the linear algebra solver
26959599516SKenneth E. Jansenc
27059599516SKenneth E. Jansen      rtmp = res(:,1:4)
27159599516SKenneth E. Jansen      call usrNew ( usr,        eqnType,          aperm,
27259599516SKenneth E. Jansen     &              atemp,      rtmp,             solinc,
27359599516SKenneth E. Jansen     &              flowDiag,   sclrDiag,         lesP,
27459599516SKenneth E. Jansen     &              lesQ,       iBC,              BC,
27559599516SKenneth E. Jansen     &              iper,       ilwork,           numpe,
27659599516SKenneth E. Jansen     &              nshg,       nshl,             nPermDims,
27759599516SKenneth E. Jansen     &              nTmpDims,   rowp,             colm,
27859599516SKenneth E. Jansen     &              lhsK,       lhsP,             rdtmp,
27959599516SKenneth E. Jansen     &              nnz_tot,    CGsol )
28059599516SKenneth E. Jansen
28159599516SKenneth E. Jansen          call myfLesSolve( lesId, usr )
28259599516SKenneth E. Jansen      endif
28359599516SKenneth E. Jansen
28459599516SKenneth E. Jansen      call cpu_time(cpusec(3))
28559599516SKenneth E. Jansen
28659599516SKenneth E. Jansen      ramg_time(1) = ramg_time(1) + cpusec(3)-cpusec(1)
28759599516SKenneth E. Jansen      ramg_time(7) = ramg_time(7) + cpusec(6)-cpusec(1)
28859599516SKenneth E. Jansen
28959599516SKenneth E. Jansen      ! ramg_time: 1 : local total
29059599516SKenneth E. Jansen      !            4 : local VG-cycle
29159599516SKenneth E. Jansen      !            7 : local setup time
29259599516SKenneth E. Jansen      !           11 : Ap-product level 1
29359599516SKenneth E. Jansen      !           12 : Ap-product level >1
29459599516SKenneth E. Jansen      !           13 : Prolongation/Restriction
29559599516SKenneth E. Jansen      !           20 : local boundary MLS time
29659599516SKenneth E. Jansen
29759599516SKenneth E. Jansen      if (myrank.eq.master) then
29859599516SKenneth E. Jansen      write(*,*)
29959599516SKenneth E. Jansen      endif
30059599516SKenneth E. Jansen      call ramg_print_time(" == AMG == Total ACUSIM Solver:",
30159599516SKenneth E. Jansen     &                    ramg_time(1))
30259599516SKenneth E. Jansen      call ramg_print_time(" == AMG == Setup: ",ramg_time(7))
30359599516SKenneth E. Jansen      call ramg_print_time(" == AMG == Prec'd cycle: ",ramg_time(4))
30459599516SKenneth E. Jansen      call ramg_print_time(" == AMG == Ap product(level=1): ",
30559599516SKenneth E. Jansen     &                      ramg_time(11))
30659599516SKenneth E. Jansen      call ramg_print_time(" == AMG == Ap product(level>=2): ",
30759599516SKenneth E. Jansen     &                      ramg_time(12))
30859599516SKenneth E. Jansen      call ramg_print_time(" == AMG == Pro/Restr ",
30959599516SKenneth E. Jansen     &                      ramg_time(13))
31059599516SKenneth E. Jansen      call ramg_print_time(" == AMG == Boundary Ap (GS only)",
31159599516SKenneth E. Jansen     &                      ramg_time(20))
31259599516SKenneth E. Jansen      if (myrank.eq.master) then
31359599516SKenneth E. Jansen      write(*,*)
31459599516SKenneth E. Jansen      endif
31559599516SKenneth E. Jansen
31659599516SKenneth E. Jansen#endif
31759599516SKenneth E. Jansen
31859599516SKenneth E. Jansen      ! End Time profiling output
31959599516SKenneth E. Jansen
32059599516SKenneth E. Jansen      call getSol ( usr, solinc )
32159599516SKenneth E. Jansen
32259599516SKenneth E. Jansen      if (numpe > 1) then
32359599516SKenneth E. Jansen         call commu ( solinc, ilwork, nflow, 'out')
32459599516SKenneth E. Jansen      endif
32579f1763eSKenneth E. Jansen      ENDIF ! end of leslib flow solve
326bd36043dSBen Matthews#endif
32759599516SKenneth E. Jansen      tlescp2 = TMRC()
32859599516SKenneth E. Jansen      impistat=0
32959599516SKenneth E. Jansen      impistat2=0
33059599516SKenneth E. Jansenc      call summary_stop()
33159599516SKenneth E. Jansen
33259599516SKenneth E. Jansen      tcorecp(1) = tcorecp(1) + telmcp2-telmcp1 ! elem. formation
33359599516SKenneth E. Jansen      tcorecp(2) = tcorecp(2) + tlescp2-tlescp1 ! linear alg. solution
33459599516SKenneth E. Jansen      call rstatic (res, y, solinc) ! output flow stats
33559599516SKenneth E. Jansenc
33659599516SKenneth E. Jansenc.... end
33759599516SKenneth E. Jansenc
33859599516SKenneth E. Jansen      return
33959599516SKenneth E. Jansen      end
34059599516SKenneth E. Jansen
34159599516SKenneth E. Jansen      subroutine SolSclr(y,          ac,         u,
34259599516SKenneth E. Jansen     &                   yold,       acold,      uold,
34359599516SKenneth E. Jansen     &                   x,          iBC,
344*da6af917SKenneth E. Jansen     &                   BC,         iper,
34559599516SKenneth E. Jansen     &                   ilwork,     shp,        shgl,
34659599516SKenneth E. Jansen     &                   shpb,       shglb,      rowp,
347*da6af917SKenneth E. Jansen     &                   colm,       solinc,
348bd36043dSBen Matthews     &                   tcorecpscal
349bd36043dSBen Matthews#ifdef HAVE_SVLS
350bd36043dSBen Matthews     &                   ,svLS_lhs,  svLS_ls,   svLS_nFaces)
351bd36043dSBen Matthews#else
352bd36043dSBen Matthews     &                   )
353bd36043dSBen Matthews#endif
35459599516SKenneth E. Jansenc
35559599516SKenneth E. Jansenc----------------------------------------------------------------------
35659599516SKenneth E. Jansenc
35759599516SKenneth E. Jansenc This is the 2nd interface routine to the Farzin's linear equation
35859599516SKenneth E. Jansenc solver library.
35959599516SKenneth E. Jansenc
36059599516SKenneth E. Jansenc input:
36159599516SKenneth E. Jansenc  y      (nshg,ndof)           : Y-variables at n+alpha_f
36259599516SKenneth E. Jansenc  ac     (nshg,ndof)           : Primvar. accel. variable n+alpha_m
36359599516SKenneth E. Jansenc  yold   (nshg,ndof)           : Y-variables at beginning of step
36459599516SKenneth E. Jansenc  x      (numnp,nsd)            : node coordinates
36559599516SKenneth E. Jansenc  iBC    (nshg)                : BC codes
36659599516SKenneth E. Jansenc  BC     (nshg,ndofBC)         : BC constraint parameters
36759599516SKenneth E. Jansenc  iper   (nshg)                : periodic nodal information
36859599516SKenneth E. Jansenc
36959599516SKenneth E. Jansenc output:
37059599516SKenneth E. Jansenc  y      (nshg,ndof)           : Y-variables at n+alpha_f
37159599516SKenneth E. Jansenc  ac     (nshg,ndof)           : Primvar. accel. variable n+alpha_m
37259599516SKenneth E. Jansenc
37359599516SKenneth E. Jansenc
37459599516SKenneth E. Jansenc The followings are preliminary steps required to use Farzin's
37559599516SKenneth E. Jansenc solver library.  New way of writing has to be used such as
37659599516SKenneth E. Jansenc
37759599516SKenneth E. Jansenc          |     E    | | dS | =  | RScal |
37859599516SKenneth E. Jansenc
37959599516SKenneth E. Jansenc----------------------------------------------------------------------
38059599516SKenneth E. Jansenc
38159599516SKenneth E. Jansen      use pointer_data
382*da6af917SKenneth E. Jansen      use solvedata
38359599516SKenneth E. Jansen
38459599516SKenneth E. Jansen      include "common.h"
38559599516SKenneth E. Jansen      include "mpif.h"
38659599516SKenneth E. Jansen      include "auxmpi.h"
387bd36043dSBen Matthews#ifdef HAVE_SVLS
388436818eeSKenneth E. Jansen        include "svLS.h"
389bd36043dSBen Matthews#endif
39059599516SKenneth E. Jansenc
39159599516SKenneth E. Jansen      real*8    y(nshg,ndof),             ac(nshg,ndof),
39259599516SKenneth E. Jansen     &          yold(nshg,ndof),          acold(nshg,ndof),
39359599516SKenneth E. Jansen     &          u(nshg,nsd),              uold(nshg,nsd),
39459599516SKenneth E. Jansen     &          x(numnp,nsd),             BC(nshg,ndofBC),
39559599516SKenneth E. Jansen     &          res(nshg,1),
39659599516SKenneth E. Jansen     &          flowDiag(nshg,4),
397*da6af917SKenneth E. Jansen     &          sclrDiag(nshg,1)
39859599516SKenneth E. Jansenc
39959599516SKenneth E. Jansen      real*8    shp(MAXTOP,maxsh,MAXQPT),
40059599516SKenneth E. Jansen     &          shgl(MAXTOP,nsd,maxsh,MAXQPT),
40159599516SKenneth E. Jansen     &          shpb(MAXTOP,maxsh,MAXQPT),
40259599516SKenneth E. Jansen     &          shglb(MAXTOP,nsd,maxsh,MAXQPT)
40359599516SKenneth E. Jansenc
40459599516SKenneth E. Jansen      integer   usr(100),                 eqnType,
40559599516SKenneth E. Jansen     &          rowp(nshg*nnz),           colm(nshg+1),
40659599516SKenneth E. Jansen     &          iBC(nshg),                ilwork(nlwork),
40759599516SKenneth E. Jansen     &          iper(nshg)
40859599516SKenneth E. Jansenc
40959599516SKenneth E. Jansen      real*8    yAlpha(nshg,ndof),        acAlpha(nshg,ndof),
41059599516SKenneth E. Jansen     &          uAlpha(nshg,nsd),
41159599516SKenneth E. Jansen     &          lesP(nshg,1),               lesQ(nshg,1),
41259599516SKenneth E. Jansen     &          solinc(nshg,1),           CGsol(nshg),
41359599516SKenneth E. Jansen     &          tcorecpscal(2)
414bd36043dSBen Matthews#ifdef HAVE_SVLS
415436818eeSKenneth E. Jansen      TYPE(svLS_lhsType), INTENT(INOUT) :: svLS_lhs
416436818eeSKenneth E. Jansen      TYPE(svLS_lsType), INTENT(INOUT) ::  svLS_ls
417bd36043dSBen Matthews      INTEGER svLS_nFaces
418bd36043dSBen Matthews#endif
419436818eeSKenneth E. Jansen      REAL*8 sumtime
420bd36043dSBen Matthews      INTEGER dof, i, j, k, l
421436818eeSKenneth E. Jansen      INTEGER, ALLOCATABLE :: incL(:)
422436818eeSKenneth E. Jansen      REAL*8, ALLOCATABLE :: faceRes(:), Res1(:,:), Val1(:,:)
42359599516SKenneth E. Jansen
42459599516SKenneth E. Jansenc
42559599516SKenneth E. Jansenc.... *******************>> Element Data Formation <<******************
42659599516SKenneth E. Jansenc
42759599516SKenneth E. Jansenc.... compute solution at n+alpha
42859599516SKenneth E. Jansenc
42959599516SKenneth E. Jansen      call itrYAlpha( uold,    yold,    acold,
43059599516SKenneth E. Jansen     &                u,       y,       ac,
43159599516SKenneth E. Jansen     &                uAlpha,  yAlpha,  acAlpha)
43259599516SKenneth E. Jansenc
43359599516SKenneth E. Jansenc.... form the LHS matrices, the residual vector (at alpha)
43459599516SKenneth E. Jansenc
43559599516SKenneth E. Jansen      impistat=2
43659599516SKenneth E. Jansen      impistat2=1
43759599516SKenneth E. Jansen      telmcp1 = TMRC()
438*da6af917SKenneth E. Jansen      jsol=nsolt+isclr
43959599516SKenneth E. Jansen      call ElmGMRSclr(yAlpha,acAlpha,    x,
44059599516SKenneth E. Jansen     &             shp,       shgl,       iBC,
44159599516SKenneth E. Jansen     &             BC,        shpb,       shglb,
44259599516SKenneth E. Jansen     &             res,       iper,       ilwork,
44359599516SKenneth E. Jansen     &             rowp,      colm,       lhsS   )
44459599516SKenneth E. Jansen      telmcp2 = TMRC()
44559599516SKenneth E. Jansen      impistat=0
44659599516SKenneth E. Jansen      impistat2=0
447436818eeSKenneth E. Jansen      statssclr(1)=0
448bd36043dSBen Matthews#ifdef HAVE_SVLS
449436818eeSKenneth E. Jansen      IF (svLSFlag .EQ. 1) THEN
450436818eeSKenneth E. Jansen
451436818eeSKenneth E. Jansenc####################################################################
452436818eeSKenneth E. Jansen!     Here calling svLS
453436818eeSKenneth E. Jansen
454436818eeSKenneth E. Jansen      ALLOCATE(faceRes(svLS_nFaces), incL(svLS_nFaces))
455ec121c45SKenneth E. Jansen      faceRes=zero
456436818eeSKenneth E. Jansen      incL = 1
457436818eeSKenneth E. Jansen      dof = 1
458436818eeSKenneth E. Jansen      IF (.NOT.ALLOCATED(Res1)) THEN
459436818eeSKenneth E. Jansen         ALLOCATE (Res1(dof,nshg), Val1(dof*dof,nnz_tot))
460436818eeSKenneth E. Jansen      END IF
461436818eeSKenneth E. Jansen
462436818eeSKenneth E. Jansen      DO i=1, nshg
463436818eeSKenneth E. Jansen         Res1(1,i) = res(i,1)
464436818eeSKenneth E. Jansen      END DO
465436818eeSKenneth E. Jansen
466436818eeSKenneth E. Jansen      DO i=1, nnz_tot
467*da6af917SKenneth E. Jansen         Val1(1,i)    = lhsS(i,jsol) ! see above jsol indexs for scalars
468436818eeSKenneth E. Jansen      END DO
469436818eeSKenneth E. Jansen
470436818eeSKenneth E. Jansen      CALL svLS_SOLVE(svLS_lhs, svLS_ls, dof, Res1, Val1, incL,
471436818eeSKenneth E. Jansen     2   faceRes)
472436818eeSKenneth E. Jansen      statssclr(1)=1.0*svLS_ls%RI%itr
473436818eeSKenneth E. Jansen      DO i=1, nshg
474436818eeSKenneth E. Jansen         solinc(i,1) = Res1(1,i)
475436818eeSKenneth E. Jansen      END DO
476bd36043dSBen Matthews      ENDIF
477bd36043dSBen Matthews#endif
478956f980fSKenneth E. Jansen#ifdef HAVE_LESLIB
47979f1763eSKenneth E. Jansen      if(leslib.eq.1) then
48059599516SKenneth E. Jansenc
48159599516SKenneth E. Jansenc.... lesSolve : main matrix solver
48259599516SKenneth E. Jansenc
48359599516SKenneth E. Jansen      lesId   = numeqns(1+nsolt+isclr)
48459599516SKenneth E. Jansen      eqnType = 2
48559599516SKenneth E. Jansenc
48659599516SKenneth E. Jansenc.... setup the linear algebra solver
48759599516SKenneth E. Jansenc
48859599516SKenneth E. Jansen
48959599516SKenneth E. Jansen      impistat=2
49059599516SKenneth E. Jansen      impistat2=1
49159599516SKenneth E. Jansen      tlescp1 = TMRC()
492*da6af917SKenneth E. Jansen      call usrNew ( usr,        eqnType,          apermS(1,1,jsol),
49359599516SKenneth E. Jansen     &              atempS,     res,              solinc,
49459599516SKenneth E. Jansen     &              flowDiag,   sclrDiag,         lesP,
49559599516SKenneth E. Jansen     &              lesQ,       iBC,              BC,
49659599516SKenneth E. Jansen     &              iper,       ilwork,           numpe,
49759599516SKenneth E. Jansen     &              nshg,       nshl,             nPermDimsS,
49859599516SKenneth E. Jansen     &              nTmpDimsS,  rowp,             colm,
49959599516SKenneth E. Jansen     &              rlhsK,      rlhsP,            lhsS,
50059599516SKenneth E. Jansen     &              nnz_tot,    CGsol )
50159599516SKenneth E. Jansenc
50259599516SKenneth E. Jansenc.... solve linear system
50359599516SKenneth E. Jansenc
50459599516SKenneth E. Jansen      call myfLesSolve ( lesId, usr )
50559599516SKenneth E. Jansen      call getSol ( usr, solinc )
50659599516SKenneth E. Jansen
50759599516SKenneth E. Jansen      if (numpe > 1) then
50859599516SKenneth E. Jansen         call commu ( solinc, ilwork, 1, 'out')
50959599516SKenneth E. Jansen      endif
51079f1763eSKenneth E. Jansen      ENDIF ! leslib conditional
511bd36043dSBen Matthews#endif
51259599516SKenneth E. Jansen      tlescp2 = TMRC()
51359599516SKenneth E. Jansen      impistat=0
51459599516SKenneth E. Jansen      impistat2=0
51559599516SKenneth E. Jansen
51659599516SKenneth E. Jansen      tcorecpscal(1) = tcorecpscal(1) + telmcp2-telmcp1 ! elem. formation
51759599516SKenneth E. Jansen      tcorecpscal(2) = tcorecpscal(2) + tlescp2-tlescp1 ! linear alg. solution
51859599516SKenneth E. Jansen
51959599516SKenneth E. Jansen      nsolsc=5+isclr
52059599516SKenneth E. Jansen      call rstaticSclr (res, y, solinc, nsolsc) ! output scalar stats
52159599516SKenneth E. Jansenc
52259599516SKenneth E. Jansenc.... end
52359599516SKenneth E. Jansenc
52459599516SKenneth E. Jansen      return
52559599516SKenneth E. Jansen      end
526