xref: /phasta/phSolver/incompressible/solfar.f (revision 956f980fcd26e82f186af0835a4a92fafc91dc7e)
159599516SKenneth E. Jansen      subroutine SolFlow(y,          ac,         u,
259599516SKenneth E. Jansen     &                   yold,       acold,      uold,
359599516SKenneth E. Jansen     &                   x,          iBC,
459599516SKenneth E. Jansen     &                   BC,         res,
559599516SKenneth E. Jansen     &                   nPermDims,  nTmpDims,  aperm,
659599516SKenneth E. Jansen     &                   atemp,      iper,
759599516SKenneth E. Jansen     &                   ilwork,     shp,        shgl,
859599516SKenneth E. Jansen     &                   shpb,       shglb,      rowp,
959599516SKenneth E. Jansen     &                   colm,       lhsK,       lhsP,
1059599516SKenneth E. Jansen     &                   solinc,     rerr,       tcorecp,
11bd36043dSBen Matthews     &                   GradV,       sumtime
12bd36043dSBen Matthews#ifdef HAVE_SVLS
13bd36043dSBen Matthews     &                   ,svLS_lhs,  svLS_ls,   svLS_nFaces)
14bd36043dSBen Matthews#else
15bd36043dSBen Matthews     &                   )
16bd36043dSBen Matthews#endif
1759599516SKenneth E. Jansenc
1859599516SKenneth E. Jansenc----------------------------------------------------------------------
1959599516SKenneth E. Jansenc
20ae8d68e4SKenneth E. Jansenc This is the 2nd interface routine to the  linear equation
2159599516SKenneth E. Jansenc solver library that uses the CGP and GMRES methods.
2259599516SKenneth E. Jansenc
2359599516SKenneth E. Jansenc input:
2459599516SKenneth E. Jansenc  y      (nshg,ndof)           : Y-variables at n+alpha_f
2559599516SKenneth E. Jansenc  ac     (nshg,ndof)           : Primvar. accel. variable n+alpha_m
2659599516SKenneth E. Jansenc  yold   (nshg,ndof)           : Y-variables at beginning of step
2759599516SKenneth E. Jansenc  acold   (nshg,ndof)          : Primvar. accel. at beginning of step
2859599516SKenneth E. Jansenc  x      (numnp,nsd)            : node coordinates
2959599516SKenneth E. Jansenc  iBC    (nshg)                : BC codes
3059599516SKenneth E. Jansenc  BC     (nshg,ndofBC)         : BC constraint parameters
3159599516SKenneth E. Jansenc  iper   (nshg)                : periodic nodal information
3259599516SKenneth E. Jansenc
3359599516SKenneth E. Jansenc output:
3459599516SKenneth E. Jansenc  res    (nshg,nflow)           : preconditioned residual
3559599516SKenneth E. Jansenc  y      (nshg,ndof)           : Y-variables at n+alpha_f
3659599516SKenneth E. Jansenc  ac     (nshg,ndof)           : Primvar. accel. variable n+alpha_m
3759599516SKenneth E. Jansenc
3859599516SKenneth E. Jansenc
3959599516SKenneth E. Jansenc The followings are preliminary steps required to use Farzin's
4059599516SKenneth E. Jansenc solver library.  New way of writing has to be used such as
4159599516SKenneth E. Jansenc
4259599516SKenneth E. Jansenc          |  K     G | | du |    | Rmom  |
4359599516SKenneth E. Jansenc          |          | |    | =  |       |
4459599516SKenneth E. Jansenc          | G^t    C | | dp |    | Rcon  |
4559599516SKenneth E. Jansenc
4659599516SKenneth E. Jansenc          |     E    | | dT | =  | Rtemp |
4759599516SKenneth E. Jansenc
4859599516SKenneth E. Jansenc     where
4959599516SKenneth E. Jansenc
5059599516SKenneth E. Jansenc      xKebe : K_ab = dRmom_a/du_b    xTe : E_ab = dRtemp_a/dT_b
5159599516SKenneth E. Jansenc
5259599516SKenneth E. Jansenc              G_ab = dRmom_a/dp_b
5359599516SKenneth E. Jansenc      xGoC  :
5459599516SKenneth E. Jansenc              C_ab = dRcon_a/dp_b
5559599516SKenneth E. Jansenc
5659599516SKenneth E. Jansenc              resf = Rmon Rcon       rest = Rtemp
5759599516SKenneth E. Jansenc
5859599516SKenneth E. Jansenc
5959599516SKenneth E. Jansenc Zdenek Johan,  Winter 1991.  (Fortran 90)
6059599516SKenneth E. Jansenc Juin Kim, Summer 1998. (Incompressible flow solver interface)
6159599516SKenneth E. Jansenc Alberto Figueroa.  CMM-FSI
6259599516SKenneth E. Jansenc----------------------------------------------------------------------
6359599516SKenneth E. Jansenc
6459599516SKenneth E. Jansen      use pointer_data
6559599516SKenneth E. Jansen#ifdef AMG
6659599516SKenneth E. Jansen      use ramg_data
6759599516SKenneth E. Jansen#endif
6859599516SKenneth E. Jansen
6959599516SKenneth E. Jansen      include "common.h"
7059599516SKenneth E. Jansen      include "mpif.h"
7159599516SKenneth E. Jansen      include "auxmpi.h"
72bd36043dSBen Matthews#ifdef HAVE_SVLS
73ae8d68e4SKenneth E. Jansen        include "svLS.h"
74bd36043dSBen Matthews#endif
7559599516SKenneth E. Jansenc
76ae8d68e4SKenneth E. JansenC
77ae8d68e4SKenneth E. JansenC     Argument variables
78ae8d68e4SKenneth E. JansenC
79ae8d68e4SKenneth E. Jansen      INTEGER            npermdims
80ae8d68e4SKenneth E. Jansen      INTEGER             ntmpdims
81ae8d68e4SKenneth E. JansenC
82ae8d68e4SKenneth E. JansenC     Local variables
83ae8d68e4SKenneth E. JansenC
84ae8d68e4SKenneth E. Jansen      INTEGER              lesid
85ae8d68e4SKenneth E. JansenC
86ae8d68e4SKenneth E. Jansen      REAL*8                rdtmp
87ae8d68e4SKenneth E. JansenC
88bd36043dSBen Matthews#ifdef HAVE_SVLS
89efb88323SKenneth E. Jansen      TYPE(svLS_lhsType), INTENT(INOUT) :: svLS_lhs
90efb88323SKenneth E. Jansen      TYPE(svLS_lsType), INTENT(INOUT) ::  svLS_ls
91bd36043dSBen Matthews#endif
92ae8d68e4SKenneth E. Jansen
9359599516SKenneth E. Jansen      real*8    y(nshg,ndof),             ac(nshg,ndof),
9459599516SKenneth E. Jansen     &          yold(nshg,ndof),          acold(nshg,ndof),
9559599516SKenneth E. Jansen     &          u(nshg,nsd),              uold(nshg,nsd),
9659599516SKenneth E. Jansen     &          x(numnp,nsd),             BC(nshg,ndofBC),
9759599516SKenneth E. Jansen     &          res(nshg,nflow),          tmpres(nshg,nflow),
9859599516SKenneth E. Jansen     &          flowDiag(nshg,4),
9959599516SKenneth E. Jansen     &          aperm(nshg,nPermDims),    atemp(nshg,nTmpDims),
10059599516SKenneth E. Jansen     &          sclrDiag(nshg,1),
10159599516SKenneth E. Jansen     &          lhsK(9,nnz_tot),          lhsP(4,nnz_tot),
10259599516SKenneth E. Jansen     &          GradV(nshg,nsdsq)
10359599516SKenneth E. Jansenc
10459599516SKenneth E. Jansen      real*8    shp(MAXTOP,maxsh,MAXQPT),
10559599516SKenneth E. Jansen     &          shgl(MAXTOP,nsd,maxsh,MAXQPT),
10659599516SKenneth E. Jansen     &          shpb(MAXTOP,maxsh,MAXQPT),
10759599516SKenneth E. Jansen     &          shglb(MAXTOP,nsd,maxsh,MAXQPT)
10859599516SKenneth E. Jansenc
10959599516SKenneth E. Jansen      integer   usr(100),                 eqnType,temp,
11059599516SKenneth E. Jansen     &          rowp(nshg*nnz),           colm(nshg+1),
11159599516SKenneth E. Jansen     &          iBC(nshg),                ilwork(nlwork),
11259599516SKenneth E. Jansen     &          iper(nshg)
11359599516SKenneth E. Jansenc
11459599516SKenneth E. Jansen      real*8    yAlpha(nshg,ndof),        acAlpha(nshg,ndof),
11559599516SKenneth E. Jansen     &          uAlpha(nshg,nsd),
11659599516SKenneth E. Jansen     &          lesP(nshg,4),             lesQ(nshg,4),
11759599516SKenneth E. Jansen     &          solinc(nshg,ndof),        CGsol(nshg)
11859599516SKenneth E. Jansen
11959599516SKenneth E. Jansen      real*8     tcorecp(2)
12059599516SKenneth E. Jansen
12159599516SKenneth E. Jansen      real*8    rerr(nshg,10),            rtmp(nshg,4),rtemp
12259599516SKenneth E. Jansen
12359599516SKenneth E. Jansen      real*8    msum(4),mval(4),cpusec(10)
124ae8d68e4SKenneth E. Jansen      REAL*8 sumtime
125bd36043dSBen Matthews#ifdef HAVE_SVLS
126bd36043dSBen Matthews      INTEGER svLS_nFaces
127bd36043dSBen Matthews#endif
128bd36043dSBen Matthews      INTEGER dof, i, j, k, l
129ae8d68e4SKenneth E. Jansen      INTEGER, ALLOCATABLE :: incL(:)
130ae8d68e4SKenneth E. Jansen      REAL*8, ALLOCATABLE :: faceRes(:), Res4(:,:), Val4(:,:)
13159599516SKenneth E. Jansen
13259599516SKenneth E. Jansenc
13359599516SKenneth E. Jansenc.... *******************>> Element Data Formation <<******************
13459599516SKenneth E. Jansenc
13559599516SKenneth E. Jansenc
13659599516SKenneth E. Jansenc.... set the parameters for flux and surface tension calculations
13759599516SKenneth E. Jansenc
13859599516SKenneth E. Jansenc
13959599516SKenneth E. Jansen        temp = npro
14059599516SKenneth E. Jansen
14159599516SKenneth E. Jansen
14259599516SKenneth E. Jansen        idflx = 0
14359599516SKenneth E. Jansen        if(idiff >= 1 )  idflx= (nflow-1) * nsd
14459599516SKenneth E. Jansen        if (isurf == 1) idflx=nflow*nsd
14559599516SKenneth E. Jansenc.... compute solution at n+alpha
14659599516SKenneth E. Jansenc
14759599516SKenneth E. Jansen      call itrYAlpha( uold,    yold,    acold,
14859599516SKenneth E. Jansen     &                u,       y,       ac,
14959599516SKenneth E. Jansen     &                uAlpha,  yAlpha,  acAlpha)
15059599516SKenneth E. Jansen
15159599516SKenneth E. Jansenc
15259599516SKenneth E. Jansenc.... form the LHS matrices, the residual vector (at alpha)
15359599516SKenneth E. Jansenc
15459599516SKenneth E. Jansenc      call summary_start()
15559599516SKenneth E. Jansen      impistat=1
15659599516SKenneth E. Jansen      impistat2=1
15759599516SKenneth E. Jansen      telmcp1 = TMRC()
15859599516SKenneth E. Jansen      call ElmGMR (uAlpha,    yAlpha,     acAlpha,    x,
15959599516SKenneth E. Jansen     &             shp,       shgl,       iBC,
16059599516SKenneth E. Jansen     &             BC,        shpb,       shglb,
16159599516SKenneth E. Jansen     &             res,       iper,       ilwork,
16259599516SKenneth E. Jansen     &             rowp,      colm,       lhsK,
16359599516SKenneth E. Jansen     &             lhsP,      rerr,       GradV   )
16459599516SKenneth E. Jansen      telmcp2 = TMRC()
16559599516SKenneth E. Jansen      impistat=0
16659599516SKenneth E. Jansen      impistat2=0
16759599516SKenneth E. Jansenc      call summary_stop()
16859599516SKenneth E. Jansen
16959599516SKenneth E. Jansen
17059599516SKenneth E. Jansen            tmpres(:,:) = res(:,:)
17159599516SKenneth E. Jansen            iblk = 1
172bd36043dSBen Matthews#ifdef HAVE_SVLS
173ae8d68e4SKenneth E. Jansen      IF (svLSFlag .EQ. 1) THEN
17459599516SKenneth E. Jansen
175ae8d68e4SKenneth E. Jansenc####################################################################
176ae8d68e4SKenneth E. Jansen!     Here calling svLS
177ae8d68e4SKenneth E. Jansen
178ae8d68e4SKenneth E. Jansen      ALLOCATE(faceRes(svLS_nFaces), incL(svLS_nFaces))
179f4e2c78fSKenneth 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
180ae8d68e4SKenneth E. Jansen      incL = 1
181ae8d68e4SKenneth E. Jansen      dof = 4
182ae8d68e4SKenneth E. Jansen      IF (.NOT.ALLOCATED(Res4)) THEN
183ae8d68e4SKenneth E. Jansen         ALLOCATE (Res4(dof,nshg), Val4(dof*dof,nnz_tot))
184ae8d68e4SKenneth E. Jansen      END IF
185ae8d68e4SKenneth E. Jansen
186ae8d68e4SKenneth E. Jansen      DO i=1, nshg
187ae8d68e4SKenneth E. Jansen         Res4(1:dof,i) = res(i,1:dof)
188ae8d68e4SKenneth E. Jansen      END DO
189ae8d68e4SKenneth E. Jansen
190ae8d68e4SKenneth E. Jansen      DO i=1, nnz_tot
191ae8d68e4SKenneth E. Jansen         Val4(1:3,i)   = lhsK(1:3,i)
192ae8d68e4SKenneth E. Jansen         Val4(5:7,i)   = lhsK(4:6,i)
193ae8d68e4SKenneth E. Jansen         Val4(9:11,i)  = lhsK(7:9,i)
194ae8d68e4SKenneth E. Jansen         Val4(13:15,i) = lhsP(1:3,i)
195ae8d68e4SKenneth E. Jansen         Val4(16,i)    = lhsP(4,i)
196ae8d68e4SKenneth E. Jansen      END DO
197ae8d68e4SKenneth E. Jansen
198ae8d68e4SKenneth E. Jansen      !Val4(4:12:4,:) = -lhsP(1:3,:)^t
199ae8d68e4SKenneth E. Jansen      DO i=1, nshg
200ae8d68e4SKenneth E. Jansen         Do j=colm(i), colm(i+1) - 1
201ae8d68e4SKenneth E. Jansen            k = rowp(j)
202ae8d68e4SKenneth E. Jansen            DO l=colm(k), colm(k+1) - 1
203ae8d68e4SKenneth E. Jansen               IF (rowp(l) .EQ. i) THEN
204ae8d68e4SKenneth E. Jansen                  Val4(4:12:4,l) = -lhsP(1:3,j)
205ae8d68e4SKenneth E. Jansen                  EXIT
206ae8d68e4SKenneth E. Jansen               END IF
207ae8d68e4SKenneth E. Jansen            END DO
208ae8d68e4SKenneth E. Jansen         END DO
209ae8d68e4SKenneth E. Jansen      END DO
210ae8d68e4SKenneth E. Jansen      CALL svLS_SOLVE(svLS_lhs, svLS_ls, dof, Res4, Val4, incL,
211ae8d68e4SKenneth E. Jansen     2   faceRes)
212ae8d68e4SKenneth E. Jansen
213ec121c45SKenneth E. Jansen      if(myrank.eq.master) write(*,*) 'svLS outer iterations', svLS_ls%RI%itr
214efb88323SKenneth E. Jansen      statsflow(1)=1.0*svLS_ls%GM%itr
215efb88323SKenneth E. Jansen      statsflow(4)=1.0*svLS_ls%CG%itr
216ae8d68e4SKenneth E. Jansen      DO i=1, nshg
217ae8d68e4SKenneth E. Jansen         solinc(i,1:dof) = Res4(1:dof,i)
218ae8d68e4SKenneth E. Jansen      END DO
219bd36043dSBen Matthews      ENDIF
220bd36043dSBen Matthews#endif
22179f1763eSKenneth E. Jansen
222*956f980fSKenneth E. Jansen#ifdef HAVE_LESLIB
22379f1763eSKenneth E. Jansen      if(leslib.eq.1) then
224ae8d68e4SKenneth E. Jansenc
22559599516SKenneth E. Jansenc.... lesSolve : main matrix solver
22659599516SKenneth E. Jansenc
22759599516SKenneth E. Jansen      lesId   = numeqns(1)
22859599516SKenneth E. Jansen      eqnType = 1
22959599516SKenneth E. Jansen
23059599516SKenneth E. Jansenc      call summary_start()
23159599516SKenneth E. Jansen      impistat=1
23259599516SKenneth E. Jansen      impistat2=1
23359599516SKenneth E. Jansen      tlescp1 = TMRC()
23459599516SKenneth E. Jansen#ifdef AMG
23559599516SKenneth E. Jansen      ! Initial Time Profiling
23659599516SKenneth E. Jansen      call cpu_time(cpusec(1))
23759599516SKenneth E. Jansen      if (irun_amg_prec.gt.0) then
23859599516SKenneth E. Jansen          call ramg_control(colm,rowp,lhsK,lhsP,
23959599516SKenneth E. Jansen     &         ilwork,BC,iBC,iper)
24059599516SKenneth E. Jansen      end if
24159599516SKenneth E. Jansen
24259599516SKenneth E. Jansen      call cpu_time(cpusec(6))
24359599516SKenneth E. Jansen      if (irun_amg_prec.gt.0) then
24459599516SKenneth E. Jansen      ramg_flag = 1
24559599516SKenneth E. Jansen      if (irun_amg_prec.eq.2) then ! Some setup work (mode a)
24659599516SKenneth E. Jansen        ramg_window = 1.0
24759599516SKenneth E. Jansen        ramg_redo = 0
24859599516SKenneth E. Jansen      endif
24959599516SKenneth E. Jansen      do while (ramg_flag.le.irun_amg_prec)
25059599516SKenneth E. Jansen      ! if smart solve, possible run solve twice
25159599516SKenneth E. Jansen      ! restart only if meets plateau
25259599516SKenneth E. Jansen#endif
25359599516SKenneth E. Jansen
25459599516SKenneth E. Jansenc
25559599516SKenneth E. Jansenc.... setup the linear algebra solver
25659599516SKenneth E. Jansenc
25759599516SKenneth E. Jansen      rtmp = res(:,1:4)
25859599516SKenneth E. Jansen      call usrNew ( usr,        eqnType,          aperm,
25959599516SKenneth E. Jansen     &              atemp,      rtmp,             solinc,
26059599516SKenneth E. Jansen     &              flowDiag,   sclrDiag,         lesP,
26159599516SKenneth E. Jansen     &              lesQ,       iBC,              BC,
26259599516SKenneth E. Jansen     &              iper,       ilwork,           numpe,
26359599516SKenneth E. Jansen     &              nshg,       nshl,             nPermDims,
26459599516SKenneth E. Jansen     &              nTmpDims,   rowp,             colm,
26559599516SKenneth E. Jansen     &              lhsK,       lhsP,             rdtmp,
26659599516SKenneth E. Jansen     &              nnz_tot,    CGsol )
26759599516SKenneth E. Jansenc
26859599516SKenneth E. Jansenc.... solve linear system
26959599516SKenneth E. Jansenc
27059599516SKenneth E. Jansen      call myfLesSolve ( lesId, usr )
27159599516SKenneth E. Jansen#ifdef AMG
27259599516SKenneth E. Jansen      ramg_flag = ramg_flag + 2 ! Default no second run in mode a
27359599516SKenneth E. Jansen      if (irun_amg_prec.eq.3) then
27459599516SKenneth E. Jansen          if (maxIters.gt.int(statsflow(4))) then
27559599516SKenneth E. Jansen          ramg_flag = ramg_flag + 1 ! Default no second run in mode b
27659599516SKenneth E. Jansen          endif
27759599516SKenneth E. Jansen      endif
27859599516SKenneth E. Jansen      enddo
27959599516SKenneth E. Jansen      else
28059599516SKenneth E. Jansenc
28159599516SKenneth E. Jansenc.... setup the linear algebra solver
28259599516SKenneth E. Jansenc
28359599516SKenneth E. Jansen      rtmp = res(:,1:4)
28459599516SKenneth E. Jansen      call usrNew ( usr,        eqnType,          aperm,
28559599516SKenneth E. Jansen     &              atemp,      rtmp,             solinc,
28659599516SKenneth E. Jansen     &              flowDiag,   sclrDiag,         lesP,
28759599516SKenneth E. Jansen     &              lesQ,       iBC,              BC,
28859599516SKenneth E. Jansen     &              iper,       ilwork,           numpe,
28959599516SKenneth E. Jansen     &              nshg,       nshl,             nPermDims,
29059599516SKenneth E. Jansen     &              nTmpDims,   rowp,             colm,
29159599516SKenneth E. Jansen     &              lhsK,       lhsP,             rdtmp,
29259599516SKenneth E. Jansen     &              nnz_tot,    CGsol )
29359599516SKenneth E. Jansen
29459599516SKenneth E. Jansen          call myfLesSolve( lesId, usr )
29559599516SKenneth E. Jansen      endif
29659599516SKenneth E. Jansen
29759599516SKenneth E. Jansen      call cpu_time(cpusec(3))
29859599516SKenneth E. Jansen
29959599516SKenneth E. Jansen      ramg_time(1) = ramg_time(1) + cpusec(3)-cpusec(1)
30059599516SKenneth E. Jansen      ramg_time(7) = ramg_time(7) + cpusec(6)-cpusec(1)
30159599516SKenneth E. Jansen
30259599516SKenneth E. Jansen      ! ramg_time: 1 : local total
30359599516SKenneth E. Jansen      !            4 : local VG-cycle
30459599516SKenneth E. Jansen      !            7 : local setup time
30559599516SKenneth E. Jansen      !           11 : Ap-product level 1
30659599516SKenneth E. Jansen      !           12 : Ap-product level >1
30759599516SKenneth E. Jansen      !           13 : Prolongation/Restriction
30859599516SKenneth E. Jansen      !           20 : local boundary MLS time
30959599516SKenneth E. Jansen
31059599516SKenneth E. Jansen      if (myrank.eq.master) then
31159599516SKenneth E. Jansen      write(*,*)
31259599516SKenneth E. Jansen      endif
31359599516SKenneth E. Jansen      call ramg_print_time(" == AMG == Total ACUSIM Solver:",
31459599516SKenneth E. Jansen     &                    ramg_time(1))
31559599516SKenneth E. Jansen      call ramg_print_time(" == AMG == Setup: ",ramg_time(7))
31659599516SKenneth E. Jansen      call ramg_print_time(" == AMG == Prec'd cycle: ",ramg_time(4))
31759599516SKenneth E. Jansen      call ramg_print_time(" == AMG == Ap product(level=1): ",
31859599516SKenneth E. Jansen     &                      ramg_time(11))
31959599516SKenneth E. Jansen      call ramg_print_time(" == AMG == Ap product(level>=2): ",
32059599516SKenneth E. Jansen     &                      ramg_time(12))
32159599516SKenneth E. Jansen      call ramg_print_time(" == AMG == Pro/Restr ",
32259599516SKenneth E. Jansen     &                      ramg_time(13))
32359599516SKenneth E. Jansen      call ramg_print_time(" == AMG == Boundary Ap (GS only)",
32459599516SKenneth E. Jansen     &                      ramg_time(20))
32559599516SKenneth E. Jansen      if (myrank.eq.master) then
32659599516SKenneth E. Jansen      write(*,*)
32759599516SKenneth E. Jansen      endif
32859599516SKenneth E. Jansen
32959599516SKenneth E. Jansen#endif
33059599516SKenneth E. Jansen
33159599516SKenneth E. Jansen      ! End Time profiling output
33259599516SKenneth E. Jansen
33359599516SKenneth E. Jansen      call getSol ( usr, solinc )
33459599516SKenneth E. Jansen
33559599516SKenneth E. Jansen      if (numpe > 1) then
33659599516SKenneth E. Jansen         call commu ( solinc, ilwork, nflow, 'out')
33759599516SKenneth E. Jansen      endif
33879f1763eSKenneth E. Jansen      ENDIF ! end of leslib flow solve
339bd36043dSBen Matthews#endif
34059599516SKenneth E. Jansen      tlescp2 = TMRC()
34159599516SKenneth E. Jansen      impistat=0
34259599516SKenneth E. Jansen      impistat2=0
34359599516SKenneth E. Jansenc      call summary_stop()
34459599516SKenneth E. Jansen
34559599516SKenneth E. Jansen      tcorecp(1) = tcorecp(1) + telmcp2-telmcp1 ! elem. formation
34659599516SKenneth E. Jansen      tcorecp(2) = tcorecp(2) + tlescp2-tlescp1 ! linear alg. solution
34759599516SKenneth E. Jansen      call rstatic (res, y, solinc) ! output flow stats
34859599516SKenneth E. Jansenc
34959599516SKenneth E. Jansenc.... end
35059599516SKenneth E. Jansenc
35159599516SKenneth E. Jansen      return
35259599516SKenneth E. Jansen      end
35359599516SKenneth E. Jansen
35459599516SKenneth E. Jansen      subroutine SolSclr(y,          ac,         u,
35559599516SKenneth E. Jansen     &                   yold,       acold,      uold,
35659599516SKenneth E. Jansen     &                   x,          iBC,
35759599516SKenneth E. Jansen     &                   BC,         nPermDimsS,  nTmpDimsS,
35859599516SKenneth E. Jansen     &                   apermS,     atempS,     iper,
35959599516SKenneth E. Jansen     &                   ilwork,     shp,        shgl,
36059599516SKenneth E. Jansen     &                   shpb,       shglb,      rowp,
36159599516SKenneth E. Jansen     &                   colm,       lhsS,       solinc,
362bd36043dSBen Matthews     &                   tcorecpscal
363bd36043dSBen Matthews#ifdef HAVE_SVLS
364bd36043dSBen Matthews     &                   ,svLS_lhs,  svLS_ls,   svLS_nFaces)
365bd36043dSBen Matthews#else
366bd36043dSBen Matthews     &                   )
367bd36043dSBen Matthews#endif
36859599516SKenneth E. Jansenc
36959599516SKenneth E. Jansenc----------------------------------------------------------------------
37059599516SKenneth E. Jansenc
37159599516SKenneth E. Jansenc This is the 2nd interface routine to the Farzin's linear equation
37259599516SKenneth E. Jansenc solver library.
37359599516SKenneth E. Jansenc
37459599516SKenneth E. Jansenc input:
37559599516SKenneth E. Jansenc  y      (nshg,ndof)           : Y-variables at n+alpha_f
37659599516SKenneth E. Jansenc  ac     (nshg,ndof)           : Primvar. accel. variable n+alpha_m
37759599516SKenneth E. Jansenc  yold   (nshg,ndof)           : Y-variables at beginning of step
37859599516SKenneth E. Jansenc  x      (numnp,nsd)            : node coordinates
37959599516SKenneth E. Jansenc  iBC    (nshg)                : BC codes
38059599516SKenneth E. Jansenc  BC     (nshg,ndofBC)         : BC constraint parameters
38159599516SKenneth E. Jansenc  iper   (nshg)                : periodic nodal information
38259599516SKenneth E. Jansenc
38359599516SKenneth E. Jansenc output:
38459599516SKenneth E. Jansenc  y      (nshg,ndof)           : Y-variables at n+alpha_f
38559599516SKenneth E. Jansenc  ac     (nshg,ndof)           : Primvar. accel. variable n+alpha_m
38659599516SKenneth E. Jansenc
38759599516SKenneth E. Jansenc
38859599516SKenneth E. Jansenc The followings are preliminary steps required to use Farzin's
38959599516SKenneth E. Jansenc solver library.  New way of writing has to be used such as
39059599516SKenneth E. Jansenc
39159599516SKenneth E. Jansenc          |     E    | | dS | =  | RScal |
39259599516SKenneth E. Jansenc
39359599516SKenneth E. Jansenc----------------------------------------------------------------------
39459599516SKenneth E. Jansenc
39559599516SKenneth E. Jansen      use pointer_data
39659599516SKenneth E. Jansen
39759599516SKenneth E. Jansen      include "common.h"
39859599516SKenneth E. Jansen      include "mpif.h"
39959599516SKenneth E. Jansen      include "auxmpi.h"
400bd36043dSBen Matthews#ifdef HAVE_SVLS
401436818eeSKenneth E. Jansen        include "svLS.h"
402bd36043dSBen Matthews#endif
40359599516SKenneth E. Jansenc
40459599516SKenneth E. Jansen      real*8    y(nshg,ndof),             ac(nshg,ndof),
40559599516SKenneth E. Jansen     &          yold(nshg,ndof),          acold(nshg,ndof),
40659599516SKenneth E. Jansen     &          u(nshg,nsd),              uold(nshg,nsd),
40759599516SKenneth E. Jansen     &          x(numnp,nsd),             BC(nshg,ndofBC),
40859599516SKenneth E. Jansen     &          res(nshg,1),
40959599516SKenneth E. Jansen     &          flowDiag(nshg,4),
41059599516SKenneth E. Jansen     &          sclrDiag(nshg,1),           lhsS(nnz_tot),
41159599516SKenneth E. Jansen     &          apermS(nshg,nPermDimsS),  atempS(nshg,nTmpDimsS)
41259599516SKenneth E. Jansen
41359599516SKenneth E. Jansenc
41459599516SKenneth E. Jansen      real*8    shp(MAXTOP,maxsh,MAXQPT),
41559599516SKenneth E. Jansen     &          shgl(MAXTOP,nsd,maxsh,MAXQPT),
41659599516SKenneth E. Jansen     &          shpb(MAXTOP,maxsh,MAXQPT),
41759599516SKenneth E. Jansen     &          shglb(MAXTOP,nsd,maxsh,MAXQPT)
41859599516SKenneth E. Jansenc
41959599516SKenneth E. Jansen      integer   usr(100),                 eqnType,
42059599516SKenneth E. Jansen     &          rowp(nshg*nnz),           colm(nshg+1),
42159599516SKenneth E. Jansen     &          iBC(nshg),                ilwork(nlwork),
42259599516SKenneth E. Jansen     &          iper(nshg)
42359599516SKenneth E. Jansenc
42459599516SKenneth E. Jansen      real*8    yAlpha(nshg,ndof),        acAlpha(nshg,ndof),
42559599516SKenneth E. Jansen     &          uAlpha(nshg,nsd),
42659599516SKenneth E. Jansen     &          lesP(nshg,1),               lesQ(nshg,1),
42759599516SKenneth E. Jansen     &          solinc(nshg,1),           CGsol(nshg),
42859599516SKenneth E. Jansen     &          tcorecpscal(2)
429bd36043dSBen Matthews#ifdef HAVE_SVLS
430436818eeSKenneth E. Jansen      TYPE(svLS_lhsType), INTENT(INOUT) :: svLS_lhs
431436818eeSKenneth E. Jansen      TYPE(svLS_lsType), INTENT(INOUT) ::  svLS_ls
432bd36043dSBen Matthews      INTEGER svLS_nFaces
433bd36043dSBen Matthews#endif
434436818eeSKenneth E. Jansen      REAL*8 sumtime
435bd36043dSBen Matthews      INTEGER dof, i, j, k, l
436436818eeSKenneth E. Jansen      INTEGER, ALLOCATABLE :: incL(:)
437436818eeSKenneth E. Jansen      REAL*8, ALLOCATABLE :: faceRes(:), Res1(:,:), Val1(:,:)
43859599516SKenneth E. Jansen
43959599516SKenneth E. Jansenc
44059599516SKenneth E. Jansenc.... *******************>> Element Data Formation <<******************
44159599516SKenneth E. Jansenc
44259599516SKenneth E. Jansenc.... compute solution at n+alpha
44359599516SKenneth E. Jansenc
44459599516SKenneth E. Jansen      call itrYAlpha( uold,    yold,    acold,
44559599516SKenneth E. Jansen     &                u,       y,       ac,
44659599516SKenneth E. Jansen     &                uAlpha,  yAlpha,  acAlpha)
44759599516SKenneth E. Jansenc
44859599516SKenneth E. Jansenc.... form the LHS matrices, the residual vector (at alpha)
44959599516SKenneth E. Jansenc
45059599516SKenneth E. Jansen      impistat=2
45159599516SKenneth E. Jansen      impistat2=1
45259599516SKenneth E. Jansen      telmcp1 = TMRC()
45359599516SKenneth E. Jansen      call ElmGMRSclr(yAlpha,acAlpha,    x,
45459599516SKenneth E. Jansen     &             shp,       shgl,       iBC,
45559599516SKenneth E. Jansen     &             BC,        shpb,       shglb,
45659599516SKenneth E. Jansen     &             res,       iper,       ilwork,
45759599516SKenneth E. Jansen     &             rowp,      colm,       lhsS   )
45859599516SKenneth E. Jansen      telmcp2 = TMRC()
45959599516SKenneth E. Jansen      impistat=0
46059599516SKenneth E. Jansen      impistat2=0
461436818eeSKenneth E. Jansen      statssclr(1)=0
462bd36043dSBen Matthews#ifdef HAVE_SVLS
463436818eeSKenneth E. Jansen      IF (svLSFlag .EQ. 1) THEN
464436818eeSKenneth E. Jansen
465436818eeSKenneth E. Jansenc####################################################################
466436818eeSKenneth E. Jansen!     Here calling svLS
467436818eeSKenneth E. Jansen
468436818eeSKenneth E. Jansen      ALLOCATE(faceRes(svLS_nFaces), incL(svLS_nFaces))
469ec121c45SKenneth E. Jansen      faceRes=zero
470436818eeSKenneth E. Jansen      incL = 1
471436818eeSKenneth E. Jansen      dof = 1
472436818eeSKenneth E. Jansen      IF (.NOT.ALLOCATED(Res1)) THEN
473436818eeSKenneth E. Jansen         ALLOCATE (Res1(dof,nshg), Val1(dof*dof,nnz_tot))
474436818eeSKenneth E. Jansen      END IF
475436818eeSKenneth E. Jansen
476436818eeSKenneth E. Jansen      DO i=1, nshg
477436818eeSKenneth E. Jansen         Res1(1,i) = res(i,1)
478436818eeSKenneth E. Jansen      END DO
479436818eeSKenneth E. Jansen
480436818eeSKenneth E. Jansen      DO i=1, nnz_tot
481436818eeSKenneth E. Jansen         Val1(1,i)    = lhsS(i)
482436818eeSKenneth E. Jansen      END DO
483436818eeSKenneth E. Jansen
484436818eeSKenneth E. Jansen      CALL svLS_SOLVE(svLS_lhs, svLS_ls, dof, Res1, Val1, incL,
485436818eeSKenneth E. Jansen     2   faceRes)
486436818eeSKenneth E. Jansen      statssclr(1)=1.0*svLS_ls%RI%itr
487436818eeSKenneth E. Jansen      DO i=1, nshg
488436818eeSKenneth E. Jansen         solinc(i,1) = Res1(1,i)
489436818eeSKenneth E. Jansen      END DO
490bd36043dSBen Matthews      ENDIF
491bd36043dSBen Matthews#endif
492*956f980fSKenneth E. Jansen#ifdef HAVE_LESLIB
49379f1763eSKenneth E. Jansen      if(leslib.eq.1) then
49459599516SKenneth E. Jansenc
49559599516SKenneth E. Jansenc.... lesSolve : main matrix solver
49659599516SKenneth E. Jansenc
49759599516SKenneth E. Jansen      lesId   = numeqns(1+nsolt+isclr)
49859599516SKenneth E. Jansen      eqnType = 2
49959599516SKenneth E. Jansenc
50059599516SKenneth E. Jansenc.... setup the linear algebra solver
50159599516SKenneth E. Jansenc
50259599516SKenneth E. Jansen
50359599516SKenneth E. Jansen      impistat=2
50459599516SKenneth E. Jansen      impistat2=1
50559599516SKenneth E. Jansen      tlescp1 = TMRC()
50659599516SKenneth E. Jansen      call usrNew ( usr,        eqnType,          apermS,
50759599516SKenneth E. Jansen     &              atempS,     res,              solinc,
50859599516SKenneth E. Jansen     &              flowDiag,   sclrDiag,         lesP,
50959599516SKenneth E. Jansen     &              lesQ,       iBC,              BC,
51059599516SKenneth E. Jansen     &              iper,       ilwork,           numpe,
51159599516SKenneth E. Jansen     &              nshg,       nshl,             nPermDimsS,
51259599516SKenneth E. Jansen     &              nTmpDimsS,  rowp,             colm,
51359599516SKenneth E. Jansen     &              rlhsK,      rlhsP,            lhsS,
51459599516SKenneth E. Jansen     &              nnz_tot,    CGsol )
51559599516SKenneth E. Jansenc
51659599516SKenneth E. Jansenc.... solve linear system
51759599516SKenneth E. Jansenc
51859599516SKenneth E. Jansen      call myfLesSolve ( lesId, usr )
51959599516SKenneth E. Jansen      call getSol ( usr, solinc )
52059599516SKenneth E. Jansen
52159599516SKenneth E. Jansen      if (numpe > 1) then
52259599516SKenneth E. Jansen         call commu ( solinc, ilwork, 1, 'out')
52359599516SKenneth E. Jansen      endif
52479f1763eSKenneth E. Jansen      ENDIF ! leslib conditional
525bd36043dSBen Matthews#endif
52659599516SKenneth E. Jansen      tlescp2 = TMRC()
52759599516SKenneth E. Jansen      impistat=0
52859599516SKenneth E. Jansen      impistat2=0
52959599516SKenneth E. Jansen
53059599516SKenneth E. Jansen      tcorecpscal(1) = tcorecpscal(1) + telmcp2-telmcp1 ! elem. formation
53159599516SKenneth E. Jansen      tcorecpscal(2) = tcorecpscal(2) + tlescp2-tlescp1 ! linear alg. solution
53259599516SKenneth E. Jansen
53359599516SKenneth E. Jansen      nsolsc=5+isclr
53459599516SKenneth E. Jansen      call rstaticSclr (res, y, solinc, nsolsc) ! output scalar stats
53559599516SKenneth E. Jansenc
53659599516SKenneth E. Jansenc.... end
53759599516SKenneth E. Jansenc
53859599516SKenneth E. Jansen      return
53959599516SKenneth E. Jansen      end
54059599516SKenneth E. Jansen
54159599516SKenneth E. Jansen
54259599516SKenneth E. Jansen
54359599516SKenneth E. Jansen
54459599516SKenneth E. Jansen
545