1cb512458SBarry Smith #ifndef lint 2*2eb8c8abSBarry Smith static char vcid[] = "$Id: mpiaij.c,v 1.19 1995/03/26 04:43:10 bsmith 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); 669e25ed09SBarry Smith MEMSET(aij->colmap,0,aij->N*sizeof(int)); 67abc0e9e4SLois Curfman McInnes for ( i=0; i<n; i++ ) aij->colmap[aij->garray[i]] = i+1; 689e25ed09SBarry Smith return 0; 699e25ed09SBarry Smith } 709e25ed09SBarry Smith 7144a69424SLois Curfman McInnes static int MatInsertValues_MPIAIJ(Mat mat,int m,int *idxm,int n, 721eb62cbbSBarry Smith int *idxn,Scalar *v,InsertMode addv) 738a729477SBarry Smith { 7444a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 751eb62cbbSBarry Smith int ierr,i,j, rstart = aij->rstart, rend = aij->rend; 761eb62cbbSBarry Smith int cstart = aij->cstart, cend = aij->cend,row,col; 778a729477SBarry Smith 781eb62cbbSBarry Smith if (aij->insertmode != NotSetValues && aij->insertmode != addv) { 791eb62cbbSBarry Smith SETERR(1,"You cannot mix inserts and adds"); 808a729477SBarry Smith } 811eb62cbbSBarry Smith aij->insertmode = addv; 828a729477SBarry Smith for ( i=0; i<m; i++ ) { 83da3a660dSBarry Smith if (idxm[i] < 0) SETERR(1,"Negative row index"); 84da3a660dSBarry Smith if (idxm[i] >= aij->M) SETERR(1,"Row index too large"); 851eb62cbbSBarry Smith if (idxm[i] >= rstart && idxm[i] < rend) { 861eb62cbbSBarry Smith row = idxm[i] - rstart; 871eb62cbbSBarry Smith for ( j=0; j<n; j++ ) { 88da3a660dSBarry Smith if (idxn[j] < 0) SETERR(1,"Negative column index"); 89da3a660dSBarry Smith if (idxn[j] >= aij->N) SETERR(1,"Column index too large"); 901eb62cbbSBarry Smith if (idxn[j] >= cstart && idxn[j] < cend){ 911eb62cbbSBarry Smith col = idxn[j] - cstart; 921eb62cbbSBarry Smith ierr = MatSetValues(aij->A,1,&row,1,&col,v+i*n+j,addv);CHKERR(ierr); 931eb62cbbSBarry Smith } 941eb62cbbSBarry Smith else { 95d6dfbf8fSBarry Smith if (aij->assembled) { 969e25ed09SBarry Smith if (!aij->colmap) {ierr = CreateColmap(mat); CHKERR(ierr);} 979e25ed09SBarry Smith col = aij->colmap[idxn[j]] - 1; 989e25ed09SBarry Smith if (col < 0) { 999e25ed09SBarry Smith SETERR(1,"Cannot insert new off diagonal block nonzero in\ 1009e25ed09SBarry Smith already\ 101d6dfbf8fSBarry Smith assembled matrix. Contact petsc-maint@mcs.anl.gov\ 102d6dfbf8fSBarry Smith if your need this feature"); 103d6dfbf8fSBarry Smith } 1049e25ed09SBarry Smith } 1059e25ed09SBarry Smith else col = idxn[j]; 1061eb62cbbSBarry Smith ierr = MatSetValues(aij->B,1,&row,1,&col,v+i*n+j,addv);CHKERR(ierr); 1071eb62cbbSBarry Smith } 1081eb62cbbSBarry Smith } 1091eb62cbbSBarry Smith } 1101eb62cbbSBarry Smith else { 1111eb62cbbSBarry Smith ierr = StashValues(&aij->stash,idxm[i],n,idxn,v+i*n,addv);CHKERR(ierr); 1121eb62cbbSBarry Smith } 1138a729477SBarry Smith } 1148a729477SBarry Smith return 0; 1158a729477SBarry Smith } 1168a729477SBarry Smith 1178a729477SBarry Smith /* 1181eb62cbbSBarry Smith the assembly code is alot like the code for vectors, we should 1191eb62cbbSBarry Smith sometime derive a single assembly code that can be used for 1201eb62cbbSBarry Smith either case. 1218a729477SBarry Smith */ 1228a729477SBarry Smith 12344a69424SLois Curfman McInnes static int MatBeginAssemble_MPIAIJ(Mat mat) 1248a729477SBarry Smith { 12544a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 126d6dfbf8fSBarry Smith MPI_Comm comm = mat->comm; 1276abc6512SBarry Smith int numtids = aij->numtids, *owners = aij->rowners; 1281eb62cbbSBarry Smith int mytid = aij->mytid; 1291eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 1306abc6512SBarry Smith int *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work; 1311eb62cbbSBarry Smith int tag = 50, *owner,*starts,count; 1321eb62cbbSBarry Smith InsertMode addv; 1331eb62cbbSBarry Smith Scalar *rvalues,*svalues; 1341eb62cbbSBarry Smith 1351eb62cbbSBarry Smith /* make sure all processors are either in INSERTMODE or ADDMODE */ 13628988994SBarry Smith MPI_Allreduce((void *) &aij->insertmode,(void *) &addv,1,MPI_INT, 1371eb62cbbSBarry Smith MPI_BOR,comm); 1381eb62cbbSBarry Smith if (addv == (AddValues|InsertValues)) { 1391eb62cbbSBarry Smith SETERR(1,"Some processors have inserted while others have added"); 1401eb62cbbSBarry Smith } 1411eb62cbbSBarry Smith aij->insertmode = addv; /* in case this processor had no cache */ 1421eb62cbbSBarry Smith 1431eb62cbbSBarry Smith /* first count number of contributors to each processor */ 1441eb62cbbSBarry Smith nprocs = (int *) MALLOC( 2*numtids*sizeof(int) ); CHKPTR(nprocs); 1451eb62cbbSBarry Smith MEMSET(nprocs,0,2*numtids*sizeof(int)); procs = nprocs + numtids; 1461eb62cbbSBarry Smith owner = (int *) MALLOC( (aij->stash.n+1)*sizeof(int) ); CHKPTR(owner); 1471eb62cbbSBarry Smith for ( i=0; i<aij->stash.n; i++ ) { 1481eb62cbbSBarry Smith idx = aij->stash.idx[i]; 1491eb62cbbSBarry Smith for ( j=0; j<numtids; j++ ) { 1501eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 1511eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; break; 1528a729477SBarry Smith } 1538a729477SBarry Smith } 1548a729477SBarry Smith } 1551eb62cbbSBarry Smith nsends = 0; for ( i=0; i<numtids; i++ ) { nsends += procs[i];} 1561eb62cbbSBarry Smith 1571eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 1581eb62cbbSBarry Smith work = (int *) MALLOC( numtids*sizeof(int) ); CHKPTR(work); 1591eb62cbbSBarry Smith MPI_Allreduce((void *) procs,(void *) work,numtids,MPI_INT,MPI_SUM,comm); 1601eb62cbbSBarry Smith nreceives = work[mytid]; 1611eb62cbbSBarry Smith MPI_Allreduce((void *) nprocs,(void *) work,numtids,MPI_INT,MPI_MAX,comm); 1621eb62cbbSBarry Smith nmax = work[mytid]; 1631eb62cbbSBarry Smith FREE(work); 1641eb62cbbSBarry Smith 1651eb62cbbSBarry Smith /* post receives: 1661eb62cbbSBarry Smith 1) each message will consist of ordered pairs 1671eb62cbbSBarry Smith (global index,value) we store the global index as a double 168d6dfbf8fSBarry Smith to simplify the message passing. 1691eb62cbbSBarry Smith 2) since we don't know how long each individual message is we 1701eb62cbbSBarry Smith allocate the largest needed buffer for each receive. Potentially 1711eb62cbbSBarry Smith this is a lot of wasted space. 1721eb62cbbSBarry Smith 1731eb62cbbSBarry Smith 1741eb62cbbSBarry Smith This could be done better. 1751eb62cbbSBarry Smith */ 17628988994SBarry Smith rvalues = (Scalar *) MALLOC(3*(nreceives+1)*(nmax+1)*sizeof(Scalar)); 1771eb62cbbSBarry Smith CHKPTR(rvalues); 1781eb62cbbSBarry Smith recv_waits = (MPI_Request *) MALLOC((nreceives+1)*sizeof(MPI_Request)); 1791eb62cbbSBarry Smith CHKPTR(recv_waits); 1801eb62cbbSBarry Smith for ( i=0; i<nreceives; i++ ) { 1811eb62cbbSBarry Smith MPI_Irecv((void *)(rvalues+3*nmax*i),3*nmax,MPI_SCALAR,MPI_ANY_SOURCE,tag, 1821eb62cbbSBarry Smith comm,recv_waits+i); 1831eb62cbbSBarry Smith } 1841eb62cbbSBarry Smith 1851eb62cbbSBarry Smith /* do sends: 1861eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 1871eb62cbbSBarry Smith the ith processor 1881eb62cbbSBarry Smith */ 1891eb62cbbSBarry Smith svalues = (Scalar *) MALLOC( 3*(aij->stash.n+1)*sizeof(Scalar) ); 1901eb62cbbSBarry Smith CHKPTR(svalues); 1911eb62cbbSBarry Smith send_waits = (MPI_Request *) MALLOC( (nsends+1)*sizeof(MPI_Request)); 1921eb62cbbSBarry Smith CHKPTR(send_waits); 1931eb62cbbSBarry Smith starts = (int *) MALLOC( numtids*sizeof(int) ); CHKPTR(starts); 1941eb62cbbSBarry Smith starts[0] = 0; 1951eb62cbbSBarry Smith for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 1961eb62cbbSBarry Smith for ( i=0; i<aij->stash.n; i++ ) { 1971eb62cbbSBarry Smith svalues[3*starts[owner[i]]] = (Scalar) aij->stash.idx[i]; 1981eb62cbbSBarry Smith svalues[3*starts[owner[i]]+1] = (Scalar) aij->stash.idy[i]; 1991eb62cbbSBarry Smith svalues[3*(starts[owner[i]]++)+2] = aij->stash.array[i]; 2001eb62cbbSBarry Smith } 2011eb62cbbSBarry Smith FREE(owner); 2021eb62cbbSBarry Smith starts[0] = 0; 2031eb62cbbSBarry Smith for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 2041eb62cbbSBarry Smith count = 0; 2051eb62cbbSBarry Smith for ( i=0; i<numtids; i++ ) { 2061eb62cbbSBarry Smith if (procs[i]) { 2071eb62cbbSBarry Smith MPI_Isend((void*)(svalues+3*starts[i]),3*nprocs[i],MPI_SCALAR,i,tag, 2081eb62cbbSBarry Smith comm,send_waits+count++); 2091eb62cbbSBarry Smith } 2101eb62cbbSBarry Smith } 2111eb62cbbSBarry Smith FREE(starts); FREE(nprocs); 2121eb62cbbSBarry Smith 2131eb62cbbSBarry Smith /* Free cache space */ 2141eb62cbbSBarry Smith aij->stash.nmax = aij->stash.n = 0; 2151eb62cbbSBarry Smith if (aij->stash.array){ FREE(aij->stash.array); aij->stash.array = 0;} 2161eb62cbbSBarry Smith 2171eb62cbbSBarry Smith aij->svalues = svalues; aij->rvalues = rvalues; 2181eb62cbbSBarry Smith aij->nsends = nsends; aij->nrecvs = nreceives; 2191eb62cbbSBarry Smith aij->send_waits = send_waits; aij->recv_waits = recv_waits; 2201eb62cbbSBarry Smith aij->rmax = nmax; 2211eb62cbbSBarry Smith 2221eb62cbbSBarry Smith return 0; 2231eb62cbbSBarry Smith } 22444a69424SLois Curfman McInnes extern int MatSetUpMultiply_MPIAIJ(Mat); 2251eb62cbbSBarry Smith 22644a69424SLois Curfman McInnes static int MatEndAssemble_MPIAIJ(Mat mat) 2271eb62cbbSBarry Smith { 2281eb62cbbSBarry Smith int ierr; 22944a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 2301eb62cbbSBarry Smith 2311eb62cbbSBarry Smith MPI_Status *send_status,recv_status; 2326abc6512SBarry Smith int imdex,nrecvs = aij->nrecvs, count = nrecvs, i, n; 2331eb62cbbSBarry Smith int row,col; 2341eb62cbbSBarry Smith Scalar *values,val; 2351eb62cbbSBarry Smith InsertMode addv = aij->insertmode; 2361eb62cbbSBarry Smith 2371eb62cbbSBarry Smith /* wait on receives */ 2381eb62cbbSBarry Smith while (count) { 239d6dfbf8fSBarry Smith MPI_Waitany(nrecvs,aij->recv_waits,&imdex,&recv_status); 2401eb62cbbSBarry Smith /* unpack receives into our local space */ 241d6dfbf8fSBarry Smith values = aij->rvalues + 3*imdex*aij->rmax; 2421eb62cbbSBarry Smith MPI_Get_count(&recv_status,MPI_SCALAR,&n); 2431eb62cbbSBarry Smith n = n/3; 2441eb62cbbSBarry Smith for ( i=0; i<n; i++ ) { 2451eb62cbbSBarry Smith row = (int) PETSCREAL(values[3*i]) - aij->rstart; 2461eb62cbbSBarry Smith col = (int) PETSCREAL(values[3*i+1]); 2471eb62cbbSBarry Smith val = values[3*i+2]; 2481eb62cbbSBarry Smith if (col >= aij->cstart && col < aij->cend) { 2491eb62cbbSBarry Smith col -= aij->cstart; 2501eb62cbbSBarry Smith MatSetValues(aij->A,1,&row,1,&col,&val,addv); 2511eb62cbbSBarry Smith } 2521eb62cbbSBarry Smith else { 253d6dfbf8fSBarry Smith if (aij->assembled) { 2549e25ed09SBarry Smith if (!aij->colmap) {ierr = CreateColmap(mat); CHKERR(ierr);} 2559e25ed09SBarry Smith col = aij->colmap[col] - 1; 2569e25ed09SBarry Smith if (col < 0) { 2579e25ed09SBarry Smith SETERR(1,"Cannot insert new off diagonal block nonzero in\ 2589e25ed09SBarry Smith already\ 259d6dfbf8fSBarry Smith assembled matrix. Contact petsc-maint@mcs.anl.gov\ 260d6dfbf8fSBarry Smith if your need this feature"); 261d6dfbf8fSBarry Smith } 2629e25ed09SBarry Smith } 2631eb62cbbSBarry Smith MatSetValues(aij->B,1,&row,1,&col,&val,addv); 2641eb62cbbSBarry Smith } 2651eb62cbbSBarry Smith } 2661eb62cbbSBarry Smith count--; 2671eb62cbbSBarry Smith } 2681eb62cbbSBarry Smith FREE(aij->recv_waits); FREE(aij->rvalues); 2691eb62cbbSBarry Smith 2701eb62cbbSBarry Smith /* wait on sends */ 2711eb62cbbSBarry Smith if (aij->nsends) { 2721eb62cbbSBarry Smith send_status = (MPI_Status *) MALLOC( aij->nsends*sizeof(MPI_Status) ); 2731eb62cbbSBarry Smith CHKPTR(send_status); 2741eb62cbbSBarry Smith MPI_Waitall(aij->nsends,aij->send_waits,send_status); 2751eb62cbbSBarry Smith FREE(send_status); 2761eb62cbbSBarry Smith } 2771eb62cbbSBarry Smith FREE(aij->send_waits); FREE(aij->svalues); 2781eb62cbbSBarry Smith 2791eb62cbbSBarry Smith aij->insertmode = NotSetValues; 2801eb62cbbSBarry Smith ierr = MatBeginAssembly(aij->A); CHKERR(ierr); 2811eb62cbbSBarry Smith ierr = MatEndAssembly(aij->A); CHKERR(ierr); 2821eb62cbbSBarry Smith 2839e25ed09SBarry Smith if (!aij->assembled) { 284d3c9fed9SLois Curfman McInnes ierr = MatSetUpMultiply_MPIAIJ(mat); CHKERR(ierr); 2859e25ed09SBarry Smith } 2861eb62cbbSBarry Smith ierr = MatBeginAssembly(aij->B); CHKERR(ierr); 2871eb62cbbSBarry Smith ierr = MatEndAssembly(aij->B); CHKERR(ierr); 288d6dfbf8fSBarry Smith 289d6dfbf8fSBarry Smith aij->assembled = 1; 2908a729477SBarry Smith return 0; 2918a729477SBarry Smith } 2928a729477SBarry Smith 29344a69424SLois Curfman McInnes static int MatZero_MPIAIJ(Mat A) 2941eb62cbbSBarry Smith { 29544a69424SLois Curfman McInnes Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data; 2961eb62cbbSBarry Smith 2971eb62cbbSBarry Smith MatZeroEntries(l->A); MatZeroEntries(l->B); 2981eb62cbbSBarry Smith return 0; 2991eb62cbbSBarry Smith } 3001eb62cbbSBarry Smith 3011eb62cbbSBarry Smith /* again this uses the same basic stratagy as in the assembly and 3021eb62cbbSBarry Smith scatter create routines, we should try to do it systemamatically 3031eb62cbbSBarry Smith if we can figure out the proper level of generality. */ 3041eb62cbbSBarry Smith 3051eb62cbbSBarry Smith /* the code does not do the diagonal entries correctly unless the 3061eb62cbbSBarry Smith matrix is square and the column and row owerships are identical. 3071eb62cbbSBarry Smith This is a BUG. The only way to fix it seems to be to access 3081eb62cbbSBarry Smith aij->A and aij->B directly and not through the MatZeroRows() 3091eb62cbbSBarry Smith routine. 3101eb62cbbSBarry Smith */ 31144a69424SLois Curfman McInnes static int MatZeroRows_MPIAIJ(Mat A,IS is,Scalar *diag) 3121eb62cbbSBarry Smith { 31344a69424SLois Curfman McInnes Mat_MPIAIJ *l = (Mat_MPIAIJ *) A->data; 3141eb62cbbSBarry Smith int i,ierr,N, *rows,*owners = l->rowners,numtids = l->numtids; 3156abc6512SBarry Smith int *procs,*nprocs,j,found,idx,nsends,*work; 3161eb62cbbSBarry Smith int nmax,*svalues,*starts,*owner,nrecvs,mytid = l->mytid; 3176abc6512SBarry Smith int *rvalues,tag = 67,count,base,slen,n,*source; 318d6dfbf8fSBarry Smith int *lens,imdex,*lrows,*values; 319d6dfbf8fSBarry Smith MPI_Comm comm = A->comm; 3201eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 3211eb62cbbSBarry Smith MPI_Status recv_status,*send_status; 3221eb62cbbSBarry Smith IS istmp; 3231eb62cbbSBarry Smith 32444a69424SLois Curfman McInnes if (!l->assembled) SETERR(1,"MatZeroRows_MPIAIJ: must assemble matrix first"); 3251eb62cbbSBarry Smith ierr = ISGetLocalSize(is,&N); CHKERR(ierr); 3261eb62cbbSBarry Smith ierr = ISGetIndices(is,&rows); CHKERR(ierr); 3271eb62cbbSBarry Smith 3281eb62cbbSBarry Smith /* first count number of contributors to each processor */ 3291eb62cbbSBarry Smith nprocs = (int *) MALLOC( 2*numtids*sizeof(int) ); CHKPTR(nprocs); 3301eb62cbbSBarry Smith MEMSET(nprocs,0,2*numtids*sizeof(int)); procs = nprocs + numtids; 3311eb62cbbSBarry Smith owner = (int *) MALLOC((N+1)*sizeof(int)); CHKPTR(owner); /* see note*/ 3321eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 3331eb62cbbSBarry Smith idx = rows[i]; 3341eb62cbbSBarry Smith found = 0; 3351eb62cbbSBarry Smith for ( j=0; j<numtids; j++ ) { 3361eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 3371eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 3381eb62cbbSBarry Smith } 3391eb62cbbSBarry Smith } 340d6dfbf8fSBarry Smith if (!found) SETERR(1,"Imdex out of range"); 3411eb62cbbSBarry Smith } 3421eb62cbbSBarry Smith nsends = 0; for ( i=0; i<numtids; i++ ) { nsends += procs[i];} 3431eb62cbbSBarry Smith 3441eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 3451eb62cbbSBarry Smith work = (int *) MALLOC( numtids*sizeof(int) ); CHKPTR(work); 3461eb62cbbSBarry Smith MPI_Allreduce((void *) procs,(void *) work,numtids,MPI_INT,MPI_SUM,comm); 3471eb62cbbSBarry Smith nrecvs = work[mytid]; 3481eb62cbbSBarry Smith MPI_Allreduce((void *) nprocs,(void *) work,numtids,MPI_INT,MPI_MAX,comm); 3491eb62cbbSBarry Smith nmax = work[mytid]; 3501eb62cbbSBarry Smith FREE(work); 3511eb62cbbSBarry Smith 3521eb62cbbSBarry Smith /* post receives: */ 35328988994SBarry Smith rvalues = (int *) MALLOC((nrecvs+1)*(nmax+1)*sizeof(int)); /*see note */ 3541eb62cbbSBarry Smith CHKPTR(rvalues); 3551eb62cbbSBarry Smith recv_waits = (MPI_Request *) MALLOC((nrecvs+1)*sizeof(MPI_Request)); 3561eb62cbbSBarry Smith CHKPTR(recv_waits); 3571eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 3581eb62cbbSBarry Smith MPI_Irecv((void *)(rvalues+nmax*i),nmax,MPI_INT,MPI_ANY_SOURCE,tag, 3591eb62cbbSBarry Smith comm,recv_waits+i); 3601eb62cbbSBarry Smith } 3611eb62cbbSBarry Smith 3621eb62cbbSBarry Smith /* do sends: 3631eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 3641eb62cbbSBarry Smith the ith processor 3651eb62cbbSBarry Smith */ 3661eb62cbbSBarry Smith svalues = (int *) MALLOC( (N+1)*sizeof(int) ); CHKPTR(svalues); 3671eb62cbbSBarry Smith send_waits = (MPI_Request *) MALLOC( (nsends+1)*sizeof(MPI_Request)); 3681eb62cbbSBarry Smith CHKPTR(send_waits); 3691eb62cbbSBarry Smith starts = (int *) MALLOC( (numtids+1)*sizeof(int) ); CHKPTR(starts); 3701eb62cbbSBarry Smith starts[0] = 0; 3711eb62cbbSBarry Smith for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3721eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 3731eb62cbbSBarry Smith svalues[starts[owner[i]]++] = rows[i]; 3741eb62cbbSBarry Smith } 3751eb62cbbSBarry Smith ISRestoreIndices(is,&rows); 3761eb62cbbSBarry Smith 3771eb62cbbSBarry Smith starts[0] = 0; 3781eb62cbbSBarry Smith for ( i=1; i<numtids+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3791eb62cbbSBarry Smith count = 0; 3801eb62cbbSBarry Smith for ( i=0; i<numtids; i++ ) { 3811eb62cbbSBarry Smith if (procs[i]) { 3821eb62cbbSBarry Smith MPI_Isend((void*)(svalues+starts[i]),nprocs[i],MPI_INT,i,tag, 3831eb62cbbSBarry Smith comm,send_waits+count++); 3841eb62cbbSBarry Smith } 3851eb62cbbSBarry Smith } 3861eb62cbbSBarry Smith FREE(starts); 3871eb62cbbSBarry Smith 3881eb62cbbSBarry Smith base = owners[mytid]; 3891eb62cbbSBarry Smith 3901eb62cbbSBarry Smith /* wait on receives */ 3911eb62cbbSBarry Smith lens = (int *) MALLOC( 2*(nrecvs+1)*sizeof(int) ); CHKPTR(lens); 3921eb62cbbSBarry Smith source = lens + nrecvs; 3931eb62cbbSBarry Smith count = nrecvs; slen = 0; 3941eb62cbbSBarry Smith while (count) { 395d6dfbf8fSBarry Smith MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status); 3961eb62cbbSBarry Smith /* unpack receives into our local space */ 3971eb62cbbSBarry Smith MPI_Get_count(&recv_status,MPI_INT,&n); 398d6dfbf8fSBarry Smith source[imdex] = recv_status.MPI_SOURCE; 399d6dfbf8fSBarry Smith lens[imdex] = n; 4001eb62cbbSBarry Smith slen += n; 4011eb62cbbSBarry Smith count--; 4021eb62cbbSBarry Smith } 4031eb62cbbSBarry Smith FREE(recv_waits); 4041eb62cbbSBarry Smith 4051eb62cbbSBarry Smith /* move the data into the send scatter */ 4061eb62cbbSBarry Smith lrows = (int *) MALLOC( slen*sizeof(int) ); CHKPTR(lrows); 4071eb62cbbSBarry Smith count = 0; 4081eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 4091eb62cbbSBarry Smith values = rvalues + i*nmax; 4101eb62cbbSBarry Smith for ( j=0; j<lens[i]; j++ ) { 4111eb62cbbSBarry Smith lrows[count++] = values[j] - base; 4121eb62cbbSBarry Smith } 4131eb62cbbSBarry Smith } 4141eb62cbbSBarry Smith FREE(rvalues); FREE(lens); 4151eb62cbbSBarry Smith FREE(owner); FREE(nprocs); 4161eb62cbbSBarry Smith 4171eb62cbbSBarry Smith /* actually zap the local rows */ 4181eb62cbbSBarry Smith ierr = ISCreateSequential(slen,lrows,&istmp); CHKERR(ierr); FREE(lrows); 4191eb62cbbSBarry Smith ierr = MatZeroRows(l->A,istmp,diag); CHKERR(ierr); 4201eb62cbbSBarry Smith ierr = MatZeroRows(l->B,istmp,0); CHKERR(ierr); 4211eb62cbbSBarry Smith ierr = ISDestroy(istmp); CHKERR(ierr); 4221eb62cbbSBarry Smith 4231eb62cbbSBarry Smith /* wait on sends */ 4241eb62cbbSBarry Smith if (nsends) { 4251eb62cbbSBarry Smith send_status = (MPI_Status *) MALLOC( nsends*sizeof(MPI_Status) ); 4261eb62cbbSBarry Smith CHKPTR(send_status); 4271eb62cbbSBarry Smith MPI_Waitall(nsends,send_waits,send_status); 4281eb62cbbSBarry Smith FREE(send_status); 4291eb62cbbSBarry Smith } 4301eb62cbbSBarry Smith FREE(send_waits); FREE(svalues); 4311eb62cbbSBarry Smith 4321eb62cbbSBarry Smith 4331eb62cbbSBarry Smith return 0; 4341eb62cbbSBarry Smith } 4351eb62cbbSBarry Smith 43644a69424SLois Curfman McInnes static int MatMult_MPIAIJ(Mat aijin,Vec xx,Vec yy) 4371eb62cbbSBarry Smith { 43844a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data; 4391eb62cbbSBarry Smith int ierr; 44044a69424SLois Curfman McInnes if (!aij->assembled) SETERR(1,"MatMult_MPIAIJ: must assemble matrix first"); 441d6dfbf8fSBarry Smith ierr = VecScatterBegin(xx,0,aij->lvec,0,InsertValues,ScatterAll,aij->Mvctx); 4421eb62cbbSBarry Smith CHKERR(ierr); 4431eb62cbbSBarry Smith ierr = MatMult(aij->A,xx,yy); CHKERR(ierr); 444d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,aij->lvec,0,InsertValues,ScatterAll,aij->Mvctx); 445d6dfbf8fSBarry Smith CHKERR(ierr); 4461eb62cbbSBarry Smith ierr = MatMultAdd(aij->B,aij->lvec,yy,yy); CHKERR(ierr); 4471eb62cbbSBarry Smith return 0; 4481eb62cbbSBarry Smith } 4491eb62cbbSBarry Smith 45044a69424SLois Curfman McInnes static int MatMultAdd_MPIAIJ(Mat aijin,Vec xx,Vec yy,Vec zz) 451da3a660dSBarry Smith { 45244a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data; 453da3a660dSBarry Smith int ierr; 45444a69424SLois Curfman McInnes if (!aij->assembled) SETERR(1,"MatMult_MPIAIJ: must assemble matrix first"); 455da3a660dSBarry Smith ierr = VecScatterBegin(xx,0,aij->lvec,0,InsertValues,ScatterAll,aij->Mvctx); 456da3a660dSBarry Smith CHKERR(ierr); 457da3a660dSBarry Smith ierr = MatMultAdd(aij->A,xx,yy,zz); CHKERR(ierr); 458da3a660dSBarry Smith ierr = VecScatterEnd(xx,0,aij->lvec,0,InsertValues,ScatterAll,aij->Mvctx); 459da3a660dSBarry Smith CHKERR(ierr); 460da3a660dSBarry Smith ierr = MatMultAdd(aij->B,aij->lvec,zz,zz); CHKERR(ierr); 461da3a660dSBarry Smith return 0; 462da3a660dSBarry Smith } 463da3a660dSBarry Smith 46444a69424SLois Curfman McInnes static int MatMultTrans_MPIAIJ(Mat aijin,Vec xx,Vec yy) 465da3a660dSBarry Smith { 46644a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data; 467da3a660dSBarry Smith int ierr; 468da3a660dSBarry Smith 46944a69424SLois Curfman McInnes if (!aij->assembled) 47044a69424SLois Curfman McInnes SETERR(1,"MatMulTrans_MPIAIJ: must assemble matrix first"); 471da3a660dSBarry Smith /* do nondiagonal part */ 472da3a660dSBarry Smith ierr = MatMultTrans(aij->B,xx,aij->lvec); CHKERR(ierr); 473da3a660dSBarry Smith /* send it on its way */ 474da3a660dSBarry Smith ierr = VecScatterBegin(aij->lvec,0,yy,0,AddValues, 475da3a660dSBarry Smith ScatterAll|ScatterReverse,aij->Mvctx); CHKERR(ierr); 476da3a660dSBarry Smith /* do local part */ 477da3a660dSBarry Smith ierr = MatMultTrans(aij->A,xx,yy); CHKERR(ierr); 478da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 479da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 480da3a660dSBarry Smith /* but is not perhaps always true. */ 481da3a660dSBarry Smith ierr = VecScatterEnd(aij->lvec,0,yy,0,AddValues,ScatterAll|ScatterReverse, 482da3a660dSBarry Smith aij->Mvctx); CHKERR(ierr); 483da3a660dSBarry Smith return 0; 484da3a660dSBarry Smith } 485da3a660dSBarry Smith 48644a69424SLois Curfman McInnes static int MatMultTransAdd_MPIAIJ(Mat aijin,Vec xx,Vec yy,Vec zz) 487da3a660dSBarry Smith { 48844a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data; 489da3a660dSBarry Smith int ierr; 490da3a660dSBarry Smith 49144a69424SLois Curfman McInnes if (!aij->assembled) 49244a69424SLois Curfman McInnes SETERR(1,"MatMulTransAdd_MPIAIJ: must assemble matrix first"); 493da3a660dSBarry Smith /* do nondiagonal part */ 494da3a660dSBarry Smith ierr = MatMultTrans(aij->B,xx,aij->lvec); CHKERR(ierr); 495da3a660dSBarry Smith /* send it on its way */ 496da3a660dSBarry Smith ierr = VecScatterBegin(aij->lvec,0,zz,0,AddValues, 497da3a660dSBarry Smith ScatterAll|ScatterReverse,aij->Mvctx); CHKERR(ierr); 498da3a660dSBarry Smith /* do local part */ 499da3a660dSBarry Smith ierr = MatMultTransAdd(aij->A,xx,yy,zz); CHKERR(ierr); 500da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 501da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 502da3a660dSBarry Smith /* but is not perhaps always true. */ 503da3a660dSBarry Smith ierr = VecScatterEnd(aij->lvec,0,zz,0,AddValues,ScatterAll|ScatterReverse, 504da3a660dSBarry Smith aij->Mvctx); CHKERR(ierr); 505da3a660dSBarry Smith return 0; 506da3a660dSBarry Smith } 507da3a660dSBarry Smith 5081eb62cbbSBarry Smith /* 5091eb62cbbSBarry Smith This only works correctly for square matrices where the subblock A->A is the 5101eb62cbbSBarry Smith diagonal block 5111eb62cbbSBarry Smith */ 51244a69424SLois Curfman McInnes static int MatGetDiag_MPIAIJ(Mat Ain,Vec v) 5131eb62cbbSBarry Smith { 51444a69424SLois Curfman McInnes Mat_MPIAIJ *A = (Mat_MPIAIJ *) Ain->data; 51544a69424SLois Curfman McInnes if (!A->assembled) SETERR(1,"MatGetDiag_MPIAIJ: must assemble matrix first"); 5161eb62cbbSBarry Smith return MatGetDiagonal(A->A,v); 5171eb62cbbSBarry Smith } 5181eb62cbbSBarry Smith 51944a69424SLois Curfman McInnes static int MatDestroy_MPIAIJ(PetscObject obj) 5201eb62cbbSBarry Smith { 5211eb62cbbSBarry Smith Mat mat = (Mat) obj; 52244a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 5231eb62cbbSBarry Smith int ierr; 524a5a9c739SBarry Smith #if defined(PETSC_LOG) 525a5a9c739SBarry Smith PLogObjectState(obj,"Rows %d Cols %d",aij->M,aij->N); 526a5a9c739SBarry Smith #endif 5271eb62cbbSBarry Smith FREE(aij->rowners); 5281eb62cbbSBarry Smith ierr = MatDestroy(aij->A); CHKERR(ierr); 5291eb62cbbSBarry Smith ierr = MatDestroy(aij->B); CHKERR(ierr); 5309e25ed09SBarry Smith if (aij->colmap) FREE(aij->colmap); 5319e25ed09SBarry Smith if (aij->garray) FREE(aij->garray); 5321eb62cbbSBarry Smith if (aij->lvec) VecDestroy(aij->lvec); 5331eb62cbbSBarry Smith if (aij->Mvctx) VecScatterCtxDestroy(aij->Mvctx); 534a5a9c739SBarry Smith FREE(aij); 535a5a9c739SBarry Smith PLogObjectDestroy(mat); 536a5a9c739SBarry Smith PETSCHEADERDESTROY(mat); 5371eb62cbbSBarry Smith return 0; 5381eb62cbbSBarry Smith } 5391eb62cbbSBarry Smith 54044a69424SLois Curfman McInnes static int MatView_MPIAIJ(PetscObject obj,Viewer viewer) 5411eb62cbbSBarry Smith { 5421eb62cbbSBarry Smith Mat mat = (Mat) obj; 54344a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data; 5441eb62cbbSBarry Smith int ierr; 545d13ab20cSBarry Smith PetscObject vobj = (PetscObject) viewer; 5461eb62cbbSBarry Smith 54744a69424SLois Curfman McInnes if (!aij->assembled) SETERR(1,"MatView_MPIAIJ: must assemble matrix first"); 548d13ab20cSBarry Smith if (vobj->cookie == VIEWER_COOKIE) { 549d13ab20cSBarry Smith FILE *fd = ViewerFileGetPointer(viewer); 550d13ab20cSBarry Smith if (vobj->type == FILE_VIEWER) { 551d6dfbf8fSBarry Smith MPE_Seq_begin(mat->comm,1); 552d13ab20cSBarry Smith fprintf(fd,"[%d] rows %d starts %d ends %d cols %d starts %d ends %d\n", 5531eb62cbbSBarry Smith aij->mytid,aij->m,aij->rstart,aij->rend,aij->n,aij->cstart, 5541eb62cbbSBarry Smith aij->cend); 555da69df5fSBarry Smith ierr = MatView(aij->A,viewer); CHKERR(ierr); 556da69df5fSBarry Smith ierr = MatView(aij->B,viewer); CHKERR(ierr); 557d13ab20cSBarry Smith fflush(fd); 558d6dfbf8fSBarry Smith MPE_Seq_end(mat->comm,1); 559d13ab20cSBarry Smith } 56095373324SBarry Smith else if (vobj->type == FILES_VIEWER) { 56195373324SBarry Smith int numtids = aij->numtids, mytid = aij->mytid; 56295373324SBarry Smith if (numtids == 1) { 56395373324SBarry Smith ierr = MatView(aij->A,viewer); CHKERR(ierr); 56495373324SBarry Smith } 56595373324SBarry Smith else { 56695373324SBarry Smith /* assemble the entire matrix onto first processor. */ 56795373324SBarry Smith Mat A; 56895373324SBarry Smith Mat_AIJ *Aaij; 569*2eb8c8abSBarry Smith int M = aij->M, N = aij->N,m,*ai,*aj,row,*cols,i,*ct; 57095373324SBarry Smith Scalar *a; 57195373324SBarry Smith if (!mytid) { 57295373324SBarry Smith ierr = MatCreateMPIAIJ(mat->comm,M,N,M,N,0,0,0,0,&A); 57395373324SBarry Smith } 57495373324SBarry Smith else { 57595373324SBarry Smith ierr = MatCreateMPIAIJ(mat->comm,0,0,M,N,0,0,0,0,&A); 57695373324SBarry Smith } 57795373324SBarry Smith CHKERR(ierr); 57895373324SBarry Smith 57995373324SBarry Smith /* copy over the A part */ 58095373324SBarry Smith Aaij = (Mat_AIJ*) aij->A->data; 581*2eb8c8abSBarry Smith m = Aaij->m; ai = Aaij->i; aj = Aaij->j; a = Aaij->a; 58295373324SBarry Smith row = aij->rstart; 58395373324SBarry Smith for ( i=0; i<ai[m]; i++ ) {aj[i] += aij->cstart - 1;} 58495373324SBarry Smith for ( i=0; i<m; i++ ) { 58595373324SBarry Smith ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],aj,a,InsertValues); 58695373324SBarry Smith CHKERR(ierr); 58795373324SBarry Smith row++; a += ai[i+1]-ai[i]; aj += ai[i+1]-ai[i]; 58895373324SBarry Smith } 58995373324SBarry Smith aj = Aaij->j; 59095373324SBarry Smith for ( i=0; i<ai[m]; i++ ) {aj[i] -= aij->cstart - 1;} 59195373324SBarry Smith 59295373324SBarry Smith /* copy over the B part */ 59395373324SBarry Smith Aaij = (Mat_AIJ*) aij->B->data; 594*2eb8c8abSBarry Smith m = Aaij->m; ai = Aaij->i; aj = Aaij->j; a = Aaij->a; 59595373324SBarry Smith row = aij->rstart; 59695373324SBarry Smith ct = cols = (int *) MALLOC( (ai[m]+1)*sizeof(int) ); CHKPTR(cols); 59795373324SBarry Smith for ( i=0; i<ai[m]; i++ ) {cols[i] = aij->garray[aj[i]-1];} 59895373324SBarry Smith for ( i=0; i<m; i++ ) { 59995373324SBarry Smith ierr = MatSetValues(A,1,&row,ai[i+1]-ai[i],cols,a,InsertValues); 60095373324SBarry Smith CHKERR(ierr); 60195373324SBarry Smith row++; a += ai[i+1]-ai[i]; cols += ai[i+1]-ai[i]; 60295373324SBarry Smith } 60395373324SBarry Smith FREE(ct); 60495373324SBarry Smith 60595373324SBarry Smith ierr = MatBeginAssembly(A); CHKERR(ierr); 60695373324SBarry Smith ierr = MatEndAssembly(A); CHKERR(ierr); 60795373324SBarry Smith if (!mytid) { 60895373324SBarry Smith ierr = MatView(((Mat_MPIAIJ*)(A->data))->A,viewer); CHKERR(ierr); 60995373324SBarry Smith } 61095373324SBarry Smith ierr = MatDestroy(A); CHKERR(ierr); 61195373324SBarry Smith } 61295373324SBarry Smith } 613d13ab20cSBarry Smith } 6141eb62cbbSBarry Smith return 0; 6151eb62cbbSBarry Smith } 6161eb62cbbSBarry Smith 617d3c9fed9SLois Curfman McInnes extern int MatMarkDiag_AIJ(Mat_AIJ *); 6181eb62cbbSBarry Smith /* 6191eb62cbbSBarry Smith This has to provide several versions. 6201eb62cbbSBarry Smith 6211eb62cbbSBarry Smith 1) per sequential 6221eb62cbbSBarry Smith 2) a) use only local smoothing updating outer values only once. 6231eb62cbbSBarry Smith b) local smoothing updating outer values each inner iteration 624d6dfbf8fSBarry Smith 3) color updating out values betwen colors. 6251eb62cbbSBarry Smith */ 62644a69424SLois Curfman McInnes static int MatRelax_MPIAIJ(Mat matin,Vec bb,double omega,int flag,double shift, 6278a729477SBarry Smith int its,Vec xx) 6288a729477SBarry Smith { 62944a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 630d6dfbf8fSBarry Smith Mat AA = mat->A, BB = mat->B; 63144a69424SLois Curfman McInnes Mat_AIJ *A = (Mat_AIJ *) AA->data, *B = (Mat_AIJ *)BB->data; 6326abc6512SBarry Smith Scalar zero = 0.0,*b,*x,*xs,*ls,d,*v,sum,scale,*t,*ts; 6336abc6512SBarry Smith int ierr,*idx, *diag; 6346abc6512SBarry Smith int n = mat->n, m = mat->m, i; 635da3a660dSBarry Smith Vec tt; 6368a729477SBarry Smith 63744a69424SLois Curfman McInnes if (!mat->assembled) SETERR(1,"MatRelax_MPIAIJ: must assemble matrix first"); 6381eb62cbbSBarry Smith 639d6dfbf8fSBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); VecGetArray(mat->lvec,&ls); 640d6dfbf8fSBarry Smith xs = x -1; /* shift by one for index start of 1 */ 641da3a660dSBarry Smith ls--; 642d3c9fed9SLois Curfman McInnes if (!A->diag) {if ((ierr = MatMarkDiag_AIJ(A))) return ierr;} 643d6dfbf8fSBarry Smith diag = A->diag; 644acb40c82SBarry Smith if (flag == SOR_APPLY_UPPER || flag == SOR_APPLY_LOWER) { 645acb40c82SBarry Smith SETERR(1,"That option not yet support for parallel AIJ matrices"); 646acb40c82SBarry Smith } 647da3a660dSBarry Smith if (flag & SOR_EISENSTAT) { 648da3a660dSBarry Smith /* Let A = L + U + D; where L is lower trianglar, 649da3a660dSBarry Smith U is upper triangular, E is diagonal; This routine applies 650da3a660dSBarry Smith 651da3a660dSBarry Smith (L + E)^{-1} A (U + E)^{-1} 652da3a660dSBarry Smith 653da3a660dSBarry Smith to a vector efficiently using Eisenstat's trick. This is for 654da3a660dSBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 655da3a660dSBarry Smith is the relaxation factor. 656da3a660dSBarry Smith */ 657da3a660dSBarry Smith ierr = VecCreate(xx,&tt); CHKERR(ierr); 658da3a660dSBarry Smith VecGetArray(tt,&t); 659da3a660dSBarry Smith scale = (2.0/omega) - 1.0; 660da3a660dSBarry Smith /* x = (E + U)^{-1} b */ 661da3a660dSBarry Smith VecSet(&zero,mat->lvec); 662da3a660dSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineUp, 663da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 664da3a660dSBarry Smith for ( i=m-1; i>-1; i-- ) { 665da3a660dSBarry Smith n = A->i[i+1] - diag[i] - 1; 666da3a660dSBarry Smith idx = A->j + diag[i]; 667da3a660dSBarry Smith v = A->a + diag[i]; 668da3a660dSBarry Smith sum = b[i]; 669da3a660dSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 670da3a660dSBarry Smith d = shift + A->a[diag[i]-1]; 671da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 672da3a660dSBarry Smith idx = B->j + B->i[i] - 1; 673da3a660dSBarry Smith v = B->a + B->i[i] - 1; 674da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 675da3a660dSBarry Smith x[i] = omega*(sum/d); 676da3a660dSBarry Smith } 677da3a660dSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineUp, 678da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 679da3a660dSBarry Smith 680da3a660dSBarry Smith /* t = b - (2*E - D)x */ 681da3a660dSBarry Smith v = A->a; 682da3a660dSBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ - 1])*x[i]; } 683da3a660dSBarry Smith 684da3a660dSBarry Smith /* t = (E + L)^{-1}t */ 685da3a660dSBarry Smith ts = t - 1; /* shifted by one for index start of a or mat->j*/ 686da3a660dSBarry Smith diag = A->diag; 687da3a660dSBarry Smith VecSet(&zero,mat->lvec); 688da3a660dSBarry Smith ierr = VecPipelineBegin(tt,0,mat->lvec,0,InsertValues,PipelineDown, 689da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 690da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 691da3a660dSBarry Smith n = diag[i] - A->i[i]; 692da3a660dSBarry Smith idx = A->j + A->i[i] - 1; 693da3a660dSBarry Smith v = A->a + A->i[i] - 1; 694da3a660dSBarry Smith sum = t[i]; 695da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 696da3a660dSBarry Smith d = shift + A->a[diag[i]-1]; 697da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 698da3a660dSBarry Smith idx = B->j + B->i[i] - 1; 699da3a660dSBarry Smith v = B->a + B->i[i] - 1; 700da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 701da3a660dSBarry Smith t[i] = omega*(sum/d); 702da3a660dSBarry Smith } 703da3a660dSBarry Smith ierr = VecPipelineEnd(tt,0,mat->lvec,0,InsertValues,PipelineDown, 704da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 705da3a660dSBarry Smith /* x = x + t */ 706da3a660dSBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 707da3a660dSBarry Smith VecDestroy(tt); 708da3a660dSBarry Smith return 0; 709da3a660dSBarry Smith } 710da3a660dSBarry Smith 7111eb62cbbSBarry Smith 712d6dfbf8fSBarry Smith if ((flag & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP){ 713da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 714da3a660dSBarry Smith VecSet(&zero,mat->lvec); VecSet(&zero,xx); 715da3a660dSBarry Smith } 716da3a660dSBarry Smith else { 717d6dfbf8fSBarry Smith ierr=VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterUp,mat->Mvctx); 718d6dfbf8fSBarry Smith CHKERR(ierr); 719d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterUp,mat->Mvctx); 720d6dfbf8fSBarry Smith CHKERR(ierr); 721da3a660dSBarry Smith } 722d6dfbf8fSBarry Smith while (its--) { 723d6dfbf8fSBarry Smith /* go down through the rows */ 724d6dfbf8fSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineDown, 725d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 726d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 727d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 728d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 729d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 730d6dfbf8fSBarry Smith sum = b[i]; 731d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 732d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 733d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 734d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 735d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 736d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 737d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 738d6dfbf8fSBarry Smith } 739d6dfbf8fSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineDown, 740d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 741d6dfbf8fSBarry Smith /* come up through the rows */ 742d6dfbf8fSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineUp, 743d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 744d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 745d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 746d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 747d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 748d6dfbf8fSBarry Smith sum = b[i]; 749d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 750d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 751d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 752d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 753d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 754d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 755d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 756d6dfbf8fSBarry Smith } 757d6dfbf8fSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineUp, 758d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 759d6dfbf8fSBarry Smith } 760d6dfbf8fSBarry Smith } 761d6dfbf8fSBarry Smith else if (flag & SOR_FORWARD_SWEEP){ 762da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 763da3a660dSBarry Smith VecSet(&zero,mat->lvec); 764da3a660dSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineDown, 765da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 766da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 767da3a660dSBarry Smith n = diag[i] - A->i[i]; 768da3a660dSBarry Smith idx = A->j + A->i[i] - 1; 769da3a660dSBarry Smith v = A->a + A->i[i] - 1; 770da3a660dSBarry Smith sum = b[i]; 771da3a660dSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 772da3a660dSBarry Smith d = shift + A->a[diag[i]-1]; 773da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 774da3a660dSBarry Smith idx = B->j + B->i[i] - 1; 775da3a660dSBarry Smith v = B->a + B->i[i] - 1; 776da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 777da3a660dSBarry Smith x[i] = omega*(sum/d); 778da3a660dSBarry Smith } 779da3a660dSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineDown, 780da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 781da3a660dSBarry Smith its--; 782da3a660dSBarry Smith } 783d6dfbf8fSBarry Smith while (its--) { 784d6dfbf8fSBarry Smith ierr=VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterUp,mat->Mvctx); 785d6dfbf8fSBarry Smith CHKERR(ierr); 786d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterUp,mat->Mvctx); 787d6dfbf8fSBarry Smith CHKERR(ierr); 788d6dfbf8fSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineDown, 789d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 790d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 791d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 792d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 793d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 794d6dfbf8fSBarry Smith sum = b[i]; 795d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 796d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 797d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 798d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 799d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 800d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 801d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 802d6dfbf8fSBarry Smith } 803d6dfbf8fSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineDown, 804d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 805d6dfbf8fSBarry Smith } 806d6dfbf8fSBarry Smith } 807d6dfbf8fSBarry Smith else if (flag & SOR_BACKWARD_SWEEP){ 808da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 809da3a660dSBarry Smith VecSet(&zero,mat->lvec); 810da3a660dSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineUp, 811da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 812da3a660dSBarry Smith for ( i=m-1; i>-1; i-- ) { 813da3a660dSBarry Smith n = A->i[i+1] - diag[i] - 1; 814da3a660dSBarry Smith idx = A->j + diag[i]; 815da3a660dSBarry Smith v = A->a + diag[i]; 816da3a660dSBarry Smith sum = b[i]; 817da3a660dSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 818da3a660dSBarry Smith d = shift + A->a[diag[i]-1]; 819da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 820da3a660dSBarry Smith idx = B->j + B->i[i] - 1; 821da3a660dSBarry Smith v = B->a + B->i[i] - 1; 822da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 823da3a660dSBarry Smith x[i] = omega*(sum/d); 824da3a660dSBarry Smith } 825da3a660dSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineUp, 826da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 827da3a660dSBarry Smith its--; 828da3a660dSBarry Smith } 829d6dfbf8fSBarry Smith while (its--) { 830d6dfbf8fSBarry Smith ierr = VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterDown, 831d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 832d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterDown, 833d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 834d6dfbf8fSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineUp, 835d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 836d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 837d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 838d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 839d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 840d6dfbf8fSBarry Smith sum = b[i]; 841d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 842d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 843d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 844d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 845d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 846d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 847d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 848d6dfbf8fSBarry Smith } 849d6dfbf8fSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineUp, 850d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 851d6dfbf8fSBarry Smith } 852d6dfbf8fSBarry Smith } 853d6dfbf8fSBarry Smith else if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){ 854da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 855da3a660dSBarry Smith return MatRelax(mat->A,bb,omega,flag,shift,its,xx); 856da3a660dSBarry Smith } 857d6dfbf8fSBarry Smith ierr=VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterAll,mat->Mvctx); 858d6dfbf8fSBarry Smith CHKERR(ierr); 859d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterAll,mat->Mvctx); 860d6dfbf8fSBarry Smith CHKERR(ierr); 861d6dfbf8fSBarry Smith while (its--) { 862d6dfbf8fSBarry Smith /* go down through the rows */ 863d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 864d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 865d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 866d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 867d6dfbf8fSBarry Smith sum = b[i]; 868d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 869d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 870d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 871d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 872d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 873d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 874d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 875d6dfbf8fSBarry Smith } 876d6dfbf8fSBarry Smith /* come up through the rows */ 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 } 892d6dfbf8fSBarry Smith else if (flag & SOR_LOCAL_FORWARD_SWEEP){ 893da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 894da3a660dSBarry Smith return MatRelax(mat->A,bb,omega,flag,shift,its,xx); 895da3a660dSBarry Smith } 896d6dfbf8fSBarry Smith ierr=VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterAll,mat->Mvctx); 897d6dfbf8fSBarry Smith CHKERR(ierr); 898d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterAll,mat->Mvctx); 899d6dfbf8fSBarry Smith CHKERR(ierr); 900d6dfbf8fSBarry Smith while (its--) { 901d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 902d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 903d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 904d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 905d6dfbf8fSBarry Smith sum = b[i]; 906d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 907d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 908d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 909d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 910d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 911d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 912d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 913d6dfbf8fSBarry Smith } 914d6dfbf8fSBarry Smith } 915d6dfbf8fSBarry Smith } 916d6dfbf8fSBarry Smith else if (flag & SOR_LOCAL_BACKWARD_SWEEP){ 917da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 918da3a660dSBarry Smith return MatRelax(mat->A,bb,omega,flag,shift,its,xx); 919da3a660dSBarry Smith } 920d6dfbf8fSBarry Smith ierr = VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterAll, 921d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 922d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterAll, 923d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 924d6dfbf8fSBarry Smith while (its--) { 925d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 926d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 927d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 928d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 929d6dfbf8fSBarry Smith sum = b[i]; 930d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 931d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 932d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 933d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 934d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 935d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 936d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 937d6dfbf8fSBarry Smith } 938d6dfbf8fSBarry Smith } 939d6dfbf8fSBarry Smith } 9408a729477SBarry Smith return 0; 9418a729477SBarry Smith } 94244a69424SLois Curfman McInnes static int MatInsOpt_MPIAIJ(Mat aijin,int op) 943c74985f6SBarry Smith { 94444a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) aijin->data; 945c74985f6SBarry Smith 946c74985f6SBarry Smith if (op == NO_NEW_NONZERO_LOCATIONS) { 947c74985f6SBarry Smith MatSetOption(aij->A,op); 948c74985f6SBarry Smith MatSetOption(aij->B,op); 949c74985f6SBarry Smith } 950c74985f6SBarry Smith else if (op == YES_NEW_NONZERO_LOCATIONS) { 951c74985f6SBarry Smith MatSetOption(aij->A,op); 952c74985f6SBarry Smith MatSetOption(aij->B,op); 953c74985f6SBarry Smith } 954c74985f6SBarry Smith else if (op == COLUMN_ORIENTED) SETERR(1,"Column oriented not supported"); 955c74985f6SBarry Smith return 0; 956c74985f6SBarry Smith } 957c74985f6SBarry Smith 95844a69424SLois Curfman McInnes static int MatSize_MPIAIJ(Mat matin,int *m,int *n) 959c74985f6SBarry Smith { 96044a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 961c74985f6SBarry Smith *m = mat->M; *n = mat->N; 962c74985f6SBarry Smith return 0; 963c74985f6SBarry Smith } 964c74985f6SBarry Smith 96544a69424SLois Curfman McInnes static int MatLocalSize_MPIAIJ(Mat matin,int *m,int *n) 966c74985f6SBarry Smith { 96744a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 968c74985f6SBarry Smith *m = mat->m; *n = mat->n; 969c74985f6SBarry Smith return 0; 970c74985f6SBarry Smith } 971c74985f6SBarry Smith 97244a69424SLois Curfman McInnes static int MatRange_MPIAIJ(Mat matin,int *m,int *n) 973c74985f6SBarry Smith { 97444a69424SLois Curfman McInnes Mat_MPIAIJ *mat = (Mat_MPIAIJ *) matin->data; 975c74985f6SBarry Smith *m = mat->rstart; *n = mat->rend; 976c74985f6SBarry Smith return 0; 977c74985f6SBarry Smith } 978c74985f6SBarry Smith 97939e00950SLois Curfman McInnes static int MatGetRow_MPIAIJ(Mat matin,int row,int *nz,int **idx,Scalar **v) 98039e00950SLois Curfman McInnes { 98139e00950SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ *) matin->data; 98239e00950SLois Curfman McInnes Scalar *vworkA, *vworkB; 983abc0e9e4SLois Curfman McInnes int ierr, *cworkA, *cworkB, lrow, cstart = aij->cstart; 98439e00950SLois Curfman McInnes int nztot, nzA, nzB, i, rstart = aij->rstart, rend = aij->rend; 98539e00950SLois Curfman McInnes 98639e00950SLois Curfman McInnes if (!aij->assembled) 98739e00950SLois Curfman McInnes SETERR(1,"MatGetRow_MPIAIJ: Must assemble matrix first."); 98839e00950SLois Curfman McInnes if (row < rstart || row >= rend) 989abc0e9e4SLois Curfman McInnes SETERR(1,"MatGetRow_MPIAIJ: Currently you can get only local rows.") 990abc0e9e4SLois Curfman McInnes lrow = row - rstart; 991abc0e9e4SLois Curfman McInnes ierr = MatGetRow(aij->A,lrow,&nzA,&cworkA,&vworkA); CHKERR(ierr); 99239e00950SLois Curfman McInnes for (i=0; i<nzA; i++) cworkA[i] += cstart; 993abc0e9e4SLois Curfman McInnes ierr = MatGetRow(aij->B,lrow,&nzB,&cworkB,&vworkB); CHKERR(ierr); 994abc0e9e4SLois Curfman McInnes for (i=0; i<nzB; i++) cworkB[i] = aij->garray[cworkB[i]]; 99539e00950SLois Curfman McInnes 996*2eb8c8abSBarry Smith if ((nztot = nzA + nzB)) { 99739e00950SLois Curfman McInnes *idx = (int *) MALLOC( (nztot)*sizeof(int) ); CHKPTR(*idx); 99839e00950SLois Curfman McInnes *v = (Scalar *) MALLOC( (nztot)*sizeof(Scalar) ); CHKPTR(*v); 99939e00950SLois Curfman McInnes for ( i=0; i<nzA; i++ ) { 100039e00950SLois Curfman McInnes (*idx)[i] = cworkA[i]; 100139e00950SLois Curfman McInnes (*v)[i] = vworkA[i]; 100239e00950SLois Curfman McInnes } 100339e00950SLois Curfman McInnes for ( i=0; i<nzB; i++ ) { 100439e00950SLois Curfman McInnes (*idx)[i+nzA] = cworkB[i]; 1005abc0e9e4SLois Curfman McInnes (*v)[i+nzA] = vworkB[i]; 100639e00950SLois Curfman McInnes } 100739e00950SLois Curfman McInnes } 100839e00950SLois Curfman McInnes else {*idx = 0; *v=0;} 100939e00950SLois Curfman McInnes *nz = nztot; 1010abc0e9e4SLois Curfman McInnes ierr = MatRestoreRow(aij->A,lrow,&nzA,&cworkA,&vworkA); CHKERR(ierr); 1011abc0e9e4SLois Curfman McInnes ierr = MatRestoreRow(aij->B,lrow,&nzB,&cworkB,&vworkB); CHKERR(ierr); 101239e00950SLois Curfman McInnes return 0; 101339e00950SLois Curfman McInnes } 101439e00950SLois Curfman McInnes 101539e00950SLois Curfman McInnes static int MatRestoreRow_MPIAIJ(Mat mat,int row,int *nz,int **idx,Scalar **v) 101639e00950SLois Curfman McInnes { 101739e00950SLois Curfman McInnes if (*idx) FREE(*idx); 101839e00950SLois Curfman McInnes if (*v) FREE(*v); 101939e00950SLois Curfman McInnes return 0; 102039e00950SLois Curfman McInnes } 102139e00950SLois Curfman McInnes 102244a69424SLois Curfman McInnes static int MatCopy_MPIAIJ(Mat,Mat *); 102344a69424SLois Curfman McInnes extern int MatConvert_MPIAIJ(Mat,MATTYPE,Mat *); 1024d6dfbf8fSBarry Smith 10258a729477SBarry Smith /* -------------------------------------------------------------------*/ 102644a69424SLois Curfman McInnes static struct _MatOps MatOps = {MatInsertValues_MPIAIJ, 102739e00950SLois Curfman McInnes MatGetRow_MPIAIJ,MatRestoreRow_MPIAIJ, 102844a69424SLois Curfman McInnes MatMult_MPIAIJ,MatMultAdd_MPIAIJ, 102944a69424SLois Curfman McInnes MatMultTrans_MPIAIJ,MatMultTransAdd_MPIAIJ, 10301eb62cbbSBarry Smith 0,0,0,0, 10318a729477SBarry Smith 0,0, 103244a69424SLois Curfman McInnes MatRelax_MPIAIJ, 10338a729477SBarry Smith 0, 10341eb62cbbSBarry Smith 0,0,0, 103544a69424SLois Curfman McInnes MatCopy_MPIAIJ, 103644a69424SLois Curfman McInnes MatGetDiag_MPIAIJ,0,0, 103744a69424SLois Curfman McInnes MatBeginAssemble_MPIAIJ,MatEndAssemble_MPIAIJ, 10381eb62cbbSBarry Smith 0, 103944a69424SLois Curfman McInnes MatInsOpt_MPIAIJ,MatZero_MPIAIJ,MatZeroRows_MPIAIJ,0, 1040c74985f6SBarry Smith 0,0,0,0, 104144a69424SLois Curfman McInnes MatSize_MPIAIJ,MatLocalSize_MPIAIJ,MatRange_MPIAIJ, 104277915d1cSLois Curfman McInnes 0,0, 104344a69424SLois Curfman McInnes 0,MatConvert_MPIAIJ }; 10448a729477SBarry Smith 10458a729477SBarry Smith /*@ 10468a729477SBarry Smith 10471eb62cbbSBarry Smith MatCreateMPIAIJ - Creates a sparse parallel matrix 10481eb62cbbSBarry Smith in AIJ format. 10498a729477SBarry Smith 10508a729477SBarry Smith Input Parameters: 10511eb62cbbSBarry Smith . comm - MPI communicator 10521eb62cbbSBarry Smith . m,n - number of local rows and columns (or -1 to have calculated) 10531eb62cbbSBarry Smith . M,N - global rows and columns (or -1 to have calculated) 10541eb62cbbSBarry Smith . d_nz - total number nonzeros in diagonal portion of matrix 10551eb62cbbSBarry Smith . d_nzz - number of nonzeros per row in diagonal portion of matrix or null 10568a729477SBarry Smith . You must leave room for the diagonal entry even if it is zero. 10571eb62cbbSBarry Smith . o_nz - total number nonzeros in off-diagonal portion of matrix 10581eb62cbbSBarry Smith . o_nzz - number of nonzeros per row in off-diagonal portion of matrix 10591eb62cbbSBarry Smith . or null. You must have at least one nonzero per row. 10608a729477SBarry Smith 10618a729477SBarry Smith Output parameters: 10628a729477SBarry Smith . newmat - the matrix 10638a729477SBarry Smith 10641eb62cbbSBarry Smith Keywords: matrix, aij, compressed row, sparse, parallel 10658a729477SBarry Smith @*/ 10661eb62cbbSBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N, 10671eb62cbbSBarry Smith int d_nz,int *d_nnz, int o_nz,int *o_nnz,Mat *newmat) 10688a729477SBarry Smith { 10698a729477SBarry Smith Mat mat; 107044a69424SLois Curfman McInnes Mat_MPIAIJ *aij; 10716abc6512SBarry Smith int ierr, i,sum[2],work[2]; 10728a729477SBarry Smith *newmat = 0; 107344a69424SLois Curfman McInnes PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATMPIAIJ,comm); 1074a5a9c739SBarry Smith PLogObjectCreate(mat); 107544a69424SLois Curfman McInnes mat->data = (void *) (aij = NEW(Mat_MPIAIJ)); CHKPTR(aij); 10768a729477SBarry Smith mat->ops = &MatOps; 107744a69424SLois Curfman McInnes mat->destroy = MatDestroy_MPIAIJ; 107844a69424SLois Curfman McInnes mat->view = MatView_MPIAIJ; 10798a729477SBarry Smith mat->factor = 0; 10808a729477SBarry Smith mat->row = 0; 10818a729477SBarry Smith mat->col = 0; 1082d6dfbf8fSBarry Smith 1083d6dfbf8fSBarry Smith mat->comm = comm; 10841eb62cbbSBarry Smith aij->insertmode = NotSetValues; 10851eb62cbbSBarry Smith MPI_Comm_rank(comm,&aij->mytid); 10861eb62cbbSBarry Smith MPI_Comm_size(comm,&aij->numtids); 10871eb62cbbSBarry Smith 10881eb62cbbSBarry Smith if (M == -1 || N == -1) { 10891eb62cbbSBarry Smith work[0] = m; work[1] = n; 1090d6dfbf8fSBarry Smith MPI_Allreduce((void *) work,(void *) sum,2,MPI_INT,MPI_SUM,comm ); 10911eb62cbbSBarry Smith if (M == -1) M = sum[0]; 10921eb62cbbSBarry Smith if (N == -1) N = sum[1]; 10931eb62cbbSBarry Smith } 10941eb62cbbSBarry Smith if (m == -1) {m = M/aij->numtids + ((M % aij->numtids) > aij->mytid);} 10951eb62cbbSBarry Smith if (n == -1) {n = N/aij->numtids + ((N % aij->numtids) > aij->mytid);} 10968a729477SBarry Smith aij->m = m; 10978a729477SBarry Smith aij->n = n; 10981eb62cbbSBarry Smith aij->N = N; 10991eb62cbbSBarry Smith aij->M = M; 11001eb62cbbSBarry Smith 11011eb62cbbSBarry Smith /* build local table of row and column ownerships */ 11021eb62cbbSBarry Smith aij->rowners = (int *) MALLOC(2*(aij->numtids+2)*sizeof(int)); 11031eb62cbbSBarry Smith CHKPTR(aij->rowners); 11041eb62cbbSBarry Smith aij->cowners = aij->rowners + aij->numtids + 1; 11051eb62cbbSBarry Smith MPI_Allgather(&m,1,MPI_INT,aij->rowners+1,1,MPI_INT,comm); 11061eb62cbbSBarry Smith aij->rowners[0] = 0; 11071eb62cbbSBarry Smith for ( i=2; i<=aij->numtids; i++ ) { 11081eb62cbbSBarry Smith aij->rowners[i] += aij->rowners[i-1]; 11098a729477SBarry Smith } 11101eb62cbbSBarry Smith aij->rstart = aij->rowners[aij->mytid]; 11111eb62cbbSBarry Smith aij->rend = aij->rowners[aij->mytid+1]; 11121eb62cbbSBarry Smith MPI_Allgather(&n,1,MPI_INT,aij->cowners+1,1,MPI_INT,comm); 11131eb62cbbSBarry Smith aij->cowners[0] = 0; 11141eb62cbbSBarry Smith for ( i=2; i<=aij->numtids; i++ ) { 11151eb62cbbSBarry Smith aij->cowners[i] += aij->cowners[i-1]; 11168a729477SBarry Smith } 11171eb62cbbSBarry Smith aij->cstart = aij->cowners[aij->mytid]; 11181eb62cbbSBarry Smith aij->cend = aij->cowners[aij->mytid+1]; 11198a729477SBarry Smith 11208a729477SBarry Smith 11211eb62cbbSBarry Smith ierr = MatCreateSequentialAIJ(m,n,d_nz,d_nnz,&aij->A); CHKERR(ierr); 1122a5a9c739SBarry Smith PLogObjectParent(mat,aij->A); 11231eb62cbbSBarry Smith ierr = MatCreateSequentialAIJ(m,N,o_nz,o_nnz,&aij->B); CHKERR(ierr); 1124a5a9c739SBarry Smith PLogObjectParent(mat,aij->B); 11258a729477SBarry Smith 11261eb62cbbSBarry Smith /* build cache for off array entries formed */ 11271eb62cbbSBarry Smith aij->stash.nmax = CHUNCKSIZE; /* completely arbratray number */ 11281eb62cbbSBarry Smith aij->stash.n = 0; 11291eb62cbbSBarry Smith aij->stash.array = (Scalar *) MALLOC( aij->stash.nmax*(2*sizeof(int) + 11301eb62cbbSBarry Smith sizeof(Scalar))); CHKPTR(aij->stash.array); 11311eb62cbbSBarry Smith aij->stash.idx = (int *) (aij->stash.array + aij->stash.nmax); 11321eb62cbbSBarry Smith aij->stash.idy = (int *) (aij->stash.idx + aij->stash.nmax); 11339e25ed09SBarry Smith aij->colmap = 0; 11349e25ed09SBarry Smith aij->garray = 0; 11358a729477SBarry Smith 11361eb62cbbSBarry Smith /* stuff used for matrix vector multiply */ 11371eb62cbbSBarry Smith aij->lvec = 0; 11381eb62cbbSBarry Smith aij->Mvctx = 0; 1139d6dfbf8fSBarry Smith aij->assembled = 0; 11408a729477SBarry Smith 1141d6dfbf8fSBarry Smith *newmat = mat; 1142d6dfbf8fSBarry Smith return 0; 1143d6dfbf8fSBarry Smith } 1144c74985f6SBarry Smith 114544a69424SLois Curfman McInnes static int MatCopy_MPIAIJ(Mat matin,Mat *newmat) 1146d6dfbf8fSBarry Smith { 1147d6dfbf8fSBarry Smith Mat mat; 114844a69424SLois Curfman McInnes Mat_MPIAIJ *aij,*oldmat = (Mat_MPIAIJ *) matin->data; 11496abc6512SBarry Smith int ierr; 1150d6dfbf8fSBarry Smith *newmat = 0; 1151d6dfbf8fSBarry Smith 1152d6dfbf8fSBarry Smith if (!oldmat->assembled) SETERR(1,"Cannot copy unassembled matrix"); 115344a69424SLois Curfman McInnes PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATMPIAIJ,matin->comm); 1154a5a9c739SBarry Smith PLogObjectCreate(mat); 115544a69424SLois Curfman McInnes mat->data = (void *) (aij = NEW(Mat_MPIAIJ)); CHKPTR(aij); 1156d6dfbf8fSBarry Smith mat->ops = &MatOps; 115744a69424SLois Curfman McInnes mat->destroy = MatDestroy_MPIAIJ; 115844a69424SLois Curfman McInnes mat->view = MatView_MPIAIJ; 1159d6dfbf8fSBarry Smith mat->factor = matin->factor; 1160d6dfbf8fSBarry Smith mat->row = 0; 1161d6dfbf8fSBarry Smith mat->col = 0; 1162d6dfbf8fSBarry Smith 1163d6dfbf8fSBarry Smith aij->m = oldmat->m; 1164d6dfbf8fSBarry Smith aij->n = oldmat->n; 1165d6dfbf8fSBarry Smith aij->M = oldmat->M; 1166d6dfbf8fSBarry Smith aij->N = oldmat->N; 1167d6dfbf8fSBarry Smith 1168d6dfbf8fSBarry Smith aij->assembled = 1; 1169d6dfbf8fSBarry Smith aij->rstart = oldmat->rstart; 1170d6dfbf8fSBarry Smith aij->rend = oldmat->rend; 1171d6dfbf8fSBarry Smith aij->cstart = oldmat->cstart; 1172d6dfbf8fSBarry Smith aij->cend = oldmat->cend; 1173d6dfbf8fSBarry Smith aij->numtids = oldmat->numtids; 1174d6dfbf8fSBarry Smith aij->mytid = oldmat->mytid; 1175d6dfbf8fSBarry Smith aij->insertmode = NotSetValues; 1176d6dfbf8fSBarry Smith 1177d6dfbf8fSBarry Smith aij->rowners = (int *) MALLOC( (aij->numtids+1)*sizeof(int) ); 1178d6dfbf8fSBarry Smith CHKPTR(aij->rowners); 1179d6dfbf8fSBarry Smith MEMCPY(aij->rowners,oldmat->rowners,(aij->numtids+1)*sizeof(int)); 1180d6dfbf8fSBarry Smith aij->stash.nmax = 0; 1181d6dfbf8fSBarry Smith aij->stash.n = 0; 1182d6dfbf8fSBarry Smith aij->stash.array= 0; 11839e25ed09SBarry Smith aij->colmap = 0; 11849e25ed09SBarry Smith aij->garray = 0; 1185d6dfbf8fSBarry Smith mat->comm = matin->comm; 1186d6dfbf8fSBarry Smith 1187d6dfbf8fSBarry Smith ierr = VecCreate(oldmat->lvec,&aij->lvec); CHKERR(ierr); 1188a5a9c739SBarry Smith PLogObjectParent(mat,aij->lvec); 1189d6dfbf8fSBarry Smith ierr = VecScatterCtxCopy(oldmat->Mvctx,&aij->Mvctx); CHKERR(ierr); 1190a5a9c739SBarry Smith PLogObjectParent(mat,aij->Mvctx); 1191d6dfbf8fSBarry Smith ierr = MatCopy(oldmat->A,&aij->A); CHKERR(ierr); 1192a5a9c739SBarry Smith PLogObjectParent(mat,aij->A); 1193d6dfbf8fSBarry Smith ierr = MatCopy(oldmat->B,&aij->B); CHKERR(ierr); 1194a5a9c739SBarry Smith PLogObjectParent(mat,aij->B); 11958a729477SBarry Smith *newmat = mat; 11968a729477SBarry Smith return 0; 11978a729477SBarry Smith } 1198