18965ea79SLois Curfman McInnes #ifndef lint 2*70f55243SBarry Smith static char vcid[] = "$Id: mpidense.c,v 1.44 1996/08/06 04:02:14 bsmith Exp bsmith $"; 38965ea79SLois Curfman McInnes #endif 48965ea79SLois Curfman McInnes 5ed3cc1f0SBarry Smith /* 6ed3cc1f0SBarry Smith Basic functions for basic parallel dense matrices. 7ed3cc1f0SBarry Smith */ 8ed3cc1f0SBarry Smith 9*70f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h" 10f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 118965ea79SLois Curfman McInnes 1239ddd567SLois Curfman McInnes static int MatSetValues_MPIDense(Mat mat,int m,int *idxm,int n, 138965ea79SLois Curfman McInnes int *idxn,Scalar *v,InsertMode addv) 148965ea79SLois Curfman McInnes { 1539b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense *) mat->data; 1639b7565bSBarry Smith int ierr, i, j, rstart = A->rstart, rend = A->rend, row; 1739b7565bSBarry Smith int roworiented = A->roworiented; 188965ea79SLois Curfman McInnes 1939b7565bSBarry Smith if (A->insertmode != NOT_SET_VALUES && A->insertmode != addv) { 2039ddd567SLois Curfman McInnes SETERRQ(1,"MatSetValues_MPIDense:Cannot mix inserts and adds"); 218965ea79SLois Curfman McInnes } 2239b7565bSBarry Smith A->insertmode = addv; 238965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 2439ddd567SLois Curfman McInnes if (idxm[i] < 0) SETERRQ(1,"MatSetValues_MPIDense:Negative row"); 2539b7565bSBarry Smith if (idxm[i] >= A->M) SETERRQ(1,"MatSetValues_MPIDense:Row too large"); 268965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 278965ea79SLois Curfman McInnes row = idxm[i] - rstart; 2839b7565bSBarry Smith if (roworiented) { 2939b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv); CHKERRQ(ierr); 3039b7565bSBarry Smith } 3139b7565bSBarry Smith else { 328965ea79SLois Curfman McInnes for ( j=0; j<n; j++ ) { 3339ddd567SLois Curfman McInnes if (idxn[j] < 0) SETERRQ(1,"MatSetValues_MPIDense:Negative column"); 3439b7565bSBarry Smith if (idxn[j] >= A->N) SETERRQ(1,"MatSetValues_MPIDense:Column too large"); 3539b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv); CHKERRQ(ierr); 3639b7565bSBarry Smith } 378965ea79SLois Curfman McInnes } 388965ea79SLois Curfman McInnes } 398965ea79SLois Curfman McInnes else { 4039b7565bSBarry Smith if (roworiented) { 4139b7565bSBarry Smith ierr = StashValues_Private(&A->stash,idxm[i],n,idxn,v+i*n,addv); CHKERRQ(ierr); 4239b7565bSBarry Smith } 4339b7565bSBarry Smith else { /* must stash each seperately */ 4439b7565bSBarry Smith row = idxm[i]; 4539b7565bSBarry Smith for ( j=0; j<n; j++ ) { 4639b7565bSBarry Smith ierr = StashValues_Private(&A->stash,row,1,&idxn[j],v+i+j*m,addv); 4739b7565bSBarry Smith CHKERRQ(ierr); 4839b7565bSBarry Smith } 4939b7565bSBarry Smith } 50b49de8d1SLois Curfman McInnes } 51b49de8d1SLois Curfman McInnes } 52b49de8d1SLois Curfman McInnes return 0; 53b49de8d1SLois Curfman McInnes } 54b49de8d1SLois Curfman McInnes 55b49de8d1SLois Curfman McInnes static int MatGetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v) 56b49de8d1SLois Curfman McInnes { 57b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 58b49de8d1SLois Curfman McInnes int ierr, i, j, rstart = mdn->rstart, rend = mdn->rend, row; 59b49de8d1SLois Curfman McInnes 60b49de8d1SLois Curfman McInnes for ( i=0; i<m; i++ ) { 61b49de8d1SLois Curfman McInnes if (idxm[i] < 0) SETERRQ(1,"MatGetValues_MPIDense:Negative row"); 62b49de8d1SLois Curfman McInnes if (idxm[i] >= mdn->M) SETERRQ(1,"MatGetValues_MPIDense:Row too large"); 63b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 64b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 65b49de8d1SLois Curfman McInnes for ( j=0; j<n; j++ ) { 66b49de8d1SLois Curfman McInnes if (idxn[j] < 0) SETERRQ(1,"MatGetValues_MPIDense:Negative column"); 67b49de8d1SLois Curfman McInnes if (idxn[j] >= mdn->N) 68b49de8d1SLois Curfman McInnes SETERRQ(1,"MatGetValues_MPIDense:Column too large"); 69b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j); CHKERRQ(ierr); 70b49de8d1SLois Curfman McInnes } 71b49de8d1SLois Curfman McInnes } 72b49de8d1SLois Curfman McInnes else { 73b49de8d1SLois Curfman McInnes SETERRQ(1,"MatGetValues_MPIDense:Only local values currently supported"); 748965ea79SLois Curfman McInnes } 758965ea79SLois Curfman McInnes } 768965ea79SLois Curfman McInnes return 0; 778965ea79SLois Curfman McInnes } 788965ea79SLois Curfman McInnes 79ff14e315SSatish Balay static int MatGetArray_MPIDense(Mat A,Scalar **array) 80ff14e315SSatish Balay { 81ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense *) A->data; 82ff14e315SSatish Balay int ierr; 83ff14e315SSatish Balay 84ff14e315SSatish Balay ierr = MatGetArray(a->A,array); CHKERRQ(ierr); 85ff14e315SSatish Balay return 0; 86ff14e315SSatish Balay } 87ff14e315SSatish Balay 88ff14e315SSatish Balay static int MatRestoreArray_MPIDense(Mat A,Scalar **array) 89ff14e315SSatish Balay { 90ff14e315SSatish Balay return 0; 91ff14e315SSatish Balay } 92ff14e315SSatish Balay 9339ddd567SLois Curfman McInnes static int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 948965ea79SLois Curfman McInnes { 9539ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 968965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 9739ddd567SLois Curfman McInnes int size = mdn->size, *owners = mdn->rowners, rank = mdn->rank; 988965ea79SLois Curfman McInnes int *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work; 9939ddd567SLois Curfman McInnes int tag = mat->tag, *owner,*starts,count,ierr; 1008965ea79SLois Curfman McInnes InsertMode addv; 10139ddd567SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 1028965ea79SLois Curfman McInnes Scalar *rvalues,*svalues; 1038965ea79SLois Curfman McInnes 1048965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 1056d84be18SBarry Smith MPI_Allreduce(&mdn->insertmode,&addv,1,MPI_INT,MPI_BOR,comm); 10639ddd567SLois Curfman McInnes if (addv == (ADD_VALUES|INSERT_VALUES)) { SETERRQ(1, 10739ddd567SLois Curfman McInnes "MatAssemblyBegin_MPIDense:Cannot mix adds/inserts on different procs"); 1088965ea79SLois Curfman McInnes } 10939ddd567SLois Curfman McInnes mdn->insertmode = addv; /* in case this processor had no cache */ 1108965ea79SLois Curfman McInnes 1118965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 1120452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 113cddf8d76SBarry Smith PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 1140452661fSBarry Smith owner = (int *) PetscMalloc( (mdn->stash.n+1)*sizeof(int) ); CHKPTRQ(owner); 11539ddd567SLois Curfman McInnes for ( i=0; i<mdn->stash.n; i++ ) { 11639ddd567SLois Curfman McInnes idx = mdn->stash.idx[i]; 1178965ea79SLois Curfman McInnes for ( j=0; j<size; j++ ) { 1188965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 1198965ea79SLois Curfman McInnes nprocs[j]++; procs[j] = 1; owner[i] = j; break; 1208965ea79SLois Curfman McInnes } 1218965ea79SLois Curfman McInnes } 1228965ea79SLois Curfman McInnes } 1238965ea79SLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 1248965ea79SLois Curfman McInnes 1258965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 1260452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 1276d84be18SBarry Smith MPI_Allreduce(procs,work,size,MPI_INT,MPI_SUM,comm); 1288965ea79SLois Curfman McInnes nreceives = work[rank]; 1296d84be18SBarry Smith MPI_Allreduce(nprocs,work,size,MPI_INT,MPI_MAX,comm); 1308965ea79SLois Curfman McInnes nmax = work[rank]; 1310452661fSBarry Smith PetscFree(work); 1328965ea79SLois Curfman McInnes 1338965ea79SLois Curfman McInnes /* post receives: 1348965ea79SLois Curfman McInnes 1) each message will consist of ordered pairs 1358965ea79SLois Curfman McInnes (global index,value) we store the global index as a double 1368965ea79SLois Curfman McInnes to simplify the message passing. 1378965ea79SLois Curfman McInnes 2) since we don't know how long each individual message is we 1388965ea79SLois Curfman McInnes allocate the largest needed buffer for each receive. Potentially 1398965ea79SLois Curfman McInnes this is a lot of wasted space. 1408965ea79SLois Curfman McInnes 1418965ea79SLois Curfman McInnes This could be done better. 1428965ea79SLois Curfman McInnes */ 1430452661fSBarry Smith rvalues = (Scalar *) PetscMalloc(3*(nreceives+1)*(nmax+1)*sizeof(Scalar)); 1448965ea79SLois Curfman McInnes CHKPTRQ(rvalues); 1450452661fSBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nreceives+1)*sizeof(MPI_Request)); 1468965ea79SLois Curfman McInnes CHKPTRQ(recv_waits); 1478965ea79SLois Curfman McInnes for ( i=0; i<nreceives; i++ ) { 1486d84be18SBarry Smith MPI_Irecv(rvalues+3*nmax*i,3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag, 1498965ea79SLois Curfman McInnes comm,recv_waits+i); 1508965ea79SLois Curfman McInnes } 1518965ea79SLois Curfman McInnes 1528965ea79SLois Curfman McInnes /* do sends: 1538965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 1548965ea79SLois Curfman McInnes the ith processor 1558965ea79SLois Curfman McInnes */ 1560452661fSBarry Smith svalues = (Scalar *) PetscMalloc( 3*(mdn->stash.n+1)*sizeof(Scalar) ); 15739ddd567SLois Curfman McInnes CHKPTRQ(svalues); 1580452661fSBarry Smith send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request)); 1598965ea79SLois Curfman McInnes CHKPTRQ(send_waits); 1600452661fSBarry Smith starts = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(starts); 1618965ea79SLois Curfman McInnes starts[0] = 0; 1628965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 16339ddd567SLois Curfman McInnes for ( i=0; i<mdn->stash.n; i++ ) { 16439ddd567SLois Curfman McInnes svalues[3*starts[owner[i]]] = (Scalar) mdn->stash.idx[i]; 16539ddd567SLois Curfman McInnes svalues[3*starts[owner[i]]+1] = (Scalar) mdn->stash.idy[i]; 16639ddd567SLois Curfman McInnes svalues[3*(starts[owner[i]]++)+2] = mdn->stash.array[i]; 1678965ea79SLois Curfman McInnes } 1680452661fSBarry Smith PetscFree(owner); 1698965ea79SLois Curfman McInnes starts[0] = 0; 1708965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 1718965ea79SLois Curfman McInnes count = 0; 1728965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 1738965ea79SLois Curfman McInnes if (procs[i]) { 1746d84be18SBarry Smith MPI_Isend(svalues+3*starts[i],3*nprocs[i],MPIU_SCALAR,i,tag, 1758965ea79SLois Curfman McInnes comm,send_waits+count++); 1768965ea79SLois Curfman McInnes } 1778965ea79SLois Curfman McInnes } 1780452661fSBarry Smith PetscFree(starts); PetscFree(nprocs); 1798965ea79SLois Curfman McInnes 1808965ea79SLois Curfman McInnes /* Free cache space */ 18139ddd567SLois Curfman McInnes ierr = StashDestroy_Private(&mdn->stash); CHKERRQ(ierr); 1828965ea79SLois Curfman McInnes 18339ddd567SLois Curfman McInnes mdn->svalues = svalues; mdn->rvalues = rvalues; 18439ddd567SLois Curfman McInnes mdn->nsends = nsends; mdn->nrecvs = nreceives; 18539ddd567SLois Curfman McInnes mdn->send_waits = send_waits; mdn->recv_waits = recv_waits; 18639ddd567SLois Curfman McInnes mdn->rmax = nmax; 1878965ea79SLois Curfman McInnes 1888965ea79SLois Curfman McInnes return 0; 1898965ea79SLois Curfman McInnes } 19039ddd567SLois Curfman McInnes extern int MatSetUpMultiply_MPIDense(Mat); 1918965ea79SLois Curfman McInnes 19239ddd567SLois Curfman McInnes static int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 1938965ea79SLois Curfman McInnes { 19439ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 1958965ea79SLois Curfman McInnes MPI_Status *send_status,recv_status; 19639ddd567SLois Curfman McInnes int imdex, nrecvs=mdn->nrecvs, count=nrecvs, i, n, ierr, row, col; 1978965ea79SLois Curfman McInnes Scalar *values,val; 19839ddd567SLois Curfman McInnes InsertMode addv = mdn->insertmode; 1998965ea79SLois Curfman McInnes 2008965ea79SLois Curfman McInnes /* wait on receives */ 2018965ea79SLois Curfman McInnes while (count) { 20239ddd567SLois Curfman McInnes MPI_Waitany(nrecvs,mdn->recv_waits,&imdex,&recv_status); 2038965ea79SLois Curfman McInnes /* unpack receives into our local space */ 20439ddd567SLois Curfman McInnes values = mdn->rvalues + 3*imdex*mdn->rmax; 2058965ea79SLois Curfman McInnes MPI_Get_count(&recv_status,MPIU_SCALAR,&n); 2068965ea79SLois Curfman McInnes n = n/3; 2078965ea79SLois Curfman McInnes for ( i=0; i<n; i++ ) { 208227d817aSBarry Smith row = (int) PetscReal(values[3*i]) - mdn->rstart; 209227d817aSBarry Smith col = (int) PetscReal(values[3*i+1]); 2108965ea79SLois Curfman McInnes val = values[3*i+2]; 21139ddd567SLois Curfman McInnes if (col >= 0 && col < mdn->N) { 21239ddd567SLois Curfman McInnes MatSetValues(mdn->A,1,&row,1,&col,&val,addv); 2138965ea79SLois Curfman McInnes } 21439ddd567SLois Curfman McInnes else {SETERRQ(1,"MatAssemblyEnd_MPIDense:Invalid column");} 2158965ea79SLois Curfman McInnes } 2168965ea79SLois Curfman McInnes count--; 2178965ea79SLois Curfman McInnes } 2180452661fSBarry Smith PetscFree(mdn->recv_waits); PetscFree(mdn->rvalues); 2198965ea79SLois Curfman McInnes 2208965ea79SLois Curfman McInnes /* wait on sends */ 22139ddd567SLois Curfman McInnes if (mdn->nsends) { 2220452661fSBarry Smith send_status = (MPI_Status *) PetscMalloc( mdn->nsends*sizeof(MPI_Status) ); 2238965ea79SLois Curfman McInnes CHKPTRQ(send_status); 22439ddd567SLois Curfman McInnes MPI_Waitall(mdn->nsends,mdn->send_waits,send_status); 2250452661fSBarry Smith PetscFree(send_status); 2268965ea79SLois Curfman McInnes } 2270452661fSBarry Smith PetscFree(mdn->send_waits); PetscFree(mdn->svalues); 2288965ea79SLois Curfman McInnes 22939ddd567SLois Curfman McInnes mdn->insertmode = NOT_SET_VALUES; 23039ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode); CHKERRQ(ierr); 23139ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode); CHKERRQ(ierr); 2328965ea79SLois Curfman McInnes 2336d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 23439ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat); CHKERRQ(ierr); 2358965ea79SLois Curfman McInnes } 2368965ea79SLois Curfman McInnes return 0; 2378965ea79SLois Curfman McInnes } 2388965ea79SLois Curfman McInnes 23939ddd567SLois Curfman McInnes static int MatZeroEntries_MPIDense(Mat A) 2408965ea79SLois Curfman McInnes { 24139ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense *) A->data; 24239ddd567SLois Curfman McInnes return MatZeroEntries(l->A); 2438965ea79SLois Curfman McInnes } 2448965ea79SLois Curfman McInnes 2458965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 2468965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 2478965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 2483501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 2498965ea79SLois Curfman McInnes routine. 2508965ea79SLois Curfman McInnes */ 25139ddd567SLois Curfman McInnes static int MatZeroRows_MPIDense(Mat A,IS is,Scalar *diag) 2528965ea79SLois Curfman McInnes { 25339ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense *) A->data; 2548965ea79SLois Curfman McInnes int i,ierr,N, *rows,*owners = l->rowners,size = l->size; 2558965ea79SLois Curfman McInnes int *procs,*nprocs,j,found,idx,nsends,*work; 2568965ea79SLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 2578965ea79SLois Curfman McInnes int *rvalues,tag = A->tag,count,base,slen,n,*source; 2588965ea79SLois Curfman McInnes int *lens,imdex,*lrows,*values; 2598965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 2608965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 2618965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 2628965ea79SLois Curfman McInnes IS istmp; 2638965ea79SLois Curfman McInnes 26477c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 2658965ea79SLois Curfman McInnes ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 2668965ea79SLois Curfman McInnes 2678965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 2680452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 269cddf8d76SBarry Smith PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 2700452661fSBarry Smith owner = (int *) PetscMalloc((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/ 2718965ea79SLois Curfman McInnes for ( i=0; i<N; i++ ) { 2728965ea79SLois Curfman McInnes idx = rows[i]; 2738965ea79SLois Curfman McInnes found = 0; 2748965ea79SLois Curfman McInnes for ( j=0; j<size; j++ ) { 2758965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 2768965ea79SLois Curfman McInnes nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 2778965ea79SLois Curfman McInnes } 2788965ea79SLois Curfman McInnes } 27939ddd567SLois Curfman McInnes if (!found) SETERRQ(1,"MatZeroRows_MPIDense:Index out of range"); 2808965ea79SLois Curfman McInnes } 2818965ea79SLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 2828965ea79SLois Curfman McInnes 2838965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 2840452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 2858965ea79SLois Curfman McInnes MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm); 2868965ea79SLois Curfman McInnes nrecvs = work[rank]; 2878965ea79SLois Curfman McInnes MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm); 2888965ea79SLois Curfman McInnes nmax = work[rank]; 2890452661fSBarry Smith PetscFree(work); 2908965ea79SLois Curfman McInnes 2918965ea79SLois Curfman McInnes /* post receives: */ 2920452661fSBarry Smith rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int)); /*see note */ 2938965ea79SLois Curfman McInnes CHKPTRQ(rvalues); 2940452661fSBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request)); 2958965ea79SLois Curfman McInnes CHKPTRQ(recv_waits); 2968965ea79SLois Curfman McInnes for ( i=0; i<nrecvs; i++ ) { 2978965ea79SLois Curfman McInnes MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i); 2988965ea79SLois Curfman McInnes } 2998965ea79SLois Curfman McInnes 3008965ea79SLois Curfman McInnes /* do sends: 3018965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 3028965ea79SLois Curfman McInnes the ith processor 3038965ea79SLois Curfman McInnes */ 3040452661fSBarry Smith svalues = (int *) PetscMalloc( (N+1)*sizeof(int) ); CHKPTRQ(svalues); 3050452661fSBarry Smith send_waits = (MPI_Request *) PetscMalloc( (nsends+1)*sizeof(MPI_Request)); 3068965ea79SLois Curfman McInnes CHKPTRQ(send_waits); 3070452661fSBarry Smith starts = (int *) PetscMalloc( (size+1)*sizeof(int) ); CHKPTRQ(starts); 3088965ea79SLois Curfman McInnes starts[0] = 0; 3098965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3108965ea79SLois Curfman McInnes for ( i=0; i<N; i++ ) { 3118965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 3128965ea79SLois Curfman McInnes } 3138965ea79SLois Curfman McInnes ISRestoreIndices(is,&rows); 3148965ea79SLois Curfman McInnes 3158965ea79SLois Curfman McInnes starts[0] = 0; 3168965ea79SLois Curfman McInnes for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3178965ea79SLois Curfman McInnes count = 0; 3188965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 3198965ea79SLois Curfman McInnes if (procs[i]) { 3208965ea79SLois Curfman McInnes MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++); 3218965ea79SLois Curfman McInnes } 3228965ea79SLois Curfman McInnes } 3230452661fSBarry Smith PetscFree(starts); 3248965ea79SLois Curfman McInnes 3258965ea79SLois Curfman McInnes base = owners[rank]; 3268965ea79SLois Curfman McInnes 3278965ea79SLois Curfman McInnes /* wait on receives */ 3280452661fSBarry Smith lens = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens); 3298965ea79SLois Curfman McInnes source = lens + nrecvs; 3308965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 3318965ea79SLois Curfman McInnes while (count) { 3328965ea79SLois Curfman McInnes MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status); 3338965ea79SLois Curfman McInnes /* unpack receives into our local space */ 3348965ea79SLois Curfman McInnes MPI_Get_count(&recv_status,MPI_INT,&n); 3358965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 3368965ea79SLois Curfman McInnes lens[imdex] = n; 3378965ea79SLois Curfman McInnes slen += n; 3388965ea79SLois Curfman McInnes count--; 3398965ea79SLois Curfman McInnes } 3400452661fSBarry Smith PetscFree(recv_waits); 3418965ea79SLois Curfman McInnes 3428965ea79SLois Curfman McInnes /* move the data into the send scatter */ 3430452661fSBarry Smith lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) ); CHKPTRQ(lrows); 3448965ea79SLois Curfman McInnes count = 0; 3458965ea79SLois Curfman McInnes for ( i=0; i<nrecvs; i++ ) { 3468965ea79SLois Curfman McInnes values = rvalues + i*nmax; 3478965ea79SLois Curfman McInnes for ( j=0; j<lens[i]; j++ ) { 3488965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 3498965ea79SLois Curfman McInnes } 3508965ea79SLois Curfman McInnes } 3510452661fSBarry Smith PetscFree(rvalues); PetscFree(lens); 3520452661fSBarry Smith PetscFree(owner); PetscFree(nprocs); 3538965ea79SLois Curfman McInnes 3548965ea79SLois Curfman McInnes /* actually zap the local rows */ 3558965ea79SLois Curfman McInnes ierr = ISCreateSeq(MPI_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 3568965ea79SLois Curfman McInnes PLogObjectParent(A,istmp); 3570452661fSBarry Smith PetscFree(lrows); 3588965ea79SLois Curfman McInnes ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr); 3598965ea79SLois Curfman McInnes ierr = ISDestroy(istmp); CHKERRQ(ierr); 3608965ea79SLois Curfman McInnes 3618965ea79SLois Curfman McInnes /* wait on sends */ 3628965ea79SLois Curfman McInnes if (nsends) { 3630452661fSBarry Smith send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status)); 3648965ea79SLois Curfman McInnes CHKPTRQ(send_status); 3658965ea79SLois Curfman McInnes MPI_Waitall(nsends,send_waits,send_status); 3660452661fSBarry Smith PetscFree(send_status); 3678965ea79SLois Curfman McInnes } 3680452661fSBarry Smith PetscFree(send_waits); PetscFree(svalues); 3698965ea79SLois Curfman McInnes 3708965ea79SLois Curfman McInnes return 0; 3718965ea79SLois Curfman McInnes } 3728965ea79SLois Curfman McInnes 37339ddd567SLois Curfman McInnes static int MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 3748965ea79SLois Curfman McInnes { 37539ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3768965ea79SLois Curfman McInnes int ierr; 377c456f294SBarry Smith 37839ddd567SLois Curfman McInnes ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx); 37939ddd567SLois Curfman McInnes CHKERRQ(ierr); 38039ddd567SLois Curfman McInnes ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx); 38139ddd567SLois Curfman McInnes CHKERRQ(ierr); 38244cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy); CHKERRQ(ierr); 3838965ea79SLois Curfman McInnes return 0; 3848965ea79SLois Curfman McInnes } 3858965ea79SLois Curfman McInnes 38639ddd567SLois Curfman McInnes static int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 3878965ea79SLois Curfman McInnes { 38839ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3898965ea79SLois Curfman McInnes int ierr; 390c456f294SBarry Smith 391c456f294SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx);CHKERRQ(ierr); 392c456f294SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_ALL,mdn->Mvctx);CHKERRQ(ierr); 39344cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz); CHKERRQ(ierr); 3948965ea79SLois Curfman McInnes return 0; 3958965ea79SLois Curfman McInnes } 3968965ea79SLois Curfman McInnes 397096963f5SLois Curfman McInnes static int MatMultTrans_MPIDense(Mat A,Vec xx,Vec yy) 398096963f5SLois Curfman McInnes { 399096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 400096963f5SLois Curfman McInnes int ierr; 4013501a2bdSLois Curfman McInnes Scalar zero = 0.0; 402096963f5SLois Curfman McInnes 4033501a2bdSLois Curfman McInnes ierr = VecSet(&zero,yy); CHKERRQ(ierr); 40444cd7ae7SLois Curfman McInnes ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec); CHKERRQ(ierr); 405096963f5SLois Curfman McInnes ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES, 406096963f5SLois Curfman McInnes (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr); 407096963f5SLois Curfman McInnes ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES, 408096963f5SLois Curfman McInnes (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr); 409096963f5SLois Curfman McInnes return 0; 410096963f5SLois Curfman McInnes } 411096963f5SLois Curfman McInnes 412096963f5SLois Curfman McInnes static int MatMultTransAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 413096963f5SLois Curfman McInnes { 414096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 415096963f5SLois Curfman McInnes int ierr; 416096963f5SLois Curfman McInnes 4173501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz); CHKERRQ(ierr); 41844cd7ae7SLois Curfman McInnes ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec); CHKERRQ(ierr); 419096963f5SLois Curfman McInnes ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES, 420096963f5SLois Curfman McInnes (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr); 421096963f5SLois Curfman McInnes ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES, 422096963f5SLois Curfman McInnes (ScatterMode)(SCATTER_ALL|SCATTER_REVERSE),a->Mvctx); CHKERRQ(ierr); 423096963f5SLois Curfman McInnes return 0; 424096963f5SLois Curfman McInnes } 425096963f5SLois Curfman McInnes 42639ddd567SLois Curfman McInnes static int MatGetDiagonal_MPIDense(Mat A,Vec v) 4278965ea79SLois Curfman McInnes { 42839ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 429096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense *) a->A->data; 43044cd7ae7SLois Curfman McInnes int ierr, len, i, n, m = a->m, radd; 43144cd7ae7SLois Curfman McInnes Scalar *x, zero = 0.0; 432ed3cc1f0SBarry Smith 43344cd7ae7SLois Curfman McInnes VecSet(&zero,v); 434096963f5SLois Curfman McInnes ierr = VecGetArray(v,&x); CHKERRQ(ierr); 435096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n); CHKERRQ(ierr); 436096963f5SLois Curfman McInnes if (n != a->M) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec"); 43744cd7ae7SLois Curfman McInnes len = PetscMin(aloc->m,aloc->n); 4387ddc982cSLois Curfman McInnes radd = a->rstart*m; 43944cd7ae7SLois Curfman McInnes for ( i=0; i<len; i++ ) { 440096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 441096963f5SLois Curfman McInnes } 442096963f5SLois Curfman McInnes return 0; 4438965ea79SLois Curfman McInnes } 4448965ea79SLois Curfman McInnes 44539ddd567SLois Curfman McInnes static int MatDestroy_MPIDense(PetscObject obj) 4468965ea79SLois Curfman McInnes { 4478965ea79SLois Curfman McInnes Mat mat = (Mat) obj; 4483501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 4498965ea79SLois Curfman McInnes int ierr; 450ed3cc1f0SBarry Smith 4518965ea79SLois Curfman McInnes #if defined(PETSC_LOG) 4523501a2bdSLois Curfman McInnes PLogObjectState(obj,"Rows=%d, Cols=%d",mdn->M,mdn->N); 4538965ea79SLois Curfman McInnes #endif 4540452661fSBarry Smith PetscFree(mdn->rowners); 4553501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A); CHKERRQ(ierr); 4563501a2bdSLois Curfman McInnes if (mdn->lvec) VecDestroy(mdn->lvec); 4573501a2bdSLois Curfman McInnes if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 458622d7880SLois Curfman McInnes if (mdn->factor) { 459622d7880SLois Curfman McInnes if (mdn->factor->temp) PetscFree(mdn->factor->temp); 460622d7880SLois Curfman McInnes if (mdn->factor->tag) PetscFree(mdn->factor->tag); 461622d7880SLois Curfman McInnes if (mdn->factor->pivots) PetscFree(mdn->factor->pivots); 462622d7880SLois Curfman McInnes PetscFree(mdn->factor); 463622d7880SLois Curfman McInnes } 4640452661fSBarry Smith PetscFree(mdn); 4658965ea79SLois Curfman McInnes PLogObjectDestroy(mat); 4660452661fSBarry Smith PetscHeaderDestroy(mat); 4678965ea79SLois Curfman McInnes return 0; 4688965ea79SLois Curfman McInnes } 46939ddd567SLois Curfman McInnes 4708965ea79SLois Curfman McInnes #include "pinclude/pviewer.h" 4718965ea79SLois Curfman McInnes 47239ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer) 4738965ea79SLois Curfman McInnes { 47439ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 4758965ea79SLois Curfman McInnes int ierr; 47639ddd567SLois Curfman McInnes if (mdn->size == 1) { 47739ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer); CHKERRQ(ierr); 4788965ea79SLois Curfman McInnes } 47939ddd567SLois Curfman McInnes else SETERRQ(1,"MatView_MPIDense_Binary:Only uniprocessor output supported"); 4808965ea79SLois Curfman McInnes return 0; 4818965ea79SLois Curfman McInnes } 4828965ea79SLois Curfman McInnes 48339ddd567SLois Curfman McInnes static int MatView_MPIDense_ASCII(Mat mat,Viewer viewer) 4848965ea79SLois Curfman McInnes { 48539ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 48639ddd567SLois Curfman McInnes int ierr, format; 4878965ea79SLois Curfman McInnes FILE *fd; 48819bcc07fSBarry Smith ViewerType vtype; 4898965ea79SLois Curfman McInnes 49019bcc07fSBarry Smith ViewerGetType(viewer,&vtype); 49190ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 49290ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 49390ace30eSBarry Smith if (format == ASCII_FORMAT_INFO_DETAILED) { 494096963f5SLois Curfman McInnes int nz, nzalloc, mem, rank; 495096963f5SLois Curfman McInnes MPI_Comm_rank(mat->comm,&rank); 496096963f5SLois Curfman McInnes ierr = MatGetInfo(mat,MAT_LOCAL,&nz,&nzalloc,&mem); 49777c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 4983501a2bdSLois Curfman McInnes fprintf(fd," [%d] local rows %d nz %d nz alloced %d mem %d \n", 499096963f5SLois Curfman McInnes rank,mdn->m,nz,nzalloc,mem); 500096963f5SLois Curfman McInnes fflush(fd); 50177c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 5023501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer); CHKERRQ(ierr); 5033501a2bdSLois Curfman McInnes return 0; 5043501a2bdSLois Curfman McInnes } 50590ace30eSBarry Smith else if (format == ASCII_FORMAT_INFO) { 506096963f5SLois Curfman McInnes return 0; 5078965ea79SLois Curfman McInnes } 50819bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER) { 50977c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 51039ddd567SLois Curfman McInnes fprintf(fd,"[%d] rows %d starts %d ends %d cols %d\n", 51139ddd567SLois Curfman McInnes mdn->rank,mdn->m,mdn->rstart,mdn->rend,mdn->n); 51239ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer); CHKERRQ(ierr); 5138965ea79SLois Curfman McInnes fflush(fd); 51477c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 5158965ea79SLois Curfman McInnes } 5168965ea79SLois Curfman McInnes else { 51739ddd567SLois Curfman McInnes int size = mdn->size, rank = mdn->rank; 5188965ea79SLois Curfman McInnes if (size == 1) { 51939ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer); CHKERRQ(ierr); 5208965ea79SLois Curfman McInnes } 5218965ea79SLois Curfman McInnes else { 5228965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 5238965ea79SLois Curfman McInnes Mat A; 52439ddd567SLois Curfman McInnes int M = mdn->M, N = mdn->N,m,row,i, nz, *cols; 52539ddd567SLois Curfman McInnes Scalar *vals; 52639ddd567SLois Curfman McInnes Mat_SeqDense *Amdn = (Mat_SeqDense*) mdn->A->data; 5278965ea79SLois Curfman McInnes 5288965ea79SLois Curfman McInnes if (!rank) { 529b4fd4287SBarry Smith ierr = MatCreateMPIDense(mat->comm,M,M,N,N,PETSC_NULL,&A); CHKERRQ(ierr); 5308965ea79SLois Curfman McInnes } 5318965ea79SLois Curfman McInnes else { 532b4fd4287SBarry Smith ierr = MatCreateMPIDense(mat->comm,0,M,N,N,PETSC_NULL,&A); CHKERRQ(ierr); 5338965ea79SLois Curfman McInnes } 5348965ea79SLois Curfman McInnes PLogObjectParent(mat,A); 5358965ea79SLois Curfman McInnes 53639ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 53739ddd567SLois Curfman McInnes but it's quick for now */ 53839ddd567SLois Curfman McInnes row = mdn->rstart; m = Amdn->m; 5398965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 54039ddd567SLois Curfman McInnes ierr = MatGetRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr); 54139ddd567SLois Curfman McInnes ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES); CHKERRQ(ierr); 54239ddd567SLois Curfman McInnes ierr = MatRestoreRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr); 54339ddd567SLois Curfman McInnes row++; 5448965ea79SLois Curfman McInnes } 5458965ea79SLois Curfman McInnes 5466d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 5476d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 5488965ea79SLois Curfman McInnes if (!rank) { 54939ddd567SLois Curfman McInnes ierr = MatView(((Mat_MPIDense*)(A->data))->A,viewer); CHKERRQ(ierr); 5508965ea79SLois Curfman McInnes } 5518965ea79SLois Curfman McInnes ierr = MatDestroy(A); CHKERRQ(ierr); 5528965ea79SLois Curfman McInnes } 5538965ea79SLois Curfman McInnes } 5548965ea79SLois Curfman McInnes return 0; 5558965ea79SLois Curfman McInnes } 5568965ea79SLois Curfman McInnes 55739ddd567SLois Curfman McInnes static int MatView_MPIDense(PetscObject obj,Viewer viewer) 5588965ea79SLois Curfman McInnes { 5598965ea79SLois Curfman McInnes Mat mat = (Mat) obj; 56039ddd567SLois Curfman McInnes int ierr; 561bcd2baecSBarry Smith ViewerType vtype; 5628965ea79SLois Curfman McInnes 563bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 564bcd2baecSBarry Smith if (vtype == ASCII_FILE_VIEWER) { 56539ddd567SLois Curfman McInnes ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr); 5668965ea79SLois Curfman McInnes } 567bcd2baecSBarry Smith else if (vtype == ASCII_FILES_VIEWER) { 56839ddd567SLois Curfman McInnes ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr); 5698965ea79SLois Curfman McInnes } 570bcd2baecSBarry Smith else if (vtype == BINARY_FILE_VIEWER) { 57139ddd567SLois Curfman McInnes return MatView_MPIDense_Binary(mat,viewer); 5728965ea79SLois Curfman McInnes } 5738965ea79SLois Curfman McInnes return 0; 5748965ea79SLois Curfman McInnes } 5758965ea79SLois Curfman McInnes 5763501a2bdSLois Curfman McInnes static int MatGetInfo_MPIDense(Mat A,MatInfoType flag,int *nz, 5778965ea79SLois Curfman McInnes int *nzalloc,int *mem) 5788965ea79SLois Curfman McInnes { 5793501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 5803501a2bdSLois Curfman McInnes Mat mdn = mat->A; 58139ddd567SLois Curfman McInnes int ierr, isend[3], irecv[3]; 5828965ea79SLois Curfman McInnes 5833501a2bdSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,&isend[0],&isend[1],&isend[2]); CHKERRQ(ierr); 5848965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 585bcd2baecSBarry Smith if (nz) *nz = isend[0]; 586bcd2baecSBarry Smith if (nzalloc) *nzalloc = isend[1]; 587bcd2baecSBarry Smith if (mem) *mem = isend[2]; 5888965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 5893501a2bdSLois Curfman McInnes MPI_Allreduce(isend,irecv,3,MPI_INT,MPI_MAX,A->comm); 590bcd2baecSBarry Smith if (nz) *nz = irecv[0]; 591bcd2baecSBarry Smith if (nzalloc) *nzalloc = irecv[1]; 592bcd2baecSBarry Smith if (mem) *mem = irecv[2]; 5938965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 5943501a2bdSLois Curfman McInnes MPI_Allreduce(isend,irecv,3,MPI_INT,MPI_SUM,A->comm); 595bcd2baecSBarry Smith if (nz) *nz = irecv[0]; 596bcd2baecSBarry Smith if (nzalloc) *nzalloc = irecv[1]; 597bcd2baecSBarry Smith if (mem) *mem = irecv[2]; 5988965ea79SLois Curfman McInnes } 5998965ea79SLois Curfman McInnes return 0; 6008965ea79SLois Curfman McInnes } 6018965ea79SLois Curfman McInnes 6028c469469SLois Curfman McInnes /* extern int MatLUFactorSymbolic_MPIDense(Mat,IS,IS,double,Mat*); 6038aaee692SLois Curfman McInnes extern int MatLUFactorNumeric_MPIDense(Mat,Mat*); 6048aaee692SLois Curfman McInnes extern int MatLUFactor_MPIDense(Mat,IS,IS,double); 6058aaee692SLois Curfman McInnes extern int MatSolve_MPIDense(Mat,Vec,Vec); 6068c469469SLois Curfman McInnes extern int MatSolveAdd_MPIDense(Mat,Vec,Vec,Vec); 6078aaee692SLois Curfman McInnes extern int MatSolveTrans_MPIDense(Mat,Vec,Vec); 6088aaee692SLois Curfman McInnes extern int MatSolveTransAdd_MPIDense(Mat,Vec,Vec,Vec); */ 6098aaee692SLois Curfman McInnes 61039ddd567SLois Curfman McInnes static int MatSetOption_MPIDense(Mat A,MatOption op) 6118965ea79SLois Curfman McInnes { 61239ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 6138965ea79SLois Curfman McInnes 6146d4a8577SBarry Smith if (op == MAT_NO_NEW_NONZERO_LOCATIONS || 6156d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 6166d4a8577SBarry Smith op == MAT_COLUMNS_SORTED || 6176d4a8577SBarry Smith op == MAT_ROW_ORIENTED) { 6188965ea79SLois Curfman McInnes MatSetOption(a->A,op); 6198965ea79SLois Curfman McInnes } 6206d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 6216d4a8577SBarry Smith op == MAT_SYMMETRIC || 6226d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 6236d4a8577SBarry Smith op == MAT_YES_NEW_DIAGONALS) 62494a424c1SBarry Smith PLogInfo(A,"Info:MatSetOption_MPIDense:Option ignored\n"); 6256d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) 62639b7565bSBarry Smith { a->roworiented = 0; MatSetOption(a->A,op);} 6276d4a8577SBarry Smith else if (op == MAT_NO_NEW_DIAGONALS) 6286d4a8577SBarry Smith {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIDense:MAT_NO_NEW_DIAGONALS");} 6298965ea79SLois Curfman McInnes else 63039ddd567SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_MPIDense:unknown option");} 6318965ea79SLois Curfman McInnes return 0; 6328965ea79SLois Curfman McInnes } 6338965ea79SLois Curfman McInnes 6343501a2bdSLois Curfman McInnes static int MatGetSize_MPIDense(Mat A,int *m,int *n) 6358965ea79SLois Curfman McInnes { 6363501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6378965ea79SLois Curfman McInnes *m = mat->M; *n = mat->N; 6388965ea79SLois Curfman McInnes return 0; 6398965ea79SLois Curfman McInnes } 6408965ea79SLois Curfman McInnes 6413501a2bdSLois Curfman McInnes static int MatGetLocalSize_MPIDense(Mat A,int *m,int *n) 6428965ea79SLois Curfman McInnes { 6433501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6448965ea79SLois Curfman McInnes *m = mat->m; *n = mat->N; 6458965ea79SLois Curfman McInnes return 0; 6468965ea79SLois Curfman McInnes } 6478965ea79SLois Curfman McInnes 6483501a2bdSLois Curfman McInnes static int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n) 6498965ea79SLois Curfman McInnes { 6503501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6518965ea79SLois Curfman McInnes *m = mat->rstart; *n = mat->rend; 6528965ea79SLois Curfman McInnes return 0; 6538965ea79SLois Curfman McInnes } 6548965ea79SLois Curfman McInnes 6553501a2bdSLois Curfman McInnes static int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v) 6568965ea79SLois Curfman McInnes { 6573501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 65839ddd567SLois Curfman McInnes int lrow, rstart = mat->rstart, rend = mat->rend; 6598965ea79SLois Curfman McInnes 66039ddd567SLois Curfman McInnes if (row < rstart || row >= rend) SETERRQ(1,"MatGetRow_MPIDense:only local rows") 6618965ea79SLois Curfman McInnes lrow = row - rstart; 66239ddd567SLois Curfman McInnes return MatGetRow(mat->A,lrow,nz,idx,v); 6638965ea79SLois Curfman McInnes } 6648965ea79SLois Curfman McInnes 66539ddd567SLois Curfman McInnes static int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v) 6668965ea79SLois Curfman McInnes { 6670452661fSBarry Smith if (idx) PetscFree(*idx); 6680452661fSBarry Smith if (v) PetscFree(*v); 6698965ea79SLois Curfman McInnes return 0; 6708965ea79SLois Curfman McInnes } 6718965ea79SLois Curfman McInnes 672cddf8d76SBarry Smith static int MatNorm_MPIDense(Mat A,NormType type,double *norm) 673096963f5SLois Curfman McInnes { 6743501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) A->data; 6753501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*) mdn->A->data; 6763501a2bdSLois Curfman McInnes int ierr, i, j; 6773501a2bdSLois Curfman McInnes double sum = 0.0; 6783501a2bdSLois Curfman McInnes Scalar *v = mat->v; 6793501a2bdSLois Curfman McInnes 6803501a2bdSLois Curfman McInnes if (mdn->size == 1) { 6813501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,norm); CHKERRQ(ierr); 6823501a2bdSLois Curfman McInnes } else { 6833501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 6843501a2bdSLois Curfman McInnes for (i=0; i<mat->n*mat->m; i++ ) { 6853501a2bdSLois Curfman McInnes #if defined(PETSC_COMPLEX) 6863501a2bdSLois Curfman McInnes sum += real(conj(*v)*(*v)); v++; 6873501a2bdSLois Curfman McInnes #else 6883501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 6893501a2bdSLois Curfman McInnes #endif 6903501a2bdSLois Curfman McInnes } 6916d84be18SBarry Smith MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm); 6923501a2bdSLois Curfman McInnes *norm = sqrt(*norm); 6933501a2bdSLois Curfman McInnes PLogFlops(2*mat->n*mat->m); 6943501a2bdSLois Curfman McInnes } 6953501a2bdSLois Curfman McInnes else if (type == NORM_1) { 6963501a2bdSLois Curfman McInnes double *tmp, *tmp2; 6970452661fSBarry Smith tmp = (double *) PetscMalloc( 2*mdn->N*sizeof(double) ); CHKPTRQ(tmp); 6983501a2bdSLois Curfman McInnes tmp2 = tmp + mdn->N; 699cddf8d76SBarry Smith PetscMemzero(tmp,2*mdn->N*sizeof(double)); 700096963f5SLois Curfman McInnes *norm = 0.0; 7013501a2bdSLois Curfman McInnes v = mat->v; 7023501a2bdSLois Curfman McInnes for ( j=0; j<mat->n; j++ ) { 7033501a2bdSLois Curfman McInnes for ( i=0; i<mat->m; i++ ) { 70467e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 7053501a2bdSLois Curfman McInnes } 7063501a2bdSLois Curfman McInnes } 7076d84be18SBarry Smith MPI_Allreduce(tmp,tmp2,mdn->N,MPI_DOUBLE,MPI_SUM,A->comm); 7083501a2bdSLois Curfman McInnes for ( j=0; j<mdn->N; j++ ) { 7093501a2bdSLois Curfman McInnes if (tmp2[j] > *norm) *norm = tmp2[j]; 7103501a2bdSLois Curfman McInnes } 7110452661fSBarry Smith PetscFree(tmp); 7123501a2bdSLois Curfman McInnes PLogFlops(mat->n*mat->m); 7133501a2bdSLois Curfman McInnes } 7143501a2bdSLois Curfman McInnes else if (type == NORM_INFINITY) { /* max row norm */ 7153501a2bdSLois Curfman McInnes double ntemp; 7163501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp); CHKERRQ(ierr); 7176d84be18SBarry Smith MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm); 7183501a2bdSLois Curfman McInnes } 7193501a2bdSLois Curfman McInnes else { 7203501a2bdSLois Curfman McInnes SETERRQ(1,"MatNorm_MPIDense:No support for two norm"); 7213501a2bdSLois Curfman McInnes } 7223501a2bdSLois Curfman McInnes } 7233501a2bdSLois Curfman McInnes return 0; 7243501a2bdSLois Curfman McInnes } 7253501a2bdSLois Curfman McInnes 7263501a2bdSLois Curfman McInnes static int MatTranspose_MPIDense(Mat A,Mat *matout) 7273501a2bdSLois Curfman McInnes { 7283501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 7293501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense *) a->A->data; 7303501a2bdSLois Curfman McInnes Mat B; 7313501a2bdSLois Curfman McInnes int M = a->M, N = a->N, m, n, *rwork, rstart = a->rstart; 7323501a2bdSLois Curfman McInnes int j, i, ierr; 7333501a2bdSLois Curfman McInnes Scalar *v; 7343501a2bdSLois Curfman McInnes 7353638b69dSLois Curfman McInnes if (matout == PETSC_NULL && M != N) 7363501a2bdSLois Curfman McInnes SETERRQ(1,"MatTranspose_MPIDense:Supports square matrix only in-place"); 737b4fd4287SBarry Smith ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B); 738ed2daf61SLois Curfman McInnes CHKERRQ(ierr); 7393501a2bdSLois Curfman McInnes 7403501a2bdSLois Curfman McInnes m = Aloc->m; n = Aloc->n; v = Aloc->v; 7410452661fSBarry Smith rwork = (int *) PetscMalloc(n*sizeof(int)); CHKPTRQ(rwork); 7423501a2bdSLois Curfman McInnes for ( j=0; j<n; j++ ) { 7433501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 7443501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES); CHKERRQ(ierr); 7453501a2bdSLois Curfman McInnes v += m; 7463501a2bdSLois Curfman McInnes } 7470452661fSBarry Smith PetscFree(rwork); 7486d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 7496d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 7503638b69dSLois Curfman McInnes if (matout != PETSC_NULL) { 7513501a2bdSLois Curfman McInnes *matout = B; 7523501a2bdSLois Curfman McInnes } else { 7533501a2bdSLois Curfman McInnes /* This isn't really an in-place transpose, but free data struct from a */ 7540452661fSBarry Smith PetscFree(a->rowners); 7553501a2bdSLois Curfman McInnes ierr = MatDestroy(a->A); CHKERRQ(ierr); 7563501a2bdSLois Curfman McInnes if (a->lvec) VecDestroy(a->lvec); 7573501a2bdSLois Curfman McInnes if (a->Mvctx) VecScatterDestroy(a->Mvctx); 7580452661fSBarry Smith PetscFree(a); 7593501a2bdSLois Curfman McInnes PetscMemcpy(A,B,sizeof(struct _Mat)); 7600452661fSBarry Smith PetscHeaderDestroy(B); 7613501a2bdSLois Curfman McInnes } 762096963f5SLois Curfman McInnes return 0; 763096963f5SLois Curfman McInnes } 764096963f5SLois Curfman McInnes 76544cd7ae7SLois Curfman McInnes #include "pinclude/plapack.h" 76644cd7ae7SLois Curfman McInnes static int MatScale_MPIDense(Scalar *alpha,Mat inA) 76744cd7ae7SLois Curfman McInnes { 76844cd7ae7SLois Curfman McInnes Mat_MPIDense *A = (Mat_MPIDense *) inA->data; 76944cd7ae7SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->A->data; 77044cd7ae7SLois Curfman McInnes int one = 1, nz; 77144cd7ae7SLois Curfman McInnes 77244cd7ae7SLois Curfman McInnes nz = a->m*a->n; 77344cd7ae7SLois Curfman McInnes BLscal_( &nz, alpha, a->v, &one ); 77444cd7ae7SLois Curfman McInnes PLogFlops(nz); 77544cd7ae7SLois Curfman McInnes return 0; 77644cd7ae7SLois Curfman McInnes } 77744cd7ae7SLois Curfman McInnes 7783d1612f7SBarry Smith static int MatConvertSameType_MPIDense(Mat,Mat *,int); 7798965ea79SLois Curfman McInnes 7808965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 78139ddd567SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIDense, 78239ddd567SLois Curfman McInnes MatGetRow_MPIDense,MatRestoreRow_MPIDense, 78339ddd567SLois Curfman McInnes MatMult_MPIDense,MatMultAdd_MPIDense, 784096963f5SLois Curfman McInnes MatMultTrans_MPIDense,MatMultTransAdd_MPIDense, 785e7ca0642SLois Curfman McInnes /* MatSolve_MPIDense,0, */ 786e7ca0642SLois Curfman McInnes 0,0, 78739ddd567SLois Curfman McInnes 0,0, 7888c469469SLois Curfman McInnes 0,0, 7898c469469SLois Curfman McInnes /* MatLUFactor_MPIDense,0, */ 7908aaee692SLois Curfman McInnes 0,MatTranspose_MPIDense, 79139ddd567SLois Curfman McInnes MatGetInfo_MPIDense,0, 792096963f5SLois Curfman McInnes MatGetDiagonal_MPIDense,0,MatNorm_MPIDense, 79339ddd567SLois Curfman McInnes MatAssemblyBegin_MPIDense,MatAssemblyEnd_MPIDense, 7948965ea79SLois Curfman McInnes 0, 79539ddd567SLois Curfman McInnes MatSetOption_MPIDense,MatZeroEntries_MPIDense,MatZeroRows_MPIDense, 7968c469469SLois Curfman McInnes 0,0,0, 7978c469469SLois Curfman McInnes /* 0,MatLUFactorSymbolic_MPIDense,MatLUFactorNumeric_MPIDense, */ 7988aaee692SLois Curfman McInnes 0,0, 79939ddd567SLois Curfman McInnes MatGetSize_MPIDense,MatGetLocalSize_MPIDense, 80039ddd567SLois Curfman McInnes MatGetOwnershipRange_MPIDense, 801ff14e315SSatish Balay 0,0, MatGetArray_MPIDense, MatRestoreArray_MPIDense, 802905e6a2fSBarry Smith 0,MatConvertSameType_MPIDense, 803b49de8d1SLois Curfman McInnes 0,0,0,0,0, 80444cd7ae7SLois Curfman McInnes 0,0,MatGetValues_MPIDense,0,0,MatScale_MPIDense}; 8058965ea79SLois Curfman McInnes 8068965ea79SLois Curfman McInnes /*@C 80739ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 8088965ea79SLois Curfman McInnes 8098965ea79SLois Curfman McInnes Input Parameters: 8108965ea79SLois Curfman McInnes . comm - MPI communicator 8118965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 8128965ea79SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated 8138965ea79SLois Curfman McInnes if N is given) 8148965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 8158965ea79SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated 8168965ea79SLois Curfman McInnes if n is given) 817b4fd4287SBarry Smith . data - optional location of matrix data. Set data=PETSC_NULL for PETSc 818dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 8198965ea79SLois Curfman McInnes 8208965ea79SLois Curfman McInnes Output Parameter: 821477f1c0bSLois Curfman McInnes . A - the matrix 8228965ea79SLois Curfman McInnes 8238965ea79SLois Curfman McInnes Notes: 82439ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 82539ddd567SLois Curfman McInnes storage by columns. 8268965ea79SLois Curfman McInnes 82718f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 82818f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 829b4fd4287SBarry Smith set data=PETSC_NULL. 83018f449edSLois Curfman McInnes 8318965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 8328965ea79SLois Curfman McInnes (possibly both). 8338965ea79SLois Curfman McInnes 8343501a2bdSLois Curfman McInnes Currently, the only parallel dense matrix decomposition is by rows, 8353501a2bdSLois Curfman McInnes so that n=N and each submatrix owns all of the global columns. 8363501a2bdSLois Curfman McInnes 83739ddd567SLois Curfman McInnes .keywords: matrix, dense, parallel 8388965ea79SLois Curfman McInnes 83939ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 8408965ea79SLois Curfman McInnes @*/ 841477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A) 8428965ea79SLois Curfman McInnes { 8438965ea79SLois Curfman McInnes Mat mat; 84439ddd567SLois Curfman McInnes Mat_MPIDense *a; 84525cdf11fSBarry Smith int ierr, i,flg; 8468965ea79SLois Curfman McInnes 847ed2daf61SLois Curfman McInnes /* Note: For now, when data is specified above, this assumes the user correctly 848ed2daf61SLois Curfman McInnes allocates the local dense storage space. We should add error checking. */ 84918f449edSLois Curfman McInnes 850477f1c0bSLois Curfman McInnes *A = 0; 8510452661fSBarry Smith PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATMPIDENSE,comm); 8528965ea79SLois Curfman McInnes PLogObjectCreate(mat); 8530452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIDense)); CHKPTRQ(a); 8548965ea79SLois Curfman McInnes PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 85539ddd567SLois Curfman McInnes mat->destroy = MatDestroy_MPIDense; 85639ddd567SLois Curfman McInnes mat->view = MatView_MPIDense; 8578965ea79SLois Curfman McInnes mat->factor = 0; 8588965ea79SLois Curfman McInnes 859622d7880SLois Curfman McInnes a->factor = 0; 8608965ea79SLois Curfman McInnes a->insertmode = NOT_SET_VALUES; 8618965ea79SLois Curfman McInnes MPI_Comm_rank(comm,&a->rank); 8628965ea79SLois Curfman McInnes MPI_Comm_size(comm,&a->size); 8638965ea79SLois Curfman McInnes 86439ddd567SLois Curfman McInnes if (M == PETSC_DECIDE) MPI_Allreduce(&m,&M,1,MPI_INT,MPI_SUM,comm); 8658965ea79SLois Curfman McInnes if (m == PETSC_DECIDE) {m = M/a->size + ((M % a->size) > a->rank);} 86639ddd567SLois Curfman McInnes 86739ddd567SLois Curfman McInnes /* each row stores all columns */ 86839ddd567SLois Curfman McInnes if (N == PETSC_DECIDE) N = n; 869d7e8b826SBarry Smith if (n == PETSC_DECIDE) {n = N/a->size + ((N % a->size) > a->rank);} 870d7e8b826SBarry Smith /* if (n != N) SETERRQ(1,"MatCreateMPIDense:For now, only n=N is supported"); */ 871aca0ad90SLois Curfman McInnes a->N = mat->N = N; 872aca0ad90SLois Curfman McInnes a->M = mat->M = M; 873aca0ad90SLois Curfman McInnes a->m = mat->m = m; 874aca0ad90SLois Curfman McInnes a->n = mat->n = n; 8758965ea79SLois Curfman McInnes 8768965ea79SLois Curfman McInnes /* build local table of row and column ownerships */ 877d7e8b826SBarry Smith a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int)); CHKPTRQ(a->rowners); 878d7e8b826SBarry Smith a->cowners = a->rowners + a->size + 1; 879d7e8b826SBarry Smith PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _Mat)+ 88039ddd567SLois Curfman McInnes sizeof(Mat_MPIDense)); 8818965ea79SLois Curfman McInnes MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm); 8828965ea79SLois Curfman McInnes a->rowners[0] = 0; 8838965ea79SLois Curfman McInnes for ( i=2; i<=a->size; i++ ) { 8848965ea79SLois Curfman McInnes a->rowners[i] += a->rowners[i-1]; 8858965ea79SLois Curfman McInnes } 8868965ea79SLois Curfman McInnes a->rstart = a->rowners[a->rank]; 8878965ea79SLois Curfman McInnes a->rend = a->rowners[a->rank+1]; 888d7e8b826SBarry Smith MPI_Allgather(&n,1,MPI_INT,a->cowners+1,1,MPI_INT,comm); 889d7e8b826SBarry Smith a->cowners[0] = 0; 890d7e8b826SBarry Smith for ( i=2; i<=a->size; i++ ) { 891d7e8b826SBarry Smith a->cowners[i] += a->cowners[i-1]; 892d7e8b826SBarry Smith } 8938965ea79SLois Curfman McInnes 89418f449edSLois Curfman McInnes ierr = MatCreateSeqDense(MPI_COMM_SELF,m,N,data,&a->A); CHKERRQ(ierr); 8958965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 8968965ea79SLois Curfman McInnes 8978965ea79SLois Curfman McInnes /* build cache for off array entries formed */ 8988965ea79SLois Curfman McInnes ierr = StashBuild_Private(&a->stash); CHKERRQ(ierr); 8998965ea79SLois Curfman McInnes 9008965ea79SLois Curfman McInnes /* stuff used for matrix vector multiply */ 9018965ea79SLois Curfman McInnes a->lvec = 0; 9028965ea79SLois Curfman McInnes a->Mvctx = 0; 90339b7565bSBarry Smith a->roworiented = 1; 9048965ea79SLois Curfman McInnes 905477f1c0bSLois Curfman McInnes *A = mat; 90625cdf11fSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 90725cdf11fSBarry Smith if (flg) { 9088c469469SLois Curfman McInnes ierr = MatPrintHelp(mat); CHKERRQ(ierr); 9098c469469SLois Curfman McInnes } 9108965ea79SLois Curfman McInnes return 0; 9118965ea79SLois Curfman McInnes } 9128965ea79SLois Curfman McInnes 9133d1612f7SBarry Smith static int MatConvertSameType_MPIDense(Mat A,Mat *newmat,int cpvalues) 9148965ea79SLois Curfman McInnes { 9158965ea79SLois Curfman McInnes Mat mat; 9163501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense *) A->data; 91739ddd567SLois Curfman McInnes int ierr; 9182ba99913SLois Curfman McInnes FactorCtx *factor; 9198965ea79SLois Curfman McInnes 9208965ea79SLois Curfman McInnes *newmat = 0; 9210452661fSBarry Smith PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATMPIDENSE,A->comm); 9228965ea79SLois Curfman McInnes PLogObjectCreate(mat); 9230452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIDense)); CHKPTRQ(a); 9248965ea79SLois Curfman McInnes PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 92539ddd567SLois Curfman McInnes mat->destroy = MatDestroy_MPIDense; 92639ddd567SLois Curfman McInnes mat->view = MatView_MPIDense; 9273501a2bdSLois Curfman McInnes mat->factor = A->factor; 928c456f294SBarry Smith mat->assembled = PETSC_TRUE; 9298965ea79SLois Curfman McInnes 93044cd7ae7SLois Curfman McInnes a->m = mat->m = oldmat->m; 93144cd7ae7SLois Curfman McInnes a->n = mat->n = oldmat->n; 93244cd7ae7SLois Curfman McInnes a->M = mat->M = oldmat->M; 93344cd7ae7SLois Curfman McInnes a->N = mat->N = oldmat->N; 9342ba99913SLois Curfman McInnes if (oldmat->factor) { 9352ba99913SLois Curfman McInnes a->factor = (FactorCtx *) (factor = PetscNew(FactorCtx)); CHKPTRQ(factor); 9362ba99913SLois Curfman McInnes /* copy factor contents ... add this code! */ 9372ba99913SLois Curfman McInnes } else a->factor = 0; 9388965ea79SLois Curfman McInnes 9398965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 9408965ea79SLois Curfman McInnes a->rend = oldmat->rend; 9418965ea79SLois Curfman McInnes a->size = oldmat->size; 9428965ea79SLois Curfman McInnes a->rank = oldmat->rank; 9438965ea79SLois Curfman McInnes a->insertmode = NOT_SET_VALUES; 9448965ea79SLois Curfman McInnes 9450452661fSBarry Smith a->rowners = (int *) PetscMalloc((a->size+1)*sizeof(int)); CHKPTRQ(a->rowners); 94639ddd567SLois Curfman McInnes PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_MPIDense)); 9478965ea79SLois Curfman McInnes PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int)); 9488965ea79SLois Curfman McInnes ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr); 9498965ea79SLois Curfman McInnes 9508965ea79SLois Curfman McInnes ierr = VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr); 9518965ea79SLois Curfman McInnes PLogObjectParent(mat,a->lvec); 95255659b69SBarry Smith ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr); 9538965ea79SLois Curfman McInnes PLogObjectParent(mat,a->Mvctx); 9548965ea79SLois Curfman McInnes ierr = MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr); 9558965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 9568965ea79SLois Curfman McInnes *newmat = mat; 9578965ea79SLois Curfman McInnes return 0; 9588965ea79SLois Curfman McInnes } 9598965ea79SLois Curfman McInnes 96077c4ece6SBarry Smith #include "sys.h" 9618965ea79SLois Curfman McInnes 96290ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M, int N, Mat *newmat) 96390ace30eSBarry Smith { 96490ace30eSBarry Smith int *rowners, i,size,rank,m,rstart,rend,ierr,nz,j; 96590ace30eSBarry Smith Scalar *array,*vals,*vals_ptr; 96690ace30eSBarry Smith MPI_Status status; 96790ace30eSBarry Smith 96890ace30eSBarry Smith MPI_Comm_rank(comm,&rank); 96990ace30eSBarry Smith MPI_Comm_size(comm,&size); 97090ace30eSBarry Smith 97190ace30eSBarry Smith /* determine ownership of all rows */ 97290ace30eSBarry Smith m = M/size + ((M % size) > rank); 97390ace30eSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners); 97490ace30eSBarry Smith MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm); 97590ace30eSBarry Smith rowners[0] = 0; 97690ace30eSBarry Smith for ( i=2; i<=size; i++ ) { 97790ace30eSBarry Smith rowners[i] += rowners[i-1]; 97890ace30eSBarry Smith } 97990ace30eSBarry Smith rstart = rowners[rank]; 98090ace30eSBarry Smith rend = rowners[rank+1]; 98190ace30eSBarry Smith 98290ace30eSBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 98390ace30eSBarry Smith ierr = MatGetArray(*newmat,&array); CHKERRQ(ierr); 98490ace30eSBarry Smith 98590ace30eSBarry Smith if (!rank) { 98690ace30eSBarry Smith vals = (Scalar *) PetscMalloc( m*N*sizeof(Scalar) ); CHKPTRQ(vals); 98790ace30eSBarry Smith 98890ace30eSBarry Smith /* read in my part of the matrix numerical values */ 98977c4ece6SBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,BINARY_SCALAR); CHKERRQ(ierr); 99090ace30eSBarry Smith 99190ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 99290ace30eSBarry Smith vals_ptr = vals; 99390ace30eSBarry Smith for ( i=0; i<m; i++ ) { 99490ace30eSBarry Smith for ( j=0; j<N; j++ ) { 99590ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 99690ace30eSBarry Smith } 99790ace30eSBarry Smith } 99890ace30eSBarry Smith 99990ace30eSBarry Smith /* read in other processors and ship out */ 100090ace30eSBarry Smith for ( i=1; i<size; i++ ) { 100190ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 100277c4ece6SBarry Smith ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr); 100390ace30eSBarry Smith MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm); 100490ace30eSBarry Smith } 100590ace30eSBarry Smith } 100690ace30eSBarry Smith else { 100790ace30eSBarry Smith /* receive numeric values */ 100890ace30eSBarry Smith vals = (Scalar*) PetscMalloc( m*N*sizeof(Scalar) ); CHKPTRQ(vals); 100990ace30eSBarry Smith 101090ace30eSBarry Smith /* receive message of values*/ 101190ace30eSBarry Smith MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status); 101290ace30eSBarry Smith 101390ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 101490ace30eSBarry Smith vals_ptr = vals; 101590ace30eSBarry Smith for ( i=0; i<m; i++ ) { 101690ace30eSBarry Smith for ( j=0; j<N; j++ ) { 101790ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 101890ace30eSBarry Smith } 101990ace30eSBarry Smith } 102090ace30eSBarry Smith } 102190ace30eSBarry Smith PetscFree(rowners); 102290ace30eSBarry Smith PetscFree(vals); 10236d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 10246d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 102590ace30eSBarry Smith return 0; 102690ace30eSBarry Smith } 102790ace30eSBarry Smith 102890ace30eSBarry Smith 102919bcc07fSBarry Smith int MatLoad_MPIDense(Viewer viewer,MatType type,Mat *newmat) 10308965ea79SLois Curfman McInnes { 10318965ea79SLois Curfman McInnes Mat A; 10328965ea79SLois Curfman McInnes int i, nz, ierr, j,rstart, rend, fd; 10338965ea79SLois Curfman McInnes Scalar *vals,*svals; 103419bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 10358965ea79SLois Curfman McInnes MPI_Status status; 10368965ea79SLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 10378965ea79SLois Curfman McInnes int *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols; 103819bcc07fSBarry Smith int tag = ((PetscObject)viewer)->tag; 10398965ea79SLois Curfman McInnes 10408965ea79SLois Curfman McInnes MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank); 10418965ea79SLois Curfman McInnes if (!rank) { 104290ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 104377c4ece6SBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,BINARY_INT); CHKERRQ(ierr); 104439ddd567SLois Curfman McInnes if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_MPIDenseorMPIRow:not matrix object"); 10458965ea79SLois Curfman McInnes } 10468965ea79SLois Curfman McInnes 10478965ea79SLois Curfman McInnes MPI_Bcast(header+1,3,MPI_INT,0,comm); 104890ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 104990ace30eSBarry Smith 105090ace30eSBarry Smith /* 105190ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 105290ace30eSBarry Smith */ 105390ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 105490ace30eSBarry Smith return MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat); 105590ace30eSBarry Smith } 105690ace30eSBarry Smith 10578965ea79SLois Curfman McInnes /* determine ownership of all rows */ 10588965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 10590452661fSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners); 10608965ea79SLois Curfman McInnes MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm); 10618965ea79SLois Curfman McInnes rowners[0] = 0; 10628965ea79SLois Curfman McInnes for ( i=2; i<=size; i++ ) { 10638965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 10648965ea79SLois Curfman McInnes } 10658965ea79SLois Curfman McInnes rstart = rowners[rank]; 10668965ea79SLois Curfman McInnes rend = rowners[rank+1]; 10678965ea79SLois Curfman McInnes 10688965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 10690452661fSBarry Smith ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens); 10708965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 10718965ea79SLois Curfman McInnes if (!rank) { 10720452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 107377c4ece6SBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 10740452661fSBarry Smith sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts); 10758965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i]; 10768965ea79SLois Curfman McInnes MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm); 10770452661fSBarry Smith PetscFree(sndcounts); 10788965ea79SLois Curfman McInnes } 10798965ea79SLois Curfman McInnes else { 10808965ea79SLois Curfman McInnes MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm); 10818965ea79SLois Curfman McInnes } 10828965ea79SLois Curfman McInnes 10838965ea79SLois Curfman McInnes if (!rank) { 10848965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 10850452661fSBarry Smith procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz); 1086cddf8d76SBarry Smith PetscMemzero(procsnz,size*sizeof(int)); 10878965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 10888965ea79SLois Curfman McInnes for ( j=rowners[i]; j< rowners[i+1]; j++ ) { 10898965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 10908965ea79SLois Curfman McInnes } 10918965ea79SLois Curfman McInnes } 10920452661fSBarry Smith PetscFree(rowlengths); 10938965ea79SLois Curfman McInnes 10948965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 10958965ea79SLois Curfman McInnes maxnz = 0; 10968965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 10970452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 10988965ea79SLois Curfman McInnes } 10990452661fSBarry Smith cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols); 11008965ea79SLois Curfman McInnes 11018965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 11028965ea79SLois Curfman McInnes nz = procsnz[0]; 11030452661fSBarry Smith mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 110477c4ece6SBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,BINARY_INT); CHKERRQ(ierr); 11058965ea79SLois Curfman McInnes 11068965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 11078965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { 11088965ea79SLois Curfman McInnes nz = procsnz[i]; 110977c4ece6SBarry Smith ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr); 11108965ea79SLois Curfman McInnes MPI_Send(cols,nz,MPI_INT,i,tag,comm); 11118965ea79SLois Curfman McInnes } 11120452661fSBarry Smith PetscFree(cols); 11138965ea79SLois Curfman McInnes } 11148965ea79SLois Curfman McInnes else { 11158965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 11168965ea79SLois Curfman McInnes nz = 0; 11178965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 11188965ea79SLois Curfman McInnes nz += ourlens[i]; 11198965ea79SLois Curfman McInnes } 11200452661fSBarry Smith mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 11218965ea79SLois Curfman McInnes 11228965ea79SLois Curfman McInnes /* receive message of column indices*/ 11238965ea79SLois Curfman McInnes MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status); 11248965ea79SLois Curfman McInnes MPI_Get_count(&status,MPI_INT,&maxnz); 112539ddd567SLois Curfman McInnes if (maxnz != nz) SETERRQ(1,"MatLoad_MPIDenseorMPIRow:something is wrong with file"); 11268965ea79SLois Curfman McInnes } 11278965ea79SLois Curfman McInnes 11288965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 1129cddf8d76SBarry Smith PetscMemzero(offlens,m*sizeof(int)); 11308965ea79SLois Curfman McInnes jj = 0; 11318965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 11328965ea79SLois Curfman McInnes for ( j=0; j<ourlens[i]; j++ ) { 11338965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 11348965ea79SLois Curfman McInnes jj++; 11358965ea79SLois Curfman McInnes } 11368965ea79SLois Curfman McInnes } 11378965ea79SLois Curfman McInnes 11388965ea79SLois Curfman McInnes /* create our matrix */ 11398965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 11408965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 11418965ea79SLois Curfman McInnes } 1142b4fd4287SBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 11438965ea79SLois Curfman McInnes A = *newmat; 11448965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 11458965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 11468965ea79SLois Curfman McInnes } 11478965ea79SLois Curfman McInnes 11488965ea79SLois Curfman McInnes if (!rank) { 11490452661fSBarry Smith vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals); 11508965ea79SLois Curfman McInnes 11518965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 11528965ea79SLois Curfman McInnes nz = procsnz[0]; 115377c4ece6SBarry Smith ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr); 11548965ea79SLois Curfman McInnes 11558965ea79SLois Curfman McInnes /* insert into matrix */ 11568965ea79SLois Curfman McInnes jj = rstart; 11578965ea79SLois Curfman McInnes smycols = mycols; 11588965ea79SLois Curfman McInnes svals = vals; 11598965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 11608965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 11618965ea79SLois Curfman McInnes smycols += ourlens[i]; 11628965ea79SLois Curfman McInnes svals += ourlens[i]; 11638965ea79SLois Curfman McInnes jj++; 11648965ea79SLois Curfman McInnes } 11658965ea79SLois Curfman McInnes 11668965ea79SLois Curfman McInnes /* read in other processors and ship out */ 11678965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { 11688965ea79SLois Curfman McInnes nz = procsnz[i]; 116977c4ece6SBarry Smith ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr); 11708965ea79SLois Curfman McInnes MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm); 11718965ea79SLois Curfman McInnes } 11720452661fSBarry Smith PetscFree(procsnz); 11738965ea79SLois Curfman McInnes } 11748965ea79SLois Curfman McInnes else { 11758965ea79SLois Curfman McInnes /* receive numeric values */ 11760452661fSBarry Smith vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals); 11778965ea79SLois Curfman McInnes 11788965ea79SLois Curfman McInnes /* receive message of values*/ 11798965ea79SLois Curfman McInnes MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status); 11808965ea79SLois Curfman McInnes MPI_Get_count(&status,MPIU_SCALAR,&maxnz); 118139ddd567SLois Curfman McInnes if (maxnz != nz) SETERRQ(1,"MatLoad_MPIDenseorMPIRow:something is wrong with file"); 11828965ea79SLois Curfman McInnes 11838965ea79SLois Curfman McInnes /* insert into matrix */ 11848965ea79SLois Curfman McInnes jj = rstart; 11858965ea79SLois Curfman McInnes smycols = mycols; 11868965ea79SLois Curfman McInnes svals = vals; 11878965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 11888965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 11898965ea79SLois Curfman McInnes smycols += ourlens[i]; 11908965ea79SLois Curfman McInnes svals += ourlens[i]; 11918965ea79SLois Curfman McInnes jj++; 11928965ea79SLois Curfman McInnes } 11938965ea79SLois Curfman McInnes } 11940452661fSBarry Smith PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners); 11958965ea79SLois Curfman McInnes 11966d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 11976d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 11988965ea79SLois Curfman McInnes return 0; 11998965ea79SLois Curfman McInnes } 120090ace30eSBarry Smith 120190ace30eSBarry Smith 120290ace30eSBarry Smith 120390ace30eSBarry Smith 120490ace30eSBarry Smith 1205