1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*606d414cSSatish Balay static char vcid[] = "$Id: mpidense.c,v 1.116 1999/06/09 23:19:38 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]; 226*606d414cSSatish Balay ierr = PetscFree(work);CHKERRQ(ierr); 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 } 257*606d414cSSatish Balay ierr = PetscFree(starts);CHKERRQ(ierr); 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 } 274*606d414cSSatish Balay ierr = PetscFree(recv_waits);CHKERRQ(ierr); 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 } 285*606d414cSSatish Balay ierr = PetscFree(rvalues);CHKERRQ(ierr); 286*606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 287*606d414cSSatish Balay ierr = PetscFree(owner);CHKERRQ(ierr); 288*606d414cSSatish Balay ierr = PetscFree(nprocs);CHKERRQ(ierr); 2898965ea79SLois Curfman McInnes 2908965ea79SLois Curfman McInnes /* actually zap the local rows */ 291029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 2928965ea79SLois Curfman McInnes PLogObjectParent(A,istmp); 293*606d414cSSatish Balay ierr = PetscFree(lrows);CHKERRQ(ierr); 2948965ea79SLois Curfman McInnes ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr); 2958965ea79SLois Curfman McInnes ierr = ISDestroy(istmp);CHKERRQ(ierr); 2968965ea79SLois Curfman McInnes 2978965ea79SLois Curfman McInnes /* wait on sends */ 2988965ea79SLois Curfman McInnes if (nsends) { 2997056b6fcSBarry Smith send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status); 300ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 301*606d414cSSatish Balay ierr = PetscFree(send_status);CHKERRQ(ierr); 3028965ea79SLois Curfman McInnes } 303*606d414cSSatish Balay ierr = PetscFree(send_waits);CHKERRQ(ierr); 304*606d414cSSatish Balay ierr = PetscFree(svalues);CHKERRQ(ierr); 3058965ea79SLois Curfman McInnes 3063a40ed3dSBarry Smith PetscFunctionReturn(0); 3078965ea79SLois Curfman McInnes } 3088965ea79SLois Curfman McInnes 3095615d1e5SSatish Balay #undef __FUNC__ 3105615d1e5SSatish Balay #define __FUNC__ "MatMult_MPIDense" 3118f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 3128965ea79SLois Curfman McInnes { 31339ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3148965ea79SLois Curfman McInnes int ierr; 315c456f294SBarry Smith 3163a40ed3dSBarry Smith PetscFunctionBegin; 31743a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 31843a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 31944cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 3203a40ed3dSBarry Smith PetscFunctionReturn(0); 3218965ea79SLois Curfman McInnes } 3228965ea79SLois Curfman McInnes 3235615d1e5SSatish Balay #undef __FUNC__ 3245615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_MPIDense" 3258f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 3268965ea79SLois Curfman McInnes { 32739ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3288965ea79SLois Curfman McInnes int ierr; 329c456f294SBarry Smith 3303a40ed3dSBarry Smith PetscFunctionBegin; 33143a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 33243a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 33344cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 3343a40ed3dSBarry Smith PetscFunctionReturn(0); 3358965ea79SLois Curfman McInnes } 3368965ea79SLois Curfman McInnes 3375615d1e5SSatish Balay #undef __FUNC__ 3385615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_MPIDense" 3398f6be9afSLois Curfman McInnes int MatMultTrans_MPIDense(Mat A,Vec xx,Vec yy) 340096963f5SLois Curfman McInnes { 341096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 342096963f5SLois Curfman McInnes int ierr; 3433501a2bdSLois Curfman McInnes Scalar zero = 0.0; 344096963f5SLois Curfman McInnes 3453a40ed3dSBarry Smith PetscFunctionBegin; 3463501a2bdSLois Curfman McInnes ierr = VecSet(&zero,yy);CHKERRQ(ierr); 34744cd7ae7SLois Curfman McInnes ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 348537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 349537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 3503a40ed3dSBarry Smith PetscFunctionReturn(0); 351096963f5SLois Curfman McInnes } 352096963f5SLois Curfman McInnes 3535615d1e5SSatish Balay #undef __FUNC__ 3545615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_MPIDense" 3558f6be9afSLois Curfman McInnes int MatMultTransAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 356096963f5SLois Curfman McInnes { 357096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 358096963f5SLois Curfman McInnes int ierr; 359096963f5SLois Curfman McInnes 3603a40ed3dSBarry Smith PetscFunctionBegin; 3613501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz);CHKERRQ(ierr); 36244cd7ae7SLois Curfman McInnes ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 363537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 364537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 3653a40ed3dSBarry Smith PetscFunctionReturn(0); 366096963f5SLois Curfman McInnes } 367096963f5SLois Curfman McInnes 3685615d1e5SSatish Balay #undef __FUNC__ 3695615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_MPIDense" 3708f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v) 3718965ea79SLois Curfman McInnes { 37239ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 373096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense *) a->A->data; 37444cd7ae7SLois Curfman McInnes int ierr, len, i, n, m = a->m, radd; 37544cd7ae7SLois Curfman McInnes Scalar *x, zero = 0.0; 376ed3cc1f0SBarry Smith 3773a40ed3dSBarry Smith PetscFunctionBegin; 37844cd7ae7SLois Curfman McInnes VecSet(&zero,v); 379096963f5SLois Curfman McInnes ierr = VecGetArray(v,&x);CHKERRQ(ierr); 380096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n);CHKERRQ(ierr); 381a8c6a408SBarry Smith if (n != a->M) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec"); 38244cd7ae7SLois Curfman McInnes len = PetscMin(aloc->m,aloc->n); 3837ddc982cSLois Curfman McInnes radd = a->rstart*m; 38444cd7ae7SLois Curfman McInnes for ( i=0; i<len; i++ ) { 385096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 386096963f5SLois Curfman McInnes } 3879a8c540fSSatish Balay ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 3883a40ed3dSBarry Smith PetscFunctionReturn(0); 3898965ea79SLois Curfman McInnes } 3908965ea79SLois Curfman McInnes 3915615d1e5SSatish Balay #undef __FUNC__ 3925615d1e5SSatish Balay #define __FUNC__ "MatDestroy_MPIDense" 393e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat) 3948965ea79SLois Curfman McInnes { 3953501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3968965ea79SLois Curfman McInnes int ierr; 397ed3cc1f0SBarry Smith 3983a40ed3dSBarry Smith PetscFunctionBegin; 39994d884c6SBarry Smith if (--mat->refct > 0) PetscFunctionReturn(0); 40094d884c6SBarry Smith 40194d884c6SBarry Smith if (mat->mapping) { 40294d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->mapping);CHKERRQ(ierr); 40394d884c6SBarry Smith } 40494d884c6SBarry Smith if (mat->bmapping) { 40594d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->bmapping);CHKERRQ(ierr); 40694d884c6SBarry Smith } 407aa482453SBarry Smith #if defined(PETSC_USE_LOG) 408e1311b90SBarry Smith PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mdn->M,mdn->N); 4098965ea79SLois Curfman McInnes #endif 4108798bf22SSatish Balay ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 411*606d414cSSatish Balay ierr = PetscFree(mdn->rowners);CHKERRQ(ierr); 4123501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A);CHKERRQ(ierr); 4133501a2bdSLois Curfman McInnes if (mdn->lvec) VecDestroy(mdn->lvec); 4143501a2bdSLois Curfman McInnes if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 415622d7880SLois Curfman McInnes if (mdn->factor) { 416*606d414cSSatish Balay if (mdn->factor->temp) {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);} 417*606d414cSSatish Balay if (mdn->factor->tag) {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);} 418*606d414cSSatish Balay if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);} 419*606d414cSSatish Balay ierr = PetscFree(mdn->factor);CHKERRQ(ierr); 420622d7880SLois Curfman McInnes } 421*606d414cSSatish Balay ierr = PetscFree(mdn);CHKERRQ(ierr); 42261b13de0SBarry Smith if (mat->rmap) { 42361b13de0SBarry Smith ierr = MapDestroy(mat->rmap);CHKERRQ(ierr); 42461b13de0SBarry Smith } 42561b13de0SBarry Smith if (mat->cmap) { 42661b13de0SBarry Smith ierr = MapDestroy(mat->cmap);CHKERRQ(ierr); 42790f02eecSBarry Smith } 4288965ea79SLois Curfman McInnes PLogObjectDestroy(mat); 4290452661fSBarry Smith PetscHeaderDestroy(mat); 4303a40ed3dSBarry Smith PetscFunctionReturn(0); 4318965ea79SLois Curfman McInnes } 43239ddd567SLois Curfman McInnes 4335615d1e5SSatish Balay #undef __FUNC__ 4345615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense_Binary" 43539ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer) 4368965ea79SLois Curfman McInnes { 43739ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 4388965ea79SLois Curfman McInnes int ierr; 4397056b6fcSBarry Smith 4403a40ed3dSBarry Smith PetscFunctionBegin; 44139ddd567SLois Curfman McInnes if (mdn->size == 1) { 44239ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 4438965ea79SLois Curfman McInnes } 444a8c6a408SBarry Smith else SETERRQ(PETSC_ERR_SUP,0,"Only uniprocessor output supported"); 4453a40ed3dSBarry Smith PetscFunctionReturn(0); 4468965ea79SLois Curfman McInnes } 4478965ea79SLois Curfman McInnes 4485615d1e5SSatish Balay #undef __FUNC__ 4495615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense_ASCII" 45039ddd567SLois Curfman McInnes static int MatView_MPIDense_ASCII(Mat mat,Viewer viewer) 4518965ea79SLois Curfman McInnes { 45239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 45377ed5343SBarry Smith int ierr, format, size = mdn->size, rank = mdn->rank; 4548965ea79SLois Curfman McInnes FILE *fd; 45519bcc07fSBarry Smith ViewerType vtype; 4568965ea79SLois Curfman McInnes 4573a40ed3dSBarry Smith PetscFunctionBegin; 4583a40ed3dSBarry Smith ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 45990ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd);CHKERRQ(ierr); 460888f2ed8SSatish Balay ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 461639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 4624e220ebcSLois Curfman McInnes MatInfo info; 463888f2ed8SSatish Balay ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 46477c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 4654e220ebcSLois Curfman McInnes fprintf(fd," [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mdn->m, 4664e220ebcSLois Curfman McInnes (int)info.nz_used,(int)info.nz_allocated,(int)info.memory); 467096963f5SLois Curfman McInnes fflush(fd); 46877c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 4693501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 4703a40ed3dSBarry Smith PetscFunctionReturn(0); 47196f6c058SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_INFO) { 4723a40ed3dSBarry Smith PetscFunctionReturn(0); 4738965ea79SLois Curfman McInnes } 47477ed5343SBarry Smith 4758965ea79SLois Curfman McInnes if (size == 1) { 47639ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 4773a40ed3dSBarry Smith } else { 4788965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 4798965ea79SLois Curfman McInnes Mat A; 48039ddd567SLois Curfman McInnes int M = mdn->M, N = mdn->N,m,row,i, nz, *cols; 48139ddd567SLois Curfman McInnes Scalar *vals; 48239ddd567SLois Curfman McInnes Mat_SeqDense *Amdn = (Mat_SeqDense*) mdn->A->data; 4838965ea79SLois Curfman McInnes 4848965ea79SLois Curfman McInnes if (!rank) { 4850513a670SBarry Smith ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr); 4863a40ed3dSBarry Smith } else { 4870513a670SBarry Smith ierr = MatCreateMPIDense(mat->comm,0,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr); 4888965ea79SLois Curfman McInnes } 4898965ea79SLois Curfman McInnes PLogObjectParent(mat,A); 4908965ea79SLois Curfman McInnes 49139ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 49239ddd567SLois Curfman McInnes but it's quick for now */ 49339ddd567SLois Curfman McInnes row = mdn->rstart; m = Amdn->m; 4948965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 49539ddd567SLois Curfman McInnes ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 49639ddd567SLois Curfman McInnes ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 49739ddd567SLois Curfman McInnes ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 49839ddd567SLois Curfman McInnes row++; 4998965ea79SLois Curfman McInnes } 5008965ea79SLois Curfman McInnes 5016d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5026d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5038965ea79SLois Curfman McInnes if (!rank) { 50439ddd567SLois Curfman McInnes ierr = MatView(((Mat_MPIDense*)(A->data))->A,viewer);CHKERRQ(ierr); 5058965ea79SLois Curfman McInnes } 5068965ea79SLois Curfman McInnes ierr = MatDestroy(A);CHKERRQ(ierr); 5078965ea79SLois Curfman McInnes } 5083a40ed3dSBarry Smith PetscFunctionReturn(0); 5098965ea79SLois Curfman McInnes } 5108965ea79SLois Curfman McInnes 5115615d1e5SSatish Balay #undef __FUNC__ 5125615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense" 513e1311b90SBarry Smith int MatView_MPIDense(Mat mat,Viewer viewer) 5148965ea79SLois Curfman McInnes { 51539ddd567SLois Curfman McInnes int ierr; 516bcd2baecSBarry Smith ViewerType vtype; 5178965ea79SLois Curfman McInnes 518bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 5193f1db9ecSBarry Smith if (PetscTypeCompare(vtype,ASCII_VIEWER)) { 52039ddd567SLois Curfman McInnes ierr = MatView_MPIDense_ASCII(mat,viewer);CHKERRQ(ierr); 5213f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,BINARY_VIEWER)) { 5223a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 5235cd90555SBarry Smith } else { 5245cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 5258965ea79SLois Curfman McInnes } 5263a40ed3dSBarry Smith PetscFunctionReturn(0); 5278965ea79SLois Curfman McInnes } 5288965ea79SLois Curfman McInnes 5295615d1e5SSatish Balay #undef __FUNC__ 5305615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_MPIDense" 5318f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 5328965ea79SLois Curfman McInnes { 5333501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 5343501a2bdSLois Curfman McInnes Mat mdn = mat->A; 5354e220ebcSLois Curfman McInnes int ierr; 5364e220ebcSLois Curfman McInnes double isend[5], irecv[5]; 5378965ea79SLois Curfman McInnes 5383a40ed3dSBarry Smith PetscFunctionBegin; 5394e220ebcSLois Curfman McInnes info->rows_global = (double)mat->M; 5404e220ebcSLois Curfman McInnes info->columns_global = (double)mat->N; 5414e220ebcSLois Curfman McInnes info->rows_local = (double)mat->m; 5424e220ebcSLois Curfman McInnes info->columns_local = (double)mat->N; 5434e220ebcSLois Curfman McInnes info->block_size = 1.0; 5444e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 5454e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 5464e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 5478965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 5484e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 5494e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 5504e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 5514e220ebcSLois Curfman McInnes info->memory = isend[3]; 5524e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 5538965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 554f7cdd7c9SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr); 5554e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 5564e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 5574e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 5584e220ebcSLois Curfman McInnes info->memory = irecv[3]; 5594e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 5608965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 561f7cdd7c9SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 5624e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 5634e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 5644e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 5654e220ebcSLois Curfman McInnes info->memory = irecv[3]; 5664e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 5678965ea79SLois Curfman McInnes } 5684e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 5694e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 5704e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 5713a40ed3dSBarry Smith PetscFunctionReturn(0); 5728965ea79SLois Curfman McInnes } 5738965ea79SLois Curfman McInnes 5748c469469SLois Curfman McInnes /* extern int MatLUFactorSymbolic_MPIDense(Mat,IS,IS,double,Mat*); 5758aaee692SLois Curfman McInnes extern int MatLUFactorNumeric_MPIDense(Mat,Mat*); 5768aaee692SLois Curfman McInnes extern int MatLUFactor_MPIDense(Mat,IS,IS,double); 5778aaee692SLois Curfman McInnes extern int MatSolve_MPIDense(Mat,Vec,Vec); 5788c469469SLois Curfman McInnes extern int MatSolveAdd_MPIDense(Mat,Vec,Vec,Vec); 5798aaee692SLois Curfman McInnes extern int MatSolveTrans_MPIDense(Mat,Vec,Vec); 5808aaee692SLois Curfman McInnes extern int MatSolveTransAdd_MPIDense(Mat,Vec,Vec,Vec); */ 5818aaee692SLois Curfman McInnes 5825615d1e5SSatish Balay #undef __FUNC__ 5835615d1e5SSatish Balay #define __FUNC__ "MatSetOption_MPIDense" 5848f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op) 5858965ea79SLois Curfman McInnes { 58639ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 5878965ea79SLois Curfman McInnes 5883a40ed3dSBarry Smith PetscFunctionBegin; 5896d4a8577SBarry Smith if (op == MAT_NO_NEW_NONZERO_LOCATIONS || 5906d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 5914787f768SSatish Balay op == MAT_NEW_NONZERO_LOCATION_ERR || 5924787f768SSatish Balay op == MAT_NEW_NONZERO_ALLOCATION_ERR || 593219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_SORTED || 594219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED) { 595b1fbbac0SLois Curfman McInnes MatSetOption(a->A,op); 596b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROW_ORIENTED) { 597aeafbbfcSLois Curfman McInnes a->roworiented = 1; 5988965ea79SLois Curfman McInnes MatSetOption(a->A,op); 599b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROWS_SORTED || 600219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 6016d4a8577SBarry Smith op == MAT_SYMMETRIC || 6026d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 603b51ba29fSSatish Balay op == MAT_YES_NEW_DIAGONALS || 604b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 605981c4779SBarry Smith PLogInfo(A,"MatSetOption_MPIDense:Option ignored\n"); 6063a40ed3dSBarry Smith } else if (op == MAT_COLUMN_ORIENTED) { 6073a40ed3dSBarry Smith a->roworiented = 0; MatSetOption(a->A,op); 6083782ba37SSatish Balay } else if (op == MAT_IGNORE_OFF_PROC_ENTRIES) { 6093782ba37SSatish Balay a->donotstash = 1; 6103a40ed3dSBarry Smith } else if (op == MAT_NO_NEW_DIAGONALS) { 6113a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 6123a40ed3dSBarry Smith } else { 6133a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 6143a40ed3dSBarry Smith } 6153a40ed3dSBarry Smith PetscFunctionReturn(0); 6168965ea79SLois Curfman McInnes } 6178965ea79SLois Curfman McInnes 6185615d1e5SSatish Balay #undef __FUNC__ 6195615d1e5SSatish Balay #define __FUNC__ "MatGetSize_MPIDense" 6208f6be9afSLois Curfman McInnes int MatGetSize_MPIDense(Mat A,int *m,int *n) 6218965ea79SLois Curfman McInnes { 6223501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6233a40ed3dSBarry Smith 6243a40ed3dSBarry Smith PetscFunctionBegin; 6258965ea79SLois Curfman McInnes *m = mat->M; *n = mat->N; 6263a40ed3dSBarry Smith PetscFunctionReturn(0); 6278965ea79SLois Curfman McInnes } 6288965ea79SLois Curfman McInnes 6295615d1e5SSatish Balay #undef __FUNC__ 6305615d1e5SSatish Balay #define __FUNC__ "MatGetLocalSize_MPIDense" 6318f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIDense(Mat A,int *m,int *n) 6328965ea79SLois Curfman McInnes { 6333501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6343a40ed3dSBarry Smith 6353a40ed3dSBarry Smith PetscFunctionBegin; 6368965ea79SLois Curfman McInnes *m = mat->m; *n = mat->N; 6373a40ed3dSBarry Smith PetscFunctionReturn(0); 6388965ea79SLois Curfman McInnes } 6398965ea79SLois Curfman McInnes 6405615d1e5SSatish Balay #undef __FUNC__ 6415615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_MPIDense" 6428f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n) 6438965ea79SLois Curfman McInnes { 6443501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6453a40ed3dSBarry Smith 6463a40ed3dSBarry Smith PetscFunctionBegin; 6478965ea79SLois Curfman McInnes *m = mat->rstart; *n = mat->rend; 6483a40ed3dSBarry Smith PetscFunctionReturn(0); 6498965ea79SLois Curfman McInnes } 6508965ea79SLois Curfman McInnes 6515615d1e5SSatish Balay #undef __FUNC__ 6525615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIDense" 6538f6be9afSLois Curfman McInnes int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v) 6548965ea79SLois Curfman McInnes { 6553501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6563a40ed3dSBarry Smith int lrow, rstart = mat->rstart, rend = mat->rend,ierr; 6578965ea79SLois Curfman McInnes 6583a40ed3dSBarry Smith PetscFunctionBegin; 659a8c6a408SBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,0,"only local rows") 6608965ea79SLois Curfman McInnes lrow = row - rstart; 6613a40ed3dSBarry Smith ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr); 6623a40ed3dSBarry Smith PetscFunctionReturn(0); 6638965ea79SLois Curfman McInnes } 6648965ea79SLois Curfman McInnes 6655615d1e5SSatish Balay #undef __FUNC__ 6665615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_MPIDense" 6678f6be9afSLois Curfman McInnes int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v) 6688965ea79SLois Curfman McInnes { 669*606d414cSSatish Balay int ierr; 670*606d414cSSatish Balay 6713a40ed3dSBarry Smith PetscFunctionBegin; 672*606d414cSSatish Balay if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 673*606d414cSSatish Balay if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 6743a40ed3dSBarry Smith PetscFunctionReturn(0); 6758965ea79SLois Curfman McInnes } 6768965ea79SLois Curfman McInnes 6775615d1e5SSatish Balay #undef __FUNC__ 6785615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIDense" 6798f6be9afSLois Curfman McInnes int MatNorm_MPIDense(Mat A,NormType type,double *norm) 680096963f5SLois Curfman McInnes { 6813501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) A->data; 6823501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*) mdn->A->data; 6833501a2bdSLois Curfman McInnes int ierr, i, j; 6843501a2bdSLois Curfman McInnes double sum = 0.0; 6853501a2bdSLois Curfman McInnes Scalar *v = mat->v; 6863501a2bdSLois Curfman McInnes 6873a40ed3dSBarry Smith PetscFunctionBegin; 6883501a2bdSLois Curfman McInnes if (mdn->size == 1) { 6893501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,norm);CHKERRQ(ierr); 6903501a2bdSLois Curfman McInnes } else { 6913501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 6923501a2bdSLois Curfman McInnes for (i=0; i<mat->n*mat->m; i++ ) { 693aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 694e20fef11SSatish Balay sum += PetscReal(PetscConj(*v)*(*v)); v++; 6953501a2bdSLois Curfman McInnes #else 6963501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 6973501a2bdSLois Curfman McInnes #endif 6983501a2bdSLois Curfman McInnes } 699ca161407SBarry Smith ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 7003501a2bdSLois Curfman McInnes *norm = sqrt(*norm); 7013501a2bdSLois Curfman McInnes PLogFlops(2*mat->n*mat->m); 7023a40ed3dSBarry Smith } else if (type == NORM_1) { 7033501a2bdSLois Curfman McInnes double *tmp, *tmp2; 7040452661fSBarry Smith tmp = (double *) PetscMalloc( 2*mdn->N*sizeof(double) );CHKPTRQ(tmp); 7053501a2bdSLois Curfman McInnes tmp2 = tmp + mdn->N; 706549d3d68SSatish Balay ierr = PetscMemzero(tmp,2*mdn->N*sizeof(double));CHKERRQ(ierr); 707096963f5SLois Curfman McInnes *norm = 0.0; 7083501a2bdSLois Curfman McInnes v = mat->v; 7093501a2bdSLois Curfman McInnes for ( j=0; j<mat->n; j++ ) { 7103501a2bdSLois Curfman McInnes for ( i=0; i<mat->m; i++ ) { 71167e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 7123501a2bdSLois Curfman McInnes } 7133501a2bdSLois Curfman McInnes } 714ca161407SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,mdn->N,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 7153501a2bdSLois Curfman McInnes for ( j=0; j<mdn->N; j++ ) { 7163501a2bdSLois Curfman McInnes if (tmp2[j] > *norm) *norm = tmp2[j]; 7173501a2bdSLois Curfman McInnes } 718*606d414cSSatish Balay ierr = PetscFree(tmp);CHKERRQ(ierr); 7193501a2bdSLois Curfman McInnes PLogFlops(mat->n*mat->m); 7203a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 7213501a2bdSLois Curfman McInnes double ntemp; 7223501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 723ca161407SBarry Smith ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr); 7243a40ed3dSBarry Smith } else { 725a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"No support for two norm"); 7263501a2bdSLois Curfman McInnes } 7273501a2bdSLois Curfman McInnes } 7283a40ed3dSBarry Smith PetscFunctionReturn(0); 7293501a2bdSLois Curfman McInnes } 7303501a2bdSLois Curfman McInnes 7315615d1e5SSatish Balay #undef __FUNC__ 7325615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIDense" 7338f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout) 7343501a2bdSLois Curfman McInnes { 7353501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 7363501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense *) a->A->data; 7373501a2bdSLois Curfman McInnes Mat B; 7383501a2bdSLois Curfman McInnes int M = a->M, N = a->N, m, n, *rwork, rstart = a->rstart; 7393501a2bdSLois Curfman McInnes int j, i, ierr; 7403501a2bdSLois Curfman McInnes Scalar *v; 7413501a2bdSLois Curfman McInnes 7423a40ed3dSBarry Smith PetscFunctionBegin; 7437056b6fcSBarry Smith if (matout == PETSC_NULL && M != N) { 744a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Supports square matrix only in-place"); 7457056b6fcSBarry Smith } 7467056b6fcSBarry Smith ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr); 7473501a2bdSLois Curfman McInnes 7483501a2bdSLois Curfman McInnes m = Aloc->m; n = Aloc->n; v = Aloc->v; 7490452661fSBarry Smith rwork = (int *) PetscMalloc(n*sizeof(int));CHKPTRQ(rwork); 7503501a2bdSLois Curfman McInnes for ( j=0; j<n; j++ ) { 7513501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 7523501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 7533501a2bdSLois Curfman McInnes v += m; 7543501a2bdSLois Curfman McInnes } 755*606d414cSSatish Balay ierr = PetscFree(rwork);CHKERRQ(ierr); 7566d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7576d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7583638b69dSLois Curfman McInnes if (matout != PETSC_NULL) { 7593501a2bdSLois Curfman McInnes *matout = B; 7603501a2bdSLois Curfman McInnes } else { 761f830108cSBarry Smith PetscOps *Abops; 76209dc0095SBarry Smith MatOps Aops; 763f830108cSBarry Smith 7643501a2bdSLois Curfman McInnes /* This isn't really an in-place transpose, but free data struct from a */ 765*606d414cSSatish Balay ierr = PetscFree(a->rowners);CHKERRQ(ierr); 7663501a2bdSLois Curfman McInnes ierr = MatDestroy(a->A);CHKERRQ(ierr); 7673501a2bdSLois Curfman McInnes if (a->lvec) VecDestroy(a->lvec); 7683501a2bdSLois Curfman McInnes if (a->Mvctx) VecScatterDestroy(a->Mvctx); 769*606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 770f830108cSBarry Smith 771f830108cSBarry Smith /* 772f830108cSBarry Smith This is horrible, horrible code. We need to keep the 773f830108cSBarry Smith A pointers for the bops and ops but copy everything 774f830108cSBarry Smith else from C. 775f830108cSBarry Smith */ 776f830108cSBarry Smith Abops = A->bops; 777f830108cSBarry Smith Aops = A->ops; 778549d3d68SSatish Balay ierr = PetscMemcpy(A,B,sizeof(struct _p_Mat));CHKERRQ(ierr); 779f830108cSBarry Smith A->bops = Abops; 780f830108cSBarry Smith A->ops = Aops; 781f830108cSBarry Smith 7820452661fSBarry Smith PetscHeaderDestroy(B); 7833501a2bdSLois Curfman McInnes } 7843a40ed3dSBarry Smith PetscFunctionReturn(0); 785096963f5SLois Curfman McInnes } 786096963f5SLois Curfman McInnes 787eadb2fb4SBarry Smith #include "pinclude/blaslapack.h" 7885615d1e5SSatish Balay #undef __FUNC__ 7895615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIDense" 7908f6be9afSLois Curfman McInnes int MatScale_MPIDense(Scalar *alpha,Mat inA) 79144cd7ae7SLois Curfman McInnes { 79244cd7ae7SLois Curfman McInnes Mat_MPIDense *A = (Mat_MPIDense *) inA->data; 79344cd7ae7SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->A->data; 79444cd7ae7SLois Curfman McInnes int one = 1, nz; 79544cd7ae7SLois Curfman McInnes 7963a40ed3dSBarry Smith PetscFunctionBegin; 79744cd7ae7SLois Curfman McInnes nz = a->m*a->n; 79844cd7ae7SLois Curfman McInnes BLscal_( &nz, alpha, a->v, &one ); 79944cd7ae7SLois Curfman McInnes PLogFlops(nz); 8003a40ed3dSBarry Smith PetscFunctionReturn(0); 80144cd7ae7SLois Curfman McInnes } 80244cd7ae7SLois Curfman McInnes 8035609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 8047b2a1423SBarry Smith extern int MatGetSubMatrices_MPIDense(Mat,int,IS *,IS *,MatReuse,Mat **); 8058965ea79SLois Curfman McInnes 8068965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 80709dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 80809dc0095SBarry Smith MatGetRow_MPIDense, 80909dc0095SBarry Smith MatRestoreRow_MPIDense, 81009dc0095SBarry Smith MatMult_MPIDense, 81109dc0095SBarry Smith MatMultAdd_MPIDense, 81209dc0095SBarry Smith MatMultTrans_MPIDense, 81309dc0095SBarry Smith MatMultTransAdd_MPIDense, 8148965ea79SLois Curfman McInnes 0, 81509dc0095SBarry Smith 0, 81609dc0095SBarry Smith 0, 81709dc0095SBarry Smith 0, 81809dc0095SBarry Smith 0, 81909dc0095SBarry Smith 0, 82009dc0095SBarry Smith 0, 82109dc0095SBarry Smith MatTranspose_MPIDense, 82209dc0095SBarry Smith MatGetInfo_MPIDense,0, 82309dc0095SBarry Smith MatGetDiagonal_MPIDense, 82409dc0095SBarry Smith 0, 82509dc0095SBarry Smith MatNorm_MPIDense, 82609dc0095SBarry Smith MatAssemblyBegin_MPIDense, 82709dc0095SBarry Smith MatAssemblyEnd_MPIDense, 82809dc0095SBarry Smith 0, 82909dc0095SBarry Smith MatSetOption_MPIDense, 83009dc0095SBarry Smith MatZeroEntries_MPIDense, 83109dc0095SBarry Smith MatZeroRows_MPIDense, 83209dc0095SBarry Smith 0, 83309dc0095SBarry Smith 0, 83409dc0095SBarry Smith 0, 83509dc0095SBarry Smith 0, 83609dc0095SBarry Smith MatGetSize_MPIDense, 83709dc0095SBarry Smith MatGetLocalSize_MPIDense, 83839ddd567SLois Curfman McInnes MatGetOwnershipRange_MPIDense, 83909dc0095SBarry Smith 0, 84009dc0095SBarry Smith 0, 84109dc0095SBarry Smith MatGetArray_MPIDense, 84209dc0095SBarry Smith MatRestoreArray_MPIDense, 8435609ef8eSBarry Smith MatDuplicate_MPIDense, 84409dc0095SBarry Smith 0, 84509dc0095SBarry Smith 0, 84609dc0095SBarry Smith 0, 84709dc0095SBarry Smith 0, 84809dc0095SBarry Smith 0, 8492ce60cd0SSatish Balay MatGetSubMatrices_MPIDense, 85009dc0095SBarry Smith 0, 85109dc0095SBarry Smith MatGetValues_MPIDense, 85209dc0095SBarry Smith 0, 85309dc0095SBarry Smith 0, 85409dc0095SBarry Smith MatScale_MPIDense, 85509dc0095SBarry Smith 0, 85609dc0095SBarry Smith 0, 85709dc0095SBarry Smith 0, 85809dc0095SBarry Smith MatGetBlockSize_MPIDense, 85909dc0095SBarry Smith 0, 86009dc0095SBarry Smith 0, 86109dc0095SBarry Smith 0, 86209dc0095SBarry Smith 0, 86309dc0095SBarry Smith 0, 86409dc0095SBarry Smith 0, 86509dc0095SBarry Smith 0, 86609dc0095SBarry Smith 0, 86709dc0095SBarry Smith 0, 86809dc0095SBarry Smith 0, 86909dc0095SBarry Smith 0, 87009dc0095SBarry Smith 0, 87109dc0095SBarry Smith MatGetMaps_Petsc}; 8728965ea79SLois Curfman McInnes 8735615d1e5SSatish Balay #undef __FUNC__ 8745615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIDense" 8758965ea79SLois Curfman McInnes /*@C 87639ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 8778965ea79SLois Curfman McInnes 878db81eaa0SLois Curfman McInnes Collective on MPI_Comm 879db81eaa0SLois Curfman McInnes 8808965ea79SLois Curfman McInnes Input Parameters: 881db81eaa0SLois Curfman McInnes + comm - MPI communicator 8828965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 883db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 8848965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 885db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 886db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 887dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 8888965ea79SLois Curfman McInnes 8898965ea79SLois Curfman McInnes Output Parameter: 890477f1c0bSLois Curfman McInnes . A - the matrix 8918965ea79SLois Curfman McInnes 892b259b22eSLois Curfman McInnes Notes: 89339ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 89439ddd567SLois Curfman McInnes storage by columns. 8958965ea79SLois Curfman McInnes 89618f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 89718f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 898b4fd4287SBarry Smith set data=PETSC_NULL. 89918f449edSLois Curfman McInnes 9008965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 9018965ea79SLois Curfman McInnes (possibly both). 9028965ea79SLois Curfman McInnes 9033501a2bdSLois Curfman McInnes Currently, the only parallel dense matrix decomposition is by rows, 9043501a2bdSLois Curfman McInnes so that n=N and each submatrix owns all of the global columns. 9053501a2bdSLois Curfman McInnes 906027ccd11SLois Curfman McInnes Level: intermediate 907027ccd11SLois Curfman McInnes 90839ddd567SLois Curfman McInnes .keywords: matrix, dense, parallel 9098965ea79SLois Curfman McInnes 91039ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 9118965ea79SLois Curfman McInnes @*/ 912477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A) 9138965ea79SLois Curfman McInnes { 9148965ea79SLois Curfman McInnes Mat mat; 91539ddd567SLois Curfman McInnes Mat_MPIDense *a; 91625cdf11fSBarry Smith int ierr, i,flg; 9178965ea79SLois Curfman McInnes 9183a40ed3dSBarry Smith PetscFunctionBegin; 919ed2daf61SLois Curfman McInnes /* Note: For now, when data is specified above, this assumes the user correctly 920ed2daf61SLois Curfman McInnes allocates the local dense storage space. We should add error checking. */ 92118f449edSLois Curfman McInnes 922477f1c0bSLois Curfman McInnes *A = 0; 9233f1db9ecSBarry Smith PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,"Mat",comm,MatDestroy,MatView); 9248965ea79SLois Curfman McInnes PLogObjectCreate(mat); 9250452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIDense));CHKPTRQ(a); 926549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 927e1311b90SBarry Smith mat->ops->destroy = MatDestroy_MPIDense; 928e1311b90SBarry Smith mat->ops->view = MatView_MPIDense; 9298965ea79SLois Curfman McInnes mat->factor = 0; 93090f02eecSBarry Smith mat->mapping = 0; 9318965ea79SLois Curfman McInnes 932622d7880SLois Curfman McInnes a->factor = 0; 933e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 934d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&a->rank);CHKERRQ(ierr); 935d132466eSBarry Smith ierr = MPI_Comm_size(comm,&a->size);CHKERRQ(ierr); 9368965ea79SLois Curfman McInnes 93796f6c058SBarry Smith ierr = PetscSplitOwnership(comm,&m,&M);CHKERRQ(ierr); 93839ddd567SLois Curfman McInnes 939c7fcc2eaSBarry Smith /* 940c7fcc2eaSBarry Smith The computation of n is wrong below, n should represent the number of local 941c7fcc2eaSBarry Smith rows in the right (column vector) 942c7fcc2eaSBarry Smith */ 943c7fcc2eaSBarry Smith 94439ddd567SLois Curfman McInnes /* each row stores all columns */ 94539ddd567SLois Curfman McInnes if (N == PETSC_DECIDE) N = n; 946d7e8b826SBarry Smith if (n == PETSC_DECIDE) {n = N/a->size + ((N % a->size) > a->rank);} 947a8c6a408SBarry Smith /* if (n != N) SETERRQ(PETSC_ERR_SUP,0,"For now, only n=N is supported"); */ 948aca0ad90SLois Curfman McInnes a->N = mat->N = N; 949aca0ad90SLois Curfman McInnes a->M = mat->M = M; 950aca0ad90SLois Curfman McInnes a->m = mat->m = m; 951aca0ad90SLois Curfman McInnes a->n = mat->n = n; 9528965ea79SLois Curfman McInnes 953c7fcc2eaSBarry Smith /* the information in the maps duplicates the information computed below, eventually 954c7fcc2eaSBarry Smith we should remove the duplicate information that is not contained in the maps */ 955488ecbafSBarry Smith ierr = MapCreateMPI(comm,m,M,&mat->rmap);CHKERRQ(ierr); 95696f6c058SBarry Smith ierr = MapCreateMPI(comm,PETSC_DECIDE,N,&mat->cmap);CHKERRQ(ierr); 957c7fcc2eaSBarry Smith 9588965ea79SLois Curfman McInnes /* build local table of row and column ownerships */ 959d7e8b826SBarry Smith a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int));CHKPTRQ(a->rowners); 960d7e8b826SBarry Smith a->cowners = a->rowners + a->size + 1; 961f09e8eb9SSatish Balay PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 962ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 9638965ea79SLois Curfman McInnes a->rowners[0] = 0; 9648965ea79SLois Curfman McInnes for ( i=2; i<=a->size; i++ ) { 9658965ea79SLois Curfman McInnes a->rowners[i] += a->rowners[i-1]; 9668965ea79SLois Curfman McInnes } 9678965ea79SLois Curfman McInnes a->rstart = a->rowners[a->rank]; 9688965ea79SLois Curfman McInnes a->rend = a->rowners[a->rank+1]; 969ca161407SBarry Smith ierr = MPI_Allgather(&n,1,MPI_INT,a->cowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 970d7e8b826SBarry Smith a->cowners[0] = 0; 971d7e8b826SBarry Smith for ( i=2; i<=a->size; i++ ) { 972d7e8b826SBarry Smith a->cowners[i] += a->cowners[i-1]; 973d7e8b826SBarry Smith } 9748965ea79SLois Curfman McInnes 975029af93fSBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,m,N,data,&a->A);CHKERRQ(ierr); 9768965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 9778965ea79SLois Curfman McInnes 9788965ea79SLois Curfman McInnes /* build cache for off array entries formed */ 9793782ba37SSatish Balay a->donotstash = 0; 9808798bf22SSatish Balay ierr = MatStashCreate_Private(comm,1,&mat->stash);CHKERRQ(ierr); 9818965ea79SLois Curfman McInnes 9828965ea79SLois Curfman McInnes /* stuff used for matrix vector multiply */ 9838965ea79SLois Curfman McInnes a->lvec = 0; 9848965ea79SLois Curfman McInnes a->Mvctx = 0; 98539b7565bSBarry Smith a->roworiented = 1; 9868965ea79SLois Curfman McInnes 987477f1c0bSLois Curfman McInnes *A = mat; 98825cdf11fSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr); 98925cdf11fSBarry Smith if (flg) { 9908c469469SLois Curfman McInnes ierr = MatPrintHelp(mat);CHKERRQ(ierr); 9918c469469SLois Curfman McInnes } 9923a40ed3dSBarry Smith PetscFunctionReturn(0); 9938965ea79SLois Curfman McInnes } 9948965ea79SLois Curfman McInnes 9955615d1e5SSatish Balay #undef __FUNC__ 9965609ef8eSBarry Smith #define __FUNC__ "MatDuplicate_MPIDense" 9975609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 9988965ea79SLois Curfman McInnes { 9998965ea79SLois Curfman McInnes Mat mat; 10003501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense *) A->data; 100139ddd567SLois Curfman McInnes int ierr; 10022ba99913SLois Curfman McInnes FactorCtx *factor; 10038965ea79SLois Curfman McInnes 10043a40ed3dSBarry Smith PetscFunctionBegin; 10058965ea79SLois Curfman McInnes *newmat = 0; 10063f1db9ecSBarry Smith PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,"Mat",A->comm,MatDestroy,MatView); 10078965ea79SLois Curfman McInnes PLogObjectCreate(mat); 10080452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIDense));CHKPTRQ(a); 1009549d3d68SSatish Balay ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1010e1311b90SBarry Smith mat->ops->destroy = MatDestroy_MPIDense; 1011e1311b90SBarry Smith mat->ops->view = MatView_MPIDense; 10123501a2bdSLois Curfman McInnes mat->factor = A->factor; 1013c456f294SBarry Smith mat->assembled = PETSC_TRUE; 10148965ea79SLois Curfman McInnes 101544cd7ae7SLois Curfman McInnes a->m = mat->m = oldmat->m; 101644cd7ae7SLois Curfman McInnes a->n = mat->n = oldmat->n; 101744cd7ae7SLois Curfman McInnes a->M = mat->M = oldmat->M; 101844cd7ae7SLois Curfman McInnes a->N = mat->N = oldmat->N; 10192ba99913SLois Curfman McInnes if (oldmat->factor) { 10202ba99913SLois Curfman McInnes a->factor = (FactorCtx *) (factor = PetscNew(FactorCtx));CHKPTRQ(factor); 10212ba99913SLois Curfman McInnes /* copy factor contents ... add this code! */ 10222ba99913SLois Curfman McInnes } else a->factor = 0; 10238965ea79SLois Curfman McInnes 10248965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 10258965ea79SLois Curfman McInnes a->rend = oldmat->rend; 10268965ea79SLois Curfman McInnes a->size = oldmat->size; 10278965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1028e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 10293782ba37SSatish Balay a->donotstash = oldmat->donotstash; 10300452661fSBarry Smith a->rowners = (int *) PetscMalloc((a->size+1)*sizeof(int));CHKPTRQ(a->rowners); 1031f09e8eb9SSatish Balay PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1032549d3d68SSatish Balay ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr); 10338798bf22SSatish Balay ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr); 10348965ea79SLois Curfman McInnes 10358965ea79SLois Curfman McInnes ierr = VecDuplicate(oldmat->lvec,&a->lvec);CHKERRQ(ierr); 10368965ea79SLois Curfman McInnes PLogObjectParent(mat,a->lvec); 103755659b69SBarry Smith ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx);CHKERRQ(ierr); 10388965ea79SLois Curfman McInnes PLogObjectParent(mat,a->Mvctx); 10395609ef8eSBarry Smith ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 10408965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 10418965ea79SLois Curfman McInnes *newmat = mat; 10423a40ed3dSBarry Smith PetscFunctionReturn(0); 10438965ea79SLois Curfman McInnes } 10448965ea79SLois Curfman McInnes 104577c4ece6SBarry Smith #include "sys.h" 10468965ea79SLois Curfman McInnes 10475615d1e5SSatish Balay #undef __FUNC__ 10485615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense_DenseInFile" 104990ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M, int N, Mat *newmat) 105090ace30eSBarry Smith { 105140011551SBarry Smith int *rowners, i,size,rank,m,ierr,nz,j; 105290ace30eSBarry Smith Scalar *array,*vals,*vals_ptr; 105390ace30eSBarry Smith MPI_Status status; 105490ace30eSBarry Smith 10553a40ed3dSBarry Smith PetscFunctionBegin; 1056d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1057d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 105890ace30eSBarry Smith 105990ace30eSBarry Smith /* determine ownership of all rows */ 106090ace30eSBarry Smith m = M/size + ((M % size) > rank); 106190ace30eSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners); 1062ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 106390ace30eSBarry Smith rowners[0] = 0; 106490ace30eSBarry Smith for ( i=2; i<=size; i++ ) { 106590ace30eSBarry Smith rowners[i] += rowners[i-1]; 106690ace30eSBarry Smith } 106790ace30eSBarry Smith 106890ace30eSBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 106990ace30eSBarry Smith ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr); 107090ace30eSBarry Smith 107190ace30eSBarry Smith if (!rank) { 107290ace30eSBarry Smith vals = (Scalar *) PetscMalloc( m*N*sizeof(Scalar) );CHKPTRQ(vals); 107390ace30eSBarry Smith 107490ace30eSBarry Smith /* read in my part of the matrix numerical values */ 10750752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 107690ace30eSBarry Smith 107790ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 107890ace30eSBarry Smith vals_ptr = vals; 107990ace30eSBarry Smith for ( i=0; i<m; i++ ) { 108090ace30eSBarry Smith for ( j=0; j<N; j++ ) { 108190ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 108290ace30eSBarry Smith } 108390ace30eSBarry Smith } 108490ace30eSBarry Smith 108590ace30eSBarry Smith /* read in other processors and ship out */ 108690ace30eSBarry Smith for ( i=1; i<size; i++ ) { 108790ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 10880752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1089ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 109090ace30eSBarry Smith } 10913a40ed3dSBarry Smith } else { 109290ace30eSBarry Smith /* receive numeric values */ 109390ace30eSBarry Smith vals = (Scalar*) PetscMalloc( m*N*sizeof(Scalar) );CHKPTRQ(vals); 109490ace30eSBarry Smith 109590ace30eSBarry Smith /* receive message of values*/ 1096ca161407SBarry Smith ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 109790ace30eSBarry Smith 109890ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 109990ace30eSBarry Smith vals_ptr = vals; 110090ace30eSBarry Smith for ( i=0; i<m; i++ ) { 110190ace30eSBarry Smith for ( j=0; j<N; j++ ) { 110290ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 110390ace30eSBarry Smith } 110490ace30eSBarry Smith } 110590ace30eSBarry Smith } 1106*606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 1107*606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 11086d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11096d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11103a40ed3dSBarry Smith PetscFunctionReturn(0); 111190ace30eSBarry Smith } 111290ace30eSBarry Smith 111390ace30eSBarry Smith 11145615d1e5SSatish Balay #undef __FUNC__ 11155615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense" 111619bcc07fSBarry Smith int MatLoad_MPIDense(Viewer viewer,MatType type,Mat *newmat) 11178965ea79SLois Curfman McInnes { 11188965ea79SLois Curfman McInnes Mat A; 11198965ea79SLois Curfman McInnes Scalar *vals,*svals; 112019bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 11218965ea79SLois Curfman McInnes MPI_Status status; 11228965ea79SLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 11238965ea79SLois Curfman McInnes int *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols; 112419bcc07fSBarry Smith int tag = ((PetscObject)viewer)->tag; 11253a40ed3dSBarry Smith int i, nz, ierr, j,rstart, rend, fd; 11268965ea79SLois Curfman McInnes 11273a40ed3dSBarry Smith PetscFunctionBegin; 1128d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1129d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 11308965ea79SLois Curfman McInnes if (!rank) { 113190ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 11320752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr); 1133a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object"); 11348965ea79SLois Curfman McInnes } 11358965ea79SLois Curfman McInnes 1136ca161407SBarry Smith ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr); 113790ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 113890ace30eSBarry Smith 113990ace30eSBarry Smith /* 114090ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 114190ace30eSBarry Smith */ 114290ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 11433a40ed3dSBarry Smith ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr); 11443a40ed3dSBarry Smith PetscFunctionReturn(0); 114590ace30eSBarry Smith } 114690ace30eSBarry Smith 11478965ea79SLois Curfman McInnes /* determine ownership of all rows */ 11488965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 11490452661fSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int));CHKPTRQ(rowners); 1150ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 11518965ea79SLois Curfman McInnes rowners[0] = 0; 11528965ea79SLois Curfman McInnes for ( i=2; i<=size; i++ ) { 11538965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 11548965ea79SLois Curfman McInnes } 11558965ea79SLois Curfman McInnes rstart = rowners[rank]; 11568965ea79SLois Curfman McInnes rend = rowners[rank+1]; 11578965ea79SLois Curfman McInnes 11588965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 11590452661fSBarry Smith ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) );CHKPTRQ(ourlens); 11608965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 11618965ea79SLois Curfman McInnes if (!rank) { 11620452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) );CHKPTRQ(rowlengths); 11630752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 11640452661fSBarry Smith sndcounts = (int*) PetscMalloc( size*sizeof(int) );CHKPTRQ(sndcounts); 11658965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i]; 1166ca161407SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 1167*606d414cSSatish Balay ierr = PetscFree(sndcounts);CHKERRQ(ierr); 1168ca161407SBarry Smith } else { 1169ca161407SBarry Smith ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);CHKERRQ(ierr); 11708965ea79SLois Curfman McInnes } 11718965ea79SLois Curfman McInnes 11728965ea79SLois Curfman McInnes if (!rank) { 11738965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 11740452661fSBarry Smith procsnz = (int*) PetscMalloc( size*sizeof(int) );CHKPTRQ(procsnz); 1175549d3d68SSatish Balay ierr = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr); 11768965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 11778965ea79SLois Curfman McInnes for ( j=rowners[i]; j< rowners[i+1]; j++ ) { 11788965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 11798965ea79SLois Curfman McInnes } 11808965ea79SLois Curfman McInnes } 1181*606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 11828965ea79SLois Curfman McInnes 11838965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 11848965ea79SLois Curfman McInnes maxnz = 0; 11858965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 11860452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 11878965ea79SLois Curfman McInnes } 11880452661fSBarry Smith cols = (int *) PetscMalloc( maxnz*sizeof(int) );CHKPTRQ(cols); 11898965ea79SLois Curfman McInnes 11908965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 11918965ea79SLois Curfman McInnes nz = procsnz[0]; 11920452661fSBarry Smith mycols = (int *) PetscMalloc( nz*sizeof(int) );CHKPTRQ(mycols); 11930752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 11948965ea79SLois Curfman McInnes 11958965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 11968965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { 11978965ea79SLois Curfman McInnes nz = procsnz[i]; 11980752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 1199ca161407SBarry Smith ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr); 12008965ea79SLois Curfman McInnes } 1201*606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 12023a40ed3dSBarry Smith } else { 12038965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 12048965ea79SLois Curfman McInnes nz = 0; 12058965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12068965ea79SLois Curfman McInnes nz += ourlens[i]; 12078965ea79SLois Curfman McInnes } 12080452661fSBarry Smith mycols = (int*) PetscMalloc( nz*sizeof(int) );CHKPTRQ(mycols); 12098965ea79SLois Curfman McInnes 12108965ea79SLois Curfman McInnes /* receive message of column indices*/ 1211ca161407SBarry Smith ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr); 1212ca161407SBarry Smith ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr); 1213a8c6a408SBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 12148965ea79SLois Curfman McInnes } 12158965ea79SLois Curfman McInnes 12168965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 1217549d3d68SSatish Balay ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr); 12188965ea79SLois Curfman McInnes jj = 0; 12198965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12208965ea79SLois Curfman McInnes for ( j=0; j<ourlens[i]; j++ ) { 12218965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 12228965ea79SLois Curfman McInnes jj++; 12238965ea79SLois Curfman McInnes } 12248965ea79SLois Curfman McInnes } 12258965ea79SLois Curfman McInnes 12268965ea79SLois Curfman McInnes /* create our matrix */ 12278965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12288965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 12298965ea79SLois Curfman McInnes } 1230b4fd4287SBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 12318965ea79SLois Curfman McInnes A = *newmat; 12328965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12338965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 12348965ea79SLois Curfman McInnes } 12358965ea79SLois Curfman McInnes 12368965ea79SLois Curfman McInnes if (!rank) { 12370452661fSBarry Smith vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) );CHKPTRQ(vals); 12388965ea79SLois Curfman McInnes 12398965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 12408965ea79SLois Curfman McInnes nz = procsnz[0]; 12410752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 12428965ea79SLois Curfman McInnes 12438965ea79SLois Curfman McInnes /* insert into matrix */ 12448965ea79SLois Curfman McInnes jj = rstart; 12458965ea79SLois Curfman McInnes smycols = mycols; 12468965ea79SLois Curfman McInnes svals = vals; 12478965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12488965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 12498965ea79SLois Curfman McInnes smycols += ourlens[i]; 12508965ea79SLois Curfman McInnes svals += ourlens[i]; 12518965ea79SLois Curfman McInnes jj++; 12528965ea79SLois Curfman McInnes } 12538965ea79SLois Curfman McInnes 12548965ea79SLois Curfman McInnes /* read in other processors and ship out */ 12558965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { 12568965ea79SLois Curfman McInnes nz = procsnz[i]; 12570752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1258ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 12598965ea79SLois Curfman McInnes } 1260*606d414cSSatish Balay ierr = PetscFree(procsnz);CHKERRQ(ierr); 12613a40ed3dSBarry Smith } else { 12628965ea79SLois Curfman McInnes /* receive numeric values */ 12630452661fSBarry Smith vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) );CHKPTRQ(vals); 12648965ea79SLois Curfman McInnes 12658965ea79SLois Curfman McInnes /* receive message of values*/ 1266ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1267ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 1268a8c6a408SBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 12698965ea79SLois Curfman McInnes 12708965ea79SLois Curfman McInnes /* insert into matrix */ 12718965ea79SLois Curfman McInnes jj = rstart; 12728965ea79SLois Curfman McInnes smycols = mycols; 12738965ea79SLois Curfman McInnes svals = vals; 12748965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12758965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 12768965ea79SLois Curfman McInnes smycols += ourlens[i]; 12778965ea79SLois Curfman McInnes svals += ourlens[i]; 12788965ea79SLois Curfman McInnes jj++; 12798965ea79SLois Curfman McInnes } 12808965ea79SLois Curfman McInnes } 1281*606d414cSSatish Balay ierr = PetscFree(ourlens);CHKERRQ(ierr); 1282*606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 1283*606d414cSSatish Balay ierr = PetscFree(mycols);CHKERRQ(ierr); 1284*606d414cSSatish Balay ierr = PetscFree(rowners);CHKERRQ(ierr); 12858965ea79SLois Curfman McInnes 12866d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12876d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12883a40ed3dSBarry Smith PetscFunctionReturn(0); 12898965ea79SLois Curfman McInnes } 129090ace30eSBarry Smith 129190ace30eSBarry Smith 129290ace30eSBarry Smith 129390ace30eSBarry Smith 129490ace30eSBarry Smith 1295