1cb512458SBarry Smith #ifndef lint 2*28988994SBarry Smith static char vcid[] = "$Id: mpiaij.c,v 1.7 1995/03/06 04:04:51 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 551eb62cbbSBarry Smith static int MatiAIJInsertValues(Mat mat,int m,int *idxm,int n, 561eb62cbbSBarry Smith int *idxn,Scalar *v,InsertMode addv) 578a729477SBarry Smith { 581eb62cbbSBarry Smith Matimpiaij *aij = (Matimpiaij *) mat->data; 591eb62cbbSBarry Smith int ierr,i,j, rstart = aij->rstart, rend = aij->rend; 601eb62cbbSBarry Smith int cstart = aij->cstart, cend = aij->cend,row,col; 618a729477SBarry Smith 621eb62cbbSBarry Smith if (aij->insertmode != NotSetValues && aij->insertmode != addv) { 631eb62cbbSBarry Smith SETERR(1,"You cannot mix inserts and adds"); 648a729477SBarry Smith } 651eb62cbbSBarry Smith aij->insertmode = addv; 668a729477SBarry Smith for ( i=0; i<m; i++ ) { 67da3a660dSBarry Smith if (idxm[i] < 0) SETERR(1,"Negative row index"); 68da3a660dSBarry Smith if (idxm[i] >= aij->M) SETERR(1,"Row index too large"); 691eb62cbbSBarry Smith if (idxm[i] >= rstart && idxm[i] < rend) { 701eb62cbbSBarry Smith row = idxm[i] - rstart; 711eb62cbbSBarry Smith for ( j=0; j<n; j++ ) { 72da3a660dSBarry Smith if (idxn[j] < 0) SETERR(1,"Negative column index"); 73da3a660dSBarry Smith if (idxn[j] >= aij->N) SETERR(1,"Column index too large"); 741eb62cbbSBarry Smith if (idxn[j] >= cstart && idxn[j] < cend){ 751eb62cbbSBarry Smith col = idxn[j] - cstart; 761eb62cbbSBarry Smith ierr = MatSetValues(aij->A,1,&row,1,&col,v+i*n+j,addv);CHKERR(ierr); 771eb62cbbSBarry Smith } 781eb62cbbSBarry Smith else { 79d6dfbf8fSBarry Smith if (aij->assembled) { 80d6dfbf8fSBarry Smith SETERR(1,"Cannot insert off diagonal block in already\ 81d6dfbf8fSBarry Smith assembled matrix. Contact petsc-maint@mcs.anl.gov\ 82d6dfbf8fSBarry Smith if your need this feature"); 83d6dfbf8fSBarry Smith } 841eb62cbbSBarry Smith col = idxn[j]; 851eb62cbbSBarry Smith ierr = MatSetValues(aij->B,1,&row,1,&col,v+i*n+j,addv);CHKERR(ierr); 861eb62cbbSBarry Smith } 871eb62cbbSBarry Smith } 881eb62cbbSBarry Smith } 891eb62cbbSBarry Smith else { 901eb62cbbSBarry Smith ierr = StashValues(&aij->stash,idxm[i],n,idxn,v+i*n,addv);CHKERR(ierr); 911eb62cbbSBarry Smith } 928a729477SBarry Smith } 938a729477SBarry Smith return 0; 948a729477SBarry Smith } 958a729477SBarry Smith 968a729477SBarry Smith /* 971eb62cbbSBarry Smith the assembly code is alot like the code for vectors, we should 981eb62cbbSBarry Smith sometime derive a single assembly code that can be used for 991eb62cbbSBarry Smith either case. 1008a729477SBarry Smith */ 1018a729477SBarry Smith 1021eb62cbbSBarry Smith static int MatiAIJBeginAssemble(Mat mat) 1038a729477SBarry Smith { 1041eb62cbbSBarry Smith Matimpiaij *aij = (Matimpiaij *) mat->data; 105d6dfbf8fSBarry Smith MPI_Comm comm = mat->comm; 1066abc6512SBarry Smith int numtids = aij->numtids, *owners = aij->rowners; 1071eb62cbbSBarry Smith int mytid = aij->mytid; 1081eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 1096abc6512SBarry Smith int *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work; 1101eb62cbbSBarry Smith int tag = 50, *owner,*starts,count; 1111eb62cbbSBarry Smith InsertMode addv; 1121eb62cbbSBarry Smith Scalar *rvalues,*svalues; 1131eb62cbbSBarry Smith 1141eb62cbbSBarry Smith /* make sure all processors are either in INSERTMODE or ADDMODE */ 115*28988994SBarry Smith MPI_Allreduce((void *) &aij->insertmode,(void *) &addv,1,MPI_INT, 1161eb62cbbSBarry Smith MPI_BOR,comm); 1171eb62cbbSBarry Smith if (addv == (AddValues|InsertValues)) { 1181eb62cbbSBarry Smith SETERR(1,"Some processors have inserted while others have added"); 1191eb62cbbSBarry Smith } 1201eb62cbbSBarry Smith aij->insertmode = addv; /* in case this processor had no cache */ 1211eb62cbbSBarry Smith 1221eb62cbbSBarry Smith /* first count number of contributors to each processor */ 1231eb62cbbSBarry Smith nprocs = (int *) MALLOC( 2*numtids*sizeof(int) ); CHKPTR(nprocs); 1241eb62cbbSBarry Smith MEMSET(nprocs,0,2*numtids*sizeof(int)); procs = nprocs + numtids; 1251eb62cbbSBarry Smith owner = (int *) MALLOC( (aij->stash.n+1)*sizeof(int) ); CHKPTR(owner); 1261eb62cbbSBarry Smith for ( i=0; i<aij->stash.n; i++ ) { 1271eb62cbbSBarry Smith idx = aij->stash.idx[i]; 1281eb62cbbSBarry Smith for ( j=0; j<numtids; j++ ) { 1291eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 1301eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; break; 1318a729477SBarry Smith } 1328a729477SBarry Smith } 1338a729477SBarry Smith } 1341eb62cbbSBarry Smith nsends = 0; for ( i=0; i<numtids; i++ ) { nsends += procs[i];} 1351eb62cbbSBarry Smith 1361eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 1371eb62cbbSBarry Smith work = (int *) MALLOC( numtids*sizeof(int) ); CHKPTR(work); 1381eb62cbbSBarry Smith MPI_Allreduce((void *) procs,(void *) work,numtids,MPI_INT,MPI_SUM,comm); 1391eb62cbbSBarry Smith nreceives = work[mytid]; 1401eb62cbbSBarry Smith MPI_Allreduce((void *) nprocs,(void *) work,numtids,MPI_INT,MPI_MAX,comm); 1411eb62cbbSBarry Smith nmax = work[mytid]; 1421eb62cbbSBarry Smith FREE(work); 1431eb62cbbSBarry Smith 1441eb62cbbSBarry Smith /* post receives: 1451eb62cbbSBarry Smith 1) each message will consist of ordered pairs 1461eb62cbbSBarry Smith (global index,value) we store the global index as a double 147d6dfbf8fSBarry Smith to simplify the message passing. 1481eb62cbbSBarry Smith 2) since we don't know how long each individual message is we 1491eb62cbbSBarry Smith allocate the largest needed buffer for each receive. Potentially 1501eb62cbbSBarry Smith this is a lot of wasted space. 1511eb62cbbSBarry Smith 1521eb62cbbSBarry Smith 1531eb62cbbSBarry Smith This could be done better. 1541eb62cbbSBarry Smith */ 155*28988994SBarry Smith rvalues = (Scalar *) MALLOC(3*(nreceives+1)*(nmax+1)*sizeof(Scalar)); 1561eb62cbbSBarry Smith CHKPTR(rvalues); 1571eb62cbbSBarry Smith recv_waits = (MPI_Request *) MALLOC((nreceives+1)*sizeof(MPI_Request)); 1581eb62cbbSBarry Smith CHKPTR(recv_waits); 1591eb62cbbSBarry Smith for ( i=0; i<nreceives; i++ ) { 1601eb62cbbSBarry Smith MPI_Irecv((void *)(rvalues+3*nmax*i),3*nmax,MPI_SCALAR,MPI_ANY_SOURCE,tag, 1611eb62cbbSBarry Smith comm,recv_waits+i); 1621eb62cbbSBarry Smith } 1631eb62cbbSBarry Smith 1641eb62cbbSBarry Smith /* do sends: 1651eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 1661eb62cbbSBarry Smith the ith processor 1671eb62cbbSBarry Smith */ 1681eb62cbbSBarry Smith svalues = (Scalar *) MALLOC( 3*(aij->stash.n+1)*sizeof(Scalar) ); 1691eb62cbbSBarry Smith CHKPTR(svalues); 1701eb62cbbSBarry Smith send_waits = (MPI_Request *) MALLOC( (nsends+1)*sizeof(MPI_Request)); 1711eb62cbbSBarry Smith CHKPTR(send_waits); 1721eb62cbbSBarry Smith starts = (int *) MALLOC( numtids*sizeof(int) ); CHKPTR(starts); 1731eb62cbbSBarry Smith starts[0] = 0; 1741eb62cbbSBarry Smith for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 1751eb62cbbSBarry Smith for ( i=0; i<aij->stash.n; i++ ) { 1761eb62cbbSBarry Smith svalues[3*starts[owner[i]]] = (Scalar) aij->stash.idx[i]; 1771eb62cbbSBarry Smith svalues[3*starts[owner[i]]+1] = (Scalar) aij->stash.idy[i]; 1781eb62cbbSBarry Smith svalues[3*(starts[owner[i]]++)+2] = aij->stash.array[i]; 1791eb62cbbSBarry Smith } 1801eb62cbbSBarry Smith FREE(owner); 1811eb62cbbSBarry Smith starts[0] = 0; 1821eb62cbbSBarry Smith for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 1831eb62cbbSBarry Smith count = 0; 1841eb62cbbSBarry Smith for ( i=0; i<numtids; i++ ) { 1851eb62cbbSBarry Smith if (procs[i]) { 1861eb62cbbSBarry Smith MPI_Isend((void*)(svalues+3*starts[i]),3*nprocs[i],MPI_SCALAR,i,tag, 1871eb62cbbSBarry Smith comm,send_waits+count++); 1881eb62cbbSBarry Smith } 1891eb62cbbSBarry Smith } 1901eb62cbbSBarry Smith FREE(starts); FREE(nprocs); 1911eb62cbbSBarry Smith 1921eb62cbbSBarry Smith /* Free cache space */ 1931eb62cbbSBarry Smith aij->stash.nmax = aij->stash.n = 0; 1941eb62cbbSBarry Smith if (aij->stash.array){ FREE(aij->stash.array); aij->stash.array = 0;} 1951eb62cbbSBarry Smith 1961eb62cbbSBarry Smith aij->svalues = svalues; aij->rvalues = rvalues; 1971eb62cbbSBarry Smith aij->nsends = nsends; aij->nrecvs = nreceives; 1981eb62cbbSBarry Smith aij->send_waits = send_waits; aij->recv_waits = recv_waits; 1991eb62cbbSBarry Smith aij->rmax = nmax; 2001eb62cbbSBarry Smith 2011eb62cbbSBarry Smith return 0; 2021eb62cbbSBarry Smith } 2031eb62cbbSBarry Smith extern int MPIAIJSetUpMultiply(Mat); 2041eb62cbbSBarry Smith 2051eb62cbbSBarry Smith static int MatiAIJEndAssemble(Mat mat) 2061eb62cbbSBarry Smith { 2071eb62cbbSBarry Smith int ierr; 2081eb62cbbSBarry Smith Matimpiaij *aij = (Matimpiaij *) mat->data; 2091eb62cbbSBarry Smith 2101eb62cbbSBarry Smith MPI_Status *send_status,recv_status; 2116abc6512SBarry Smith int imdex,nrecvs = aij->nrecvs, count = nrecvs, i, n; 2121eb62cbbSBarry Smith int row,col; 2131eb62cbbSBarry Smith Scalar *values,val; 2141eb62cbbSBarry Smith InsertMode addv = aij->insertmode; 2151eb62cbbSBarry Smith 2161eb62cbbSBarry Smith /* wait on receives */ 2171eb62cbbSBarry Smith while (count) { 218d6dfbf8fSBarry Smith MPI_Waitany(nrecvs,aij->recv_waits,&imdex,&recv_status); 2191eb62cbbSBarry Smith /* unpack receives into our local space */ 220d6dfbf8fSBarry Smith values = aij->rvalues + 3*imdex*aij->rmax; 2211eb62cbbSBarry Smith MPI_Get_count(&recv_status,MPI_SCALAR,&n); 2221eb62cbbSBarry Smith n = n/3; 2231eb62cbbSBarry Smith for ( i=0; i<n; i++ ) { 2241eb62cbbSBarry Smith row = (int) PETSCREAL(values[3*i]) - aij->rstart; 2251eb62cbbSBarry Smith col = (int) PETSCREAL(values[3*i+1]); 2261eb62cbbSBarry Smith val = values[3*i+2]; 2271eb62cbbSBarry Smith if (col >= aij->cstart && col < aij->cend) { 2281eb62cbbSBarry Smith col -= aij->cstart; 2291eb62cbbSBarry Smith MatSetValues(aij->A,1,&row,1,&col,&val,addv); 2301eb62cbbSBarry Smith } 2311eb62cbbSBarry Smith else { 232d6dfbf8fSBarry Smith if (aij->assembled) { 233d6dfbf8fSBarry Smith SETERR(1,"Cannot insert off diagonal block in already\ 234d6dfbf8fSBarry Smith assembled matrix. Contact petsc-maint@mcs.anl.gov\ 235d6dfbf8fSBarry Smith if your need this feature"); 236d6dfbf8fSBarry Smith } 2371eb62cbbSBarry Smith MatSetValues(aij->B,1,&row,1,&col,&val,addv); 2381eb62cbbSBarry Smith } 2391eb62cbbSBarry Smith } 2401eb62cbbSBarry Smith count--; 2411eb62cbbSBarry Smith } 2421eb62cbbSBarry Smith FREE(aij->recv_waits); FREE(aij->rvalues); 2431eb62cbbSBarry Smith 2441eb62cbbSBarry Smith /* wait on sends */ 2451eb62cbbSBarry Smith if (aij->nsends) { 2461eb62cbbSBarry Smith send_status = (MPI_Status *) MALLOC( aij->nsends*sizeof(MPI_Status) ); 2471eb62cbbSBarry Smith CHKPTR(send_status); 2481eb62cbbSBarry Smith MPI_Waitall(aij->nsends,aij->send_waits,send_status); 2491eb62cbbSBarry Smith FREE(send_status); 2501eb62cbbSBarry Smith } 2511eb62cbbSBarry Smith FREE(aij->send_waits); FREE(aij->svalues); 2521eb62cbbSBarry Smith 2531eb62cbbSBarry Smith aij->insertmode = NotSetValues; 2541eb62cbbSBarry Smith ierr = MatBeginAssembly(aij->A); CHKERR(ierr); 2551eb62cbbSBarry Smith ierr = MatEndAssembly(aij->A); CHKERR(ierr); 2561eb62cbbSBarry Smith 257d6dfbf8fSBarry Smith 258d6dfbf8fSBarry Smith if (aij->assembled) { 259d6dfbf8fSBarry Smith /* this is because the present implimentation doesn't support 260d6dfbf8fSBarry Smith inserting values off the diagonal block.*/ 261d6dfbf8fSBarry Smith return 0; 262d6dfbf8fSBarry Smith } 263d6dfbf8fSBarry Smith 2641eb62cbbSBarry Smith ierr = MPIAIJSetUpMultiply(mat); CHKERR(ierr); 2651eb62cbbSBarry Smith ierr = MatBeginAssembly(aij->B); CHKERR(ierr); 2661eb62cbbSBarry Smith ierr = MatEndAssembly(aij->B); CHKERR(ierr); 267d6dfbf8fSBarry Smith 268d6dfbf8fSBarry Smith aij->assembled = 1; 2698a729477SBarry Smith return 0; 2708a729477SBarry Smith } 2718a729477SBarry Smith 2721eb62cbbSBarry Smith static int MatiZero(Mat A) 2731eb62cbbSBarry Smith { 2741eb62cbbSBarry Smith Matimpiaij *l = (Matimpiaij *) A->data; 2751eb62cbbSBarry Smith 2761eb62cbbSBarry Smith MatZeroEntries(l->A); MatZeroEntries(l->B); 2771eb62cbbSBarry Smith return 0; 2781eb62cbbSBarry Smith } 2791eb62cbbSBarry Smith 2801eb62cbbSBarry Smith /* again this uses the same basic stratagy as in the assembly and 2811eb62cbbSBarry Smith scatter create routines, we should try to do it systemamatically 2821eb62cbbSBarry Smith if we can figure out the proper level of generality. */ 2831eb62cbbSBarry Smith 2841eb62cbbSBarry Smith /* the code does not do the diagonal entries correctly unless the 2851eb62cbbSBarry Smith matrix is square and the column and row owerships are identical. 2861eb62cbbSBarry Smith This is a BUG. The only way to fix it seems to be to access 2871eb62cbbSBarry Smith aij->A and aij->B directly and not through the MatZeroRows() 2881eb62cbbSBarry Smith routine. 2891eb62cbbSBarry Smith */ 2901eb62cbbSBarry Smith static int MatiZerorows(Mat A,IS is,Scalar *diag) 2911eb62cbbSBarry Smith { 2921eb62cbbSBarry Smith Matimpiaij *l = (Matimpiaij *) A->data; 2931eb62cbbSBarry Smith int i,ierr,N, *rows,*owners = l->rowners,numtids = l->numtids; 2946abc6512SBarry Smith int *procs,*nprocs,j,found,idx,nsends,*work; 2951eb62cbbSBarry Smith int nmax,*svalues,*starts,*owner,nrecvs,mytid = l->mytid; 2966abc6512SBarry Smith int *rvalues,tag = 67,count,base,slen,n,*source; 297d6dfbf8fSBarry Smith int *lens,imdex,*lrows,*values; 298d6dfbf8fSBarry Smith MPI_Comm comm = A->comm; 2991eb62cbbSBarry Smith MPI_Request *send_waits,*recv_waits; 3001eb62cbbSBarry Smith MPI_Status recv_status,*send_status; 3011eb62cbbSBarry Smith IS istmp; 3021eb62cbbSBarry Smith 303da3a660dSBarry Smith if (!l->assembled) SETERR(1,"MatiZerorows: must assmble matrix first"); 3041eb62cbbSBarry Smith ierr = ISGetLocalSize(is,&N); CHKERR(ierr); 3051eb62cbbSBarry Smith ierr = ISGetIndices(is,&rows); CHKERR(ierr); 3061eb62cbbSBarry Smith 3071eb62cbbSBarry Smith /* first count number of contributors to each processor */ 3081eb62cbbSBarry Smith nprocs = (int *) MALLOC( 2*numtids*sizeof(int) ); CHKPTR(nprocs); 3091eb62cbbSBarry Smith MEMSET(nprocs,0,2*numtids*sizeof(int)); procs = nprocs + numtids; 3101eb62cbbSBarry Smith owner = (int *) MALLOC((N+1)*sizeof(int)); CHKPTR(owner); /* see note*/ 3111eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 3121eb62cbbSBarry Smith idx = rows[i]; 3131eb62cbbSBarry Smith found = 0; 3141eb62cbbSBarry Smith for ( j=0; j<numtids; j++ ) { 3151eb62cbbSBarry Smith if (idx >= owners[j] && idx < owners[j+1]) { 3161eb62cbbSBarry Smith nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 3171eb62cbbSBarry Smith } 3181eb62cbbSBarry Smith } 319d6dfbf8fSBarry Smith if (!found) SETERR(1,"Imdex out of range"); 3201eb62cbbSBarry Smith } 3211eb62cbbSBarry Smith nsends = 0; for ( i=0; i<numtids; i++ ) { nsends += procs[i];} 3221eb62cbbSBarry Smith 3231eb62cbbSBarry Smith /* inform other processors of number of messages and max length*/ 3241eb62cbbSBarry Smith work = (int *) MALLOC( numtids*sizeof(int) ); CHKPTR(work); 3251eb62cbbSBarry Smith MPI_Allreduce((void *) procs,(void *) work,numtids,MPI_INT,MPI_SUM,comm); 3261eb62cbbSBarry Smith nrecvs = work[mytid]; 3271eb62cbbSBarry Smith MPI_Allreduce((void *) nprocs,(void *) work,numtids,MPI_INT,MPI_MAX,comm); 3281eb62cbbSBarry Smith nmax = work[mytid]; 3291eb62cbbSBarry Smith FREE(work); 3301eb62cbbSBarry Smith 3311eb62cbbSBarry Smith /* post receives: */ 332*28988994SBarry Smith rvalues = (int *) MALLOC((nrecvs+1)*(nmax+1)*sizeof(int)); /*see note */ 3331eb62cbbSBarry Smith CHKPTR(rvalues); 3341eb62cbbSBarry Smith recv_waits = (MPI_Request *) MALLOC((nrecvs+1)*sizeof(MPI_Request)); 3351eb62cbbSBarry Smith CHKPTR(recv_waits); 3361eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 3371eb62cbbSBarry Smith MPI_Irecv((void *)(rvalues+nmax*i),nmax,MPI_INT,MPI_ANY_SOURCE,tag, 3381eb62cbbSBarry Smith comm,recv_waits+i); 3391eb62cbbSBarry Smith } 3401eb62cbbSBarry Smith 3411eb62cbbSBarry Smith /* do sends: 3421eb62cbbSBarry Smith 1) starts[i] gives the starting index in svalues for stuff going to 3431eb62cbbSBarry Smith the ith processor 3441eb62cbbSBarry Smith */ 3451eb62cbbSBarry Smith svalues = (int *) MALLOC( (N+1)*sizeof(int) ); CHKPTR(svalues); 3461eb62cbbSBarry Smith send_waits = (MPI_Request *) MALLOC( (nsends+1)*sizeof(MPI_Request)); 3471eb62cbbSBarry Smith CHKPTR(send_waits); 3481eb62cbbSBarry Smith starts = (int *) MALLOC( (numtids+1)*sizeof(int) ); CHKPTR(starts); 3491eb62cbbSBarry Smith starts[0] = 0; 3501eb62cbbSBarry Smith for ( i=1; i<numtids; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3511eb62cbbSBarry Smith for ( i=0; i<N; i++ ) { 3521eb62cbbSBarry Smith svalues[starts[owner[i]]++] = rows[i]; 3531eb62cbbSBarry Smith } 3541eb62cbbSBarry Smith ISRestoreIndices(is,&rows); 3551eb62cbbSBarry Smith 3561eb62cbbSBarry Smith starts[0] = 0; 3571eb62cbbSBarry Smith for ( i=1; i<numtids+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3581eb62cbbSBarry Smith count = 0; 3591eb62cbbSBarry Smith for ( i=0; i<numtids; i++ ) { 3601eb62cbbSBarry Smith if (procs[i]) { 3611eb62cbbSBarry Smith MPI_Isend((void*)(svalues+starts[i]),nprocs[i],MPI_INT,i,tag, 3621eb62cbbSBarry Smith comm,send_waits+count++); 3631eb62cbbSBarry Smith } 3641eb62cbbSBarry Smith } 3651eb62cbbSBarry Smith FREE(starts); 3661eb62cbbSBarry Smith 3671eb62cbbSBarry Smith base = owners[mytid]; 3681eb62cbbSBarry Smith 3691eb62cbbSBarry Smith /* wait on receives */ 3701eb62cbbSBarry Smith lens = (int *) MALLOC( 2*(nrecvs+1)*sizeof(int) ); CHKPTR(lens); 3711eb62cbbSBarry Smith source = lens + nrecvs; 3721eb62cbbSBarry Smith count = nrecvs; slen = 0; 3731eb62cbbSBarry Smith while (count) { 374d6dfbf8fSBarry Smith MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status); 3751eb62cbbSBarry Smith /* unpack receives into our local space */ 3761eb62cbbSBarry Smith MPI_Get_count(&recv_status,MPI_INT,&n); 377d6dfbf8fSBarry Smith source[imdex] = recv_status.MPI_SOURCE; 378d6dfbf8fSBarry Smith lens[imdex] = n; 3791eb62cbbSBarry Smith slen += n; 3801eb62cbbSBarry Smith count--; 3811eb62cbbSBarry Smith } 3821eb62cbbSBarry Smith FREE(recv_waits); 3831eb62cbbSBarry Smith 3841eb62cbbSBarry Smith /* move the data into the send scatter */ 3851eb62cbbSBarry Smith lrows = (int *) MALLOC( slen*sizeof(int) ); CHKPTR(lrows); 3861eb62cbbSBarry Smith count = 0; 3871eb62cbbSBarry Smith for ( i=0; i<nrecvs; i++ ) { 3881eb62cbbSBarry Smith values = rvalues + i*nmax; 3891eb62cbbSBarry Smith for ( j=0; j<lens[i]; j++ ) { 3901eb62cbbSBarry Smith lrows[count++] = values[j] - base; 3911eb62cbbSBarry Smith } 3921eb62cbbSBarry Smith } 3931eb62cbbSBarry Smith FREE(rvalues); FREE(lens); 3941eb62cbbSBarry Smith FREE(owner); FREE(nprocs); 3951eb62cbbSBarry Smith 3961eb62cbbSBarry Smith /* actually zap the local rows */ 3971eb62cbbSBarry Smith ierr = ISCreateSequential(slen,lrows,&istmp); CHKERR(ierr); FREE(lrows); 3981eb62cbbSBarry Smith ierr = MatZeroRows(l->A,istmp,diag); CHKERR(ierr); 3991eb62cbbSBarry Smith ierr = MatZeroRows(l->B,istmp,0); CHKERR(ierr); 4001eb62cbbSBarry Smith ierr = ISDestroy(istmp); CHKERR(ierr); 4011eb62cbbSBarry Smith 4021eb62cbbSBarry Smith /* wait on sends */ 4031eb62cbbSBarry Smith if (nsends) { 4041eb62cbbSBarry Smith send_status = (MPI_Status *) MALLOC( nsends*sizeof(MPI_Status) ); 4051eb62cbbSBarry Smith CHKPTR(send_status); 4061eb62cbbSBarry Smith MPI_Waitall(nsends,send_waits,send_status); 4071eb62cbbSBarry Smith FREE(send_status); 4081eb62cbbSBarry Smith } 4091eb62cbbSBarry Smith FREE(send_waits); FREE(svalues); 4101eb62cbbSBarry Smith 4111eb62cbbSBarry Smith 4121eb62cbbSBarry Smith return 0; 4131eb62cbbSBarry Smith } 4141eb62cbbSBarry Smith 4151eb62cbbSBarry Smith static int MatiAIJMult(Mat aijin,Vec xx,Vec yy) 4161eb62cbbSBarry Smith { 4171eb62cbbSBarry Smith Matimpiaij *aij = (Matimpiaij *) aijin->data; 4181eb62cbbSBarry Smith int ierr; 419da3a660dSBarry Smith if (!aij->assembled) SETERR(1,"MatiAIJMul: must assmble matrix first"); 420d6dfbf8fSBarry Smith ierr = VecScatterBegin(xx,0,aij->lvec,0,InsertValues,ScatterAll,aij->Mvctx); 4211eb62cbbSBarry Smith CHKERR(ierr); 4221eb62cbbSBarry Smith ierr = MatMult(aij->A,xx,yy); CHKERR(ierr); 423d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,aij->lvec,0,InsertValues,ScatterAll,aij->Mvctx); 424d6dfbf8fSBarry Smith CHKERR(ierr); 4251eb62cbbSBarry Smith ierr = MatMultAdd(aij->B,aij->lvec,yy,yy); CHKERR(ierr); 4261eb62cbbSBarry Smith return 0; 4271eb62cbbSBarry Smith } 4281eb62cbbSBarry Smith 429da3a660dSBarry Smith static int MatiAIJMultadd(Mat aijin,Vec xx,Vec yy,Vec zz) 430da3a660dSBarry Smith { 431da3a660dSBarry Smith Matimpiaij *aij = (Matimpiaij *) aijin->data; 432da3a660dSBarry Smith int ierr; 433da3a660dSBarry Smith if (!aij->assembled) SETERR(1,"MatiAIJMul: must assmble matrix first"); 434da3a660dSBarry Smith ierr = VecScatterBegin(xx,0,aij->lvec,0,InsertValues,ScatterAll,aij->Mvctx); 435da3a660dSBarry Smith CHKERR(ierr); 436da3a660dSBarry Smith ierr = MatMultAdd(aij->A,xx,yy,zz); CHKERR(ierr); 437da3a660dSBarry Smith ierr = VecScatterEnd(xx,0,aij->lvec,0,InsertValues,ScatterAll,aij->Mvctx); 438da3a660dSBarry Smith CHKERR(ierr); 439da3a660dSBarry Smith ierr = MatMultAdd(aij->B,aij->lvec,zz,zz); CHKERR(ierr); 440da3a660dSBarry Smith return 0; 441da3a660dSBarry Smith } 442da3a660dSBarry Smith 443da3a660dSBarry Smith static int MatiAIJMultTrans(Mat aijin,Vec xx,Vec yy) 444da3a660dSBarry Smith { 445da3a660dSBarry Smith Matimpiaij *aij = (Matimpiaij *) aijin->data; 446da3a660dSBarry Smith int ierr; 447da3a660dSBarry Smith 448da3a660dSBarry Smith if (!aij->assembled) SETERR(1,"MatiAIJMulTrans: must assmble matrix first"); 449da3a660dSBarry Smith /* do nondiagonal part */ 450da3a660dSBarry Smith ierr = MatMultTrans(aij->B,xx,aij->lvec); CHKERR(ierr); 451da3a660dSBarry Smith /* send it on its way */ 452da3a660dSBarry Smith ierr = VecScatterBegin(aij->lvec,0,yy,0,AddValues, 453da3a660dSBarry Smith ScatterAll|ScatterReverse,aij->Mvctx); CHKERR(ierr); 454da3a660dSBarry Smith /* do local part */ 455da3a660dSBarry Smith ierr = MatMultTrans(aij->A,xx,yy); CHKERR(ierr); 456da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 457da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 458da3a660dSBarry Smith /* but is not perhaps always true. */ 459da3a660dSBarry Smith ierr = VecScatterEnd(aij->lvec,0,yy,0,AddValues,ScatterAll|ScatterReverse, 460da3a660dSBarry Smith aij->Mvctx); CHKERR(ierr); 461da3a660dSBarry Smith return 0; 462da3a660dSBarry Smith } 463da3a660dSBarry Smith 464da3a660dSBarry Smith static int MatiAIJMultTransadd(Mat aijin,Vec xx,Vec yy,Vec zz) 465da3a660dSBarry Smith { 466da3a660dSBarry Smith Matimpiaij *aij = (Matimpiaij *) aijin->data; 467da3a660dSBarry Smith int ierr; 468da3a660dSBarry Smith 469da3a660dSBarry Smith if (!aij->assembled) SETERR(1,"MatiAIJMulTrans: must assmble matrix first"); 470da3a660dSBarry Smith /* do nondiagonal part */ 471da3a660dSBarry Smith ierr = MatMultTrans(aij->B,xx,aij->lvec); CHKERR(ierr); 472da3a660dSBarry Smith /* send it on its way */ 473da3a660dSBarry Smith ierr = VecScatterBegin(aij->lvec,0,zz,0,AddValues, 474da3a660dSBarry Smith ScatterAll|ScatterReverse,aij->Mvctx); CHKERR(ierr); 475da3a660dSBarry Smith /* do local part */ 476da3a660dSBarry Smith ierr = MatMultTransAdd(aij->A,xx,yy,zz); CHKERR(ierr); 477da3a660dSBarry Smith /* receive remote parts: note this assumes the values are not actually */ 478da3a660dSBarry Smith /* inserted in yy until the next line, which is true for my implementation*/ 479da3a660dSBarry Smith /* but is not perhaps always true. */ 480da3a660dSBarry Smith ierr = VecScatterEnd(aij->lvec,0,zz,0,AddValues,ScatterAll|ScatterReverse, 481da3a660dSBarry Smith aij->Mvctx); CHKERR(ierr); 482da3a660dSBarry Smith return 0; 483da3a660dSBarry Smith } 484da3a660dSBarry Smith 4851eb62cbbSBarry Smith /* 4861eb62cbbSBarry Smith This only works correctly for square matrices where the subblock A->A is the 4871eb62cbbSBarry Smith diagonal block 4881eb62cbbSBarry Smith */ 4891eb62cbbSBarry Smith static int MatiAIJgetdiag(Mat Ain,Vec v) 4901eb62cbbSBarry Smith { 4911eb62cbbSBarry Smith Matimpiaij *A = (Matimpiaij *) Ain->data; 492da3a660dSBarry Smith if (!A->assembled) SETERR(1,"MatiAIJgetdiag: must assmble matrix first"); 4931eb62cbbSBarry Smith return MatGetDiagonal(A->A,v); 4941eb62cbbSBarry Smith } 4951eb62cbbSBarry Smith 4961eb62cbbSBarry Smith static int MatiAIJdestroy(PetscObject obj) 4971eb62cbbSBarry Smith { 4981eb62cbbSBarry Smith Mat mat = (Mat) obj; 4991eb62cbbSBarry Smith Matimpiaij *aij = (Matimpiaij *) mat->data; 5001eb62cbbSBarry Smith int ierr; 5011eb62cbbSBarry Smith FREE(aij->rowners); 5021eb62cbbSBarry Smith ierr = MatDestroy(aij->A); CHKERR(ierr); 5031eb62cbbSBarry Smith ierr = MatDestroy(aij->B); CHKERR(ierr); 5041eb62cbbSBarry Smith FREE(aij); FREE(mat); 5051eb62cbbSBarry Smith if (aij->lvec) VecDestroy(aij->lvec); 5061eb62cbbSBarry Smith if (aij->Mvctx) VecScatterCtxDestroy(aij->Mvctx); 5071eb62cbbSBarry Smith return 0; 5081eb62cbbSBarry Smith } 5091eb62cbbSBarry Smith 5101eb62cbbSBarry Smith static int MatiView(PetscObject obj,Viewer viewer) 5111eb62cbbSBarry Smith { 5121eb62cbbSBarry Smith Mat mat = (Mat) obj; 5131eb62cbbSBarry Smith Matimpiaij *aij = (Matimpiaij *) mat->data; 5141eb62cbbSBarry Smith int ierr; 5151eb62cbbSBarry Smith 516da3a660dSBarry Smith if (!aij->assembled) SETERR(1,"MatiAIJMulTrans: must assmble matrix first"); 517d6dfbf8fSBarry Smith MPE_Seq_begin(mat->comm,1); 5181eb62cbbSBarry Smith printf("[%d] rows %d starts %d ends %d cols %d starts %d ends %d\n", 5191eb62cbbSBarry Smith aij->mytid,aij->m,aij->rstart,aij->rend,aij->n,aij->cstart, 5201eb62cbbSBarry Smith aij->cend); 521*28988994SBarry Smith ierr = MatView(aij->A,0); CHKERR(ierr); 522*28988994SBarry Smith ierr = MatView(aij->B,0); CHKERR(ierr); 523*28988994SBarry Smith fflush(stdout); 524d6dfbf8fSBarry Smith MPE_Seq_end(mat->comm,1); 5251eb62cbbSBarry Smith return 0; 5261eb62cbbSBarry Smith } 5271eb62cbbSBarry Smith 528d6dfbf8fSBarry Smith extern int MatiAIJmarkdiag(Matiaij *); 5291eb62cbbSBarry Smith /* 5301eb62cbbSBarry Smith This has to provide several versions. 5311eb62cbbSBarry Smith 5321eb62cbbSBarry Smith 1) per sequential 5331eb62cbbSBarry Smith 2) a) use only local smoothing updating outer values only once. 5341eb62cbbSBarry Smith b) local smoothing updating outer values each inner iteration 535d6dfbf8fSBarry Smith 3) color updating out values betwen colors. 5361eb62cbbSBarry Smith */ 537d6dfbf8fSBarry Smith static int MatiAIJrelax(Mat matin,Vec bb,double omega,int flag,double shift, 5388a729477SBarry Smith int its,Vec xx) 5398a729477SBarry Smith { 5401eb62cbbSBarry Smith Matimpiaij *mat = (Matimpiaij *) matin->data; 541d6dfbf8fSBarry Smith Mat AA = mat->A, BB = mat->B; 542d6dfbf8fSBarry Smith Matiaij *A = (Matiaij *) AA->data, *B = (Matiaij *)BB->data; 5436abc6512SBarry Smith Scalar zero = 0.0,*b,*x,*xs,*ls,d,*v,sum,scale,*t,*ts; 5446abc6512SBarry Smith int ierr,*idx, *diag; 5456abc6512SBarry Smith int n = mat->n, m = mat->m, i; 546da3a660dSBarry Smith Vec tt; 5478a729477SBarry Smith 548da3a660dSBarry Smith if (!mat->assembled) SETERR(1,"MatiAIJRelax: must assmble matrix first"); 5491eb62cbbSBarry Smith 550d6dfbf8fSBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); VecGetArray(mat->lvec,&ls); 551d6dfbf8fSBarry Smith xs = x -1; /* shift by one for index start of 1 */ 552da3a660dSBarry Smith ls--; 5536abc6512SBarry Smith if (!A->diag) {if ((ierr = MatiAIJmarkdiag(A))) return ierr;} 554d6dfbf8fSBarry Smith diag = A->diag; 555da3a660dSBarry Smith if (flag & SOR_EISENSTAT) { 556da3a660dSBarry Smith /* Let A = L + U + D; where L is lower trianglar, 557da3a660dSBarry Smith U is upper triangular, E is diagonal; This routine applies 558da3a660dSBarry Smith 559da3a660dSBarry Smith (L + E)^{-1} A (U + E)^{-1} 560da3a660dSBarry Smith 561da3a660dSBarry Smith to a vector efficiently using Eisenstat's trick. This is for 562da3a660dSBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 563da3a660dSBarry Smith is the relaxation factor. 564da3a660dSBarry Smith */ 565da3a660dSBarry Smith ierr = VecCreate(xx,&tt); CHKERR(ierr); 566da3a660dSBarry Smith VecGetArray(tt,&t); 567da3a660dSBarry Smith scale = (2.0/omega) - 1.0; 568da3a660dSBarry Smith /* x = (E + U)^{-1} b */ 569da3a660dSBarry Smith VecSet(&zero,mat->lvec); 570da3a660dSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineUp, 571da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 572da3a660dSBarry Smith for ( i=m-1; i>-1; i-- ) { 573da3a660dSBarry Smith n = A->i[i+1] - diag[i] - 1; 574da3a660dSBarry Smith idx = A->j + diag[i]; 575da3a660dSBarry Smith v = A->a + diag[i]; 576da3a660dSBarry Smith sum = b[i]; 577da3a660dSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 578da3a660dSBarry Smith d = shift + A->a[diag[i]-1]; 579da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 580da3a660dSBarry Smith idx = B->j + B->i[i] - 1; 581da3a660dSBarry Smith v = B->a + B->i[i] - 1; 582da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 583da3a660dSBarry Smith x[i] = omega*(sum/d); 584da3a660dSBarry Smith } 585da3a660dSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineUp, 586da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 587da3a660dSBarry Smith 588da3a660dSBarry Smith /* t = b - (2*E - D)x */ 589da3a660dSBarry Smith v = A->a; 590da3a660dSBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ - 1])*x[i]; } 591da3a660dSBarry Smith 592da3a660dSBarry Smith /* t = (E + L)^{-1}t */ 593da3a660dSBarry Smith ts = t - 1; /* shifted by one for index start of a or mat->j*/ 594da3a660dSBarry Smith diag = A->diag; 595da3a660dSBarry Smith VecSet(&zero,mat->lvec); 596da3a660dSBarry Smith ierr = VecPipelineBegin(tt,0,mat->lvec,0,InsertValues,PipelineDown, 597da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 598da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 599da3a660dSBarry Smith n = diag[i] - A->i[i]; 600da3a660dSBarry Smith idx = A->j + A->i[i] - 1; 601da3a660dSBarry Smith v = A->a + A->i[i] - 1; 602da3a660dSBarry Smith sum = t[i]; 603da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 604da3a660dSBarry Smith d = shift + A->a[diag[i]-1]; 605da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 606da3a660dSBarry Smith idx = B->j + B->i[i] - 1; 607da3a660dSBarry Smith v = B->a + B->i[i] - 1; 608da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 609da3a660dSBarry Smith t[i] = omega*(sum/d); 610da3a660dSBarry Smith } 611da3a660dSBarry Smith ierr = VecPipelineEnd(tt,0,mat->lvec,0,InsertValues,PipelineDown, 612da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 613da3a660dSBarry Smith /* x = x + t */ 614da3a660dSBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 615da3a660dSBarry Smith VecDestroy(tt); 616da3a660dSBarry Smith return 0; 617da3a660dSBarry Smith } 618da3a660dSBarry Smith 6191eb62cbbSBarry Smith 620d6dfbf8fSBarry Smith if ((flag & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP){ 621da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 622da3a660dSBarry Smith VecSet(&zero,mat->lvec); VecSet(&zero,xx); 623da3a660dSBarry Smith } 624da3a660dSBarry Smith else { 625d6dfbf8fSBarry Smith ierr=VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterUp,mat->Mvctx); 626d6dfbf8fSBarry Smith CHKERR(ierr); 627d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterUp,mat->Mvctx); 628d6dfbf8fSBarry Smith CHKERR(ierr); 629da3a660dSBarry Smith } 630d6dfbf8fSBarry Smith while (its--) { 631d6dfbf8fSBarry Smith /* go down through the rows */ 632d6dfbf8fSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineDown, 633d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 634d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 635d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 636d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 637d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 638d6dfbf8fSBarry Smith sum = b[i]; 639d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 640d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 641d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 642d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 643d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 644d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 645d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 646d6dfbf8fSBarry Smith } 647d6dfbf8fSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineDown, 648d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 649d6dfbf8fSBarry Smith /* come up through the rows */ 650d6dfbf8fSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineUp, 651d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 652d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 653d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 654d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 655d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 656d6dfbf8fSBarry Smith sum = b[i]; 657d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 658d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 659d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 660d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 661d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 662d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 663d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 664d6dfbf8fSBarry Smith } 665d6dfbf8fSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineUp, 666d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 667d6dfbf8fSBarry Smith } 668d6dfbf8fSBarry Smith } 669d6dfbf8fSBarry Smith else if (flag & SOR_FORWARD_SWEEP){ 670da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 671da3a660dSBarry Smith VecSet(&zero,mat->lvec); 672da3a660dSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineDown, 673da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 674da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 675da3a660dSBarry Smith n = diag[i] - A->i[i]; 676da3a660dSBarry Smith idx = A->j + A->i[i] - 1; 677da3a660dSBarry Smith v = A->a + A->i[i] - 1; 678da3a660dSBarry Smith sum = b[i]; 679da3a660dSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 680da3a660dSBarry Smith d = shift + A->a[diag[i]-1]; 681da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 682da3a660dSBarry Smith idx = B->j + B->i[i] - 1; 683da3a660dSBarry Smith v = B->a + B->i[i] - 1; 684da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 685da3a660dSBarry Smith x[i] = omega*(sum/d); 686da3a660dSBarry Smith } 687da3a660dSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineDown, 688da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 689da3a660dSBarry Smith its--; 690da3a660dSBarry Smith } 691d6dfbf8fSBarry Smith while (its--) { 692d6dfbf8fSBarry Smith ierr=VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterUp,mat->Mvctx); 693d6dfbf8fSBarry Smith CHKERR(ierr); 694d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterUp,mat->Mvctx); 695d6dfbf8fSBarry Smith CHKERR(ierr); 696d6dfbf8fSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineDown, 697d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 698d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 699d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 700d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 701d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 702d6dfbf8fSBarry Smith sum = b[i]; 703d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 704d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 705d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 706d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 707d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 708d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 709d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 710d6dfbf8fSBarry Smith } 711d6dfbf8fSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineDown, 712d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 713d6dfbf8fSBarry Smith } 714d6dfbf8fSBarry Smith } 715d6dfbf8fSBarry Smith else if (flag & SOR_BACKWARD_SWEEP){ 716da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 717da3a660dSBarry Smith VecSet(&zero,mat->lvec); 718da3a660dSBarry Smith ierr = VecPipelineBegin(xx,0,mat->lvec,0,InsertValues,PipelineUp, 719da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 720da3a660dSBarry Smith for ( i=m-1; i>-1; i-- ) { 721da3a660dSBarry Smith n = A->i[i+1] - diag[i] - 1; 722da3a660dSBarry Smith idx = A->j + diag[i]; 723da3a660dSBarry Smith v = A->a + diag[i]; 724da3a660dSBarry Smith sum = b[i]; 725da3a660dSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 726da3a660dSBarry Smith d = shift + A->a[diag[i]-1]; 727da3a660dSBarry Smith n = B->i[i+1] - B->i[i]; 728da3a660dSBarry Smith idx = B->j + B->i[i] - 1; 729da3a660dSBarry Smith v = B->a + B->i[i] - 1; 730da3a660dSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 731da3a660dSBarry Smith x[i] = omega*(sum/d); 732da3a660dSBarry Smith } 733da3a660dSBarry Smith ierr = VecPipelineEnd(xx,0,mat->lvec,0,InsertValues,PipelineUp, 734da3a660dSBarry Smith mat->Mvctx); CHKERR(ierr); 735da3a660dSBarry Smith its--; 736da3a660dSBarry Smith } 737d6dfbf8fSBarry Smith while (its--) { 738d6dfbf8fSBarry Smith ierr = VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterDown, 739d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 740d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterDown, 741d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 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_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){ 762da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 763da3a660dSBarry Smith return MatRelax(mat->A,bb,omega,flag,shift,its,xx); 764da3a660dSBarry Smith } 765d6dfbf8fSBarry Smith ierr=VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterAll,mat->Mvctx); 766d6dfbf8fSBarry Smith CHKERR(ierr); 767d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterAll,mat->Mvctx); 768d6dfbf8fSBarry Smith CHKERR(ierr); 769d6dfbf8fSBarry Smith while (its--) { 770d6dfbf8fSBarry Smith /* go down through the rows */ 771d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 772d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 773d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 774d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 775d6dfbf8fSBarry Smith sum = b[i]; 776d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 777d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 778d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 779d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 780d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 781d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 782d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 783d6dfbf8fSBarry Smith } 784d6dfbf8fSBarry Smith /* come up through the rows */ 785d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 786d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 787d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 788d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 789d6dfbf8fSBarry Smith sum = b[i]; 790d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 791d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 792d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 793d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 794d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 795d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 796d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 797d6dfbf8fSBarry Smith } 798d6dfbf8fSBarry Smith } 799d6dfbf8fSBarry Smith } 800d6dfbf8fSBarry Smith else if (flag & SOR_LOCAL_FORWARD_SWEEP){ 801da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 802da3a660dSBarry Smith return MatRelax(mat->A,bb,omega,flag,shift,its,xx); 803da3a660dSBarry Smith } 804d6dfbf8fSBarry Smith ierr=VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterAll,mat->Mvctx); 805d6dfbf8fSBarry Smith CHKERR(ierr); 806d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterAll,mat->Mvctx); 807d6dfbf8fSBarry Smith CHKERR(ierr); 808d6dfbf8fSBarry Smith while (its--) { 809d6dfbf8fSBarry Smith for ( i=0; i<m; i++ ) { 810d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 811d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 812d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 813d6dfbf8fSBarry Smith sum = b[i]; 814d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 815d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 816d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 817d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 818d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 819d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 820d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 821d6dfbf8fSBarry Smith } 822d6dfbf8fSBarry Smith } 823d6dfbf8fSBarry Smith } 824d6dfbf8fSBarry Smith else if (flag & SOR_LOCAL_BACKWARD_SWEEP){ 825da3a660dSBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 826da3a660dSBarry Smith return MatRelax(mat->A,bb,omega,flag,shift,its,xx); 827da3a660dSBarry Smith } 828d6dfbf8fSBarry Smith ierr = VecScatterBegin(xx,0,mat->lvec,0,InsertValues,ScatterAll, 829d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 830d6dfbf8fSBarry Smith ierr = VecScatterEnd(xx,0,mat->lvec,0,InsertValues,ScatterAll, 831d6dfbf8fSBarry Smith mat->Mvctx); CHKERR(ierr); 832d6dfbf8fSBarry Smith while (its--) { 833d6dfbf8fSBarry Smith for ( i=m-1; i>-1; i-- ) { 834d6dfbf8fSBarry Smith n = A->i[i+1] - A->i[i]; 835d6dfbf8fSBarry Smith idx = A->j + A->i[i] - 1; 836d6dfbf8fSBarry Smith v = A->a + A->i[i] - 1; 837d6dfbf8fSBarry Smith sum = b[i]; 838d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 839d6dfbf8fSBarry Smith d = shift + A->a[diag[i]-1]; 840d6dfbf8fSBarry Smith n = B->i[i+1] - B->i[i]; 841d6dfbf8fSBarry Smith idx = B->j + B->i[i] - 1; 842d6dfbf8fSBarry Smith v = B->a + B->i[i] - 1; 843d6dfbf8fSBarry Smith SPARSEDENSEMDOT(sum,ls,v,idx,n); 844d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 845d6dfbf8fSBarry Smith } 846d6dfbf8fSBarry Smith } 847d6dfbf8fSBarry Smith } 8488a729477SBarry Smith return 0; 8498a729477SBarry Smith } 850c74985f6SBarry Smith static int MatiAIJinsopt(Mat aijin,int op) 851c74985f6SBarry Smith { 852c74985f6SBarry Smith Matimpiaij *aij = (Matimpiaij *) aijin->data; 853c74985f6SBarry Smith 854c74985f6SBarry Smith if (op == NO_NEW_NONZERO_LOCATIONS) { 855c74985f6SBarry Smith MatSetOption(aij->A,op); 856c74985f6SBarry Smith MatSetOption(aij->B,op); 857c74985f6SBarry Smith } 858c74985f6SBarry Smith else if (op == YES_NEW_NONZERO_LOCATIONS) { 859c74985f6SBarry Smith MatSetOption(aij->A,op); 860c74985f6SBarry Smith MatSetOption(aij->B,op); 861c74985f6SBarry Smith } 862c74985f6SBarry Smith else if (op == COLUMN_ORIENTED) SETERR(1,"Column oriented not supported"); 863c74985f6SBarry Smith return 0; 864c74985f6SBarry Smith } 865c74985f6SBarry Smith 866c74985f6SBarry Smith static int MatiAIJsize(Mat matin,int *m,int *n) 867c74985f6SBarry Smith { 868c74985f6SBarry Smith Matimpiaij *mat = (Matimpiaij *) matin->data; 869c74985f6SBarry Smith *m = mat->M; *n = mat->N; 870c74985f6SBarry Smith return 0; 871c74985f6SBarry Smith } 872c74985f6SBarry Smith 873c74985f6SBarry Smith static int MatiAIJlocalsize(Mat matin,int *m,int *n) 874c74985f6SBarry Smith { 875c74985f6SBarry Smith Matimpiaij *mat = (Matimpiaij *) matin->data; 876c74985f6SBarry Smith *m = mat->m; *n = mat->n; 877c74985f6SBarry Smith return 0; 878c74985f6SBarry Smith } 879c74985f6SBarry Smith 880c74985f6SBarry Smith static int MatiAIJrange(Mat matin,int *m,int *n) 881c74985f6SBarry Smith { 882c74985f6SBarry Smith Matimpiaij *mat = (Matimpiaij *) matin->data; 883c74985f6SBarry Smith *m = mat->rstart; *n = mat->rend; 884c74985f6SBarry Smith return 0; 885c74985f6SBarry Smith } 886c74985f6SBarry Smith 887d6dfbf8fSBarry Smith static int MatiCopy(Mat,Mat *); 888d6dfbf8fSBarry Smith 8898a729477SBarry Smith /* -------------------------------------------------------------------*/ 8901eb62cbbSBarry Smith static struct _MatOps MatOps = {MatiAIJInsertValues, 8918a729477SBarry Smith 0, 0, 892da3a660dSBarry Smith MatiAIJMult,MatiAIJMultadd,MatiAIJMultTrans,MatiAIJMultTransadd, 8931eb62cbbSBarry Smith 0,0,0,0, 8948a729477SBarry Smith 0,0, 8958a729477SBarry Smith MatiAIJrelax, 8968a729477SBarry Smith 0, 8971eb62cbbSBarry Smith 0,0,0, 898d6dfbf8fSBarry Smith MatiCopy, 8998a729477SBarry Smith MatiAIJgetdiag,0,0, 9001eb62cbbSBarry Smith MatiAIJBeginAssemble,MatiAIJEndAssemble, 9011eb62cbbSBarry Smith 0, 902c74985f6SBarry Smith MatiAIJinsopt,MatiZero,MatiZerorows,0, 903c74985f6SBarry Smith 0,0,0,0, 904c74985f6SBarry Smith MatiAIJsize,MatiAIJlocalsize,MatiAIJrange }; 9058a729477SBarry Smith 9068a729477SBarry Smith 9078a729477SBarry Smith 9088a729477SBarry Smith /*@ 9098a729477SBarry Smith 9101eb62cbbSBarry Smith MatCreateMPIAIJ - Creates a sparse parallel matrix 9111eb62cbbSBarry Smith in AIJ format. 9128a729477SBarry Smith 9138a729477SBarry Smith Input Parameters: 9141eb62cbbSBarry Smith . comm - MPI communicator 9151eb62cbbSBarry Smith . m,n - number of local rows and columns (or -1 to have calculated) 9161eb62cbbSBarry Smith . M,N - global rows and columns (or -1 to have calculated) 9171eb62cbbSBarry Smith . d_nz - total number nonzeros in diagonal portion of matrix 9181eb62cbbSBarry Smith . d_nzz - number of nonzeros per row in diagonal portion of matrix or null 9198a729477SBarry Smith . You must leave room for the diagonal entry even if it is zero. 9201eb62cbbSBarry Smith . o_nz - total number nonzeros in off-diagonal portion of matrix 9211eb62cbbSBarry Smith . o_nzz - number of nonzeros per row in off-diagonal portion of matrix 9221eb62cbbSBarry Smith . or null. You must have at least one nonzero per row. 9238a729477SBarry Smith 9248a729477SBarry Smith Output parameters: 9258a729477SBarry Smith . newmat - the matrix 9268a729477SBarry Smith 9271eb62cbbSBarry Smith Keywords: matrix, aij, compressed row, sparse, parallel 9288a729477SBarry Smith @*/ 9291eb62cbbSBarry Smith int MatCreateMPIAIJ(MPI_Comm comm,int m,int n,int M,int N, 9301eb62cbbSBarry Smith int d_nz,int *d_nnz, int o_nz,int *o_nnz,Mat *newmat) 9318a729477SBarry Smith { 9328a729477SBarry Smith Mat mat; 9331eb62cbbSBarry Smith Matimpiaij *aij; 9346abc6512SBarry Smith int ierr, i,sum[2],work[2]; 9358a729477SBarry Smith *newmat = 0; 9368a729477SBarry Smith CREATEHEADER(mat,_Mat); 9371eb62cbbSBarry Smith mat->data = (void *) (aij = NEW(Matimpiaij)); CHKPTR(aij); 9388a729477SBarry Smith mat->cookie = MAT_COOKIE; 9398a729477SBarry Smith mat->ops = &MatOps; 9408a729477SBarry Smith mat->destroy = MatiAIJdestroy; 9411eb62cbbSBarry Smith mat->view = MatiView; 9421eb62cbbSBarry Smith mat->type = MATAIJMPI; 9438a729477SBarry Smith mat->factor = 0; 9448a729477SBarry Smith mat->row = 0; 9458a729477SBarry Smith mat->col = 0; 946d6dfbf8fSBarry Smith 947d6dfbf8fSBarry Smith mat->comm = comm; 9481eb62cbbSBarry Smith aij->insertmode = NotSetValues; 9491eb62cbbSBarry Smith MPI_Comm_rank(comm,&aij->mytid); 9501eb62cbbSBarry Smith MPI_Comm_size(comm,&aij->numtids); 9511eb62cbbSBarry Smith 9521eb62cbbSBarry Smith if (M == -1 || N == -1) { 9531eb62cbbSBarry Smith work[0] = m; work[1] = n; 954d6dfbf8fSBarry Smith MPI_Allreduce((void *) work,(void *) sum,2,MPI_INT,MPI_SUM,comm ); 9551eb62cbbSBarry Smith if (M == -1) M = sum[0]; 9561eb62cbbSBarry Smith if (N == -1) N = sum[1]; 9571eb62cbbSBarry Smith } 9581eb62cbbSBarry Smith if (m == -1) {m = M/aij->numtids + ((M % aij->numtids) > aij->mytid);} 9591eb62cbbSBarry Smith if (n == -1) {n = N/aij->numtids + ((N % aij->numtids) > aij->mytid);} 9608a729477SBarry Smith aij->m = m; 9618a729477SBarry Smith aij->n = n; 9621eb62cbbSBarry Smith aij->N = N; 9631eb62cbbSBarry Smith aij->M = M; 9641eb62cbbSBarry Smith 9651eb62cbbSBarry Smith /* build local table of row and column ownerships */ 9661eb62cbbSBarry Smith aij->rowners = (int *) MALLOC(2*(aij->numtids+2)*sizeof(int)); 9671eb62cbbSBarry Smith CHKPTR(aij->rowners); 9681eb62cbbSBarry Smith aij->cowners = aij->rowners + aij->numtids + 1; 9691eb62cbbSBarry Smith MPI_Allgather(&m,1,MPI_INT,aij->rowners+1,1,MPI_INT,comm); 9701eb62cbbSBarry Smith aij->rowners[0] = 0; 9711eb62cbbSBarry Smith for ( i=2; i<=aij->numtids; i++ ) { 9721eb62cbbSBarry Smith aij->rowners[i] += aij->rowners[i-1]; 9738a729477SBarry Smith } 9741eb62cbbSBarry Smith aij->rstart = aij->rowners[aij->mytid]; 9751eb62cbbSBarry Smith aij->rend = aij->rowners[aij->mytid+1]; 9761eb62cbbSBarry Smith MPI_Allgather(&n,1,MPI_INT,aij->cowners+1,1,MPI_INT,comm); 9771eb62cbbSBarry Smith aij->cowners[0] = 0; 9781eb62cbbSBarry Smith for ( i=2; i<=aij->numtids; i++ ) { 9791eb62cbbSBarry Smith aij->cowners[i] += aij->cowners[i-1]; 9808a729477SBarry Smith } 9811eb62cbbSBarry Smith aij->cstart = aij->cowners[aij->mytid]; 9821eb62cbbSBarry Smith aij->cend = aij->cowners[aij->mytid+1]; 9838a729477SBarry Smith 9848a729477SBarry Smith 9851eb62cbbSBarry Smith ierr = MatCreateSequentialAIJ(m,n,d_nz,d_nnz,&aij->A); CHKERR(ierr); 9861eb62cbbSBarry Smith ierr = MatCreateSequentialAIJ(m,N,o_nz,o_nnz,&aij->B); CHKERR(ierr); 9878a729477SBarry Smith 9881eb62cbbSBarry Smith /* build cache for off array entries formed */ 9891eb62cbbSBarry Smith aij->stash.nmax = CHUNCKSIZE; /* completely arbratray number */ 9901eb62cbbSBarry Smith aij->stash.n = 0; 9911eb62cbbSBarry Smith aij->stash.array = (Scalar *) MALLOC( aij->stash.nmax*(2*sizeof(int) + 9921eb62cbbSBarry Smith sizeof(Scalar))); CHKPTR(aij->stash.array); 9931eb62cbbSBarry Smith aij->stash.idx = (int *) (aij->stash.array + aij->stash.nmax); 9941eb62cbbSBarry Smith aij->stash.idy = (int *) (aij->stash.idx + aij->stash.nmax); 9958a729477SBarry Smith 9961eb62cbbSBarry Smith /* stuff used for matrix vector multiply */ 9971eb62cbbSBarry Smith aij->lvec = 0; 9981eb62cbbSBarry Smith aij->Mvctx = 0; 999d6dfbf8fSBarry Smith aij->assembled = 0; 10008a729477SBarry Smith 1001d6dfbf8fSBarry Smith *newmat = mat; 1002d6dfbf8fSBarry Smith return 0; 1003d6dfbf8fSBarry Smith } 1004c74985f6SBarry Smith 1005d6dfbf8fSBarry Smith static int MatiCopy(Mat matin,Mat *newmat) 1006d6dfbf8fSBarry Smith { 1007d6dfbf8fSBarry Smith Mat mat; 1008d6dfbf8fSBarry Smith Matimpiaij *aij,*oldmat = (Matimpiaij *) matin->data; 10096abc6512SBarry Smith int ierr; 1010d6dfbf8fSBarry Smith *newmat = 0; 1011d6dfbf8fSBarry Smith 1012d6dfbf8fSBarry Smith if (!oldmat->assembled) SETERR(1,"Cannot copy unassembled matrix"); 1013d6dfbf8fSBarry Smith CREATEHEADER(mat,_Mat); 1014d6dfbf8fSBarry Smith mat->data = (void *) (aij = NEW(Matimpiaij)); CHKPTR(aij); 1015d6dfbf8fSBarry Smith mat->cookie = MAT_COOKIE; 1016d6dfbf8fSBarry Smith mat->ops = &MatOps; 1017d6dfbf8fSBarry Smith mat->destroy = MatiAIJdestroy; 1018d6dfbf8fSBarry Smith mat->view = MatiView; 1019*28988994SBarry Smith mat->type = MATAIJMPI; 1020d6dfbf8fSBarry Smith mat->factor = matin->factor; 1021d6dfbf8fSBarry Smith mat->row = 0; 1022d6dfbf8fSBarry Smith mat->col = 0; 1023d6dfbf8fSBarry Smith 1024d6dfbf8fSBarry Smith aij->m = oldmat->m; 1025d6dfbf8fSBarry Smith aij->n = oldmat->n; 1026d6dfbf8fSBarry Smith aij->M = oldmat->M; 1027d6dfbf8fSBarry Smith aij->N = oldmat->N; 1028d6dfbf8fSBarry Smith 1029d6dfbf8fSBarry Smith aij->assembled = 1; 1030d6dfbf8fSBarry Smith aij->rstart = oldmat->rstart; 1031d6dfbf8fSBarry Smith aij->rend = oldmat->rend; 1032d6dfbf8fSBarry Smith aij->cstart = oldmat->cstart; 1033d6dfbf8fSBarry Smith aij->cend = oldmat->cend; 1034d6dfbf8fSBarry Smith aij->numtids = oldmat->numtids; 1035d6dfbf8fSBarry Smith aij->mytid = oldmat->mytid; 1036d6dfbf8fSBarry Smith aij->insertmode = NotSetValues; 1037d6dfbf8fSBarry Smith 1038d6dfbf8fSBarry Smith aij->rowners = (int *) MALLOC( (aij->numtids+1)*sizeof(int) ); 1039d6dfbf8fSBarry Smith CHKPTR(aij->rowners); 1040d6dfbf8fSBarry Smith MEMCPY(aij->rowners,oldmat->rowners,(aij->numtids+1)*sizeof(int)); 1041d6dfbf8fSBarry Smith aij->stash.nmax = 0; 1042d6dfbf8fSBarry Smith aij->stash.n = 0; 1043d6dfbf8fSBarry Smith aij->stash.array= 0; 1044d6dfbf8fSBarry Smith mat->comm = matin->comm; 1045d6dfbf8fSBarry Smith 1046d6dfbf8fSBarry Smith ierr = VecCreate(oldmat->lvec,&aij->lvec); CHKERR(ierr); 1047d6dfbf8fSBarry Smith ierr = VecScatterCtxCopy(oldmat->Mvctx,&aij->Mvctx); CHKERR(ierr); 1048d6dfbf8fSBarry Smith ierr = MatCopy(oldmat->A,&aij->A); CHKERR(ierr); 1049d6dfbf8fSBarry Smith ierr = MatCopy(oldmat->B,&aij->B); CHKERR(ierr); 10508a729477SBarry Smith *newmat = mat; 10518a729477SBarry Smith return 0; 10528a729477SBarry Smith } 1053