xref: /phasta/phSolver/compressible/solgmrpetsc.c (revision d06028c1ed09954bf598b0145fe1ac8aecabad11)
117860365SKenneth E. Jansen #ifdef HAVE_PETSC
210167291SKenneth E. Jansen #include <stdio.h>
310167291SKenneth E. Jansen 
410167291SKenneth E. Jansen #include "petscsys.h"
510167291SKenneth E. Jansen #include "petscvec.h"
610167291SKenneth E. Jansen #include "petscmat.h"
710167291SKenneth E. Jansen #include "petscpc.h"
810167291SKenneth E. Jansen #include "petscksp.h"
910167291SKenneth E. Jansen 
1010167291SKenneth E. Jansen #include "assert.h"
1110167291SKenneth E. Jansen #include "common_c.h"
1210167291SKenneth E. Jansen 
1310167291SKenneth E. Jansen #include <FCMangle.h>
1410167291SKenneth E. Jansen #define SolGMRp FortranCInterface_GLOBAL_(solgmrp,SOLGMRP)
1510167291SKenneth E. Jansen #define SolGMRpSclr FortranCInterface_GLOBAL_(solgmrpsclr,SOLGMRPSCLR)
1610167291SKenneth E. Jansen #define ElmGMRPETSc FortranCInterface_GLOBAL_(elmgmrpetsc, ELMGMRPETSC)
1710167291SKenneth E. Jansen #define ElmGMRPETScSclr FortranCInterface_GLOBAL_(elmgmrpetscsclr, ELMGMRPETSCSCLR)
182801f607SKenneth E. Jansen #define rstatp FortranCInterface_GLOBAL_(rstatp, RSTATP)
192801f607SKenneth E. Jansen #define rstatpSclr FortranCInterface_GLOBAL_(rstatpsclr, RSTATPSCLR)
2010167291SKenneth E. Jansen #define min(a,b) (((a) < (b)) ? (a) : (b))
2110167291SKenneth E. Jansen #define max(a,b) (((a) > (b)) ? (a) : (b))
2210167291SKenneth E. Jansen #define get_time FortranCInterface_GLOBAL_(get_time,GET_TIME)
2310167291SKenneth E. Jansen #define get_max_time_diff FortranCInterface_GLOBAL_(get_max_time_diff,GET_MAX_TIME_DIFF)
2410167291SKenneth E. Jansen 
2510167291SKenneth E. Jansen typedef long long int gcorp_t;
2610167291SKenneth E. Jansen 
2710167291SKenneth E. Jansen void get_time(uint64_t* rv, uint64_t* cycle);
2810167291SKenneth E. Jansen void get_max_time_diff(uint64_t* first, uint64_t* last, uint64_t* c_first, uint64_t* c_last, char* lbl);
2910167291SKenneth E. Jansen 
3010167291SKenneth E. Jansen //#include "auxmpi.h"
3110167291SKenneth E. Jansen //      static PetscOffset poff;
3210167291SKenneth E. Jansen 
3310167291SKenneth E. Jansen       static Mat lhsP;
3410167291SKenneth E. Jansen       static PC pc;
3510167291SKenneth E. Jansen       static KSP ksp;
360be30ed5SKenneth E. Jansen       static Vec DyP, resP, DyPLocal;
3710167291SKenneth E. Jansen       static PetscErrorCode ierr;
3810167291SKenneth E. Jansen       static PetscInt PetscOne, PetscRow, PetscCol, LocalRow, LocalCol;
398ae99c59SKenneth E. Jansen       static IS LocalIndexSet;
4010167291SKenneth E. Jansen       static ISLocalToGlobalMapping VectorMapping;
4110167291SKenneth E. Jansen       static  VecScatter scatter7;
4210167291SKenneth E. Jansen       static int firstpetsccall = 1;
43*6d194905SKenneth E. Jansen       static  PetscInt maxitsHist;
44*6d194905SKenneth E. Jansen       static  PetscReal* resHist;
4510167291SKenneth E. Jansen 
4610167291SKenneth E. Jansen       static Mat lhsPs;
4710167291SKenneth E. Jansen       static PC pcs;
4810167291SKenneth E. Jansen       static KSP ksps;
490be30ed5SKenneth E. Jansen       static Vec DyPs, resPs, DyPLocals;
508ae99c59SKenneth E. Jansen       static IS LocalIndexSets;
5110167291SKenneth E. Jansen       static ISLocalToGlobalMapping VectorMappings;
5210167291SKenneth E. Jansen       static VecScatter scatter7s;
5310167291SKenneth E. Jansen       static int firstpetsccalls = 1;
5410167291SKenneth E. Jansen 
55175e1b6bSKenneth E. Jansen       static int rankdump=-1;  // 8121 was the problem rank with 3.5.3
562801f607SKenneth E. Jansen       PetscReal resNrm;
57175e1b6bSKenneth E. Jansen 
58175e1b6bSKenneth E. Jansen 
SolGMRp(double * y,double * ac,double * yold,double * x,int * iBC,double * BC,int * col,int * row,double * lhsk,double * res,double * BDiag,int * iper,int * ilwork,double * shp,double * shgl,double * shpb,double * shglb,double * Dy,double * rerr,long long int * fncorp)5910167291SKenneth E. Jansen void     SolGMRp(double* y,         double* ac,        double* yold,
6010167291SKenneth E. Jansen      	double* x,         int* iBC,       double* BC,
6110167291SKenneth E. Jansen      	int* col,       int* row,       double* lhsk,
6210167291SKenneth E. Jansen      	double* res,       double* BDiag,
6310167291SKenneth E. Jansen      	int* iper,
6410167291SKenneth E. Jansen      	int* ilwork,    double* shp,       double* shgl,      double* shpb,
6510167291SKenneth E. Jansen      	double* shglb,     double* Dy,        double* rerr, long long int* fncorp)
6610167291SKenneth E. Jansen {
6710167291SKenneth E. Jansen //
6810167291SKenneth E. Jansen // ----------------------------------------------------------------------
6910167291SKenneth E. Jansen //
7010167291SKenneth E. Jansen //  This is the preconditioned GMRES driver routine.
7110167291SKenneth E. Jansen //
7210167291SKenneth E. Jansen // input:
7310167291SKenneth E. Jansen //  y      (nshg,ndof)           : Y-variables at n+alpha_v
7410167291SKenneth E. Jansen //  ac     (nshg,ndof)           : Primvar. accel. variable n+alpha_m
7510167291SKenneth E. Jansen //  yold   (nshg,ndof)           : Y-variables at beginning of step
7610167291SKenneth E. Jansen //  acold  (nshg,ndof)           : Primvar. accel. variable at begng step
7710167291SKenneth E. Jansen //  x      (numnp,nsd)            : node coordinates
7810167291SKenneth E. Jansen //  iBC    (nshg)                : BC codes
7910167291SKenneth E. Jansen //  BC     (nshg,ndofBC)         : BC constraint parameters
8010167291SKenneth E. Jansen //  shp(b) (nen,maxsh,melCat)     : element shape functions (boundary)
8110167291SKenneth E. Jansen //  shgl(b)(nsd,nen,maxsh,melCat) : local gradients of shape functions
8210167291SKenneth E. Jansen //
8310167291SKenneth E. Jansen // output:
8410167291SKenneth E. Jansen //  res    (nshg,nflow)           : preconditioned residual
8510167291SKenneth E. Jansen //  BDiag  (nshg,nflow,nflow)      : block-diagonal preconditioner
8610167291SKenneth E. Jansen //
8710167291SKenneth E. Jansen //
8810167291SKenneth E. Jansen // Zdenek Johan,  Winter 1991.  (Fortran 90)
8910167291SKenneth E. Jansen // ----------------------------------------------------------------------
9010167291SKenneth E. Jansen //
9110167291SKenneth E. Jansen 
9210167291SKenneth E. Jansen 
9310167291SKenneth E. Jansen // Get variables from common_c.h
9410167291SKenneth E. Jansen       int nshg, nflow, nsd, iownnodes;
9510167291SKenneth E. Jansen       nshg  = conpar.nshg;
9610167291SKenneth E. Jansen       nflow = conpar.nflow;
9710167291SKenneth E. Jansen       nsd = NSD;
9810167291SKenneth E. Jansen       iownnodes = conpar.iownnodes;
9910167291SKenneth E. Jansen 
10010167291SKenneth E. Jansen       gcorp_t nshgt;
10110167291SKenneth E. Jansen       gcorp_t mbeg;
10210167291SKenneth E. Jansen       gcorp_t mend;
10310167291SKenneth E. Jansen       nshgt = (gcorp_t) newdim.nshgt; //Fix nshgt in common_c.h
10410167291SKenneth E. Jansen       mbeg = (gcorp_t) newdim.minowned;
10510167291SKenneth E. Jansen       mend = (gcorp_t) newdim.maxowned;
10610167291SKenneth E. Jansen 
10710167291SKenneth E. Jansen 
10810167291SKenneth E. Jansen       int node, element, var, eqn;
10910167291SKenneth E. Jansen       double valtoinsert;
11010167291SKenneth E. Jansen       int nenl, iel, lelCat, lcsyst, iorder, nshl;
11110167291SKenneth E. Jansen       int mattyp, ndofl, nsymdl, npro, ngauss, nppro;
11210167291SKenneth E. Jansen       double DyFlat[nshg*nflow];
11310167291SKenneth E. Jansen       double DyFlatPhasta[nshg*nflow];
11410167291SKenneth E. Jansen       double rmes[nshg*nflow];
11510167291SKenneth E. Jansen // DEBUG
11610167291SKenneth E. Jansen       int i,j,k,l,m;
11710167291SKenneth E. Jansen 
11810167291SKenneth E. Jansen       // FIXME: PetscScalar
11910167291SKenneth E. Jansen       double  real_rtol, real_abstol, real_dtol;
12010167291SKenneth E. Jansen // /DEBUG
12110167291SKenneth E. Jansen       //double parray[1]; // Should be a PetscScalar
12210167291SKenneth E. Jansen       double *parray; // Should be a PetscScalar
12310167291SKenneth E. Jansen       PetscInt petsc_bs, petsc_m, petsc_M,petsc_PA;
12410167291SKenneth E. Jansen       PetscInt petsc_n;
12510167291SKenneth E. Jansen       PetscOne = 1;
12610167291SKenneth E. Jansen       uint64_t duration[4];
12710167291SKenneth E. Jansen 
12810167291SKenneth E. Jansen       gcorp_t glbNZ;
12910167291SKenneth E. Jansen 
13010167291SKenneth E. Jansen       if(firstpetsccall == 1) {
1318ae99c59SKenneth E. Jansen //Everthing in this conditional block should be moved to a function to estimate the size of PETSc's matrix which improves time on the first matrix fill
1328ae99c59SKenneth E. Jansen //
133f4d0b58bSKenneth E. Jansen         PetscInt* idiagnz= (PetscInt*) malloc(sizeof(PetscInt)*iownnodes);
134f4d0b58bSKenneth E. Jansen         PetscInt* iodiagnz= (PetscInt*) malloc(sizeof(PetscInt)*iownnodes);
13510167291SKenneth E. Jansen         for(i=0;i<iownnodes;i++) {
13610167291SKenneth E. Jansen              idiagnz[i]=0;
13710167291SKenneth E. Jansen              iodiagnz[i]=0;
13810167291SKenneth E. Jansen         }
13910167291SKenneth E. Jansen         i=0;
14010167291SKenneth E. Jansen         for(k=0;k<nshg;k++) {
14110167291SKenneth E. Jansen           if((fncorp[k] < mbeg) || (fncorp[k] >mend)){
14210167291SKenneth E. Jansen // this node is not owned by this rank so we skip
14310167291SKenneth E. Jansen           } else {
14410167291SKenneth E. Jansen            for(j=col[i]-1;j<col[i+1]-1;j++) {
14579f1763eSKenneth E. Jansen //          assert(row[j]<=nshg);
14679f1763eSKenneth E. Jansen //          assert(fncorp[row[j]-1]<=nshgt);
14779f1763eSKenneth E. Jansen              glbNZ=fncorp[row[j]-1];
14810167291SKenneth E. Jansen              if((glbNZ < mbeg) || (glbNZ > mend)) {
14910167291SKenneth E. Jansen                 iodiagnz[i]++;
15010167291SKenneth E. Jansen              } else {
15110167291SKenneth E. Jansen                 idiagnz[i]++;
15210167291SKenneth E. Jansen              }
15310167291SKenneth E. Jansen            }
15410167291SKenneth E. Jansen            i++;
15510167291SKenneth E. Jansen           }
15610167291SKenneth E. Jansen         }
15710167291SKenneth E. Jansen         gcorp_t mind=1000;
15810167291SKenneth E. Jansen         gcorp_t mino=1000;
15910167291SKenneth E. Jansen         gcorp_t maxd=0;
16010167291SKenneth E. Jansen         gcorp_t maxo=0;
16110167291SKenneth E. Jansen         for(i=0;i<iownnodes;i++) {
16210167291SKenneth E. Jansen            mind=min(mind,idiagnz[i]);
16310167291SKenneth E. Jansen            mino=min(mino,iodiagnz[i]);
16410167291SKenneth E. Jansen            maxd=max(maxd,idiagnz[i]);
16510167291SKenneth E. Jansen            maxo=max(maxo,iodiagnz[i]);
16610167291SKenneth E. Jansen //           iodiagnz[i]=max(iodiagnz[i],10);
16710167291SKenneth E. Jansen //           idiagnz[i]=max(idiagnz[i],10);
16810167291SKenneth E. Jansen //           iodiagnz[i]=2*iodiagnz[i];   //  estimate a bit higher for off-part interactions
16910167291SKenneth E. Jansen //           idiagnz[i]=2*idiagnz[i];   //  estimate a bit higher for off-part interactions
17010167291SKenneth E. Jansen         }
17110167291SKenneth E. Jansen // the above was pretty good but below is faster and not too much more memory...of course once you do this
17210167291SKenneth E. Jansen // could just use the constant fill parameters in create but keep it alive for potential later optimization
17310167291SKenneth E. Jansen 
17410167291SKenneth E. Jansen         for(i=0;i<iownnodes;i++) {
17510167291SKenneth E. Jansen            iodiagnz[i]=1.3*maxd;
17610167291SKenneth E. Jansen            idiagnz[i]=1.3*maxd;
17710167291SKenneth E. Jansen         }
17810167291SKenneth E. Jansen 
17910167291SKenneth E. Jansen 
18010167291SKenneth E. Jansen 
18110167291SKenneth E. Jansen         if(workfc.numpe < 200){
18210167291SKenneth E. Jansen           printf("myrank,i,iownnodes,nshg %d %d %d %d \n",workfc.myrank,i,iownnodes,nshg);
18310167291SKenneth E. Jansen           printf("myrank,mind,maxd,mino,maxo %d %d %d %d %d \n",workfc.myrank,mind,maxd,mino,maxo);
18410167291SKenneth E. Jansen         }
18510167291SKenneth E. Jansen         // Print debug info
18610167291SKenneth E. Jansen         if(nshgt < 200){
18710167291SKenneth E. Jansen           int irank;
18810167291SKenneth E. Jansen           for(irank=0;irank<workfc.numpe;irank++) {
18910167291SKenneth E. Jansen             if(irank == workfc.myrank){
19010167291SKenneth E. Jansen               printf("mbeg,mend,myrank,idiagnz, iodiagnz %d %d %d \n",mbeg,mend,workfc.myrank);
19110167291SKenneth E. Jansen               for(i=0;i<iownnodes;i++) {
19210167291SKenneth E. Jansen                 printf("%d %ld %ld \n",workfc.myrank,idiagnz[i],iodiagnz[i]);
19310167291SKenneth E. Jansen               }
19410167291SKenneth E. Jansen             }
19510167291SKenneth E. Jansen            MPI_Barrier(MPI_COMM_WORLD);
19610167291SKenneth E. Jansen           }
19710167291SKenneth E. Jansen         }
19810167291SKenneth E. Jansen         petsc_bs = (PetscInt) nflow;
19910167291SKenneth E. Jansen         petsc_m  = (PetscInt) nflow* (PetscInt) iownnodes;
20010167291SKenneth E. Jansen         petsc_M  = (PetscInt) nshgt * (PetscInt) nflow;
20110167291SKenneth E. Jansen         petsc_PA  = (PetscInt) 40;
20210167291SKenneth E. Jansen         ierr = MatCreateBAIJ(PETSC_COMM_WORLD, petsc_bs, petsc_m, petsc_m, petsc_M, petsc_M,
20310167291SKenneth E. Jansen                              0, idiagnz, 0, iodiagnz, &lhsP);
20410167291SKenneth E. Jansen 
20510167291SKenneth E. Jansen         ierr = MatSetOption(lhsP, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
20610167291SKenneth E. Jansen 
20710167291SKenneth E. Jansen // Next is Jed Brown's improvement to imprint Assembly to make that stage scalable after the first call
2088ae99c59SKenneth E. Jansen #ifdef JEDBROWN
20910167291SKenneth E. Jansen         ierr = MatSetOption(lhsP, MAT_SUBSET_OFF_PROC_ENTRIES, PETSC_TRUE);
210dcce770dSKenneth E. Jansen #endif
21110167291SKenneth E. Jansen         ierr = MatSetUp(lhsP);
21210167291SKenneth E. Jansen 
21310167291SKenneth E. Jansen       PetscInt myMatStart, myMatEnd;
21410167291SKenneth E. Jansen       ierr = MatGetOwnershipRange(lhsP, &myMatStart, &myMatEnd);
215175e1b6bSKenneth E. Jansen //debug
216175e1b6bSKenneth E. Jansen       if(workfc.myrank == rankdump) printf("Flow myrank,myMatStart,myMatEnd %d,%d,%d, \n", workfc.myrank,myMatStart,myMatEnd);
21710167291SKenneth E. Jansen       }
21810167291SKenneth E. Jansen       if(genpar.lhs ==1) ierr = MatZeroEntries(lhsP);
21910167291SKenneth E. Jansen 
2208ae99c59SKenneth E. Jansen //
2218ae99c59SKenneth E. Jansen // .... *******************>> Element Data Formation <<******************
2228ae99c59SKenneth E. Jansen //
2238ae99c59SKenneth E. Jansen //
2248ae99c59SKenneth E. Jansen // .... set the parameters for flux and surface tension calculations
2258ae99c59SKenneth E. Jansen //
2268ae99c59SKenneth E. Jansen //
2278ae99c59SKenneth E. Jansen       genpar.idflx = 0 ;
2288ae99c59SKenneth E. Jansen       if(genpar.idiff >= 1)  genpar.idflx= genpar.idflx + (nflow-1) * nsd;
2298ae99c59SKenneth E. Jansen       if (genpar.isurf == 1) genpar.idflx=genpar.idflx + nsd;
23010167291SKenneth E. Jansen 
2318ae99c59SKenneth E. Jansen       get_time(duration, (duration+1));
23210167291SKenneth E. Jansen 
23310167291SKenneth E. Jansen       ElmGMRPETSc(y,             ac,            x,
23410167291SKenneth E. Jansen                   shp,           shgl,          iBC,
23510167291SKenneth E. Jansen                   BC,            shpb,
23610167291SKenneth E. Jansen                   shglb,         res,
23710167291SKenneth E. Jansen                   rmes,                   iper,
23810167291SKenneth E. Jansen                   ilwork,        rerr ,         &lhsP);
2398ae99c59SKenneth E. Jansen 
24010167291SKenneth E. Jansen       get_time((duration+2), (duration+3));
24110167291SKenneth E. Jansen       get_max_time_diff((duration), (duration+2),
24210167291SKenneth E. Jansen                         (duration+1), (duration+3),
24310167291SKenneth E. Jansen                         "ElmGMRPETSc \0");  // char(0))
24410167291SKenneth E. Jansen        if(firstpetsccall == 1) {
24510167291SKenneth E. Jansen       // Setup IndexSet. For now, we mimic vector insertion procedure
24610167291SKenneth E. Jansen       // Since we always reference by global indexes this doesn't matter
24710167291SKenneth E. Jansen       // except for cache performance)
24810167291SKenneth E. Jansen       // TODO: Better arrangment?
2498636783dSKenneth E. Jansen       PetscInt* indexsetary = malloc(sizeof(PetscInt)*nflow*nshg);
2508636783dSKenneth E. Jansen       PetscInt nodetoinsert;
25110167291SKenneth E. Jansen         nodetoinsert = 0;
25210167291SKenneth E. Jansen         k=0;
2538ae99c59SKenneth E. Jansen //debug
254175e1b6bSKenneth E. Jansen          if(workfc.myrank == rankdump) {
255175e1b6bSKenneth E. Jansen              printf("myrank,i,iownnodes,nshg %d %d %d %d \n",workfc.myrank,i,iownnodes,nshg);
256175e1b6bSKenneth E. Jansen              printf("myrank,mbeg,mend %d %d %d \n",workfc.myrank,mbeg,mend);
257175e1b6bSKenneth E. Jansen            }
25810167291SKenneth E. Jansen         if(workfc.numpe > 1) {
25910167291SKenneth E. Jansen           for (i=0; i<nshg ; i++) {
26010167291SKenneth E. Jansen             nodetoinsert = fncorp[i]-1;
261175e1b6bSKenneth E. Jansen //debug
262175e1b6bSKenneth E. Jansen          if(workfc.myrank == rankdump) {
263175e1b6bSKenneth E. Jansen              printf("myrank,i,nodetoinsert %d %d %d \n",workfc.myrank,i,nodetoinsert);
264175e1b6bSKenneth E. Jansen          }
265175e1b6bSKenneth E. Jansen 
26610167291SKenneth E. Jansen //            assert(fncorp[i]>=0);
26710167291SKenneth E. Jansen             for (j=1; j<=nflow; j++) {
26810167291SKenneth E. Jansen               indexsetary[k] = nodetoinsert*nflow+(j-1);
26910167291SKenneth E. Jansen               assert(indexsetary[k]>=0);
27010167291SKenneth E. Jansen //              assert(fncorp[i]>=0);
27110167291SKenneth E. Jansen               k = k+1;
27210167291SKenneth E. Jansen             }
27310167291SKenneth E. Jansen           }
27410167291SKenneth E. Jansen         }
27510167291SKenneth E. Jansen         else {
27610167291SKenneth E. Jansen           for (i=0; i<nshg ; i++) {
27710167291SKenneth E. Jansen             nodetoinsert = i;
27810167291SKenneth E. Jansen             for (j=1; j<=nflow; j++) {
27910167291SKenneth E. Jansen               indexsetary[k] = nodetoinsert*nflow+(j-1);
28010167291SKenneth E. Jansen               k = k+1;
28110167291SKenneth E. Jansen             }
28210167291SKenneth E. Jansen           }
28310167291SKenneth E. Jansen         }
28410167291SKenneth E. Jansen 
28510167291SKenneth E. Jansen //  Create Vector Index Maps
28610167291SKenneth E. Jansen         petsc_n  = (PetscInt) nshg * (PetscInt) nflow;
28710167291SKenneth E. Jansen         ierr = ISCreateGeneral(PETSC_COMM_SELF, petsc_n, indexsetary,
28810167291SKenneth E. Jansen      	       PETSC_COPY_VALUES, &LocalIndexSet);
2898636783dSKenneth E. Jansen       free(indexsetary);
29010167291SKenneth E. Jansen       }
29110167291SKenneth E. Jansen       if(genpar.lhs == 1) {
29210167291SKenneth E. Jansen       get_time((duration), (duration+1));
29310167291SKenneth E. Jansen       ierr = MatAssemblyBegin(lhsP, MAT_FINAL_ASSEMBLY);
29410167291SKenneth E. Jansen       ierr = MatAssemblyEnd(lhsP, MAT_FINAL_ASSEMBLY);
29510167291SKenneth E. Jansen       get_time((duration+2),(duration+3));
29610167291SKenneth E. Jansen       get_max_time_diff((duration), (duration+2),
29710167291SKenneth E. Jansen                         (duration+1), (duration+3),
29810167291SKenneth E. Jansen                         "MatAssembly \0"); // char(0))
29910167291SKenneth E. Jansen       get_time(duration, (duration+1));
30010167291SKenneth E. Jansen       }
30110167291SKenneth E. Jansen       if(firstpetsccall == 1) {
30210167291SKenneth E. Jansen         ierr = MatGetLocalSize(lhsP, &LocalRow, &LocalCol);
30310167291SKenneth E. Jansen         ierr = VecCreateMPI(PETSC_COMM_WORLD, LocalRow, petsc_M, &resP);
30410167291SKenneth E. Jansen         ierr = VecCreateMPI(PETSC_COMM_WORLD, LocalRow, petsc_M, &DyP);
30510167291SKenneth E. Jansen       }
30610167291SKenneth E. Jansen       ierr = VecZeroEntries(resP);
30710167291SKenneth E. Jansen       if(firstpetsccall == 1) {
30810167291SKenneth E. Jansen         ierr = VecCreateSeq(PETSC_COMM_SELF, petsc_n, &DyPLocal);
30910167291SKenneth E. Jansen       }
31010167291SKenneth E. Jansen 
31110167291SKenneth E. Jansen       PetscRow=0;
31210167291SKenneth E. Jansen       k = 0;
31310167291SKenneth E. Jansen       int index;
31410167291SKenneth E. Jansen       for (i=0; i<nshg; i++ ){
31510167291SKenneth E. Jansen         for (j = 1; j<=nflow; j++){
31610167291SKenneth E. Jansen           index = i + (j-1)*nshg;
31710167291SKenneth E. Jansen           valtoinsert = res[index];
31810167291SKenneth E. Jansen           if(workfc.numpe > 1) {
31910167291SKenneth E. Jansen             PetscRow = (fncorp[i]-1)*nflow+(j-1);
32010167291SKenneth E. Jansen           }
32110167291SKenneth E. Jansen           else {
32210167291SKenneth E. Jansen             PetscRow = i*nflow+(j-1);
32310167291SKenneth E. Jansen           }
32410167291SKenneth E. Jansen           assert(fncorp[i]<=nshgt);
32510167291SKenneth E. Jansen           assert(fncorp[i]>0);
32610167291SKenneth E. Jansen           assert(PetscRow>=0);
32710167291SKenneth E. Jansen           assert(PetscRow<=nshgt*nflow);
3282801f607SKenneth E. Jansen           ierr =  VecSetValue(resP, PetscRow, valtoinsert, ADD_VALUES);
32910167291SKenneth E. Jansen         }
33010167291SKenneth E. Jansen       }
33110167291SKenneth E. Jansen       ierr = VecAssemblyBegin(resP);
33210167291SKenneth E. Jansen       ierr = VecAssemblyEnd(resP);
3332801f607SKenneth E. Jansen       ierr = VecNorm(resP,NORM_2,&resNrm);
33410167291SKenneth E. Jansen       get_time((duration+2), (duration+3));
33510167291SKenneth E. Jansen       get_max_time_diff((duration), (duration+2),
33610167291SKenneth E. Jansen                         (duration+1), (duration+3),
33710167291SKenneth E. Jansen                         "VectorWorkPre \0"); // char(0))
33810167291SKenneth E. Jansen 
33910167291SKenneth E. Jansen       get_time((duration),(duration+1));
34010167291SKenneth E. Jansen 
34110167291SKenneth E. Jansen       if(firstpetsccall == 1) {
34210167291SKenneth E. Jansen         ierr = KSPCreate(PETSC_COMM_WORLD, &ksp);
34310167291SKenneth E. Jansen         ierr = KSPSetOperators(ksp, lhsP, lhsP);
34410167291SKenneth E. Jansen         ierr = KSPGetPC(ksp, &pc);
34510167291SKenneth E. Jansen         ierr = PCSetType(pc, PCPBJACOBI);
34610167291SKenneth E. Jansen         PetscInt maxits;
34710167291SKenneth E. Jansen         maxits = (PetscInt)  solpar.nGMRES * (PetscInt) solpar.Kspace;
34810167291SKenneth E. Jansen         ierr = KSPSetTolerances(ksp, timdat.etol, PETSC_DEFAULT, PETSC_DEFAULT, maxits);
34910167291SKenneth E. Jansen         ierr = KSPSetFromOptions(ksp);
350*6d194905SKenneth E. Jansen         maxitsHist=1000;
351*6d194905SKenneth E. Jansen         resHist= (PetscReal*) malloc(sizeof(PetscReal)*maxitsHist);
352*6d194905SKenneth E. Jansen         ierr = KSPSetResidualHistory(ksp,resHist,maxitsHist, PETSC_TRUE);
35310167291SKenneth E. Jansen       }
35410167291SKenneth E. Jansen       ierr = KSPSolve(ksp, resP, DyP);
355*6d194905SKenneth E. Jansen       ierr = KSPGetResidualHistory(ksp,&resHist,&maxitsHist);
356*6d194905SKenneth E. Jansen       PetscReal resNrmP=resHist[0];
35710167291SKenneth E. Jansen       get_time((duration+2),(duration+3));
35810167291SKenneth E. Jansen       get_max_time_diff((duration), (duration+2),
35910167291SKenneth E. Jansen                         (duration+1), (duration+3),
36010167291SKenneth E. Jansen                         "KSPSolve \0"); // char(0))
36110167291SKenneth E. Jansen       get_time((duration),(duration+1));
36210167291SKenneth E. Jansen       if(firstpetsccall == 1) {
3638ae99c59SKenneth E. Jansen       ierr = VecScatterCreate(DyP, LocalIndexSet, DyPLocal, NULL, &scatter7);
36410167291SKenneth E. Jansen       }
36510167291SKenneth E. Jansen       ierr = VecScatterBegin(scatter7, DyP, DyPLocal, INSERT_VALUES, SCATTER_FORWARD);
36610167291SKenneth E. Jansen       ierr = VecScatterEnd(scatter7, DyP, DyPLocal, INSERT_VALUES, SCATTER_FORWARD);
36710167291SKenneth E. Jansen       ierr = VecGetArray(DyPLocal, &parray);
36810167291SKenneth E. Jansen       PetscRow = 0;
36910167291SKenneth E. Jansen       for ( node = 0; node< nshg; node++) {
37010167291SKenneth E. Jansen         for (var = 1; var<= nflow; var++) {
37110167291SKenneth E. Jansen           index = node + (var-1)*nshg;
37210167291SKenneth E. Jansen           Dy[index] = parray[PetscRow];
37310167291SKenneth E. Jansen           PetscRow = PetscRow+1;
37410167291SKenneth E. Jansen         }
37510167291SKenneth E. Jansen       }
37610167291SKenneth E. Jansen       ierr = VecRestoreArray(DyPLocal, &parray);
37710167291SKenneth E. Jansen 
37810167291SKenneth E. Jansen       firstpetsccall = 0;
3798636783dSKenneth E. Jansen 
38010167291SKenneth E. Jansen // .... output the statistics
38110167291SKenneth E. Jansen //
38210167291SKenneth E. Jansen       itrpar.iKs=0; // see rstat()
38310167291SKenneth E. Jansen       PetscInt its;
38410167291SKenneth E. Jansen       ierr = KSPGetIterationNumber(ksp, &its);
38510167291SKenneth E. Jansen       itrpar.iKs = (int) its;
386*6d194905SKenneth E. Jansen       /*
387*6d194905SKenneth E. Jansen       PetscReal scale=1.0/sqrt(1.0*nshgt);
388*6d194905SKenneth E. Jansen       if(workfc.myrank ==0) {
389*6d194905SKenneth E. Jansen         printf("node resNrmP rosqrtNshgt\n");
390*6d194905SKenneth E. Jansen         for ( node = 0; node<its; node++) {
391*6d194905SKenneth E. Jansen           printf(" %d %f %f \n",node,resHist[node],scale*resHist[node]);
392*6d194905SKenneth E. Jansen         }
393*6d194905SKenneth E. Jansen       }
394*6d194905SKenneth E. Jansen */
39510167291SKenneth E. Jansen       get_time((duration+2),(duration+3));
39610167291SKenneth E. Jansen       get_max_time_diff((duration), (duration+2),
39710167291SKenneth E. Jansen                         (duration+1), (duration+3),
39810167291SKenneth E. Jansen                         "solWork \0"); // char(0))
39910167291SKenneth E. Jansen       itrpar.ntotGM += (int) its;
400*6d194905SKenneth E. Jansen       rstatp (&resNrm,&resNrmP);
40110167291SKenneth E. Jansen //
40210167291SKenneth E. Jansen // .... end
40310167291SKenneth E. Jansen //
40410167291SKenneth E. Jansen }
40510167291SKenneth E. Jansen 
SolGMRpSclr(double * y,double * ac,double * x,double * elDw,int * iBC,double * BC,int * col,int * row,int * iper,int * ilwork,double * shp,double * shgl,double * shpb,double * shglb,double * res,double * Dy,long long int * fncorp)40610167291SKenneth E. Jansen void     SolGMRpSclr(double* y,         double* ac,
40710167291SKenneth E. Jansen      	double* x,      double* elDw,   int* iBC,       double* BC,
40810167291SKenneth E. Jansen      	int* col,       int* row,
40910167291SKenneth E. Jansen      	int* iper,
41010167291SKenneth E. Jansen      	int* ilwork,    double* shp,       double* shgl,      double* shpb,
41110167291SKenneth E. Jansen      	double* shglb,  double* res,   double* Dy,     long long int* fncorp)
41210167291SKenneth E. Jansen {
41310167291SKenneth E. Jansen //
41410167291SKenneth E. Jansen // ----------------------------------------------------------------------
41510167291SKenneth E. Jansen //
41610167291SKenneth E. Jansen //  This is the preconditioned GMRES driver routine.
41710167291SKenneth E. Jansen //
41810167291SKenneth E. Jansen // input:
41910167291SKenneth E. Jansen //  y      (nshg,ndof)           : Y-variables at n+alpha_v
42010167291SKenneth E. Jansen //  ac     (nshg,ndof)           : Primvar. accel. variable n+alpha_m
42110167291SKenneth E. Jansen //  yold   (nshg,ndof)           : Y-variables at beginning of step
42210167291SKenneth E. Jansen //  acold  (nshg,ndof)           : Primvar. accel. variable at begng step
42310167291SKenneth E. Jansen //  x      (numnp,nsd)            : node coordinates
42410167291SKenneth E. Jansen //  iBC    (nshg)                : BC codes
42510167291SKenneth E. Jansen //  BC     (nshg,ndofBC)         : BC constraint parameters
42610167291SKenneth E. Jansen //  shp(b) (nen,maxsh,melCat)     : element shape functions (boundary)
42710167291SKenneth E. Jansen //  shgl(b)(nsd,nen,maxsh,melCat) : local gradients of shape functions
42810167291SKenneth E. Jansen //
42910167291SKenneth E. Jansen // ----------------------------------------------------------------------
43010167291SKenneth E. Jansen //
43110167291SKenneth E. Jansen 
43210167291SKenneth E. Jansen 
43310167291SKenneth E. Jansen // Get variables from common_c.h
43410167291SKenneth E. Jansen       int nshg, nflow, nsd, iownnodes;
43510167291SKenneth E. Jansen       nshg  = conpar.nshg;
43610167291SKenneth E. Jansen       nsd = NSD;
43710167291SKenneth E. Jansen       iownnodes = conpar.iownnodes;
43810167291SKenneth E. Jansen 
43910167291SKenneth E. Jansen       gcorp_t nshgt;
44010167291SKenneth E. Jansen       gcorp_t mbeg;
44110167291SKenneth E. Jansen       gcorp_t mend;
44210167291SKenneth E. Jansen       nshgt = (gcorp_t) newdim.nshgt; //Fix nshgt in common_c.h
44310167291SKenneth E. Jansen       mbeg = (gcorp_t) newdim.minowned;
44410167291SKenneth E. Jansen       mend = (gcorp_t) newdim.maxowned;
44510167291SKenneth E. Jansen 
44610167291SKenneth E. Jansen 
44710167291SKenneth E. Jansen       int node, element, var, eqn;
44810167291SKenneth E. Jansen       double valtoinsert;
44910167291SKenneth E. Jansen       int nenl, iel, lelCat, lcsyst, iorder, nshl;
45010167291SKenneth E. Jansen       int mattyp, ndofl, nsymdl, npro, ngauss, nppro;
45110167291SKenneth E. Jansen       double DyFlats[nshg];
45210167291SKenneth E. Jansen       double DyFlatPhastas[nshg];
45310167291SKenneth E. Jansen       double rmes[nshg];
45410167291SKenneth E. Jansen // DEBUG
45510167291SKenneth E. Jansen       int i,j,k,l,m;
45610167291SKenneth E. Jansen 
45710167291SKenneth E. Jansen       double  real_rtol, real_abstol, real_dtol;
45810167291SKenneth E. Jansen       double *parray;
45910167291SKenneth E. Jansen       PetscInt petsc_bs, petsc_m, petsc_M,petsc_PA;
46010167291SKenneth E. Jansen       PetscInt petsc_n;
46110167291SKenneth E. Jansen       PetscOne = 1;
46210167291SKenneth E. Jansen       uint64_t duration[4];
46310167291SKenneth E. Jansen 
46410167291SKenneth E. Jansen 
46510167291SKenneth E. Jansen //
46610167291SKenneth E. Jansen //
46710167291SKenneth E. Jansen // .... *******************>> Element Data Formation <<******************
46810167291SKenneth E. Jansen //
46910167291SKenneth E. Jansen //
47010167291SKenneth E. Jansen // .... set the parameters for flux and surface tension calculations
47110167291SKenneth E. Jansen //
47210167291SKenneth E. Jansen //
47310167291SKenneth E. Jansen       genpar.idflx = 0 ;
47410167291SKenneth E. Jansen       if(genpar.idiff >= 1)  genpar.idflx= genpar.idflx + (nflow-1) * nsd;
47510167291SKenneth E. Jansen       if (genpar.isurf == 1) genpar.idflx=genpar.idflx + nsd;
47610167291SKenneth E. Jansen 
47710167291SKenneth E. Jansen 
47810167291SKenneth E. Jansen 
47910167291SKenneth E. Jansen       gcorp_t glbNZ;
48010167291SKenneth E. Jansen 
48110167291SKenneth E. Jansen       if(firstpetsccalls == 1) {
482f4d0b58bSKenneth E. Jansen         PetscInt* idiagnz= (PetscInt*) malloc(sizeof(PetscInt)*iownnodes);
483f4d0b58bSKenneth E. Jansen         PetscInt* iodiagnz= (PetscInt*) malloc(sizeof(PetscInt)*iownnodes);
48410167291SKenneth E. Jansen         for(i=0;i<iownnodes;i++) {
48510167291SKenneth E. Jansen              idiagnz[i]=0;
48610167291SKenneth E. Jansen              iodiagnz[i]=0;
48710167291SKenneth E. Jansen         }
48810167291SKenneth E. Jansen         i=0;
48910167291SKenneth E. Jansen         for(k=0;k<nshg;k++) {
49010167291SKenneth E. Jansen           if((fncorp[k] < mbeg) || (fncorp[k] >mend)){
49110167291SKenneth E. Jansen // this node is not owned by this rank so we skip
49210167291SKenneth E. Jansen           } else {
49310167291SKenneth E. Jansen            for(j=col[i]-1;j<col[i+1]-1;j++) {
494f4d0b58bSKenneth E. Jansen              glbNZ=fncorp[row[j]-1];
49510167291SKenneth E. Jansen              if((glbNZ < mbeg) || (glbNZ > mend)) {
49610167291SKenneth E. Jansen                 iodiagnz[i]++;
49710167291SKenneth E. Jansen              } else {
49810167291SKenneth E. Jansen                 idiagnz[i]++;
49910167291SKenneth E. Jansen              }
50010167291SKenneth E. Jansen            }
50110167291SKenneth E. Jansen            i++;
50210167291SKenneth E. Jansen           }
50310167291SKenneth E. Jansen         }
50410167291SKenneth E. Jansen         gcorp_t mind=1000;
50510167291SKenneth E. Jansen         gcorp_t mino=1000;
50610167291SKenneth E. Jansen         gcorp_t maxd=0;
50710167291SKenneth E. Jansen         gcorp_t maxo=0;
50810167291SKenneth E. Jansen         for(i=0;i<iownnodes;i++) {
50910167291SKenneth E. Jansen            mind=min(mind,idiagnz[i]);
51010167291SKenneth E. Jansen            mino=min(mino,iodiagnz[i]);
51110167291SKenneth E. Jansen            maxd=max(maxd,idiagnz[i]);
51210167291SKenneth E. Jansen            maxo=max(maxo,iodiagnz[i]);
51310167291SKenneth E. Jansen         }
51410167291SKenneth E. Jansen 
51510167291SKenneth E. Jansen         for(i=0;i<iownnodes;i++) {
51610167291SKenneth E. Jansen            iodiagnz[i]=1.3*maxd;
51710167291SKenneth E. Jansen            idiagnz[i]=1.3*maxd;
51810167291SKenneth E. Jansen         }
51910167291SKenneth E. Jansen 
52010167291SKenneth E. Jansen 
52110167291SKenneth E. Jansen 
522f4d0b58bSKenneth E. Jansen //       }
523f4d0b58bSKenneth E. Jansen //       if(firstpetsccalls == 1) {
52410167291SKenneth E. Jansen 
52510167291SKenneth E. Jansen         petsc_m  = (PetscInt) iownnodes;
52610167291SKenneth E. Jansen         petsc_M  = (PetscInt) nshgt;
52710167291SKenneth E. Jansen         petsc_PA  = (PetscInt) 40;
52810167291SKenneth E. Jansen 
52910167291SKenneth E. Jansen         ierr = MatCreateAIJ(PETSC_COMM_WORLD, petsc_m, petsc_m, petsc_M, petsc_M,
53010167291SKenneth E. Jansen                             0, idiagnz, 0, iodiagnz, &lhsPs);
53110167291SKenneth E. Jansen 
53210167291SKenneth E. Jansen         ierr = MatSetOption(lhsPs, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
53310167291SKenneth E. Jansen 
53410167291SKenneth E. Jansen // Next is Jed Brown's improvement to imprint Assembly to make that stage scalable after the first call
535175e1b6bSKenneth E. Jansen #ifdef HIDEJEDBROWN
53610167291SKenneth E. Jansen         ierr = MatSetOption(lhsPs, MAT_SUBSET_OFF_PROC_ENTRIES, PETSC_TRUE);
537dcce770dSKenneth E. Jansen #endif
53810167291SKenneth E. Jansen         ierr = MatSetUp(lhsPs);
53910167291SKenneth E. Jansen 
54010167291SKenneth E. Jansen       PetscInt myMatStart, myMatEnd;
541175e1b6bSKenneth E. Jansen       ierr = MatGetOwnershipRange(lhsPs, &myMatStart, &myMatEnd);
542175e1b6bSKenneth E. Jansen       if(workfc.myrank ==rankdump) printf("Sclr myrank,myMatStart,myMatEnd %d,%d,%d, \n", workfc.myrank,myMatStart,myMatEnd);
54310167291SKenneth E. Jansen       }
54410167291SKenneth E. Jansen //      MPI_Barrier(MPI_COMM_WORLD);
54510167291SKenneth E. Jansen  //     if(workfc.myrank ==0) printf("Before MatZeroEntries  \n");
54610167291SKenneth E. Jansen       if(genpar.lhs == 1) ierr = MatZeroEntries(lhsPs);
54710167291SKenneth E. Jansen 
54810167291SKenneth E. Jansen       get_time(duration, (duration+1));
54910167291SKenneth E. Jansen 
55010167291SKenneth E. Jansen //      MPI_Barrier(MPI_COMM_WORLD);
55110167291SKenneth E. Jansen  //     if(workfc.myrank ==0) printf("Before elmgmr  \n");
55210167291SKenneth E. Jansen //      for (i=0; i<nshg ; i++) {
55310167291SKenneth E. Jansen //            assert(fncorp[i]>0);
55410167291SKenneth E. Jansen //      }
55510167291SKenneth E. Jansen 
55610167291SKenneth E. Jansen       ElmGMRPETScSclr(y,             ac,            x, elDw,
55710167291SKenneth E. Jansen                   shp,           shgl,          iBC,
55810167291SKenneth E. Jansen                   BC,            shpb,
55910167291SKenneth E. Jansen                   shglb,         res,
56010167291SKenneth E. Jansen                   rmes,                   iper,
56110167291SKenneth E. Jansen                   ilwork,        &lhsPs);
56210167291SKenneth E. Jansen        if(firstpetsccalls == 1) {
56310167291SKenneth E. Jansen       // Setup IndexSet. For now, we mimic vector insertion procedure
56410167291SKenneth E. Jansen       // Since we always reference by global indexes this doesn't matter
56510167291SKenneth E. Jansen       // except for cache performance)
56610167291SKenneth E. Jansen       // TODO: Better arrangment?
5678636783dSKenneth E. Jansen 
5688636783dSKenneth E. Jansen       PetscInt* indexsetarys = malloc(sizeof(PetscInt)*nshg);
5698636783dSKenneth E. Jansen       PetscInt nodetoinsert;
5708636783dSKenneth E. Jansen 
57110167291SKenneth E. Jansen         nodetoinsert = 0;
57210167291SKenneth E. Jansen         k=0;
573175e1b6bSKenneth E. Jansen 
574175e1b6bSKenneth E. Jansen //debug
575175e1b6bSKenneth E. Jansen          if(workfc.myrank == rankdump) {
576175e1b6bSKenneth E. Jansen              printf("myrank,i,iownnodes,nshg %d %d %d %d \n",workfc.myrank,i,iownnodes,nshg);
577175e1b6bSKenneth E. Jansen              printf("myrank,mbeg,mend %d %d %d \n",workfc.myrank,mbeg,mend);
578175e1b6bSKenneth E. Jansen            }
579175e1b6bSKenneth E. Jansen 
58010167291SKenneth E. Jansen         if(workfc.numpe > 1) {
58110167291SKenneth E. Jansen           for (i=0; i<nshg ; i++) {
58210167291SKenneth E. Jansen             nodetoinsert = fncorp[i]-1;
583175e1b6bSKenneth E. Jansen //debug
584175e1b6bSKenneth E. Jansen          if(workfc.myrank == rankdump) {
585175e1b6bSKenneth E. Jansen              printf("myrank,i,nodetoinsert %d %d %d \n",workfc.myrank,i,nodetoinsert);
586175e1b6bSKenneth E. Jansen          }
587175e1b6bSKenneth E. Jansen 
58810167291SKenneth E. Jansen             indexsetarys[k] = nodetoinsert;
58910167291SKenneth E. Jansen             k = k+1;
59010167291SKenneth E. Jansen           }
59110167291SKenneth E. Jansen         }
59210167291SKenneth E. Jansen         else {
59310167291SKenneth E. Jansen           for (i=0; i<nshg ; i++) {
59410167291SKenneth E. Jansen             nodetoinsert = i;
59510167291SKenneth E. Jansen             indexsetarys[k] = nodetoinsert;
59610167291SKenneth E. Jansen             k = k+1;
59710167291SKenneth E. Jansen           }
59810167291SKenneth E. Jansen         }
59910167291SKenneth E. Jansen 
60010167291SKenneth E. Jansen //  Create Vector Index Maps
60110167291SKenneth E. Jansen         petsc_n  = (PetscInt) nshg;
60210167291SKenneth E. Jansen         ierr = ISCreateGeneral(PETSC_COMM_SELF, petsc_n, indexsetarys,
60310167291SKenneth E. Jansen      	       PETSC_COPY_VALUES, &LocalIndexSets);
6048636783dSKenneth E. Jansen       free(indexsetarys);
60510167291SKenneth E. Jansen       }
60610167291SKenneth E. Jansen       if(genpar.lhs ==1) {
60710167291SKenneth E. Jansen       ierr = MatAssemblyBegin(lhsPs, MAT_FINAL_ASSEMBLY);
60810167291SKenneth E. Jansen       ierr = MatAssemblyEnd(lhsPs, MAT_FINAL_ASSEMBLY);
60910167291SKenneth E. Jansen       }
61010167291SKenneth E. Jansen       if(firstpetsccalls == 1) {
61110167291SKenneth E. Jansen         ierr = MatGetLocalSize(lhsPs, &LocalRow, &LocalCol);
61210167291SKenneth E. Jansen         ierr = VecCreateMPI(PETSC_COMM_WORLD, LocalRow, petsc_M, &resPs);
61310167291SKenneth E. Jansen         ierr = VecCreateMPI(PETSC_COMM_WORLD, LocalRow, petsc_M, &DyPs);
61410167291SKenneth E. Jansen       }
61510167291SKenneth E. Jansen       ierr = VecZeroEntries(resPs);
61610167291SKenneth E. Jansen       if(firstpetsccalls == 1) {
61710167291SKenneth E. Jansen         ierr = VecCreateSeq(PETSC_COMM_SELF, petsc_n, &DyPLocals);
61810167291SKenneth E. Jansen       }
61910167291SKenneth E. Jansen 
62010167291SKenneth E. Jansen       PetscRow=0;
62110167291SKenneth E. Jansen       k = 0;
62210167291SKenneth E. Jansen       int index;
62310167291SKenneth E. Jansen       for (i=0; i<nshg; i++ ){
62410167291SKenneth E. Jansen           valtoinsert = res[i];
62510167291SKenneth E. Jansen           if(workfc.numpe > 1) {
62610167291SKenneth E. Jansen             PetscRow = (fncorp[i]-1);
62710167291SKenneth E. Jansen           }
62810167291SKenneth E. Jansen           else {
629c90fc7c7SKenneth E. Jansen             PetscRow = i;
63010167291SKenneth E. Jansen           }
63110167291SKenneth E. Jansen           assert(fncorp[i]<=nshgt);
63210167291SKenneth E. Jansen           assert(fncorp[i]>0);
63310167291SKenneth E. Jansen           assert(PetscRow>=0);
63410167291SKenneth E. Jansen           assert(PetscRow<=nshgt);
6352801f607SKenneth E. Jansen           ierr =  VecSetValue(resPs, PetscRow, valtoinsert, ADD_VALUES);
63610167291SKenneth E. Jansen       }
63710167291SKenneth E. Jansen       ierr = VecAssemblyBegin(resPs);
63810167291SKenneth E. Jansen       ierr = VecAssemblyEnd(resPs);
6392801f607SKenneth E. Jansen       ierr = VecNorm(resPs,NORM_2,&resNrm);
64010167291SKenneth E. Jansen 
64110167291SKenneth E. Jansen       if(firstpetsccalls == 1) {
64210167291SKenneth E. Jansen         ierr = KSPCreate(PETSC_COMM_WORLD, &ksps);
64310167291SKenneth E. Jansen         ierr = KSPSetOperators(ksps, lhsPs, lhsPs);
64410167291SKenneth E. Jansen         ierr = KSPGetPC(ksps, &pcs);
64510167291SKenneth E. Jansen         ierr = PCSetType(pcs, PCPBJACOBI);
64610167291SKenneth E. Jansen         PetscInt maxits;
64710167291SKenneth E. Jansen         maxits = (PetscInt)  solpar.nGMRES * (PetscInt) solpar.Kspace;
64810167291SKenneth E. Jansen         ierr = KSPSetTolerances(ksps, timdat.etol, PETSC_DEFAULT, PETSC_DEFAULT, maxits);
64910167291SKenneth E. Jansen         ierr = KSPSetFromOptions(ksps);
65010167291SKenneth E. Jansen       }
65110167291SKenneth E. Jansen       ierr = KSPSolve(ksps, resPs, DyPs);
65210167291SKenneth E. Jansen       if(firstpetsccalls == 1) {
6538ae99c59SKenneth E. Jansen       ierr = VecScatterCreate(DyPs, LocalIndexSets, DyPLocals, NULL, &scatter7s);
65410167291SKenneth E. Jansen       }
65510167291SKenneth E. Jansen       ierr = VecScatterBegin(scatter7s, DyPs, DyPLocals, INSERT_VALUES, SCATTER_FORWARD);
65610167291SKenneth E. Jansen       ierr = VecScatterEnd(scatter7s, DyPs, DyPLocals, INSERT_VALUES, SCATTER_FORWARD);
65710167291SKenneth E. Jansen       ierr = VecGetArray(DyPLocals, &parray);
65810167291SKenneth E. Jansen       PetscRow = 0;
65910167291SKenneth E. Jansen       for ( node = 0; node< nshg; node++) {
66010167291SKenneth E. Jansen           index = node;
66110167291SKenneth E. Jansen           Dy[index] = parray[PetscRow];
66210167291SKenneth E. Jansen           PetscRow = PetscRow+1;
66310167291SKenneth E. Jansen       }
66410167291SKenneth E. Jansen       ierr = VecRestoreArray(DyPLocals, &parray);
66510167291SKenneth E. Jansen 
66610167291SKenneth E. Jansen       firstpetsccalls = 0;
66710167291SKenneth E. Jansen 
66810167291SKenneth E. Jansen // .... output the statistics
66910167291SKenneth E. Jansen //
67010167291SKenneth E. Jansen //      itrpar.iKss=0; // see rstat()
67110167291SKenneth E. Jansen       PetscInt its;
67210167291SKenneth E. Jansen       ierr = KSPGetIterationNumber(ksps, &its);
67310167291SKenneth E. Jansen       itrpar.iKss = (int) its;
67410167291SKenneth E. Jansen       itrpar.ntotGMs += (int) its;
6752801f607SKenneth E. Jansen       rstatpSclr (&resNrm);
67610167291SKenneth E. Jansen //
67710167291SKenneth E. Jansen // .... end
67810167291SKenneth E. Jansen //
67910167291SKenneth E. Jansen }
68017860365SKenneth E. Jansen #endif
681