1cb512458SBarry Smith #ifndef lint 2*d13ab20cSBarry Smith static char vcid[] = "$Id: mpiaij.c,v 1.16 1995/03/25 01:10:01 curfman Exp bsmith $"; 3cb512458SBarry Smith #endif 48a729477SBarry Smith 51eb62cbbSBarry Smith #include "mpiaij.h" 68a729477SBarry Smith #include "vec/vecimpl.h" 7d6dfbf8fSBarry Smith #include "inline/spops.h" 88a729477SBarry Smith 91eb62cbbSBarry Smith #define CHUNCKSIZE 100 101eb62cbbSBarry Smith /* 111eb62cbbSBarry Smith This is a simple minded stash. Do a linear search to determine if 121eb62cbbSBarry Smith in stash, if not add to end. 131eb62cbbSBarry Smith */ 141eb62cbbSBarry Smith static int StashValues(Stash *stash,int row,int n, int *idxn, 151eb62cbbSBarry Smith Scalar *values,InsertMode addv) 168a729477SBarry Smith { 171eb62cbbSBarry Smith int i,j,N = stash->n,found,*n_idx, *n_idy; 181eb62cbbSBarry Smith Scalar val,*n_array; 198a729477SBarry Smith 201eb62cbbSBarry Smith for ( i=0; i<n; i++ ) { 211eb62cbbSBarry Smith found = 0; 221eb62cbbSBarry Smith val = *values++; 238a729477SBarry Smith for ( j=0; j<N; j++ ) { 241eb62cbbSBarry Smith if ( stash->idx[j] == row && stash->idy[j] == idxn[i]) { 251eb62cbbSBarry Smith /* found a match */ 261eb62cbbSBarry Smith if (addv == AddValues) stash->array[j] += val; 271eb62cbbSBarry Smith else stash->array[j] = val; 281eb62cbbSBarry Smith found = 1; 298a729477SBarry Smith break; 308a729477SBarry Smith } 318a729477SBarry Smith } 321eb62cbbSBarry Smith if (!found) { /* not found so add to end */ 331eb62cbbSBarry Smith if ( stash->n == stash->nmax ) { 341eb62cbbSBarry Smith /* allocate a larger stash */ 351eb62cbbSBarry Smith n_array = (Scalar *) MALLOC( (stash->nmax + CHUNCKSIZE)*( 361eb62cbbSBarry Smith 2*sizeof(int) + sizeof(Scalar))); 371eb62cbbSBarry Smith CHKPTR(n_array); 381eb62cbbSBarry Smith n_idx = (int *) (n_array + stash->nmax + CHUNCKSIZE); 391eb62cbbSBarry Smith n_idy = (int *) (n_idx + stash->nmax + CHUNCKSIZE); 401eb62cbbSBarry Smith MEMCPY(n_array,stash->array,stash->nmax*sizeof(Scalar)); 411eb62cbbSBarry Smith MEMCPY(n_idx,stash->idx,stash->nmax*sizeof(int)); 421eb62cbbSBarry Smith MEMCPY(n_idy,stash->idy,stash->nmax*sizeof(int)); 431eb62cbbSBarry Smith if (stash->array) FREE(stash->array); 441eb62cbbSBarry Smith stash->array = n_array; stash->idx = n_idx; stash->idy = n_idy; 451eb62cbbSBarry Smith stash->nmax += CHUNCKSIZE; 461eb62cbbSBarry Smith } 471eb62cbbSBarry Smith stash->array[stash->n] = val; 481eb62cbbSBarry Smith stash->idx[stash->n] = row; 491eb62cbbSBarry Smith stash->idy[stash->n++] = idxn[i]; 501eb62cbbSBarry Smith } 518a729477SBarry Smith } 528a729477SBarry Smith return 0; 538a729477SBarry Smith } 548a729477SBarry Smith 559e25ed09SBarry Smith /* local utility routine that creates a mapping from the global column 569e25ed09SBarry Smith number to the local number in the off-diagonal part of the local 579e25ed09SBarry Smith storage of the matrix. This is done in a non scable way since the 589e25ed09SBarry Smith length of colmap equals the global matrix length. 599e25ed09SBarry Smith */ 609e25ed09SBarry Smith static int CreateColmap(Mat mat) 619e25ed09SBarry Smith { 6244a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 6344a69424SLois Curfman McInnes Mat_AIJ *B = (Mat_AIJ*) aij->B->data; 649e25ed09SBarry Smith int n = B->n,i; 659e25ed09SBarry Smith aij->colmap = (int *) MALLOC( aij->N*sizeof(int) ); CHKPTR(aij->colmap); 6639e00950SLois Curfman McInnes aij->invcolmap = (int *) MALLOC( aij->N*sizeof(int) ); CHKPTR(aij->invcolmap); 679e25ed09SBarry Smith MEMSET(aij->colmap,0,aij->N*sizeof(int)); 6839e00950SLois Curfman McInnes MEMSET(aij->invcolmap,0,aij->N*sizeof(int)); 699e25ed09SBarry Smith for ( i=0; i<n; i++ ) { 709e25ed09SBarry Smith aij->colmap[aij->garray[i]] = i+1; 7139e00950SLois Curfman McInnes aij->invcolmap[i+1] = aij->garray[i]; 729e25ed09SBarry Smith } 739e25ed09SBarry Smith return 0; 749e25ed09SBarry Smith } 759e25ed09SBarry Smith 7644a69424SLois Curfman McInnes static int MatInsertValues_MPIAIJ(Mat mat,int m,int *idxm,int n, 771eb62cbbSBarry Smith int *idxn,Scalar *v,InsertMode addv) 788a729477SBarry Smith { 7944a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 801eb62cbbSBarry Smith int ierr,i,j, rstart = aij->rstart, rend = aij->rend; 811eb62cbbSBarry Smith int cstart = aij->cstart, cend = aij->cend,row,col; 828a729477SBarry Smith 831eb62cbbSBarry Smith if (aij->insertmode != NotSetValues && aij->insertmode != addv) { 841eb62cbbSBarry Smith SETERR(1,"You cannot mix inserts and adds"); 858a729477SBarry Smith } 861eb62cbbSBarry Smith aij->insertmode = addv; 878a729477SBarry Smith for ( i=0; i<m; i++ ) { 88da3a660dSBarry Smith if (idxm[i] < 0) SETERR(1,"Negative row index"); 89da3a660dSBarry Smith if (idxm[i] >= aij->M) SETERR(1,"Row index too large"); 901eb62cbbSBarry Smith if (idxm[i] >= rstart && idxm[i] < rend) { 911eb62cbbSBarry Smith row = idxm[i] - rstart; 921eb62cbbSBarry Smith for ( j=0; j<n; j++ ) { 93da3a660dSBarry Smith if (idxn[j] < 0) SETERR(1,"Negative column index"); 94da3a660dSBarry Smith if (idxn[j] >= aij->N) SETERR(1,"Column index too large"); 951eb62cbbSBarry Smith if (idxn[j] >= cstart && idxn[j] < cend){ 961eb62cbbSBarry Smith col = idxn[j] - cstart; 971eb62cbbSBarry Smith ierr = MatSetValues(aij->A,1,&row,1,&col,v+i*n+j,addv);CHKERR(ierr); 981eb62cbbSBarry Smith } 991eb62cbbSBarry Smith else { 100d6dfbf8fSBarry Smith if (aij->assembled) { 1019e25ed09SBarry Smith if (!aij->colmap) {ierr = CreateColmap(mat); CHKERR(ierr);} 1029e25ed09SBarry Smith col = aij->colmap[idxn[j]] - 1; 1039e25ed09SBarry Smith if (col < 0) { 1049e25ed09SBarry Smith SETERR(1,"Cannot insert new off diagonal block nonzero in\ 1059e25ed09SBarry Smith already\ 106d6dfbf8fSBarry Smith assembled matrix. Contact petsc-maint@mcs.anl.gov\ 107d6dfbf8fSBarry Smith if your need this feature"); 108d6dfbf8fSBarry Smith } 1099e25ed09SBarry Smith } 1109e25ed09SBarry Smith else col = idxn[j]; 1111eb62cbbSBarry Smith ierr = MatSetValues(aij->B,1,&row,1,&col,v+i*n+j,addv);CHKERR(ierr); 1121eb62cbbSBarry Smith } 1131eb62cbbSBarry Smith } 1141eb62cbbSBarry Smith } 1151eb62cbbSBarry Smith else { 1161eb62cbbSBarry Smith ierr = StashValues(&aij->stash,idxm[i],n,idxn,v+i*n,addv);CHKERR(ierr); 1171eb62cbbSBarry Smith } 1188a729477SBarry Smith } 1198a729477SBarry Smith return 0; 1208a729477SBarry Smith } 1218a729477SBarry Smith 1228a729477SBarry Smith /* 1231eb62cbbSBarry Smith the assembly code is alot like the code for vectors, we should 1241eb62cbbSBarry Smith sometime derive a single assembly code that can be used for 1251eb62cbbSBarry Smith either case. 1268a729477SBarry Smith */ 1278a729477SBarry Smith 12844a69424SLois Curfman McInnes static int MatBeginAssemble_MPIAIJ(Mat mat) 1298a729477SBarry Smith { 13044a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 131d6dfbf8fSBarry Smith MPI_Comm comm = mat->comm; 1326abc6512SBarry Smith int numtids = aij->numtids, *owners = aij->rowners; 1331eb62cbbSBarry Smith int mytid = aij->mytid; 1341eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 1356abc6512SBarry Smith int *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work; 1361eb62cbbSBarry Smith int tag = 50, *owner,*starts,count; 1371eb62cbbSBarry Smith InsertMode addv; 1381eb62cbbSBarry Smith Scalar *rvalues,*svalues; 1391eb62cbbSBarry Smith 1401eb62cbbSBarry Smith /* make sure all processors are either in INSERTMODE or ADDMODE */ 14128988994SBarry Smith MPI_Allreduce((void *) &aij->insertmode,(void *) &addv,1,MPI_INT, 1421eb62cbbSBarry Smith MPI_BOR,comm); 1431eb62cbbSBarry Smith if (addv == (AddValues|InsertValues)) { 1441eb62cbbSBarry Smith SETERR(1,"Some processors have inserted while others have added"); 1451eb62cbbSBarry Smith } 1461eb62cbbSBarry Smith aij->insertmode = addv; /* in case this processor had no cache */ 1471eb62cbbSBarry Smith 1481eb62cbbSBarry Smith /* first count number of contributors to each processor */ 1491eb62cbbSBarry Smith nprocs = (int *) MALLOC( 2*numtids*sizeof(int) ); CHKPTR(nprocs); 1501eb62cbbSBarry Smith MEMSET(nprocs,0,2*numtids*sizeof(int)); procs = nprocs + numtids; 1511eb62cbbSBarry Smith owner = (int *) MALLOC( (aij->stash.n+1)*sizeof(int) ); CHKPTR(owner); 1521eb62cbbSBarry Smith for ( i=0; i<aij->stash.n; i++ ) { 1531eb62cbbSBarry Smith idx = aij->stash.idx[i]; 1541eb62cbbSBarry Smith for ( j=0; j<numtids; j++ ) { 1551eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 1561eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; break; 1578a729477SBarry Smith } 1588a729477SBarry Smith } 1598a729477SBarry Smith } 1601eb62cbbSBarry Smith nsends = 0; for ( i=0; i<numtids; i++ ) { nsends += procs[i];} 1611eb62cbbSBarry Smith 1621eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 1631eb62cbbSBarry Smith work = (int *) MALLOC( numtids*sizeof(int) ); CHKPTR(work); 1641eb62cbbSBarry Smith MPI_Allreduce((void *) procs,(void *) work,numtids,MPI_INT,MPI_SUM,comm); 1651eb62cbbSBarry Smith nreceives = work[mytid]; 1661eb62cbbSBarry Smith MPI_Allreduce((void *) nprocs,(void *) work,numtids,MPI_INT,MPI_MAX,comm); 1671eb62cbbSBarry Smith nmax = work[mytid]; 1681eb62cbbSBarry Smith FREE(work); 1691eb62cbbSBarry Smith 1701eb62cbbSBarry Smith /* post receives: 1711eb62cbbSBarry Smith 1) each message will consist of ordered pairs 1721eb62cbbSBarry Smith (global index,value) we store the global index as a double 173d6dfbf8fSBarry Smith to simplify the message passing. 1741eb62cbbSBarry Smith 2) since we don't know how long each individual message is we 1751eb62cbbSBarry Smith allocate the largest needed buffer for each receive. Potentially 1761eb62cbbSBarry Smith this is a lot of wasted space. 1771eb62cbbSBarry Smith 1781eb62cbbSBarry Smith 1791eb62cbbSBarry Smith This could be done better. 1801eb62cbbSBarry Smith */ 18128988994SBarry Smith rvalues = (Scalar *) MALLOC(3*(nreceives+1)*(nmax+1)*sizeof(Scalar)); 1821eb62cbbSBarry Smith CHKPTR(rvalues); 1831eb62cbbSBarry Smith recv_waits = (MPI_Request *) MALLOC((nreceives+1)*sizeof(MPI_Request)); 1841eb62cbbSBarry Smith CHKPTR(recv_waits); 1851eb62cbbSBarry Smith for ( i=0; i<nreceives; i++ ) { 1861eb62cbbSBarry Smith MPI_Irecv((void *)(rvalues+3*nmax*i),3*nmax,MPI_SCALAR,MPI_ANY_SOURCE,tag, 1871eb62cbbSBarry Smith comm,recv_waits+i); 1881eb62cbbSBarry Smith } 1891eb62cbbSBarry Smith 1901eb62cbbSBarry Smith /* do sends: 1911eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 1921eb62cbbSBarry Smith the ith processor 1931eb62cbbSBarry Smith */ 1941eb62cbbSBarry Smith svalues = (Scalar *) MALLOC( 3*(aij->stash.n+1)*sizeof(Scalar) ); 1951eb62cbbSBarry Smith CHKPTR(svalues); 1961eb62cbbSBarry Smith send_waits = (MPI_Request *) MALLOC( (nsends+1)*sizeof(MPI_Request)); 1971eb62cbbSBarry Smith CHKPTR(send_waits); 1981eb62cbbSBarry Smith starts = (int *) MALLOC( numtids*sizeof(int) ); CHKPTR(starts); 1991eb62cbbSBarry Smith starts[0] = 0; 2001eb62cbbSBarry Smith for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 2011eb62cbbSBarry Smith for ( i=0; i<aij->stash.n; i++ ) { 2021eb62cbbSBarry Smith svalues[3*starts[owner[i]]] = (Scalar) aij->stash.idx[i]; 2031eb62cbbSBarry Smith svalues[3*starts[owner[i]]+1] = (Scalar) aij->stash.idy[i]; 2041eb62cbbSBarry Smith svalues[3*(starts[owner[i]]++)+2] = aij->stash.array[i]; 2051eb62cbbSBarry Smith } 2061eb62cbbSBarry Smith FREE(owner); 2071eb62cbbSBarry Smith starts[0] = 0; 2081eb62cbbSBarry Smith for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 2091eb62cbbSBarry Smith count = 0; 2101eb62cbbSBarry Smith for ( i=0; i<numtids; i++ ) { 2111eb62cbbSBarry Smith if (procs[i]) { 2121eb62cbbSBarry Smith MPI_Isend((void*)(svalues+3*starts[i]),3*nprocs[i],MPI_SCALAR,i,tag, 2131eb62cbbSBarry Smith comm,send_waits+count++); 2141eb62cbbSBarry Smith } 2151eb62cbbSBarry Smith } 2161eb62cbbSBarry Smith FREE(starts); FREE(nprocs); 2171eb62cbbSBarry Smith 2181eb62cbbSBarry Smith /* Free cache space */ 2191eb62cbbSBarry Smith aij->stash.nmax = aij->stash.n = 0; 2201eb62cbbSBarry Smith if (aij->stash.array){ FREE(aij->stash.array); aij->stash.array = 0;} 2211eb62cbbSBarry Smith 2221eb62cbbSBarry Smith aij->svalues = svalues; aij->rvalues = rvalues; 2231eb62cbbSBarry Smith aij->nsends = nsends; aij->nrecvs = nreceives; 2241eb62cbbSBarry Smith aij->send_waits = send_waits; aij->recv_waits = recv_waits; 2251eb62cbbSBarry Smith aij->rmax = nmax; 2261eb62cbbSBarry Smith 2271eb62cbbSBarry Smith return 0; 2281eb62cbbSBarry Smith } 22944a69424SLois Curfman McInnes extern int MatSetUpMultiply_MPIAIJ(Mat); 2301eb62cbbSBarry Smith 23144a69424SLois Curfman McInnes static int MatEndAssemble_MPIAIJ(Mat mat) 2321eb62cbbSBarry Smith { 2331eb62cbbSBarry Smith int ierr; 23444a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 2351eb62cbbSBarry Smith 2361eb62cbbSBarry Smith MPI_Status *send_status,recv_status; 2376abc6512SBarry Smith int imdex,nrecvs = aij->nrecvs, count = nrecvs, i, n; 2381eb62cbbSBarry Smith int row,col; 2391eb62cbbSBarry Smith Scalar *values,val; 2401eb62cbbSBarry Smith InsertMode addv = aij->insertmode; 2411eb62cbbSBarry Smith 2421eb62cbbSBarry Smith /* wait on receives */ 2431eb62cbbSBarry Smith while (count) { 244d6dfbf8fSBarry Smith MPI_Waitany(nrecvs,aij->recv_waits,&imdex,&recv_status); 2451eb62cbbSBarry Smith /* unpack receives into our local space */ 246d6dfbf8fSBarry Smith values = aij->rvalues + 3*imdex*aij->rmax; 2471eb62cbbSBarry Smith MPI_Get_count(&recv_status,MPI_SCALAR,&n); 2481eb62cbbSBarry Smith n = n/3; 2491eb62cbbSBarry Smith for ( i=0; i<n; i++ ) { 2501eb62cbbSBarry Smith row = (int) PETSCREAL(values[3*i]) - aij->rstart; 2511eb62cbbSBarry Smith col = (int) PETSCREAL(values[3*i+1]); 2521eb62cbbSBarry Smith val = values[3*i+2]; 2531eb62cbbSBarry Smith if (col >= aij->cstart && col < aij->cend) { 2541eb62cbbSBarry Smith col -= aij->cstart; 2551eb62cbbSBarry Smith MatSetValues(aij->A,1,&row,1,&col,&val,addv); 2561eb62cbbSBarry Smith } 2571eb62cbbSBarry Smith else { 258d6dfbf8fSBarry Smith if (aij->assembled) { 2599e25ed09SBarry Smith if (!aij->colmap) {ierr = CreateColmap(mat); CHKERR(ierr);} 2609e25ed09SBarry Smith col = aij->colmap[col] - 1; 2619e25ed09SBarry Smith if (col < 0) { 2629e25ed09SBarry Smith SETERR(1,"Cannot insert new off diagonal block nonzero in\ 2639e25ed09SBarry Smith already\ 264d6dfbf8fSBarry Smith assembled matrix. Contact petsc-maint@mcs.anl.gov\ 265d6dfbf8fSBarry Smith if your need this feature"); 266d6dfbf8fSBarry Smith } 2679e25ed09SBarry Smith } 2681eb62cbbSBarry Smith MatSetValues(aij->B,1,&row,1,&col,&val,addv); 2691eb62cbbSBarry Smith } 2701eb62cbbSBarry Smith } 2711eb62cbbSBarry Smith count--; 2721eb62cbbSBarry Smith } 2731eb62cbbSBarry Smith FREE(aij->recv_waits); FREE(aij->rvalues); 2741eb62cbbSBarry Smith 2751eb62cbbSBarry Smith /* wait on sends */ 2761eb62cbbSBarry Smith if (aij->nsends) { 2771eb62cbbSBarry Smith send_status = (MPI_Status *) MALLOC( aij->nsends*sizeof(MPI_Status) ); 2781eb62cbbSBarry Smith CHKPTR(send_status); 2791eb62cbbSBarry Smith MPI_Waitall(aij->nsends,aij->send_waits,send_status); 2801eb62cbbSBarry Smith FREE(send_status); 2811eb62cbbSBarry Smith } 2821eb62cbbSBarry Smith FREE(aij->send_waits); FREE(aij->svalues); 2831eb62cbbSBarry Smith 2841eb62cbbSBarry Smith aij->insertmode = NotSetValues; 2851eb62cbbSBarry Smith ierr = MatBeginAssembly(aij->A); CHKERR(ierr); 2861eb62cbbSBarry Smith ierr = MatEndAssembly(aij->A); CHKERR(ierr); 2871eb62cbbSBarry Smith 2889e25ed09SBarry Smith if (!aij->assembled) { 289d3c9fed9SLois Curfman McInnes ierr = MatSetUpMultiply_MPIAIJ(mat); CHKERR(ierr); 2909e25ed09SBarry Smith } 2911eb62cbbSBarry Smith ierr = MatBeginAssembly(aij->B); CHKERR(ierr); 2921eb62cbbSBarry Smith ierr = MatEndAssembly(aij->B); CHKERR(ierr); 293d6dfbf8fSBarry Smith 294d6dfbf8fSBarry Smith aij->assembled = 1; 2958a729477SBarry Smith return 0; 2968a729477SBarry Smith } 2978a729477SBarry Smith 29844a69424SLois Curfman McInnes static int MatZero_MPIAIJ(Mat A) 2991eb62cbbSBarry Smith { 30044a69424SLois Curfman McInnes Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data; 3011eb62cbbSBarry Smith 3021eb62cbbSBarry Smith MatZeroEntries(l->A); MatZeroEntries(l->B); 3031eb62cbbSBarry Smith return 0; 3041eb62cbbSBarry Smith } 3051eb62cbbSBarry Smith 3061eb62cbbSBarry Smith /* again this uses the same basic stratagy as in the assembly and 3071eb62cbbSBarry Smith scatter create routines, we should try to do it systemamatically 3081eb62cbbSBarry Smith if we can figure out the proper level of generality. */ 3091eb62cbbSBarry Smith 3101eb62cbbSBarry Smith /* the code does not do the diagonal entries correctly unless the 3111eb62cbbSBarry Smith matrix is square and the column and row owerships are identical. 3121eb62cbbSBarry Smith This is a BUG. The only way to fix it seems to be to access 3131eb62cbbSBarry Smith aij->A and aij->B directly and not through the MatZeroRows() 3141eb62cbbSBarry Smith routine. 3151eb62cbbSBarry Smith */ 31644a69424SLois Curfman McInnes static int MatZeroRows_MPIAIJ(Mat A,IS is,Scalar *diag) 3171eb62cbbSBarry Smith { 31844a69424SLois Curfman McInnes Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data; 3191eb62cbbSBarry Smith int i,ierr,N, *rows,*owners = l->rowners,numtids = l->numtids; 3206abc6512SBarry Smith int *procs,*nprocs,j,found,idx,nsends,*work; 3211eb62cbbSBarry Smith int nmax,*svalues,*starts,*owner,nrecvs,mytid = l->mytid; 3226abc6512SBarry Smith int *rvalues,tag = 67,count,base,slen,n,*source; 323d6dfbf8fSBarry Smith int *lens,imdex,*lrows,*values; 324d6dfbf8fSBarry Smith MPI_Comm comm = A->comm; 3251eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 3261eb62cbbSBarry Smith MPI_Status recv_status,*send_status; 3271eb62cbbSBarry Smith IS istmp; 3281eb62cbbSBarry Smith 32944a69424SLois Curfman McInnes if (!l->assembled) SETERR(1,"MatZeroRows_MPIAIJ: must assemble matrix first"); 3301eb62cbbSBarry Smith ierr = ISGetLocalSize(is,&N); CHKERR(ierr); 3311eb62cbbSBarry Smith ierr = ISGetIndices(is,&rows); CHKERR(ierr); 3321eb62cbbSBarry Smith 3331eb62cbbSBarry Smith /* first count number of contributors to each processor */ 3341eb62cbbSBarry Smith nprocs = (int *) MALLOC( 2*numtids*sizeof(int) ); CHKPTR(nprocs); 3351eb62cbbSBarry Smith MEMSET(nprocs,0,2*numtids*sizeof(int)); procs = nprocs + numtids; 3361eb62cbbSBarry Smith owner = (int *) MALLOC((N+1)*sizeof(int)); CHKPTR(owner); /* see note*/ 3371eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 3381eb62cbbSBarry Smith idx = rows[i]; 3391eb62cbbSBarry Smith found = 0; 3401eb62cbbSBarry Smith for ( j=0; j<numtids; j++ ) { 3411eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 3421eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 3431eb62cbbSBarry Smith } 3441eb62cbbSBarry Smith } 345d6dfbf8fSBarry Smith if (!found) SETERR(1,"Imdex out of range"); 3461eb62cbbSBarry Smith } 3471eb62cbbSBarry Smith nsends = 0; for ( i=0; i<numtids; i++ ) { nsends += procs[i];} 3481eb62cbbSBarry Smith 3491eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 3501eb62cbbSBarry Smith work = (int *) MALLOC( numtids*sizeof(int) ); CHKPTR(work); 3511eb62cbbSBarry Smith MPI_Allreduce((void *) procs,(void *) work,numtids,MPI_INT,MPI_SUM,comm); 3521eb62cbbSBarry Smith nrecvs = work[mytid]; 3531eb62cbbSBarry Smith MPI_Allreduce((void *) nprocs,(void *) work,numtids,MPI_INT,MPI_MAX,comm); 3541eb62cbbSBarry Smith nmax = work[mytid]; 3551eb62cbbSBarry Smith FREE(work); 3561eb62cbbSBarry Smith 3571eb62cbbSBarry Smith /* post receives: */ 35828988994SBarry Smith rvalues = (int *) MALLOC((nrecvs+1)*(nmax+1)*sizeof(int)); /*see note */ 3591eb62cbbSBarry Smith CHKPTR(rvalues); 3601eb62cbbSBarry Smith recv_waits = (MPI_Request *) MALLOC((nrecvs+1)*sizeof(MPI_Request)); 3611eb62cbbSBarry Smith CHKPTR(recv_waits); 3621eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 3631eb62cbbSBarry Smith MPI_Irecv((void *)(rvalues+nmax*i),nmax,MPI_INT,MPI_ANY_SOURCE,tag, 3641eb62cbbSBarry Smith comm,recv_waits+i); 3651eb62cbbSBarry Smith } 3661eb62cbbSBarry Smith 3671eb62cbbSBarry Smith /* do sends: 3681eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 3691eb62cbbSBarry Smith the ith processor 3701eb62cbbSBarry Smith */ 3711eb62cbbSBarry Smith svalues = (int *) MALLOC( (N+1)*sizeof(int) ); CHKPTR(svalues); 3721eb62cbbSBarry Smith send_waits = (MPI_Request *) MALLOC( (nsends+1)*sizeof(MPI_Request)); 3731eb62cbbSBarry Smith CHKPTR(send_waits); 3741eb62cbbSBarry Smith starts = (int *) MALLOC( (numtids+1)*sizeof(int) ); CHKPTR(starts); 3751eb62cbbSBarry Smith starts[0] = 0; 3761eb62cbbSBarry Smith for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3771eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 3781eb62cbbSBarry Smith svalues[starts[owner[i]]++] = rows[i]; 3791eb62cbbSBarry Smith } 3801eb62cbbSBarry Smith ISRestoreIndices(is,&rows); 3811eb62cbbSBarry Smith 3821eb62cbbSBarry Smith starts[0] = 0; 3831eb62cbbSBarry Smith for ( i=1; i<numtids+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3841eb62cbbSBarry Smith count = 0; 3851eb62cbbSBarry Smith for ( i=0; i<numtids; i++ ) { 3861eb62cbbSBarry Smith if (procs[i]) { 3871eb62cbbSBarry Smith MPI_Isend((void*)(svalues+starts[i]),nprocs[i],MPI_INT,i,tag, 3881eb62cbbSBarry Smith comm,send_waits+count++); 3891eb62cbbSBarry Smith } 3901eb62cbbSBarry Smith } 3911eb62cbbSBarry Smith FREE(starts); 3921eb62cbbSBarry Smith 3931eb62cbbSBarry Smith base = owners[mytid]; 3941eb62cbbSBarry Smith 3951eb62cbbSBarry Smith /* wait on receives */ 3961eb62cbbSBarry Smith lens = (int *) MALLOC( 2*(nrecvs+1)*sizeof(int) ); CHKPTR(lens); 3971eb62cbbSBarry Smith source = lens + nrecvs; 3981eb62cbbSBarry Smith count = nrecvs; slen = 0; 3991eb62cbbSBarry Smith while (count) { 400d6dfbf8fSBarry Smith MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status); 4011eb62cbbSBarry Smith /* unpack receives into our local space */ 4021eb62cbbSBarry Smith MPI_Get_count(&recv_status,MPI_INT,&n); 403d6dfbf8fSBarry Smith source[imdex] = recv_status.MPI_SOURCE; 404d6dfbf8fSBarry Smith lens[imdex] = n; 4051eb62cbbSBarry Smith slen += n; 4061eb62cbbSBarry Smith count--; 4071eb62cbbSBarry Smith } 4081eb62cbbSBarry Smith FREE(recv_waits); 4091eb62cbbSBarry Smith 4101eb62cbbSBarry Smith /* move the data into the send scatter */ 4111eb62cbbSBarry Smith lrows = (int *) MALLOC( slen*sizeof(int) ); CHKPTR(lrows); 4121eb62cbbSBarry Smith count = 0; 4131eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 4141eb62cbbSBarry Smith values = rvalues + i*nmax; 4151eb62cbbSBarry Smith for ( j=0; j<lens[i]; j++ ) { 4161eb62cbbSBarry Smith lrows[count++] = values[j] - base; 4171eb62cbbSBarry Smith } 4181eb62cbbSBarry Smith } 4191eb62cbbSBarry Smith FREE(rvalues); FREE(lens); 4201eb62cbbSBarry Smith FREE(owner); FREE(nprocs); 4211eb62cbbSBarry Smith 4221eb62cbbSBarry Smith /* actually zap the local rows */ 4231eb62cbbSBarry Smith ierr = ISCreateSequential(slen,lrows,&istmp); CHKERR(ierr); FREE(lrows); 4241eb62cbbSBarry Smith ierr = MatZeroRows(l->A,istmp,diag); CHKERR(ierr); 4251eb62cbbSBarry Smith ierr = MatZeroRows(l->B,istmp,0); CHKERR(ierr); 4261eb62cbbSBarry Smith ierr = ISDestroy(istmp); CHKERR(ierr); 4271eb62cbbSBarry Smith 4281eb62cbbSBarry Smith /* wait on sends */ 4291eb62cbbSBarry Smith if (nsends) { 4301eb62cbbSBarry Smith send_status = (MPI_Status *) MALLOC( nsends*sizeof(MPI_Status) ); 4311eb62cbbSBarry Smith CHKPTR(send_status); 4321eb62cbbSBarry Smith MPI_Waitall(nsends,send_waits,send_status); 4331eb62cbbSBarry Smith FREE(send_status); 4341eb62cbbSBarry Smith } 4351eb62cbbSBarry Smith FREE(send_waits); FREE(svalues); 4361eb62cbbSBarry Smith 4371eb62cbbSBarry Smith 4381eb62cbbSBarry Smith return 0; 4391eb62cbbSBarry Smith } 4401eb62cbbSBarry Smith 44144a69424SLois Curfman McInnes static int MatMult_MPIAIJ(Mat aijin,Vec xx,Vec yy) 4421eb62cbbSBarry Smith { 44344a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data; 4441eb62cbbSBarry Smith int ierr; 44544a69424SLois Curfman McInnes if (!aij->assembled) SETERR(1,"MatMult_MPIAIJ: must assemble matrix first"); 446d6dfbf8fSBarry Smith ierr = VecScatterBegin(xx,0,aij->lvec,0,InsertValues,ScatterAll,aij->Mvctx); 4471eb62cbbSBarry Smith CHKERR(ierr); 4481eb62cbbSBarry Smith ierr = MatMult(aij->A,xx,yy); CHKERR(ierr); 449d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,aij->lvec,0,InsertValues,ScatterAll,aij->Mvctx); 450d6dfbf8fSBarry Smith CHKERR(ierr); 4511eb62cbbSBarry Smith ierr = MatMultAdd(aij->B,aij->lvec,yy,yy); CHKERR(ierr); 4521eb62cbbSBarry Smith return 0; 4531eb62cbbSBarry Smith } 4541eb62cbbSBarry Smith 45544a69424SLois Curfman McInnes static int MatMultAdd_MPIAIJ(Mat aijin,Vec xx,Vec yy,Vec zz) 456da3a660dSBarry Smith { 45744a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data; 458da3a660dSBarry Smith int ierr; 45944a69424SLois Curfman McInnes if (!aij->assembled) SETERR(1,"MatMult_MPIAIJ: must assemble matrix first"); 460da3a660dSBarry Smith ierr = VecScatterBegin(xx,0,aij->lvec,0,InsertValues,ScatterAll,aij->Mvctx); 461da3a660dSBarry Smith CHKERR(ierr); 462da3a660dSBarry Smith ierr = MatMultAdd(aij->A,xx,yy,zz); CHKERR(ierr); 463da3a660dSBarry Smith ierr = VecScatterEnd(xx,0,aij->lvec,0,InsertValues,ScatterAll,aij->Mvctx); 464da3a660dSBarry Smith CHKERR(ierr); 465da3a660dSBarry Smith ierr = MatMultAdd(aij->B,aij->lvec,zz,zz); CHKERR(ierr); 466da3a660dSBarry Smith return 0; 467da3a660dSBarry Smith } 468da3a660dSBarry Smith 46944a69424SLois Curfman McInnes static int MatMultTrans_MPIAIJ(Mat aijin,Vec xx,Vec yy) 470da3a660dSBarry Smith { 47144a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data; 472da3a660dSBarry Smith int ierr; 473da3a660dSBarry Smith 47444a69424SLois Curfman McInnes if (!aij->assembled) 47544a69424SLois Curfman McInnes SETERR(1,"MatMulTrans_MPIAIJ: must assemble matrix first"); 476da3a660dSBarry Smith /* do nondiagonal part */ 477da3a660dSBarry Smith ierr = MatMultTrans(aij->B,xx,aij->lvec); CHKERR(ierr); 478da3a660dSBarry Smith /* send it on its way */ 479da3a660dSBarry Smith ierr = VecScatterBegin(aij->lvec,0,yy,0,AddValues, 480da3a660dSBarry Smith ScatterAll|ScatterReverse,aij->Mvctx); CHKERR(ierr); 481da3a660dSBarry Smith /* do local part */ 482da3a660dSBarry Smith ierr = MatMultTrans(aij->A,xx,yy); CHKERR(ierr); 483da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 484da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 485da3a660dSBarry Smith /* but is not perhaps always true. */ 486da3a660dSBarry Smith ierr = VecScatterEnd(aij->lvec,0,yy,0,AddValues,ScatterAll|ScatterReverse, 487da3a660dSBarry Smith aij->Mvctx); CHKERR(ierr); 488da3a660dSBarry Smith return 0; 489da3a660dSBarry Smith } 490da3a660dSBarry Smith 49144a69424SLois Curfman McInnes static int MatMultTransAdd_MPIAIJ(Mat aijin,Vec xx,Vec yy,Vec zz) 492da3a660dSBarry Smith { 49344a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data; 494da3a660dSBarry Smith int ierr; 495da3a660dSBarry Smith 49644a69424SLois Curfman McInnes if (!aij->assembled) 49744a69424SLois Curfman McInnes SETERR(1,"MatMulTransAdd_MPIAIJ: must assemble matrix first"); 498da3a660dSBarry Smith /* do nondiagonal part */ 499da3a660dSBarry Smith ierr = MatMultTrans(aij->B,xx,aij->lvec); CHKERR(ierr); 500da3a660dSBarry Smith /* send it on its way */ 501da3a660dSBarry Smith ierr = VecScatterBegin(aij->lvec,0,zz,0,AddValues, 502da3a660dSBarry Smith ScatterAll|ScatterReverse,aij->Mvctx); CHKERR(ierr); 503da3a660dSBarry Smith /* do local part */ 504da3a660dSBarry Smith ierr = MatMultTransAdd(aij->A,xx,yy,zz); CHKERR(ierr); 505da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 506da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 507da3a660dSBarry Smith /* but is not perhaps always true. */ 508da3a660dSBarry Smith ierr = VecScatterEnd(aij->lvec,0,zz,0,AddValues,ScatterAll|ScatterReverse, 509da3a660dSBarry Smith aij->Mvctx); CHKERR(ierr); 510da3a660dSBarry Smith return 0; 511da3a660dSBarry Smith } 512da3a660dSBarry Smith 5131eb62cbbSBarry Smith /* 5141eb62cbbSBarry Smith This only works correctly for square matrices where the subblock A->A is the 5151eb62cbbSBarry Smith diagonal block 5161eb62cbbSBarry Smith */ 51744a69424SLois Curfman McInnes static int MatGetDiag_MPIAIJ(Mat Ain,Vec v) 5181eb62cbbSBarry Smith { 51944a69424SLois Curfman McInnes Mat_MPIAIJ *A = (Mat_MPIAIJ *) Ain->data; 52044a69424SLois Curfman McInnes if (!A->assembled) SETERR(1,"MatGetDiag_MPIAIJ: must assemble matrix first"); 5211eb62cbbSBarry Smith return MatGetDiagonal(A->A,v); 5221eb62cbbSBarry Smith } 5231eb62cbbSBarry Smith 52444a69424SLois Curfman McInnes static int MatDestroy_MPIAIJ(PetscObject obj) 5251eb62cbbSBarry Smith { 5261eb62cbbSBarry Smith Mat mat = (Mat) obj; 52744a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 5281eb62cbbSBarry Smith int ierr; 529a5a9c739SBarry Smith #if defined(PETSC_LOG) 530a5a9c739SBarry Smith PLogObjectState(obj,"Rows %d Cols %d",aij->M,aij->N); 531a5a9c739SBarry Smith #endif 5321eb62cbbSBarry Smith FREE(aij->rowners); 5331eb62cbbSBarry Smith ierr = MatDestroy(aij->A); CHKERR(ierr); 5341eb62cbbSBarry Smith ierr = MatDestroy(aij->B); CHKERR(ierr); 5359e25ed09SBarry Smith if (aij->colmap) FREE(aij->colmap); 5369e25ed09SBarry Smith if (aij->garray) FREE(aij->garray); 5371eb62cbbSBarry Smith if (aij->lvec) VecDestroy(aij->lvec); 5381eb62cbbSBarry Smith if (aij->Mvctx) VecScatterCtxDestroy(aij->Mvctx); 539a5a9c739SBarry Smith FREE(aij); 540a5a9c739SBarry Smith PLogObjectDestroy(mat); 541a5a9c739SBarry Smith PETSCHEADERDESTROY(mat); 5421eb62cbbSBarry Smith return 0; 5431eb62cbbSBarry Smith } 5441eb62cbbSBarry Smith 54544a69424SLois Curfman McInnes static int MatView_MPIAIJ(PetscObject obj,Viewer viewer) 5461eb62cbbSBarry Smith { 5471eb62cbbSBarry Smith Mat mat = (Mat) obj; 54844a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 5491eb62cbbSBarry Smith int ierr; 550*d13ab20cSBarry Smith PetscObject vobj = (PetscObject) viewer; 5511eb62cbbSBarry Smith 55244a69424SLois Curfman McInnes if (!aij->assembled) SETERR(1,"MatView_MPIAIJ: must assemble matrix first"); 553*d13ab20cSBarry Smith if (vobj->cookie == VIEWER_COOKIE) { 554*d13ab20cSBarry Smith FILE *fd = ViewerFileGetPointer(viewer); 555*d13ab20cSBarry Smith if (vobj->type == FILE_VIEWER) { 556d6dfbf8fSBarry Smith MPE_Seq_begin(mat->comm,1); 557*d13ab20cSBarry Smith fprintf(fd,"[%d] rows %d starts %d ends %d cols %d starts %d ends %d\n", 5581eb62cbbSBarry Smith aij->mytid,aij->m,aij->rstart,aij->rend,aij->n,aij->cstart, 5591eb62cbbSBarry Smith aij->cend); 560da69df5fSBarry Smith ierr = MatView(aij->A,viewer); CHKERR(ierr); 561da69df5fSBarry Smith ierr = MatView(aij->B,viewer); CHKERR(ierr); 562*d13ab20cSBarry Smith fflush(fd); 563d6dfbf8fSBarry Smith MPE_Seq_end(mat->comm,1); 564*d13ab20cSBarry Smith } 565*d13ab20cSBarry Smith } 5661eb62cbbSBarry Smith return 0; 5671eb62cbbSBarry Smith } 5681eb62cbbSBarry Smith 569d3c9fed9SLois Curfman McInnes extern int MatMarkDiag_AIJ(Mat_AIJ *); 5701eb62cbbSBarry Smith /* 5711eb62cbbSBarry Smith This has to provide several versions. 5721eb62cbbSBarry Smith 5731eb62cbbSBarry Smith 1) per sequential 5741eb62cbbSBarry Smith 2) a) use only local smoothing updating outer values only once. 5751eb62cbbSBarry Smith b) local smoothing updating outer values each inner iteration 576d6dfbf8fSBarry Smith 3) color updating out values betwen colors. 5771eb62cbbSBarry Smith */ 57844a69424SLois Curfman McInnes static int MatRelax_MPIAIJ(Mat matin,Vec bb,double omega,int flag,double shift, 5798a729477SBarry Smith int its,Vec xx) 5808a729477SBarry Smith { 58144a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 582d6dfbf8fSBarry Smith Mat AA = mat->A, BB = mat->B; 58344a69424SLois Curfman McInnes Mat_AIJ *A = (Mat_AIJ *) AA->data, *B = (Mat_AIJ *)BB->data; 5846abc6512SBarry Smith Scalar zero = 0.0,*b,*x,*xs,*ls,d,*v,sum,scale,*t,*ts; 5856abc6512SBarry Smith int ierr,*idx, *diag; 5866abc6512SBarry Smith int n = mat->n, m = mat->m, i; 587da3a660dSBarry Smith Vec tt; 5888a729477SBarry Smith 58944a69424SLois Curfman McInnes if (!mat->assembled) SETERR(1,"MatRelax_MPIAIJ: must assemble matrix first"); 5901eb62cbbSBarry Smith 591d6dfbf8fSBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); VecGetArray(mat->lvec,&ls); 592d6dfbf8fSBarry Smith xs = x -1; /* shift by one for index start of 1 */ 593da3a660dSBarry Smith ls--; 594d3c9fed9SLois Curfman McInnes if (!A->diag) {if ((ierr = MatMarkDiag_AIJ(A))) return ierr;} 595d6dfbf8fSBarry Smith diag = A->diag; 596acb40c82SBarry Smith if (flag == SOR_APPLY_UPPER || flag == SOR_APPLY_LOWER) { 597acb40c82SBarry Smith SETERR(1,"That option not yet support for parallel AIJ matrices"); 598acb40c82SBarry Smith } 599da3a660dSBarry Smith if (flag & SOR_EISENSTAT) { 600da3a660dSBarry Smith /* Let A = L + U + D; where L is lower trianglar, 601da3a660dSBarry Smith U is upper triangular, E is diagonal; This routine applies 602da3a660dSBarry Smith 603da3a660dSBarry Smith (L + E)^{-1} A (U + E)^{-1} 604da3a660dSBarry Smith 605da3a660dSBarry Smith to a vector efficiently using Eisenstat's trick. This is for 606da3a660dSBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 607da3a660dSBarry Smith is the relaxation factor. 608da3a660dSBarry Smith */ 609da3a660dSBarry Smith ierr = VecCreate(xx,&tt); CHKERR(ierr); 610da3a660dSBarry Smith VecGetArray(tt,&t); 611da3a660dSBarry Smith scale = (2.0/omega) - 1.0; 612da3a660dSBarry Smith /* x = (E + U)^{-1} b */ 613da3a660dSBarry Smith VecSet(&zero,mat->lvec); 614da3a660dSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineUp, 615da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 616da3a660dSBarry Smith for ( i=m-1; i>-1; i-- ) { 617da3a660dSBarry Smith n = A->i[i+1] - diag[i] - 1; 618da3a660dSBarry Smith idx = A->j + diag[i]; 619da3a660dSBarry Smith v = A->a + diag[i]; 620da3a660dSBarry Smith sum = b[i]; 621da3a660dSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 622da3a660dSBarry Smith d = shift + A->a[diag[i]-1]; 623da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 624da3a660dSBarry Smith idx = B->j + B->i[i] - 1; 625da3a660dSBarry Smith v = B->a + B->i[i] - 1; 626da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 627da3a660dSBarry Smith x[i] = omega*(sum/d); 628da3a660dSBarry Smith } 629da3a660dSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineUp, 630da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 631da3a660dSBarry Smith 632da3a660dSBarry Smith /* t = b - (2*E - D)x */ 633da3a660dSBarry Smith v = A->a; 634da3a660dSBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ - 1])*x[i]; } 635da3a660dSBarry Smith 636da3a660dSBarry Smith /* t = (E + L)^{-1}t */ 637da3a660dSBarry Smith ts = t - 1; /* shifted by one for index start of a or mat->j*/ 638da3a660dSBarry Smith diag = A->diag; 639da3a660dSBarry Smith VecSet(&zero,mat->lvec); 640da3a660dSBarry Smith ierr = VecPipelineBegin(tt,0,mat->lvec,0,InsertValues,PipelineDown, 641da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 642da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 643da3a660dSBarry Smith n = diag[i] - A->i[i]; 644da3a660dSBarry Smith idx = A->j + A->i[i] - 1; 645da3a660dSBarry Smith v = A->a + A->i[i] - 1; 646da3a660dSBarry Smith sum = t[i]; 647da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 648da3a660dSBarry Smith d = shift + A->a[diag[i]-1]; 649da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 650da3a660dSBarry Smith idx = B->j + B->i[i] - 1; 651da3a660dSBarry Smith v = B->a + B->i[i] - 1; 652da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 653da3a660dSBarry Smith t[i] = omega*(sum/d); 654da3a660dSBarry Smith } 655da3a660dSBarry Smith ierr = VecPipelineEnd(tt,0,mat->lvec,0,InsertValues,PipelineDown, 656da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 657da3a660dSBarry Smith /* x = x + t */ 658da3a660dSBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 659da3a660dSBarry Smith VecDestroy(tt); 660da3a660dSBarry Smith return 0; 661da3a660dSBarry Smith } 662da3a660dSBarry Smith 6631eb62cbbSBarry Smith 664d6dfbf8fSBarry Smith if ((flag & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP){ 665da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 666da3a660dSBarry Smith VecSet(&zero,mat->lvec); VecSet(&zero,xx); 667da3a660dSBarry Smith } 668da3a660dSBarry Smith else { 669d6dfbf8fSBarry Smith ierr=VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterUp,mat->Mvctx); 670d6dfbf8fSBarry Smith CHKERR(ierr); 671d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterUp,mat->Mvctx); 672d6dfbf8fSBarry Smith CHKERR(ierr); 673da3a660dSBarry Smith } 674d6dfbf8fSBarry Smith while (its--) { 675d6dfbf8fSBarry Smith /* go down through the rows */ 676d6dfbf8fSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineDown, 677d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 678d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 679d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 680d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 681d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 682d6dfbf8fSBarry Smith sum = b[i]; 683d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 684d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 685d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 686d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 687d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 688d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 689d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 690d6dfbf8fSBarry Smith } 691d6dfbf8fSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineDown, 692d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 693d6dfbf8fSBarry Smith /* come up through the rows */ 694d6dfbf8fSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineUp, 695d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 696d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 697d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 698d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 699d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 700d6dfbf8fSBarry Smith sum = b[i]; 701d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 702d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 703d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 704d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 705d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 706d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 707d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 708d6dfbf8fSBarry Smith } 709d6dfbf8fSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineUp, 710d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 711d6dfbf8fSBarry Smith } 712d6dfbf8fSBarry Smith } 713d6dfbf8fSBarry Smith else if (flag & SOR_FORWARD_SWEEP){ 714da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 715da3a660dSBarry Smith VecSet(&zero,mat->lvec); 716da3a660dSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineDown, 717da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 718da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 719da3a660dSBarry Smith n = diag[i] - A->i[i]; 720da3a660dSBarry Smith idx = A->j + A->i[i] - 1; 721da3a660dSBarry Smith v = A->a + A->i[i] - 1; 722da3a660dSBarry Smith sum = b[i]; 723da3a660dSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 724da3a660dSBarry Smith d = shift + A->a[diag[i]-1]; 725da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 726da3a660dSBarry Smith idx = B->j + B->i[i] - 1; 727da3a660dSBarry Smith v = B->a + B->i[i] - 1; 728da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 729da3a660dSBarry Smith x[i] = omega*(sum/d); 730da3a660dSBarry Smith } 731da3a660dSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineDown, 732da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 733da3a660dSBarry Smith its--; 734da3a660dSBarry Smith } 735d6dfbf8fSBarry Smith while (its--) { 736d6dfbf8fSBarry Smith ierr=VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterUp,mat->Mvctx); 737d6dfbf8fSBarry Smith CHKERR(ierr); 738d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterUp,mat->Mvctx); 739d6dfbf8fSBarry Smith CHKERR(ierr); 740d6dfbf8fSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineDown, 741d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 742d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 743d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 744d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 745d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 746d6dfbf8fSBarry Smith sum = b[i]; 747d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 748d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 749d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 750d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 751d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 752d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 753d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 754d6dfbf8fSBarry Smith } 755d6dfbf8fSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineDown, 756d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 757d6dfbf8fSBarry Smith } 758d6dfbf8fSBarry Smith } 759d6dfbf8fSBarry Smith else if (flag & SOR_BACKWARD_SWEEP){ 760da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 761da3a660dSBarry Smith VecSet(&zero,mat->lvec); 762da3a660dSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineUp, 763da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 764da3a660dSBarry Smith for ( i=m-1; i>-1; i-- ) { 765da3a660dSBarry Smith n = A->i[i+1] - diag[i] - 1; 766da3a660dSBarry Smith idx = A->j + diag[i]; 767da3a660dSBarry Smith v = A->a + diag[i]; 768da3a660dSBarry Smith sum = b[i]; 769da3a660dSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 770da3a660dSBarry Smith d = shift + A->a[diag[i]-1]; 771da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 772da3a660dSBarry Smith idx = B->j + B->i[i] - 1; 773da3a660dSBarry Smith v = B->a + B->i[i] - 1; 774da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 775da3a660dSBarry Smith x[i] = omega*(sum/d); 776da3a660dSBarry Smith } 777da3a660dSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineUp, 778da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 779da3a660dSBarry Smith its--; 780da3a660dSBarry Smith } 781d6dfbf8fSBarry Smith while (its--) { 782d6dfbf8fSBarry Smith ierr = VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterDown, 783d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 784d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterDown, 785d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 786d6dfbf8fSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineUp, 787d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 788d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 789d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 790d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 791d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 792d6dfbf8fSBarry Smith sum = b[i]; 793d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 794d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 795d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 796d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 797d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 798d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 799d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 800d6dfbf8fSBarry Smith } 801d6dfbf8fSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineUp, 802d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 803d6dfbf8fSBarry Smith } 804d6dfbf8fSBarry Smith } 805d6dfbf8fSBarry Smith else if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){ 806da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 807da3a660dSBarry Smith return MatRelax(mat->A,bb,omega,flag,shift,its,xx); 808da3a660dSBarry Smith } 809d6dfbf8fSBarry Smith ierr=VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterAll,mat->Mvctx); 810d6dfbf8fSBarry Smith CHKERR(ierr); 811d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterAll,mat->Mvctx); 812d6dfbf8fSBarry Smith CHKERR(ierr); 813d6dfbf8fSBarry Smith while (its--) { 814d6dfbf8fSBarry Smith /* go down through the rows */ 815d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 816d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 817d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 818d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 819d6dfbf8fSBarry Smith sum = b[i]; 820d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 821d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 822d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 823d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 824d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 825d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 826d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 827d6dfbf8fSBarry Smith } 828d6dfbf8fSBarry Smith /* come up through the rows */ 829d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 830d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 831d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 832d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 833d6dfbf8fSBarry Smith sum = b[i]; 834d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 835d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 836d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 837d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 838d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 839d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 840d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 841d6dfbf8fSBarry Smith } 842d6dfbf8fSBarry Smith } 843d6dfbf8fSBarry Smith } 844d6dfbf8fSBarry Smith else if (flag & SOR_LOCAL_FORWARD_SWEEP){ 845da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 846da3a660dSBarry Smith return MatRelax(mat->A,bb,omega,flag,shift,its,xx); 847da3a660dSBarry Smith } 848d6dfbf8fSBarry Smith ierr=VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterAll,mat->Mvctx); 849d6dfbf8fSBarry Smith CHKERR(ierr); 850d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterAll,mat->Mvctx); 851d6dfbf8fSBarry Smith CHKERR(ierr); 852d6dfbf8fSBarry Smith while (its--) { 853d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 854d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 855d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 856d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 857d6dfbf8fSBarry Smith sum = b[i]; 858d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 859d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 860d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 861d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 862d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 863d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 864d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 865d6dfbf8fSBarry Smith } 866d6dfbf8fSBarry Smith } 867d6dfbf8fSBarry Smith } 868d6dfbf8fSBarry Smith else if (flag & SOR_LOCAL_BACKWARD_SWEEP){ 869da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 870da3a660dSBarry Smith return MatRelax(mat->A,bb,omega,flag,shift,its,xx); 871da3a660dSBarry Smith } 872d6dfbf8fSBarry Smith ierr = VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterAll, 873d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 874d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterAll, 875d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 876d6dfbf8fSBarry Smith while (its--) { 877d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 878d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 879d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 880d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 881d6dfbf8fSBarry Smith sum = b[i]; 882d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 883d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 884d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 885d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 886d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 887d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 888d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 889d6dfbf8fSBarry Smith } 890d6dfbf8fSBarry Smith } 891d6dfbf8fSBarry Smith } 8928a729477SBarry Smith return 0; 8938a729477SBarry Smith } 89444a69424SLois Curfman McInnes static int MatInsOpt_MPIAIJ(Mat aijin,int op) 895c74985f6SBarry Smith { 89644a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data; 897c74985f6SBarry Smith 898c74985f6SBarry Smith if (op == NO_NEW_NONZERO_LOCATIONS) { 899c74985f6SBarry Smith MatSetOption(aij->A,op); 900c74985f6SBarry Smith MatSetOption(aij->B,op); 901c74985f6SBarry Smith } 902c74985f6SBarry Smith else if (op == YES_NEW_NONZERO_LOCATIONS) { 903c74985f6SBarry Smith MatSetOption(aij->A,op); 904c74985f6SBarry Smith MatSetOption(aij->B,op); 905c74985f6SBarry Smith } 906c74985f6SBarry Smith else if (op == COLUMN_ORIENTED) SETERR(1,"Column oriented not supported"); 907c74985f6SBarry Smith return 0; 908c74985f6SBarry Smith } 909c74985f6SBarry Smith 91044a69424SLois Curfman McInnes static int MatSize_MPIAIJ(Mat matin,int *m,int *n) 911c74985f6SBarry Smith { 91244a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 913c74985f6SBarry Smith *m = mat->M; *n = mat->N; 914c74985f6SBarry Smith return 0; 915c74985f6SBarry Smith } 916c74985f6SBarry Smith 91744a69424SLois Curfman McInnes static int MatLocalSize_MPIAIJ(Mat matin,int *m,int *n) 918c74985f6SBarry Smith { 91944a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 920c74985f6SBarry Smith *m = mat->m; *n = mat->n; 921c74985f6SBarry Smith return 0; 922c74985f6SBarry Smith } 923c74985f6SBarry Smith 92444a69424SLois Curfman McInnes static int MatRange_MPIAIJ(Mat matin,int *m,int *n) 925c74985f6SBarry Smith { 92644a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 927c74985f6SBarry Smith *m = mat->rstart; *n = mat->rend; 928c74985f6SBarry Smith return 0; 929c74985f6SBarry Smith } 930c74985f6SBarry Smith 93139e00950SLois Curfman McInnes static int MatGetRow_MPIAIJ(Mat matin,int row,int *nz,int **idx,Scalar **v) 93239e00950SLois Curfman McInnes { 93339e00950SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) matin->data; 93439e00950SLois Curfman McInnes Mat A = aij->A, B = aij->B; 93539e00950SLois Curfman McInnes Mat_AIJ *Ad = (Mat_AIJ *)(A->data), *Bd = (Mat_AIJ *)(B->data); 93639e00950SLois Curfman McInnes Scalar *vworkA, *vworkB; 93739e00950SLois Curfman McInnes int ierr, *cworkA, *cworkB; 93839e00950SLois Curfman McInnes int nztot, nzA, nzB, i, rstart = aij->rstart, rend = aij->rend; 93939e00950SLois Curfman McInnes int cstart = aij->cstart; 94039e00950SLois Curfman McInnes 94139e00950SLois Curfman McInnes if (!aij->assembled) 94239e00950SLois Curfman McInnes SETERR(1,"MatGetRow_MPIAIJ: Must assemble matrix first."); 94339e00950SLois Curfman McInnes if (row < rstart || row >= rend) 94439e00950SLois Curfman McInnes SETERR(1,"MatGetRow_MPIAIJ: Currently you can get only the local rows.") 94539e00950SLois Curfman McInnes 94639e00950SLois Curfman McInnes ierr = MatGetRow(A,row,&nzA,&cworkA,&vworkA); CHKERR(ierr); 94739e00950SLois Curfman McInnes for (i=0; i<nzA; i++) cworkA[i] += cstart; 94839e00950SLois Curfman McInnes ierr = MatGetRow(B,row,&nzB,&cworkB,&vworkB); CHKERR(ierr); 94939e00950SLois Curfman McInnes for (i=0; i<nzB; i++) { 95039e00950SLois Curfman McInnes printf("i=%d, invcolmap[i]=%d, cwork=%d, vwork=%g\n", 95139e00950SLois Curfman McInnes i, aij->invcolmap[i], cworkB[i], vworkB[i] ); 95239e00950SLois Curfman McInnes /* cworkB[i] = aij->invcolmap[cworkB[i]]; */ 95339e00950SLois Curfman McInnes } 95439e00950SLois Curfman McInnes 95539e00950SLois Curfman McInnes if (nztot = nzA + nzB) { 95639e00950SLois Curfman McInnes *idx = (int *) MALLOC( (nztot)*sizeof(int) ); CHKPTR(*idx); 95739e00950SLois Curfman McInnes *v = (Scalar *) MALLOC( (nztot)*sizeof(Scalar) ); CHKPTR(*v); 95839e00950SLois Curfman McInnes for ( i=0; i<nzA; i++ ) { 95939e00950SLois Curfman McInnes (*idx)[i] = cworkA[i]; 96039e00950SLois Curfman McInnes (*v)[i] = vworkA[i]; 96139e00950SLois Curfman McInnes } 96239e00950SLois Curfman McInnes for ( i=0; i<nzB; i++ ) { 96339e00950SLois Curfman McInnes (*idx)[i+nzA] = cworkB[i]; 96439e00950SLois Curfman McInnes (*v)[i] = vworkB[i]; 96539e00950SLois Curfman McInnes } 96639e00950SLois Curfman McInnes } 96739e00950SLois Curfman McInnes else {*idx = 0; *v=0;} 96839e00950SLois Curfman McInnes *nz = nztot; 96939e00950SLois Curfman McInnes ierr = MatRestoreRow(A,row,&nzA,&cworkA,&vworkA); CHKERR(ierr); 97039e00950SLois Curfman McInnes ierr = MatRestoreRow(B,row,&nzB,&cworkB,&vworkB); CHKERR(ierr); 97139e00950SLois Curfman McInnes return 0; 97239e00950SLois Curfman McInnes } 97339e00950SLois Curfman McInnes 97439e00950SLois Curfman McInnes static int MatRestoreRow_MPIAIJ(Mat mat,int row,int *nz,int **idx,Scalar **v) 97539e00950SLois Curfman McInnes { 97639e00950SLois Curfman McInnes if (*idx) FREE(*idx); 97739e00950SLois Curfman McInnes if (*v) FREE(*v); 97839e00950SLois Curfman McInnes return 0; 97939e00950SLois Curfman McInnes } 98039e00950SLois Curfman McInnes 98144a69424SLois Curfman McInnes static int MatCopy_MPIAIJ(Mat,Mat *); 98244a69424SLois Curfman McInnes extern int MatConvert_MPIAIJ(Mat,MATTYPE,Mat *); 983d6dfbf8fSBarry Smith 9848a729477SBarry Smith /* -------------------------------------------------------------------*/ 98544a69424SLois Curfman McInnes static struct _MatOps MatOps = {MatInsertValues_MPIAIJ, 98639e00950SLois Curfman McInnes MatGetRow_MPIAIJ,MatRestoreRow_MPIAIJ, 98744a69424SLois Curfman McInnes MatMult_MPIAIJ,MatMultAdd_MPIAIJ, 98844a69424SLois Curfman McInnes MatMultTrans_MPIAIJ,MatMultTransAdd_MPIAIJ, 9891eb62cbbSBarry Smith 0,0,0,0, 9908a729477SBarry Smith 0,0, 99144a69424SLois Curfman McInnes MatRelax_MPIAIJ, 9928a729477SBarry Smith 0, 9931eb62cbbSBarry Smith 0,0,0, 99444a69424SLois Curfman McInnes MatCopy_MPIAIJ, 99544a69424SLois Curfman McInnes MatGetDiag_MPIAIJ,0,0, 99644a69424SLois Curfman McInnes MatBeginAssemble_MPIAIJ,MatEndAssemble_MPIAIJ, 9971eb62cbbSBarry Smith 0, 99844a69424SLois Curfman McInnes MatInsOpt_MPIAIJ,MatZero_MPIAIJ,MatZeroRows_MPIAIJ,0, 999c74985f6SBarry Smith 0,0,0,0, 100044a69424SLois Curfman McInnes MatSize_MPIAIJ,MatLocalSize_MPIAIJ,MatRange_MPIAIJ, 100177915d1cSLois Curfman McInnes 0,0, 100244a69424SLois Curfman McInnes 0,MatConvert_MPIAIJ }; 10038a729477SBarry Smith 10048a729477SBarry Smith /*@ 10058a729477SBarry Smith 10061eb62cbbSBarry Smith MatCreateMPIAIJ - Creates a sparse parallel matrix 10071eb62cbbSBarry Smith in AIJ format. 10088a729477SBarry Smith 10098a729477SBarry Smith Input Parameters: 10101eb62cbbSBarry Smith . comm - MPI communicator 10111eb62cbbSBarry Smith . m,n - number of local rows and columns (or -1 to have calculated) 10121eb62cbbSBarry Smith . M,N - global rows and columns (or -1 to have calculated) 10131eb62cbbSBarry Smith . d_nz - total number nonzeros in diagonal portion of matrix 10141eb62cbbSBarry Smith . d_nzz - number of nonzeros per row in diagonal portion of matrix or null 10158a729477SBarry Smith . You must leave room for the diagonal entry even if it is zero. 10161eb62cbbSBarry Smith . o_nz - total number nonzeros in off-diagonal portion of matrix 10171eb62cbbSBarry Smith . o_nzz - number of nonzeros per row in off-diagonal portion of matrix 10181eb62cbbSBarry Smith . or null. You must have at least one nonzero per row. 10198a729477SBarry Smith 10208a729477SBarry Smith Output parameters: 10218a729477SBarry Smith . newmat - the matrix 10228a729477SBarry Smith 10231eb62cbbSBarry Smith Keywords: matrix, aij, compressed row, sparse, parallel 10248a729477SBarry Smith @*/ 10251eb62cbbSBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N, 10261eb62cbbSBarry Smith int d_nz,int *d_nnz, int o_nz,int *o_nnz,Mat *newmat) 10278a729477SBarry Smith { 10288a729477SBarry Smith Mat mat; 102944a69424SLois Curfman McInnes Mat_MPIAIJ *aij; 10306abc6512SBarry Smith int ierr, i,sum[2],work[2]; 10318a729477SBarry Smith *newmat = 0; 103244a69424SLois Curfman McInnes PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATMPIAIJ,comm); 1033a5a9c739SBarry Smith PLogObjectCreate(mat); 103444a69424SLois Curfman McInnes mat->data = (void *) (aij = NEW(Mat_MPIAIJ)); CHKPTR(aij); 10358a729477SBarry Smith mat->ops = &MatOps; 103644a69424SLois Curfman McInnes mat->destroy = MatDestroy_MPIAIJ; 103744a69424SLois Curfman McInnes mat->view = MatView_MPIAIJ; 10388a729477SBarry Smith mat->factor = 0; 10398a729477SBarry Smith mat->row = 0; 10408a729477SBarry Smith mat->col = 0; 1041d6dfbf8fSBarry Smith 1042d6dfbf8fSBarry Smith mat->comm = comm; 10431eb62cbbSBarry Smith aij->insertmode = NotSetValues; 10441eb62cbbSBarry Smith MPI_Comm_rank(comm,&aij->mytid); 10451eb62cbbSBarry Smith MPI_Comm_size(comm,&aij->numtids); 10461eb62cbbSBarry Smith 10471eb62cbbSBarry Smith if (M == -1 || N == -1) { 10481eb62cbbSBarry Smith work[0] = m; work[1] = n; 1049d6dfbf8fSBarry Smith MPI_Allreduce((void *) work,(void *) sum,2,MPI_INT,MPI_SUM,comm ); 10501eb62cbbSBarry Smith if (M == -1) M = sum[0]; 10511eb62cbbSBarry Smith if (N == -1) N = sum[1]; 10521eb62cbbSBarry Smith } 10531eb62cbbSBarry Smith if (m == -1) {m = M/aij->numtids + ((M % aij->numtids) > aij->mytid);} 10541eb62cbbSBarry Smith if (n == -1) {n = N/aij->numtids + ((N % aij->numtids) > aij->mytid);} 10558a729477SBarry Smith aij->m = m; 10568a729477SBarry Smith aij->n = n; 10571eb62cbbSBarry Smith aij->N = N; 10581eb62cbbSBarry Smith aij->M = M; 10591eb62cbbSBarry Smith 10601eb62cbbSBarry Smith /* build local table of row and column ownerships */ 10611eb62cbbSBarry Smith aij->rowners = (int *) MALLOC(2*(aij->numtids+2)*sizeof(int)); 10621eb62cbbSBarry Smith CHKPTR(aij->rowners); 10631eb62cbbSBarry Smith aij->cowners = aij->rowners + aij->numtids + 1; 10641eb62cbbSBarry Smith MPI_Allgather(&m,1,MPI_INT,aij->rowners+1,1,MPI_INT,comm); 10651eb62cbbSBarry Smith aij->rowners[0] = 0; 10661eb62cbbSBarry Smith for ( i=2; i<=aij->numtids; i++ ) { 10671eb62cbbSBarry Smith aij->rowners[i] += aij->rowners[i-1]; 10688a729477SBarry Smith } 10691eb62cbbSBarry Smith aij->rstart = aij->rowners[aij->mytid]; 10701eb62cbbSBarry Smith aij->rend = aij->rowners[aij->mytid+1]; 10711eb62cbbSBarry Smith MPI_Allgather(&n,1,MPI_INT,aij->cowners+1,1,MPI_INT,comm); 10721eb62cbbSBarry Smith aij->cowners[0] = 0; 10731eb62cbbSBarry Smith for ( i=2; i<=aij->numtids; i++ ) { 10741eb62cbbSBarry Smith aij->cowners[i] += aij->cowners[i-1]; 10758a729477SBarry Smith } 10761eb62cbbSBarry Smith aij->cstart = aij->cowners[aij->mytid]; 10771eb62cbbSBarry Smith aij->cend = aij->cowners[aij->mytid+1]; 10788a729477SBarry Smith 10798a729477SBarry Smith 10801eb62cbbSBarry Smith ierr = MatCreateSequentialAIJ(m,n,d_nz,d_nnz,&aij->A); CHKERR(ierr); 1081a5a9c739SBarry Smith PLogObjectParent(mat,aij->A); 10821eb62cbbSBarry Smith ierr = MatCreateSequentialAIJ(m,N,o_nz,o_nnz,&aij->B); CHKERR(ierr); 1083a5a9c739SBarry Smith PLogObjectParent(mat,aij->B); 10848a729477SBarry Smith 10851eb62cbbSBarry Smith /* build cache for off array entries formed */ 10861eb62cbbSBarry Smith aij->stash.nmax = CHUNCKSIZE; /* completely arbratray number */ 10871eb62cbbSBarry Smith aij->stash.n = 0; 10881eb62cbbSBarry Smith aij->stash.array = (Scalar *) MALLOC( aij->stash.nmax*(2*sizeof(int) + 10891eb62cbbSBarry Smith sizeof(Scalar))); CHKPTR(aij->stash.array); 10901eb62cbbSBarry Smith aij->stash.idx = (int *) (aij->stash.array + aij->stash.nmax); 10911eb62cbbSBarry Smith aij->stash.idy = (int *) (aij->stash.idx + aij->stash.nmax); 10929e25ed09SBarry Smith aij->colmap = 0; 10939e25ed09SBarry Smith aij->garray = 0; 10948a729477SBarry Smith 10951eb62cbbSBarry Smith /* stuff used for matrix vector multiply */ 10961eb62cbbSBarry Smith aij->lvec = 0; 10971eb62cbbSBarry Smith aij->Mvctx = 0; 1098d6dfbf8fSBarry Smith aij->assembled = 0; 10998a729477SBarry Smith 1100d6dfbf8fSBarry Smith *newmat = mat; 1101d6dfbf8fSBarry Smith return 0; 1102d6dfbf8fSBarry Smith } 1103c74985f6SBarry Smith 110444a69424SLois Curfman McInnes static int MatCopy_MPIAIJ(Mat matin,Mat *newmat) 1105d6dfbf8fSBarry Smith { 1106d6dfbf8fSBarry Smith Mat mat; 110744a69424SLois Curfman McInnes Mat_MPIAIJ *aij,*oldmat = (Mat_MPIAIJ *) matin->data; 11086abc6512SBarry Smith int ierr; 1109d6dfbf8fSBarry Smith *newmat = 0; 1110d6dfbf8fSBarry Smith 1111d6dfbf8fSBarry Smith if (!oldmat->assembled) SETERR(1,"Cannot copy unassembled matrix"); 111244a69424SLois Curfman McInnes PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATMPIAIJ,matin->comm); 1113a5a9c739SBarry Smith PLogObjectCreate(mat); 111444a69424SLois Curfman McInnes mat->data = (void *) (aij = NEW(Mat_MPIAIJ)); CHKPTR(aij); 1115d6dfbf8fSBarry Smith mat->ops = &MatOps; 111644a69424SLois Curfman McInnes mat->destroy = MatDestroy_MPIAIJ; 111744a69424SLois Curfman McInnes mat->view = MatView_MPIAIJ; 1118d6dfbf8fSBarry Smith mat->factor = matin->factor; 1119d6dfbf8fSBarry Smith mat->row = 0; 1120d6dfbf8fSBarry Smith mat->col = 0; 1121d6dfbf8fSBarry Smith 1122d6dfbf8fSBarry Smith aij->m = oldmat->m; 1123d6dfbf8fSBarry Smith aij->n = oldmat->n; 1124d6dfbf8fSBarry Smith aij->M = oldmat->M; 1125d6dfbf8fSBarry Smith aij->N = oldmat->N; 1126d6dfbf8fSBarry Smith 1127d6dfbf8fSBarry Smith aij->assembled = 1; 1128d6dfbf8fSBarry Smith aij->rstart = oldmat->rstart; 1129d6dfbf8fSBarry Smith aij->rend = oldmat->rend; 1130d6dfbf8fSBarry Smith aij->cstart = oldmat->cstart; 1131d6dfbf8fSBarry Smith aij->cend = oldmat->cend; 1132d6dfbf8fSBarry Smith aij->numtids = oldmat->numtids; 1133d6dfbf8fSBarry Smith aij->mytid = oldmat->mytid; 1134d6dfbf8fSBarry Smith aij->insertmode = NotSetValues; 1135d6dfbf8fSBarry Smith 1136d6dfbf8fSBarry Smith aij->rowners = (int *) MALLOC( (aij->numtids+1)*sizeof(int) ); 1137d6dfbf8fSBarry Smith CHKPTR(aij->rowners); 1138d6dfbf8fSBarry Smith MEMCPY(aij->rowners,oldmat->rowners,(aij->numtids+1)*sizeof(int)); 1139d6dfbf8fSBarry Smith aij->stash.nmax = 0; 1140d6dfbf8fSBarry Smith aij->stash.n = 0; 1141d6dfbf8fSBarry Smith aij->stash.array= 0; 11429e25ed09SBarry Smith aij->colmap = 0; 11439e25ed09SBarry Smith aij->garray = 0; 1144d6dfbf8fSBarry Smith mat->comm = matin->comm; 1145d6dfbf8fSBarry Smith 1146d6dfbf8fSBarry Smith ierr = VecCreate(oldmat->lvec,&aij->lvec); CHKERR(ierr); 1147a5a9c739SBarry Smith PLogObjectParent(mat,aij->lvec); 1148d6dfbf8fSBarry Smith ierr = VecScatterCtxCopy(oldmat->Mvctx,&aij->Mvctx); CHKERR(ierr); 1149a5a9c739SBarry Smith PLogObjectParent(mat,aij->Mvctx); 1150d6dfbf8fSBarry Smith ierr = MatCopy(oldmat->A,&aij->A); CHKERR(ierr); 1151a5a9c739SBarry Smith PLogObjectParent(mat,aij->A); 1152d6dfbf8fSBarry Smith ierr = MatCopy(oldmat->B,&aij->B); CHKERR(ierr); 1153a5a9c739SBarry Smith PLogObjectParent(mat,aij->B); 11548a729477SBarry Smith *newmat = mat; 11558a729477SBarry Smith return 0; 11568a729477SBarry Smith } 1157