1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*9a8c540fSSatish Balay static char vcid[] = "$Id: mpidense.c,v 1.115 1999/06/08 22:55:41 balay Exp balay $"; 38965ea79SLois Curfman McInnes #endif 48965ea79SLois Curfman McInnes 5ed3cc1f0SBarry Smith /* 6ed3cc1f0SBarry Smith Basic functions for basic parallel dense matrices. 7ed3cc1f0SBarry Smith */ 8ed3cc1f0SBarry Smith 970f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h" 10f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 118965ea79SLois Curfman McInnes 127ef1d9bdSSatish Balay extern int MatSetUpMultiply_MPIDense(Mat); 137ef1d9bdSSatish Balay 145615d1e5SSatish Balay #undef __FUNC__ 155615d1e5SSatish Balay #define __FUNC__ "MatSetValues_MPIDense" 168f6be9afSLois Curfman McInnes int MatSetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv) 178965ea79SLois Curfman McInnes { 1839b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense *) mat->data; 1939b7565bSBarry Smith int ierr, i, j, rstart = A->rstart, rend = A->rend, row; 2039b7565bSBarry Smith int roworiented = A->roworiented; 218965ea79SLois Curfman McInnes 223a40ed3dSBarry Smith PetscFunctionBegin; 238965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 245ef9f2a5SBarry Smith if (idxm[i] < 0) continue; 25a8c6a408SBarry Smith if (idxm[i] >= A->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"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); 303a40ed3dSBarry Smith } else { 318965ea79SLois Curfman McInnes for ( j=0; j<n; j++ ) { 325ef9f2a5SBarry Smith if (idxn[j] < 0) continue; 33a8c6a408SBarry Smith if (idxn[j] >= A->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 3439b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 3539b7565bSBarry Smith } 368965ea79SLois Curfman McInnes } 373a40ed3dSBarry Smith } else { 383782ba37SSatish Balay if ( !A->donotstash) { 3939b7565bSBarry Smith if (roworiented) { 408798bf22SSatish Balay ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr); 41d36fbae8SSatish Balay } else { 428798bf22SSatish Balay ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr); 4339b7565bSBarry Smith } 44b49de8d1SLois Curfman McInnes } 45b49de8d1SLois Curfman McInnes } 463782ba37SSatish Balay } 473a40ed3dSBarry Smith PetscFunctionReturn(0); 48b49de8d1SLois Curfman McInnes } 49b49de8d1SLois Curfman McInnes 505615d1e5SSatish Balay #undef __FUNC__ 515615d1e5SSatish Balay #define __FUNC__ "MatGetValues_MPIDense" 528f6be9afSLois Curfman McInnes int MatGetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v) 53b49de8d1SLois Curfman McInnes { 54b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 55b49de8d1SLois Curfman McInnes int ierr, i, j, rstart = mdn->rstart, rend = mdn->rend, row; 56b49de8d1SLois Curfman McInnes 573a40ed3dSBarry Smith PetscFunctionBegin; 58b49de8d1SLois Curfman McInnes for ( i=0; i<m; i++ ) { 59a8c6a408SBarry Smith if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 60a8c6a408SBarry Smith if (idxm[i] >= mdn->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 61b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 62b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 63b49de8d1SLois Curfman McInnes for ( j=0; j<n; j++ ) { 64a8c6a408SBarry Smith if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 65a8c6a408SBarry Smith if (idxn[j] >= mdn->N) { 66a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 67a8c6a408SBarry Smith } 68b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr); 69b49de8d1SLois Curfman McInnes } 70a8c6a408SBarry Smith } else { 71a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Only local values currently supported"); 728965ea79SLois Curfman McInnes } 738965ea79SLois Curfman McInnes } 743a40ed3dSBarry Smith PetscFunctionReturn(0); 758965ea79SLois Curfman McInnes } 768965ea79SLois Curfman McInnes 775615d1e5SSatish Balay #undef __FUNC__ 785615d1e5SSatish Balay #define __FUNC__ "MatGetArray_MPIDense" 798f6be9afSLois Curfman McInnes int MatGetArray_MPIDense(Mat A,Scalar **array) 80ff14e315SSatish Balay { 81ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense *) A->data; 82ff14e315SSatish Balay int ierr; 83ff14e315SSatish Balay 843a40ed3dSBarry Smith PetscFunctionBegin; 85ff14e315SSatish Balay ierr = MatGetArray(a->A,array);CHKERRQ(ierr); 863a40ed3dSBarry Smith PetscFunctionReturn(0); 87ff14e315SSatish Balay } 88ff14e315SSatish Balay 895615d1e5SSatish Balay #undef __FUNC__ 905615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_MPIDense" 918f6be9afSLois Curfman McInnes int MatRestoreArray_MPIDense(Mat A,Scalar **array) 92ff14e315SSatish Balay { 933a40ed3dSBarry Smith PetscFunctionBegin; 943a40ed3dSBarry Smith PetscFunctionReturn(0); 95ff14e315SSatish Balay } 96ff14e315SSatish Balay 975615d1e5SSatish Balay #undef __FUNC__ 985615d1e5SSatish Balay #define __FUNC__ "MatAssemblyBegin_MPIDense" 998f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 1008965ea79SLois Curfman McInnes { 10139ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 1028965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 103d36fbae8SSatish Balay int ierr,nstash,reallocs; 1048965ea79SLois Curfman McInnes InsertMode addv; 1058965ea79SLois Curfman McInnes 1063a40ed3dSBarry Smith PetscFunctionBegin; 1078965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 108ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 1097056b6fcSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 110a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Cannot mix adds/inserts on different procs"); 1118965ea79SLois Curfman McInnes } 112e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 1138965ea79SLois Curfman McInnes 1148798bf22SSatish Balay ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr); 1158798bf22SSatish Balay ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 116d36fbae8SSatish Balay PLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %d entries, uses %d mallocs.\n", 117d36fbae8SSatish Balay nstash,reallocs); 1183a40ed3dSBarry Smith PetscFunctionReturn(0); 1198965ea79SLois Curfman McInnes } 1208965ea79SLois Curfman McInnes 1215615d1e5SSatish Balay #undef __FUNC__ 1225615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_MPIDense" 1238f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 1248965ea79SLois Curfman McInnes { 12539ddd567SLois Curfman McInnes Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 1267ef1d9bdSSatish Balay int i,n,ierr,*row,*col,flg,j,rstart,ncols; 1277ef1d9bdSSatish Balay Scalar *val; 128e0fa3b82SLois Curfman McInnes InsertMode addv=mat->insertmode; 1298965ea79SLois Curfman McInnes 1303a40ed3dSBarry Smith PetscFunctionBegin; 1318965ea79SLois Curfman McInnes /* wait on receives */ 1327ef1d9bdSSatish Balay while (1) { 1338798bf22SSatish Balay ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 1347ef1d9bdSSatish Balay if (!flg) break; 1358965ea79SLois Curfman McInnes 1367ef1d9bdSSatish Balay for ( i=0; i<n; ) { 1377ef1d9bdSSatish Balay /* Now identify the consecutive vals belonging to the same row */ 1387ef1d9bdSSatish Balay for ( j=i,rstart=row[j]; j<n; j++ ) { if (row[j] != rstart) break; } 1397ef1d9bdSSatish Balay if (j < n) ncols = j-i; 1407ef1d9bdSSatish Balay else ncols = n-i; 1417ef1d9bdSSatish Balay /* Now assemble all these values with a single function call */ 1427ef1d9bdSSatish Balay ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr); 1437ef1d9bdSSatish Balay i = j; 1448965ea79SLois Curfman McInnes } 1457ef1d9bdSSatish Balay } 1468798bf22SSatish Balay ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 1478965ea79SLois Curfman McInnes 14839ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 14939ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 1508965ea79SLois Curfman McInnes 1516d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 15239ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 1538965ea79SLois Curfman McInnes } 1543a40ed3dSBarry Smith PetscFunctionReturn(0); 1558965ea79SLois Curfman McInnes } 1568965ea79SLois Curfman McInnes 1575615d1e5SSatish Balay #undef __FUNC__ 1585615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_MPIDense" 1598f6be9afSLois Curfman McInnes int MatZeroEntries_MPIDense(Mat A) 1608965ea79SLois Curfman McInnes { 1613a40ed3dSBarry Smith int ierr; 16239ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense *) A->data; 1633a40ed3dSBarry Smith 1643a40ed3dSBarry Smith PetscFunctionBegin; 1653a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 1663a40ed3dSBarry Smith PetscFunctionReturn(0); 1678965ea79SLois Curfman McInnes } 1688965ea79SLois Curfman McInnes 1695615d1e5SSatish Balay #undef __FUNC__ 1705615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_MPIDense" 1718f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIDense(Mat A,int *bs) 1724e220ebcSLois Curfman McInnes { 1733a40ed3dSBarry Smith PetscFunctionBegin; 1744e220ebcSLois Curfman McInnes *bs = 1; 1753a40ed3dSBarry Smith PetscFunctionReturn(0); 1764e220ebcSLois Curfman McInnes } 1774e220ebcSLois Curfman McInnes 1788965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 1798965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 1808965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 1813501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 1828965ea79SLois Curfman McInnes routine. 1838965ea79SLois Curfman McInnes */ 1845615d1e5SSatish Balay #undef __FUNC__ 1855615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_MPIDense" 1868f6be9afSLois Curfman McInnes int MatZeroRows_MPIDense(Mat A,IS is,Scalar *diag) 1878965ea79SLois Curfman McInnes { 18839ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense *) A->data; 1898965ea79SLois Curfman McInnes int i,ierr,N, *rows,*owners = l->rowners,size = l->size; 1908965ea79SLois Curfman McInnes int *procs,*nprocs,j,found,idx,nsends,*work; 1918965ea79SLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 1928965ea79SLois Curfman McInnes int *rvalues,tag = A->tag,count,base,slen,n,*source; 1938965ea79SLois Curfman McInnes int *lens,imdex,*lrows,*values; 1948965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 1958965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 1968965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 1978965ea79SLois Curfman McInnes IS istmp; 1988965ea79SLois Curfman McInnes 1993a40ed3dSBarry Smith PetscFunctionBegin; 20077c4ece6SBarry Smith ierr = ISGetSize(is,&N);CHKERRQ(ierr); 2018965ea79SLois Curfman McInnes ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 2028965ea79SLois Curfman McInnes 2038965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 2040452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) );CHKPTRQ(nprocs); 205549d3d68SSatish Balay ierr = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr); 206549d3d68SSatish Balay procs = nprocs + size; 2070452661fSBarry Smith owner = (int *) PetscMalloc((N+1)*sizeof(int));CHKPTRQ(owner); /* see note*/ 2088965ea79SLois Curfman McInnes for ( i=0; i<N; i++ ) { 2098965ea79SLois Curfman McInnes idx = rows[i]; 2108965ea79SLois Curfman McInnes found = 0; 2118965ea79SLois Curfman McInnes for ( j=0; j<size; j++ ) { 2128965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 2138965ea79SLois Curfman McInnes nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 2148965ea79SLois Curfman McInnes } 2158965ea79SLois Curfman McInnes } 216a8c6a408SBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Index out of range"); 2178965ea79SLois Curfman McInnes } 2188965ea79SLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 2198965ea79SLois Curfman McInnes 2208965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 2210452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) );CHKPTRQ(work); 222ca161407SBarry Smith ierr = MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 2238965ea79SLois Curfman McInnes nrecvs = work[rank]; 224ca161407SBarry Smith ierr = MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); 2258965ea79SLois Curfman McInnes nmax = work[rank]; 2260452661fSBarry Smith PetscFree(work); 2278965ea79SLois Curfman McInnes 2288965ea79SLois Curfman McInnes /* post receives: */ 2293a40ed3dSBarry Smith rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int));CHKPTRQ(rvalues); 2303a40ed3dSBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits); 2318965ea79SLois Curfman McInnes for ( i=0; i<nrecvs; i++ ) { 232ca161407SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 2338965ea79SLois Curfman McInnes } 2348965ea79SLois Curfman McInnes 2358965ea79SLois Curfman McInnes /* do sends: 2368965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 2378965ea79SLois Curfman McInnes the ith processor 2388965ea79SLois Curfman McInnes */ 2390452661fSBarry Smith svalues = (int *) PetscMalloc( (N+1)*sizeof(int) );CHKPTRQ(svalues); 2407056b6fcSBarry Smith send_waits = (MPI_Request *) PetscMalloc((nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits); 2410452661fSBarry Smith starts = (int *) PetscMalloc( (size+1)*sizeof(int) );CHKPTRQ(starts); 2428965ea79SLois Curfman McInnes starts[0] = 0; 2438965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 2448965ea79SLois Curfman McInnes for ( i=0; i<N; i++ ) { 2458965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 2468965ea79SLois Curfman McInnes } 2478965ea79SLois Curfman McInnes ISRestoreIndices(is,&rows); 2488965ea79SLois Curfman McInnes 2498965ea79SLois Curfman McInnes starts[0] = 0; 2508965ea79SLois Curfman McInnes for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 2518965ea79SLois Curfman McInnes count = 0; 2528965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 2538965ea79SLois Curfman McInnes if (procs[i]) { 254ca161407SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 2558965ea79SLois Curfman McInnes } 2568965ea79SLois Curfman McInnes } 2570452661fSBarry Smith PetscFree(starts); 2588965ea79SLois Curfman McInnes 2598965ea79SLois Curfman McInnes base = owners[rank]; 2608965ea79SLois Curfman McInnes 2618965ea79SLois Curfman McInnes /* wait on receives */ 2620452661fSBarry Smith lens = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) );CHKPTRQ(lens); 2638965ea79SLois Curfman McInnes source = lens + nrecvs; 2648965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 2658965ea79SLois Curfman McInnes while (count) { 266ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 2678965ea79SLois Curfman McInnes /* unpack receives into our local space */ 268ca161407SBarry Smith ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr); 2698965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 2708965ea79SLois Curfman McInnes lens[imdex] = n; 2718965ea79SLois Curfman McInnes slen += n; 2728965ea79SLois Curfman McInnes count--; 2738965ea79SLois Curfman McInnes } 2740452661fSBarry Smith PetscFree(recv_waits); 2758965ea79SLois Curfman McInnes 2768965ea79SLois Curfman McInnes /* move the data into the send scatter */ 2770452661fSBarry Smith lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) );CHKPTRQ(lrows); 2788965ea79SLois Curfman McInnes count = 0; 2798965ea79SLois Curfman McInnes for ( i=0; i<nrecvs; i++ ) { 2808965ea79SLois Curfman McInnes values = rvalues + i*nmax; 2818965ea79SLois Curfman McInnes for ( j=0; j<lens[i]; j++ ) { 2828965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 2838965ea79SLois Curfman McInnes } 2848965ea79SLois Curfman McInnes } 2850452661fSBarry Smith PetscFree(rvalues); PetscFree(lens); 2860452661fSBarry Smith PetscFree(owner); PetscFree(nprocs); 2878965ea79SLois Curfman McInnes 2888965ea79SLois Curfman McInnes /* actually zap the local rows */ 289029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 2908965ea79SLois Curfman McInnes PLogObjectParent(A,istmp); 2910452661fSBarry Smith PetscFree(lrows); 2928965ea79SLois Curfman McInnes ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr); 2938965ea79SLois Curfman McInnes ierr = ISDestroy(istmp);CHKERRQ(ierr); 2948965ea79SLois Curfman McInnes 2958965ea79SLois Curfman McInnes /* wait on sends */ 2968965ea79SLois Curfman McInnes if (nsends) { 2977056b6fcSBarry Smith send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status); 298ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 2990452661fSBarry Smith PetscFree(send_status); 3008965ea79SLois Curfman McInnes } 3010452661fSBarry Smith PetscFree(send_waits); PetscFree(svalues); 3028965ea79SLois Curfman McInnes 3033a40ed3dSBarry Smith PetscFunctionReturn(0); 3048965ea79SLois Curfman McInnes } 3058965ea79SLois Curfman McInnes 3065615d1e5SSatish Balay #undef __FUNC__ 3075615d1e5SSatish Balay #define __FUNC__ "MatMult_MPIDense" 3088f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 3098965ea79SLois Curfman McInnes { 31039ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3118965ea79SLois Curfman McInnes int ierr; 312c456f294SBarry Smith 3133a40ed3dSBarry Smith PetscFunctionBegin; 31443a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 31543a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 31644cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 3173a40ed3dSBarry Smith PetscFunctionReturn(0); 3188965ea79SLois Curfman McInnes } 3198965ea79SLois Curfman McInnes 3205615d1e5SSatish Balay #undef __FUNC__ 3215615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_MPIDense" 3228f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 3238965ea79SLois Curfman McInnes { 32439ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3258965ea79SLois Curfman McInnes int ierr; 326c456f294SBarry Smith 3273a40ed3dSBarry Smith PetscFunctionBegin; 32843a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 32943a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 33044cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 3313a40ed3dSBarry Smith PetscFunctionReturn(0); 3328965ea79SLois Curfman McInnes } 3338965ea79SLois Curfman McInnes 3345615d1e5SSatish Balay #undef __FUNC__ 3355615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_MPIDense" 3368f6be9afSLois Curfman McInnes int MatMultTrans_MPIDense(Mat A,Vec xx,Vec yy) 337096963f5SLois Curfman McInnes { 338096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 339096963f5SLois Curfman McInnes int ierr; 3403501a2bdSLois Curfman McInnes Scalar zero = 0.0; 341096963f5SLois Curfman McInnes 3423a40ed3dSBarry Smith PetscFunctionBegin; 3433501a2bdSLois Curfman McInnes ierr = VecSet(&zero,yy);CHKERRQ(ierr); 34444cd7ae7SLois Curfman McInnes ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 345537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 346537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 3473a40ed3dSBarry Smith PetscFunctionReturn(0); 348096963f5SLois Curfman McInnes } 349096963f5SLois Curfman McInnes 3505615d1e5SSatish Balay #undef __FUNC__ 3515615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_MPIDense" 3528f6be9afSLois Curfman McInnes int MatMultTransAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 353096963f5SLois Curfman McInnes { 354096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 355096963f5SLois Curfman McInnes int ierr; 356096963f5SLois Curfman McInnes 3573a40ed3dSBarry Smith PetscFunctionBegin; 3583501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 35944cd7ae7SLois Curfman McInnes ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 360537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 361537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 3623a40ed3dSBarry Smith PetscFunctionReturn(0); 363096963f5SLois Curfman McInnes } 364096963f5SLois Curfman McInnes 3655615d1e5SSatish Balay #undef __FUNC__ 3665615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_MPIDense" 3678f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v) 3688965ea79SLois Curfman McInnes { 36939ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 370096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense *) a->A->data; 37144cd7ae7SLois Curfman McInnes int ierr, len, i, n, m = a->m, radd; 37244cd7ae7SLois Curfman McInnes Scalar *x, zero = 0.0; 373ed3cc1f0SBarry Smith 3743a40ed3dSBarry Smith PetscFunctionBegin; 37544cd7ae7SLois Curfman McInnes VecSet(&zero,v); 376096963f5SLois Curfman McInnes ierr = VecGetArray(v,&x);CHKERRQ(ierr); 377096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 378a8c6a408SBarry Smith if (n != a->M) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec"); 37944cd7ae7SLois Curfman McInnes len = PetscMin(aloc->m,aloc->n); 3807ddc982cSLois Curfman McInnes radd = a->rstart*m; 38144cd7ae7SLois Curfman McInnes for ( i=0; i<len; i++ ) { 382096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 383096963f5SLois Curfman McInnes } 384*9a8c540fSSatish Balay ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 3853a40ed3dSBarry Smith PetscFunctionReturn(0); 3868965ea79SLois Curfman McInnes } 3878965ea79SLois Curfman McInnes 3885615d1e5SSatish Balay #undef __FUNC__ 3895615d1e5SSatish Balay #define __FUNC__ "MatDestroy_MPIDense" 390e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat) 3918965ea79SLois Curfman McInnes { 3923501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3938965ea79SLois Curfman McInnes int ierr; 394ed3cc1f0SBarry Smith 3953a40ed3dSBarry Smith PetscFunctionBegin; 39694d884c6SBarry Smith if (--mat->refct > 0) PetscFunctionReturn(0); 39794d884c6SBarry Smith 39894d884c6SBarry Smith if (mat->mapping) { 39994d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->mapping);CHKERRQ(ierr); 40094d884c6SBarry Smith } 40194d884c6SBarry Smith if (mat->bmapping) { 40294d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->bmapping);CHKERRQ(ierr); 40394d884c6SBarry Smith } 404aa482453SBarry Smith #if defined(PETSC_USE_LOG) 405e1311b90SBarry Smith PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mdn->M,mdn->N); 4068965ea79SLois Curfman McInnes #endif 4078798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 4080452661fSBarry Smith PetscFree(mdn->rowners); 4093501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A);CHKERRQ(ierr); 4103501a2bdSLois Curfman McInnes if (mdn->lvec) VecDestroy(mdn->lvec); 4113501a2bdSLois Curfman McInnes if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 412622d7880SLois Curfman McInnes if (mdn->factor) { 413622d7880SLois Curfman McInnes if (mdn->factor->temp) PetscFree(mdn->factor->temp); 414622d7880SLois Curfman McInnes if (mdn->factor->tag) PetscFree(mdn->factor->tag); 415622d7880SLois Curfman McInnes if (mdn->factor->pivots) PetscFree(mdn->factor->pivots); 416622d7880SLois Curfman McInnes PetscFree(mdn->factor); 417622d7880SLois Curfman McInnes } 4180452661fSBarry Smith PetscFree(mdn); 41961b13de0SBarry Smith if (mat->rmap) { 42061b13de0SBarry Smith ierr = MapDestroy(mat->rmap);CHKERRQ(ierr); 42161b13de0SBarry Smith } 42261b13de0SBarry Smith if (mat->cmap) { 42361b13de0SBarry Smith ierr = MapDestroy(mat->cmap);CHKERRQ(ierr); 42490f02eecSBarry Smith } 4258965ea79SLois Curfman McInnes PLogObjectDestroy(mat); 4260452661fSBarry Smith PetscHeaderDestroy(mat); 4273a40ed3dSBarry Smith PetscFunctionReturn(0); 4288965ea79SLois Curfman McInnes } 42939ddd567SLois Curfman McInnes 4305615d1e5SSatish Balay #undef __FUNC__ 4315615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense_Binary" 43239ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer) 4338965ea79SLois Curfman McInnes { 43439ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 4358965ea79SLois Curfman McInnes int ierr; 4367056b6fcSBarry Smith 4373a40ed3dSBarry Smith PetscFunctionBegin; 43839ddd567SLois Curfman McInnes if (mdn->size == 1) { 43939ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 4408965ea79SLois Curfman McInnes } 441a8c6a408SBarry Smith else SETERRQ(PETSC_ERR_SUP,0,"Only uniprocessor output supported"); 4423a40ed3dSBarry Smith PetscFunctionReturn(0); 4438965ea79SLois Curfman McInnes } 4448965ea79SLois Curfman McInnes 4455615d1e5SSatish Balay #undef __FUNC__ 4465615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense_ASCII" 44739ddd567SLois Curfman McInnes static int MatView_MPIDense_ASCII(Mat mat,Viewer viewer) 4488965ea79SLois Curfman McInnes { 44939ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 45077ed5343SBarry Smith int ierr, format, size = mdn->size, rank = mdn->rank; 4518965ea79SLois Curfman McInnes FILE *fd; 45219bcc07fSBarry Smith ViewerType vtype; 4538965ea79SLois Curfman McInnes 4543a40ed3dSBarry Smith PetscFunctionBegin; 4553a40ed3dSBarry Smith ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 45690ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd);CHKERRQ(ierr); 457888f2ed8SSatish Balay ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 458639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 4594e220ebcSLois Curfman McInnes MatInfo info; 460888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 46177c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 4624e220ebcSLois Curfman McInnes fprintf(fd," [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mdn->m, 4634e220ebcSLois Curfman McInnes (int)info.nz_used,(int)info.nz_allocated,(int)info.memory); 464096963f5SLois Curfman McInnes fflush(fd); 46577c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 4663501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 4673a40ed3dSBarry Smith PetscFunctionReturn(0); 46896f6c058SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_INFO) { 4693a40ed3dSBarry Smith PetscFunctionReturn(0); 4708965ea79SLois Curfman McInnes } 47177ed5343SBarry Smith 4728965ea79SLois Curfman McInnes if (size == 1) { 47339ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 4743a40ed3dSBarry Smith } else { 4758965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 4768965ea79SLois Curfman McInnes Mat A; 47739ddd567SLois Curfman McInnes int M = mdn->M, N = mdn->N,m,row,i, nz, *cols; 47839ddd567SLois Curfman McInnes Scalar *vals; 47939ddd567SLois Curfman McInnes Mat_SeqDense *Amdn = (Mat_SeqDense*) mdn->A->data; 4808965ea79SLois Curfman McInnes 4818965ea79SLois Curfman McInnes if (!rank) { 4820513a670SBarry Smith ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr); 4833a40ed3dSBarry Smith } else { 4840513a670SBarry Smith ierr = MatCreateMPIDense(mat->comm,0,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr); 4858965ea79SLois Curfman McInnes } 4868965ea79SLois Curfman McInnes PLogObjectParent(mat,A); 4878965ea79SLois Curfman McInnes 48839ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 48939ddd567SLois Curfman McInnes but it's quick for now */ 49039ddd567SLois Curfman McInnes row = mdn->rstart; m = Amdn->m; 4918965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 49239ddd567SLois Curfman McInnes ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 49339ddd567SLois Curfman McInnes ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 49439ddd567SLois Curfman McInnes ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 49539ddd567SLois Curfman McInnes row++; 4968965ea79SLois Curfman McInnes } 4978965ea79SLois Curfman McInnes 4986d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4996d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5008965ea79SLois Curfman McInnes if (!rank) { 50139ddd567SLois Curfman McInnes ierr = MatView(((Mat_MPIDense*)(A->data))->A,viewer);CHKERRQ(ierr); 5028965ea79SLois Curfman McInnes } 5038965ea79SLois Curfman McInnes ierr = MatDestroy(A);CHKERRQ(ierr); 5048965ea79SLois Curfman McInnes } 5053a40ed3dSBarry Smith PetscFunctionReturn(0); 5068965ea79SLois Curfman McInnes } 5078965ea79SLois Curfman McInnes 5085615d1e5SSatish Balay #undef __FUNC__ 5095615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense" 510e1311b90SBarry Smith int MatView_MPIDense(Mat mat,Viewer viewer) 5118965ea79SLois Curfman McInnes { 51239ddd567SLois Curfman McInnes int ierr; 513bcd2baecSBarry Smith ViewerType vtype; 5148965ea79SLois Curfman McInnes 515bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 5163f1db9ecSBarry Smith if (PetscTypeCompare(vtype,ASCII_VIEWER)) { 51739ddd567SLois Curfman McInnes ierr = MatView_MPIDense_ASCII(mat,viewer);CHKERRQ(ierr); 5183f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,BINARY_VIEWER)) { 5193a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 5205cd90555SBarry Smith } else { 5215cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 5228965ea79SLois Curfman McInnes } 5233a40ed3dSBarry Smith PetscFunctionReturn(0); 5248965ea79SLois Curfman McInnes } 5258965ea79SLois Curfman McInnes 5265615d1e5SSatish Balay #undef __FUNC__ 5275615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_MPIDense" 5288f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 5298965ea79SLois Curfman McInnes { 5303501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 5313501a2bdSLois Curfman McInnes Mat mdn = mat->A; 5324e220ebcSLois Curfman McInnes int ierr; 5334e220ebcSLois Curfman McInnes double isend[5], irecv[5]; 5348965ea79SLois Curfman McInnes 5353a40ed3dSBarry Smith PetscFunctionBegin; 5364e220ebcSLois Curfman McInnes info->rows_global = (double)mat->M; 5374e220ebcSLois Curfman McInnes info->columns_global = (double)mat->N; 5384e220ebcSLois Curfman McInnes info->rows_local = (double)mat->m; 5394e220ebcSLois Curfman McInnes info->columns_local = (double)mat->N; 5404e220ebcSLois Curfman McInnes info->block_size = 1.0; 5414e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 5424e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 5434e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 5448965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 5454e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 5464e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 5474e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 5484e220ebcSLois Curfman McInnes info->memory = isend[3]; 5494e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 5508965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 551f7cdd7c9SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr); 5524e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 5534e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 5544e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 5554e220ebcSLois Curfman McInnes info->memory = irecv[3]; 5564e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 5578965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 558f7cdd7c9SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 5594e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 5604e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 5614e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 5624e220ebcSLois Curfman McInnes info->memory = irecv[3]; 5634e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 5648965ea79SLois Curfman McInnes } 5654e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 5664e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 5674e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 5683a40ed3dSBarry Smith PetscFunctionReturn(0); 5698965ea79SLois Curfman McInnes } 5708965ea79SLois Curfman McInnes 5718c469469SLois Curfman McInnes /* extern int MatLUFactorSymbolic_MPIDense(Mat,IS,IS,double,Mat*); 5728aaee692SLois Curfman McInnes extern int MatLUFactorNumeric_MPIDense(Mat,Mat*); 5738aaee692SLois Curfman McInnes extern int MatLUFactor_MPIDense(Mat,IS,IS,double); 5748aaee692SLois Curfman McInnes extern int MatSolve_MPIDense(Mat,Vec,Vec); 5758c469469SLois Curfman McInnes extern int MatSolveAdd_MPIDense(Mat,Vec,Vec,Vec); 5768aaee692SLois Curfman McInnes extern int MatSolveTrans_MPIDense(Mat,Vec,Vec); 5778aaee692SLois Curfman McInnes extern int MatSolveTransAdd_MPIDense(Mat,Vec,Vec,Vec); */ 5788aaee692SLois Curfman McInnes 5795615d1e5SSatish Balay #undef __FUNC__ 5805615d1e5SSatish Balay #define __FUNC__ "MatSetOption_MPIDense" 5818f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op) 5828965ea79SLois Curfman McInnes { 58339ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 5848965ea79SLois Curfman McInnes 5853a40ed3dSBarry Smith PetscFunctionBegin; 5866d4a8577SBarry Smith if (op == MAT_NO_NEW_NONZERO_LOCATIONS || 5876d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 5884787f768SSatish Balay op == MAT_NEW_NONZERO_LOCATION_ERR || 5894787f768SSatish Balay op == MAT_NEW_NONZERO_ALLOCATION_ERR || 590219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_SORTED || 591219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED) { 592b1fbbac0SLois Curfman McInnes MatSetOption(a->A,op); 593b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROW_ORIENTED) { 594aeafbbfcSLois Curfman McInnes a->roworiented = 1; 5958965ea79SLois Curfman McInnes MatSetOption(a->A,op); 596b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROWS_SORTED || 597219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 5986d4a8577SBarry Smith op == MAT_SYMMETRIC || 5996d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 600b51ba29fSSatish Balay op == MAT_YES_NEW_DIAGONALS || 601b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 602981c4779SBarry Smith PLogInfo(A,"MatSetOption_MPIDense:Option ignored\n"); 6033a40ed3dSBarry Smith } else if (op == MAT_COLUMN_ORIENTED) { 6043a40ed3dSBarry Smith a->roworiented = 0; MatSetOption(a->A,op); 6053782ba37SSatish Balay } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) { 6063782ba37SSatish Balay a->donotstash = 1; 6073a40ed3dSBarry Smith } else if (op == MAT_NO_NEW_DIAGONALS) { 6083a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 6093a40ed3dSBarry Smith } else { 6103a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 6113a40ed3dSBarry Smith } 6123a40ed3dSBarry Smith PetscFunctionReturn(0); 6138965ea79SLois Curfman McInnes } 6148965ea79SLois Curfman McInnes 6155615d1e5SSatish Balay #undef __FUNC__ 6165615d1e5SSatish Balay #define __FUNC__ "MatGetSize_MPIDense" 6178f6be9afSLois Curfman McInnes int MatGetSize_MPIDense(Mat A,int *m,int *n) 6188965ea79SLois Curfman McInnes { 6193501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6203a40ed3dSBarry Smith 6213a40ed3dSBarry Smith PetscFunctionBegin; 6228965ea79SLois Curfman McInnes *m = mat->M; *n = mat->N; 6233a40ed3dSBarry Smith PetscFunctionReturn(0); 6248965ea79SLois Curfman McInnes } 6258965ea79SLois Curfman McInnes 6265615d1e5SSatish Balay #undef __FUNC__ 6275615d1e5SSatish Balay #define __FUNC__ "MatGetLocalSize_MPIDense" 6288f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIDense(Mat A,int *m,int *n) 6298965ea79SLois Curfman McInnes { 6303501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6313a40ed3dSBarry Smith 6323a40ed3dSBarry Smith PetscFunctionBegin; 6338965ea79SLois Curfman McInnes *m = mat->m; *n = mat->N; 6343a40ed3dSBarry Smith PetscFunctionReturn(0); 6358965ea79SLois Curfman McInnes } 6368965ea79SLois Curfman McInnes 6375615d1e5SSatish Balay #undef __FUNC__ 6385615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_MPIDense" 6398f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n) 6408965ea79SLois Curfman McInnes { 6413501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6423a40ed3dSBarry Smith 6433a40ed3dSBarry Smith PetscFunctionBegin; 6448965ea79SLois Curfman McInnes *m = mat->rstart; *n = mat->rend; 6453a40ed3dSBarry Smith PetscFunctionReturn(0); 6468965ea79SLois Curfman McInnes } 6478965ea79SLois Curfman McInnes 6485615d1e5SSatish Balay #undef __FUNC__ 6495615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIDense" 6508f6be9afSLois Curfman McInnes int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v) 6518965ea79SLois Curfman McInnes { 6523501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6533a40ed3dSBarry Smith int lrow, rstart = mat->rstart, rend = mat->rend,ierr; 6548965ea79SLois Curfman McInnes 6553a40ed3dSBarry Smith PetscFunctionBegin; 656a8c6a408SBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,0,"only local rows") 6578965ea79SLois Curfman McInnes lrow = row - rstart; 6583a40ed3dSBarry Smith ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr); 6593a40ed3dSBarry Smith PetscFunctionReturn(0); 6608965ea79SLois Curfman McInnes } 6618965ea79SLois Curfman McInnes 6625615d1e5SSatish Balay #undef __FUNC__ 6635615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_MPIDense" 6648f6be9afSLois Curfman McInnes int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v) 6658965ea79SLois Curfman McInnes { 6663a40ed3dSBarry Smith PetscFunctionBegin; 6670452661fSBarry Smith if (idx) PetscFree(*idx); 6680452661fSBarry Smith if (v) PetscFree(*v); 6693a40ed3dSBarry Smith PetscFunctionReturn(0); 6708965ea79SLois Curfman McInnes } 6718965ea79SLois Curfman McInnes 6725615d1e5SSatish Balay #undef __FUNC__ 6735615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIDense" 6748f6be9afSLois Curfman McInnes int MatNorm_MPIDense(Mat A,NormType type,double *norm) 675096963f5SLois Curfman McInnes { 6763501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) A->data; 6773501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*) mdn->A->data; 6783501a2bdSLois Curfman McInnes int ierr, i, j; 6793501a2bdSLois Curfman McInnes double sum = 0.0; 6803501a2bdSLois Curfman McInnes Scalar *v = mat->v; 6813501a2bdSLois Curfman McInnes 6823a40ed3dSBarry Smith PetscFunctionBegin; 6833501a2bdSLois Curfman McInnes if (mdn->size == 1) { 6843501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,norm);CHKERRQ(ierr); 6853501a2bdSLois Curfman McInnes } else { 6863501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 6873501a2bdSLois Curfman McInnes for (i=0; i<mat->n*mat->m; i++ ) { 688aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 689e20fef11SSatish Balay sum += PetscReal(PetscConj(*v)*(*v)); v++; 6903501a2bdSLois Curfman McInnes #else 6913501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 6923501a2bdSLois Curfman McInnes #endif 6933501a2bdSLois Curfman McInnes } 694ca161407SBarry Smith ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 6953501a2bdSLois Curfman McInnes *norm = sqrt(*norm); 6963501a2bdSLois Curfman McInnes PLogFlops(2*mat->n*mat->m); 6973a40ed3dSBarry Smith } else if (type == NORM_1) { 6983501a2bdSLois Curfman McInnes double *tmp, *tmp2; 6990452661fSBarry Smith tmp = (double *) PetscMalloc( 2*mdn->N*sizeof(double) );CHKPTRQ(tmp); 7003501a2bdSLois Curfman McInnes tmp2 = tmp + mdn->N; 701549d3d68SSatish Balay ierr = PetscMemzero(tmp,2*mdn->N*sizeof(double));CHKERRQ(ierr); 702096963f5SLois Curfman McInnes *norm = 0.0; 7033501a2bdSLois Curfman McInnes v = mat->v; 7043501a2bdSLois Curfman McInnes for ( j=0; j<mat->n; j++ ) { 7053501a2bdSLois Curfman McInnes for ( i=0; i<mat->m; i++ ) { 70667e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 7073501a2bdSLois Curfman McInnes } 7083501a2bdSLois Curfman McInnes } 709ca161407SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,mdn->N,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 7103501a2bdSLois Curfman McInnes for ( j=0; j<mdn->N; j++ ) { 7113501a2bdSLois Curfman McInnes if (tmp2[j] > *norm) *norm = tmp2[j]; 7123501a2bdSLois Curfman McInnes } 7130452661fSBarry Smith PetscFree(tmp); 7143501a2bdSLois Curfman McInnes PLogFlops(mat->n*mat->m); 7153a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 7163501a2bdSLois Curfman McInnes double ntemp; 7173501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 718ca161407SBarry Smith ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr); 7193a40ed3dSBarry Smith } else { 720a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"No support for two norm"); 7213501a2bdSLois Curfman McInnes } 7223501a2bdSLois Curfman McInnes } 7233a40ed3dSBarry Smith PetscFunctionReturn(0); 7243501a2bdSLois Curfman McInnes } 7253501a2bdSLois Curfman McInnes 7265615d1e5SSatish Balay #undef __FUNC__ 7275615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIDense" 7288f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout) 7293501a2bdSLois Curfman McInnes { 7303501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 7313501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense *) a->A->data; 7323501a2bdSLois Curfman McInnes Mat B; 7333501a2bdSLois Curfman McInnes int M = a->M, N = a->N, m, n, *rwork, rstart = a->rstart; 7343501a2bdSLois Curfman McInnes int j, i, ierr; 7353501a2bdSLois Curfman McInnes Scalar *v; 7363501a2bdSLois Curfman McInnes 7373a40ed3dSBarry Smith PetscFunctionBegin; 7387056b6fcSBarry Smith if (matout == PETSC_NULL && M != N) { 739a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Supports square matrix only in-place"); 7407056b6fcSBarry Smith } 7417056b6fcSBarry Smith ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr); 7423501a2bdSLois Curfman McInnes 7433501a2bdSLois Curfman McInnes m = Aloc->m; n = Aloc->n; v = Aloc->v; 7440452661fSBarry Smith rwork = (int *) PetscMalloc(n*sizeof(int));CHKPTRQ(rwork); 7453501a2bdSLois Curfman McInnes for ( j=0; j<n; j++ ) { 7463501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 7473501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 7483501a2bdSLois Curfman McInnes v += m; 7493501a2bdSLois Curfman McInnes } 7500452661fSBarry Smith PetscFree(rwork); 7516d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7526d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7533638b69dSLois Curfman McInnes if (matout != PETSC_NULL) { 7543501a2bdSLois Curfman McInnes *matout = B; 7553501a2bdSLois Curfman McInnes } else { 756f830108cSBarry Smith PetscOps *Abops; 75709dc0095SBarry Smith MatOps Aops; 758f830108cSBarry Smith 7593501a2bdSLois Curfman McInnes /* This isn't really an in-place transpose, but free data struct from a */ 7600452661fSBarry Smith PetscFree(a->rowners); 7613501a2bdSLois Curfman McInnes ierr = MatDestroy(a->A);CHKERRQ(ierr); 7623501a2bdSLois Curfman McInnes if (a->lvec) VecDestroy(a->lvec); 7633501a2bdSLois Curfman McInnes if (a->Mvctx) VecScatterDestroy(a->Mvctx); 7640452661fSBarry Smith PetscFree(a); 765f830108cSBarry Smith 766f830108cSBarry Smith /* 767f830108cSBarry Smith This is horrible, horrible code. We need to keep the 768f830108cSBarry Smith A pointers for the bops and ops but copy everything 769f830108cSBarry Smith else from C. 770f830108cSBarry Smith */ 771f830108cSBarry Smith Abops = A->bops; 772f830108cSBarry Smith Aops = A->ops; 773549d3d68SSatish Balay ierr = PetscMemcpy(A,B,sizeof(struct _p_Mat));CHKERRQ(ierr); 774f830108cSBarry Smith A->bops = Abops; 775f830108cSBarry Smith A->ops = Aops; 776f830108cSBarry Smith 7770452661fSBarry Smith PetscHeaderDestroy(B); 7783501a2bdSLois Curfman McInnes } 7793a40ed3dSBarry Smith PetscFunctionReturn(0); 780096963f5SLois Curfman McInnes } 781096963f5SLois Curfman McInnes 782eadb2fb4SBarry Smith #include "pinclude/blaslapack.h" 7835615d1e5SSatish Balay #undef __FUNC__ 7845615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIDense" 7858f6be9afSLois Curfman McInnes int MatScale_MPIDense(Scalar *alpha,Mat inA) 78644cd7ae7SLois Curfman McInnes { 78744cd7ae7SLois Curfman McInnes Mat_MPIDense *A = (Mat_MPIDense *) inA->data; 78844cd7ae7SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->A->data; 78944cd7ae7SLois Curfman McInnes int one = 1, nz; 79044cd7ae7SLois Curfman McInnes 7913a40ed3dSBarry Smith PetscFunctionBegin; 79244cd7ae7SLois Curfman McInnes nz = a->m*a->n; 79344cd7ae7SLois Curfman McInnes BLscal_( &nz, alpha, a->v, &one ); 79444cd7ae7SLois Curfman McInnes PLogFlops(nz); 7953a40ed3dSBarry Smith PetscFunctionReturn(0); 79644cd7ae7SLois Curfman McInnes } 79744cd7ae7SLois Curfman McInnes 7985609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 7997b2a1423SBarry Smith extern int MatGetSubMatrices_MPIDense(Mat,int,IS *,IS *,MatReuse,Mat **); 8008965ea79SLois Curfman McInnes 8018965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 80209dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 80309dc0095SBarry Smith MatGetRow_MPIDense, 80409dc0095SBarry Smith MatRestoreRow_MPIDense, 80509dc0095SBarry Smith MatMult_MPIDense, 80609dc0095SBarry Smith MatMultAdd_MPIDense, 80709dc0095SBarry Smith MatMultTrans_MPIDense, 80809dc0095SBarry Smith MatMultTransAdd_MPIDense, 8098965ea79SLois Curfman McInnes 0, 81009dc0095SBarry Smith 0, 81109dc0095SBarry Smith 0, 81209dc0095SBarry Smith 0, 81309dc0095SBarry Smith 0, 81409dc0095SBarry Smith 0, 81509dc0095SBarry Smith 0, 81609dc0095SBarry Smith MatTranspose_MPIDense, 81709dc0095SBarry Smith MatGetInfo_MPIDense,0, 81809dc0095SBarry Smith MatGetDiagonal_MPIDense, 81909dc0095SBarry Smith 0, 82009dc0095SBarry Smith MatNorm_MPIDense, 82109dc0095SBarry Smith MatAssemblyBegin_MPIDense, 82209dc0095SBarry Smith MatAssemblyEnd_MPIDense, 82309dc0095SBarry Smith 0, 82409dc0095SBarry Smith MatSetOption_MPIDense, 82509dc0095SBarry Smith MatZeroEntries_MPIDense, 82609dc0095SBarry Smith MatZeroRows_MPIDense, 82709dc0095SBarry Smith 0, 82809dc0095SBarry Smith 0, 82909dc0095SBarry Smith 0, 83009dc0095SBarry Smith 0, 83109dc0095SBarry Smith MatGetSize_MPIDense, 83209dc0095SBarry Smith MatGetLocalSize_MPIDense, 83339ddd567SLois Curfman McInnes MatGetOwnershipRange_MPIDense, 83409dc0095SBarry Smith 0, 83509dc0095SBarry Smith 0, 83609dc0095SBarry Smith MatGetArray_MPIDense, 83709dc0095SBarry Smith MatRestoreArray_MPIDense, 8385609ef8eSBarry Smith MatDuplicate_MPIDense, 83909dc0095SBarry Smith 0, 84009dc0095SBarry Smith 0, 84109dc0095SBarry Smith 0, 84209dc0095SBarry Smith 0, 84309dc0095SBarry Smith 0, 8442ce60cd0SSatish Balay MatGetSubMatrices_MPIDense, 84509dc0095SBarry Smith 0, 84609dc0095SBarry Smith MatGetValues_MPIDense, 84709dc0095SBarry Smith 0, 84809dc0095SBarry Smith 0, 84909dc0095SBarry Smith MatScale_MPIDense, 85009dc0095SBarry Smith 0, 85109dc0095SBarry Smith 0, 85209dc0095SBarry Smith 0, 85309dc0095SBarry Smith MatGetBlockSize_MPIDense, 85409dc0095SBarry Smith 0, 85509dc0095SBarry Smith 0, 85609dc0095SBarry Smith 0, 85709dc0095SBarry Smith 0, 85809dc0095SBarry Smith 0, 85909dc0095SBarry Smith 0, 86009dc0095SBarry Smith 0, 86109dc0095SBarry Smith 0, 86209dc0095SBarry Smith 0, 86309dc0095SBarry Smith 0, 86409dc0095SBarry Smith 0, 86509dc0095SBarry Smith 0, 86609dc0095SBarry Smith MatGetMaps_Petsc}; 8678965ea79SLois Curfman McInnes 8685615d1e5SSatish Balay #undef __FUNC__ 8695615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIDense" 8708965ea79SLois Curfman McInnes /*@C 87139ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 8728965ea79SLois Curfman McInnes 873db81eaa0SLois Curfman McInnes Collective on MPI_Comm 874db81eaa0SLois Curfman McInnes 8758965ea79SLois Curfman McInnes Input Parameters: 876db81eaa0SLois Curfman McInnes + comm - MPI communicator 8778965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 878db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 8798965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 880db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 881db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 882dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 8838965ea79SLois Curfman McInnes 8848965ea79SLois Curfman McInnes Output Parameter: 885477f1c0bSLois Curfman McInnes . A - the matrix 8868965ea79SLois Curfman McInnes 887b259b22eSLois Curfman McInnes Notes: 88839ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 88939ddd567SLois Curfman McInnes storage by columns. 8908965ea79SLois Curfman McInnes 89118f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 89218f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 893b4fd4287SBarry Smith set data=PETSC_NULL. 89418f449edSLois Curfman McInnes 8958965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 8968965ea79SLois Curfman McInnes (possibly both). 8978965ea79SLois Curfman McInnes 8983501a2bdSLois Curfman McInnes Currently, the only parallel dense matrix decomposition is by rows, 8993501a2bdSLois Curfman McInnes so that n=N and each submatrix owns all of the global columns. 9003501a2bdSLois Curfman McInnes 901027ccd11SLois Curfman McInnes Level: intermediate 902027ccd11SLois Curfman McInnes 90339ddd567SLois Curfman McInnes .keywords: matrix, dense, parallel 9048965ea79SLois Curfman McInnes 90539ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 9068965ea79SLois Curfman McInnes @*/ 907477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A) 9088965ea79SLois Curfman McInnes { 9098965ea79SLois Curfman McInnes Mat mat; 91039ddd567SLois Curfman McInnes Mat_MPIDense *a; 91125cdf11fSBarry Smith int ierr, i,flg; 9128965ea79SLois Curfman McInnes 9133a40ed3dSBarry Smith PetscFunctionBegin; 914ed2daf61SLois Curfman McInnes /* Note: For now, when data is specified above, this assumes the user correctly 915ed2daf61SLois Curfman McInnes allocates the local dense storage space. We should add error checking. */ 91618f449edSLois Curfman McInnes 917477f1c0bSLois Curfman McInnes *A = 0; 9183f1db9ecSBarry Smith PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,"Mat",comm,MatDestroy,MatView); 9198965ea79SLois Curfman McInnes PLogObjectCreate(mat); 9200452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIDense));CHKPTRQ(a); 921549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 922e1311b90SBarry Smith mat->ops->destroy = MatDestroy_MPIDense; 923e1311b90SBarry Smith mat->ops->view = MatView_MPIDense; 9248965ea79SLois Curfman McInnes mat->factor = 0; 92590f02eecSBarry Smith mat->mapping = 0; 9268965ea79SLois Curfman McInnes 927622d7880SLois Curfman McInnes a->factor = 0; 928e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 929d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&a->rank);CHKERRQ(ierr); 930d132466eSBarry Smith ierr = MPI_Comm_size(comm,&a->size);CHKERRQ(ierr); 9318965ea79SLois Curfman McInnes 93296f6c058SBarry Smith ierr = PetscSplitOwnership(comm,&m,&M);CHKERRQ(ierr); 93339ddd567SLois Curfman McInnes 934c7fcc2eaSBarry Smith /* 935c7fcc2eaSBarry Smith The computation of n is wrong below, n should represent the number of local 936c7fcc2eaSBarry Smith rows in the right (column vector) 937c7fcc2eaSBarry Smith */ 938c7fcc2eaSBarry Smith 93939ddd567SLois Curfman McInnes /* each row stores all columns */ 94039ddd567SLois Curfman McInnes if (N == PETSC_DECIDE) N = n; 941d7e8b826SBarry Smith if (n == PETSC_DECIDE) {n = N/a->size + ((N % a->size) > a->rank);} 942a8c6a408SBarry Smith /* if (n != N) SETERRQ(PETSC_ERR_SUP,0,"For now, only n=N is supported"); */ 943aca0ad90SLois Curfman McInnes a->N = mat->N = N; 944aca0ad90SLois Curfman McInnes a->M = mat->M = M; 945aca0ad90SLois Curfman McInnes a->m = mat->m = m; 946aca0ad90SLois Curfman McInnes a->n = mat->n = n; 9478965ea79SLois Curfman McInnes 948c7fcc2eaSBarry Smith /* the information in the maps duplicates the information computed below, eventually 949c7fcc2eaSBarry Smith we should remove the duplicate information that is not contained in the maps */ 950488ecbafSBarry Smith ierr = MapCreateMPI(comm,m,M,&mat->rmap);CHKERRQ(ierr); 95196f6c058SBarry Smith ierr = MapCreateMPI(comm,PETSC_DECIDE,N,&mat->cmap);CHKERRQ(ierr); 952c7fcc2eaSBarry Smith 9538965ea79SLois Curfman McInnes /* build local table of row and column ownerships */ 954d7e8b826SBarry Smith a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int));CHKPTRQ(a->rowners); 955d7e8b826SBarry Smith a->cowners = a->rowners + a->size + 1; 956f09e8eb9SSatish Balay PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 957ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 9588965ea79SLois Curfman McInnes a->rowners[0] = 0; 9598965ea79SLois Curfman McInnes for ( i=2; i<=a->size; i++ ) { 9608965ea79SLois Curfman McInnes a->rowners[i] += a->rowners[i-1]; 9618965ea79SLois Curfman McInnes } 9628965ea79SLois Curfman McInnes a->rstart = a->rowners[a->rank]; 9638965ea79SLois Curfman McInnes a->rend = a->rowners[a->rank+1]; 964ca161407SBarry Smith ierr = MPI_Allgather(&n,1,MPI_INT,a->cowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 965d7e8b826SBarry Smith a->cowners[0] = 0; 966d7e8b826SBarry Smith for ( i=2; i<=a->size; i++ ) { 967d7e8b826SBarry Smith a->cowners[i] += a->cowners[i-1]; 968d7e8b826SBarry Smith } 9698965ea79SLois Curfman McInnes 970029af93fSBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,m,N,data,&a->A);CHKERRQ(ierr); 9718965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 9728965ea79SLois Curfman McInnes 9738965ea79SLois Curfman McInnes /* build cache for off array entries formed */ 9743782ba37SSatish Balay a->donotstash = 0; 9758798bf22SSatish Balay ierr = MatStashCreate_Private(comm,1,&mat->stash);CHKERRQ(ierr); 9768965ea79SLois Curfman McInnes 9778965ea79SLois Curfman McInnes /* stuff used for matrix vector multiply */ 9788965ea79SLois Curfman McInnes a->lvec = 0; 9798965ea79SLois Curfman McInnes a->Mvctx = 0; 98039b7565bSBarry Smith a->roworiented = 1; 9818965ea79SLois Curfman McInnes 982477f1c0bSLois Curfman McInnes *A = mat; 98325cdf11fSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr); 98425cdf11fSBarry Smith if (flg) { 9858c469469SLois Curfman McInnes ierr = MatPrintHelp(mat);CHKERRQ(ierr); 9868c469469SLois Curfman McInnes } 9873a40ed3dSBarry Smith PetscFunctionReturn(0); 9888965ea79SLois Curfman McInnes } 9898965ea79SLois Curfman McInnes 9905615d1e5SSatish Balay #undef __FUNC__ 9915609ef8eSBarry Smith #define __FUNC__ "MatDuplicate_MPIDense" 9925609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 9938965ea79SLois Curfman McInnes { 9948965ea79SLois Curfman McInnes Mat mat; 9953501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense *) A->data; 99639ddd567SLois Curfman McInnes int ierr; 9972ba99913SLois Curfman McInnes FactorCtx *factor; 9988965ea79SLois Curfman McInnes 9993a40ed3dSBarry Smith PetscFunctionBegin; 10008965ea79SLois Curfman McInnes *newmat = 0; 10013f1db9ecSBarry Smith PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,"Mat",A->comm,MatDestroy,MatView); 10028965ea79SLois Curfman McInnes PLogObjectCreate(mat); 10030452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIDense));CHKPTRQ(a); 1004549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1005e1311b90SBarry Smith mat->ops->destroy = MatDestroy_MPIDense; 1006e1311b90SBarry Smith mat->ops->view = MatView_MPIDense; 10073501a2bdSLois Curfman McInnes mat->factor = A->factor; 1008c456f294SBarry Smith mat->assembled = PETSC_TRUE; 10098965ea79SLois Curfman McInnes 101044cd7ae7SLois Curfman McInnes a->m = mat->m = oldmat->m; 101144cd7ae7SLois Curfman McInnes a->n = mat->n = oldmat->n; 101244cd7ae7SLois Curfman McInnes a->M = mat->M = oldmat->M; 101344cd7ae7SLois Curfman McInnes a->N = mat->N = oldmat->N; 10142ba99913SLois Curfman McInnes if (oldmat->factor) { 10152ba99913SLois Curfman McInnes a->factor = (FactorCtx *) (factor = PetscNew(FactorCtx));CHKPTRQ(factor); 10162ba99913SLois Curfman McInnes /* copy factor contents ... add this code! */ 10172ba99913SLois Curfman McInnes } else a->factor = 0; 10188965ea79SLois Curfman McInnes 10198965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 10208965ea79SLois Curfman McInnes a->rend = oldmat->rend; 10218965ea79SLois Curfman McInnes a->size = oldmat->size; 10228965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1023e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 10243782ba37SSatish Balay a->donotstash = oldmat->donotstash; 10250452661fSBarry Smith a->rowners = (int *) PetscMalloc((a->size+1)*sizeof(int));CHKPTRQ(a->rowners); 1026f09e8eb9SSatish Balay PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1027549d3d68SSatish Balay ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr); 10288798bf22SSatish Balay ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr); 10298965ea79SLois Curfman McInnes 10308965ea79SLois Curfman McInnes ierr = VecDuplicate(oldmat->lvec,&a->lvec);CHKERRQ(ierr); 10318965ea79SLois Curfman McInnes PLogObjectParent(mat,a->lvec); 103255659b69SBarry Smith ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx);CHKERRQ(ierr); 10338965ea79SLois Curfman McInnes PLogObjectParent(mat,a->Mvctx); 10345609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 10358965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 10368965ea79SLois Curfman McInnes *newmat = mat; 10373a40ed3dSBarry Smith PetscFunctionReturn(0); 10388965ea79SLois Curfman McInnes } 10398965ea79SLois Curfman McInnes 104077c4ece6SBarry Smith #include "sys.h" 10418965ea79SLois Curfman McInnes 10425615d1e5SSatish Balay #undef __FUNC__ 10435615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense_DenseInFile" 104490ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M, int N, Mat *newmat) 104590ace30eSBarry Smith { 104640011551SBarry Smith int *rowners, i,size,rank,m,ierr,nz,j; 104790ace30eSBarry Smith Scalar *array,*vals,*vals_ptr; 104890ace30eSBarry Smith MPI_Status status; 104990ace30eSBarry Smith 10503a40ed3dSBarry Smith PetscFunctionBegin; 1051d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1052d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 105390ace30eSBarry Smith 105490ace30eSBarry Smith /* determine ownership of all rows */ 105590ace30eSBarry Smith m = M/size + ((M % size) > rank); 105690ace30eSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners); 1057ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 105890ace30eSBarry Smith rowners[0] = 0; 105990ace30eSBarry Smith for ( i=2; i<=size; i++ ) { 106090ace30eSBarry Smith rowners[i] += rowners[i-1]; 106190ace30eSBarry Smith } 106290ace30eSBarry Smith 106390ace30eSBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 106490ace30eSBarry Smith ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr); 106590ace30eSBarry Smith 106690ace30eSBarry Smith if (!rank) { 106790ace30eSBarry Smith vals = (Scalar *) PetscMalloc( m*N*sizeof(Scalar) );CHKPTRQ(vals); 106890ace30eSBarry Smith 106990ace30eSBarry Smith /* read in my part of the matrix numerical values */ 10700752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 107190ace30eSBarry Smith 107290ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 107390ace30eSBarry Smith vals_ptr = vals; 107490ace30eSBarry Smith for ( i=0; i<m; i++ ) { 107590ace30eSBarry Smith for ( j=0; j<N; j++ ) { 107690ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 107790ace30eSBarry Smith } 107890ace30eSBarry Smith } 107990ace30eSBarry Smith 108090ace30eSBarry Smith /* read in other processors and ship out */ 108190ace30eSBarry Smith for ( i=1; i<size; i++ ) { 108290ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 10830752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1084ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 108590ace30eSBarry Smith } 10863a40ed3dSBarry Smith } else { 108790ace30eSBarry Smith /* receive numeric values */ 108890ace30eSBarry Smith vals = (Scalar*) PetscMalloc( m*N*sizeof(Scalar) );CHKPTRQ(vals); 108990ace30eSBarry Smith 109090ace30eSBarry Smith /* receive message of values*/ 1091ca161407SBarry Smith ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 109290ace30eSBarry Smith 109390ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 109490ace30eSBarry Smith vals_ptr = vals; 109590ace30eSBarry Smith for ( i=0; i<m; i++ ) { 109690ace30eSBarry Smith for ( j=0; j<N; j++ ) { 109790ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 109890ace30eSBarry Smith } 109990ace30eSBarry Smith } 110090ace30eSBarry Smith } 110190ace30eSBarry Smith PetscFree(rowners); 110290ace30eSBarry Smith PetscFree(vals); 11036d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11046d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11053a40ed3dSBarry Smith PetscFunctionReturn(0); 110690ace30eSBarry Smith } 110790ace30eSBarry Smith 110890ace30eSBarry Smith 11095615d1e5SSatish Balay #undef __FUNC__ 11105615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense" 111119bcc07fSBarry Smith int MatLoad_MPIDense(Viewer viewer,MatType type,Mat *newmat) 11128965ea79SLois Curfman McInnes { 11138965ea79SLois Curfman McInnes Mat A; 11148965ea79SLois Curfman McInnes Scalar *vals,*svals; 111519bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 11168965ea79SLois Curfman McInnes MPI_Status status; 11178965ea79SLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 11188965ea79SLois Curfman McInnes int *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols; 111919bcc07fSBarry Smith int tag = ((PetscObject)viewer)->tag; 11203a40ed3dSBarry Smith int i, nz, ierr, j,rstart, rend, fd; 11218965ea79SLois Curfman McInnes 11223a40ed3dSBarry Smith PetscFunctionBegin; 1123d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1124d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 11258965ea79SLois Curfman McInnes if (!rank) { 112690ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 11270752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr); 1128a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object"); 11298965ea79SLois Curfman McInnes } 11308965ea79SLois Curfman McInnes 1131ca161407SBarry Smith ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr); 113290ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 113390ace30eSBarry Smith 113490ace30eSBarry Smith /* 113590ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 113690ace30eSBarry Smith */ 113790ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 11383a40ed3dSBarry Smith ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr); 11393a40ed3dSBarry Smith PetscFunctionReturn(0); 114090ace30eSBarry Smith } 114190ace30eSBarry Smith 11428965ea79SLois Curfman McInnes /* determine ownership of all rows */ 11438965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 11440452661fSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners); 1145ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 11468965ea79SLois Curfman McInnes rowners[0] = 0; 11478965ea79SLois Curfman McInnes for ( i=2; i<=size; i++ ) { 11488965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 11498965ea79SLois Curfman McInnes } 11508965ea79SLois Curfman McInnes rstart = rowners[rank]; 11518965ea79SLois Curfman McInnes rend = rowners[rank+1]; 11528965ea79SLois Curfman McInnes 11538965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 11540452661fSBarry Smith ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) );CHKPTRQ(ourlens); 11558965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 11568965ea79SLois Curfman McInnes if (!rank) { 11570452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) );CHKPTRQ(rowlengths); 11580752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 11590452661fSBarry Smith sndcounts = (int*) PetscMalloc( size*sizeof(int) );CHKPTRQ(sndcounts); 11608965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i]; 1161ca161407SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 11620452661fSBarry Smith PetscFree(sndcounts); 1163ca161407SBarry Smith } else { 1164ca161407SBarry Smith ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);CHKERRQ(ierr); 11658965ea79SLois Curfman McInnes } 11668965ea79SLois Curfman McInnes 11678965ea79SLois Curfman McInnes if (!rank) { 11688965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 11690452661fSBarry Smith procsnz = (int*) PetscMalloc( size*sizeof(int) );CHKPTRQ(procsnz); 1170549d3d68SSatish Balay ierr = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr); 11718965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 11728965ea79SLois Curfman McInnes for ( j=rowners[i]; j< rowners[i+1]; j++ ) { 11738965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 11748965ea79SLois Curfman McInnes } 11758965ea79SLois Curfman McInnes } 11760452661fSBarry Smith PetscFree(rowlengths); 11778965ea79SLois Curfman McInnes 11788965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 11798965ea79SLois Curfman McInnes maxnz = 0; 11808965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 11810452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 11828965ea79SLois Curfman McInnes } 11830452661fSBarry Smith cols = (int *) PetscMalloc( maxnz*sizeof(int) );CHKPTRQ(cols); 11848965ea79SLois Curfman McInnes 11858965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 11868965ea79SLois Curfman McInnes nz = procsnz[0]; 11870452661fSBarry Smith mycols = (int *) PetscMalloc( nz*sizeof(int) );CHKPTRQ(mycols); 11880752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 11898965ea79SLois Curfman McInnes 11908965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 11918965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { 11928965ea79SLois Curfman McInnes nz = procsnz[i]; 11930752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 1194ca161407SBarry Smith ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr); 11958965ea79SLois Curfman McInnes } 11960452661fSBarry Smith PetscFree(cols); 11973a40ed3dSBarry Smith } else { 11988965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 11998965ea79SLois Curfman McInnes nz = 0; 12008965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12018965ea79SLois Curfman McInnes nz += ourlens[i]; 12028965ea79SLois Curfman McInnes } 12030452661fSBarry Smith mycols = (int*) PetscMalloc( nz*sizeof(int) );CHKPTRQ(mycols); 12048965ea79SLois Curfman McInnes 12058965ea79SLois Curfman McInnes /* receive message of column indices*/ 1206ca161407SBarry Smith ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr); 1207ca161407SBarry Smith ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr); 1208a8c6a408SBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 12098965ea79SLois Curfman McInnes } 12108965ea79SLois Curfman McInnes 12118965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 1212549d3d68SSatish Balay ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr); 12138965ea79SLois Curfman McInnes jj = 0; 12148965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12158965ea79SLois Curfman McInnes for ( j=0; j<ourlens[i]; j++ ) { 12168965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 12178965ea79SLois Curfman McInnes jj++; 12188965ea79SLois Curfman McInnes } 12198965ea79SLois Curfman McInnes } 12208965ea79SLois Curfman McInnes 12218965ea79SLois Curfman McInnes /* create our matrix */ 12228965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12238965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 12248965ea79SLois Curfman McInnes } 1225b4fd4287SBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 12268965ea79SLois Curfman McInnes A = *newmat; 12278965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12288965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 12298965ea79SLois Curfman McInnes } 12308965ea79SLois Curfman McInnes 12318965ea79SLois Curfman McInnes if (!rank) { 12320452661fSBarry Smith vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) );CHKPTRQ(vals); 12338965ea79SLois Curfman McInnes 12348965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 12358965ea79SLois Curfman McInnes nz = procsnz[0]; 12360752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 12378965ea79SLois Curfman McInnes 12388965ea79SLois Curfman McInnes /* insert into matrix */ 12398965ea79SLois Curfman McInnes jj = rstart; 12408965ea79SLois Curfman McInnes smycols = mycols; 12418965ea79SLois Curfman McInnes svals = vals; 12428965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12438965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 12448965ea79SLois Curfman McInnes smycols += ourlens[i]; 12458965ea79SLois Curfman McInnes svals += ourlens[i]; 12468965ea79SLois Curfman McInnes jj++; 12478965ea79SLois Curfman McInnes } 12488965ea79SLois Curfman McInnes 12498965ea79SLois Curfman McInnes /* read in other processors and ship out */ 12508965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { 12518965ea79SLois Curfman McInnes nz = procsnz[i]; 12520752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1253ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 12548965ea79SLois Curfman McInnes } 12550452661fSBarry Smith PetscFree(procsnz); 12563a40ed3dSBarry Smith } else { 12578965ea79SLois Curfman McInnes /* receive numeric values */ 12580452661fSBarry Smith vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) );CHKPTRQ(vals); 12598965ea79SLois Curfman McInnes 12608965ea79SLois Curfman McInnes /* receive message of values*/ 1261ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1262ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 1263a8c6a408SBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 12648965ea79SLois Curfman McInnes 12658965ea79SLois Curfman McInnes /* insert into matrix */ 12668965ea79SLois Curfman McInnes jj = rstart; 12678965ea79SLois Curfman McInnes smycols = mycols; 12688965ea79SLois Curfman McInnes svals = vals; 12698965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12708965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 12718965ea79SLois Curfman McInnes smycols += ourlens[i]; 12728965ea79SLois Curfman McInnes svals += ourlens[i]; 12738965ea79SLois Curfman McInnes jj++; 12748965ea79SLois Curfman McInnes } 12758965ea79SLois Curfman McInnes } 12760452661fSBarry Smith PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners); 12778965ea79SLois Curfman McInnes 12786d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12796d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12803a40ed3dSBarry Smith PetscFunctionReturn(0); 12818965ea79SLois Curfman McInnes } 128290ace30eSBarry Smith 128390ace30eSBarry Smith 128490ace30eSBarry Smith 128590ace30eSBarry Smith 128690ace30eSBarry Smith 1287