1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*db81eaa0SLois Curfman McInnes static char vcid[] = "$Id: mpidense.c,v 1.87 1998/04/15 22:50:39 curfman Exp curfman $"; 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 125615d1e5SSatish Balay #undef __FUNC__ 135615d1e5SSatish Balay #define __FUNC__ "MatSetValues_MPIDense" 148f6be9afSLois Curfman McInnes int MatSetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv) 158965ea79SLois Curfman McInnes { 1639b7565bSBarry Smith Mat_MPIDense *A = (Mat_MPIDense *) mat->data; 1739b7565bSBarry Smith int ierr, i, j, rstart = A->rstart, rend = A->rend, row; 1839b7565bSBarry Smith int roworiented = A->roworiented; 198965ea79SLois Curfman McInnes 203a40ed3dSBarry Smith PetscFunctionBegin; 218965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 22a8c6a408SBarry Smith if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 23a8c6a408SBarry Smith if (idxm[i] >= A->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 248965ea79SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 258965ea79SLois Curfman McInnes row = idxm[i] - rstart; 2639b7565bSBarry Smith if (roworiented) { 2739b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv); CHKERRQ(ierr); 283a40ed3dSBarry Smith } else { 298965ea79SLois Curfman McInnes for ( j=0; j<n; j++ ) { 30a8c6a408SBarry Smith if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 31a8c6a408SBarry Smith if (idxn[j] >= A->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 3239b7565bSBarry Smith ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv); CHKERRQ(ierr); 3339b7565bSBarry Smith } 348965ea79SLois Curfman McInnes } 353a40ed3dSBarry Smith } else { 3639b7565bSBarry Smith if (roworiented) { 3739b7565bSBarry Smith ierr = StashValues_Private(&A->stash,idxm[i],n,idxn,v+i*n,addv); CHKERRQ(ierr); 383a40ed3dSBarry Smith } else { /* must stash each seperately */ 3939b7565bSBarry Smith row = idxm[i]; 4039b7565bSBarry Smith for ( j=0; j<n; j++ ) { 417056b6fcSBarry Smith ierr = StashValues_Private(&A->stash,row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 4239b7565bSBarry Smith } 4339b7565bSBarry Smith } 44b49de8d1SLois Curfman McInnes } 45b49de8d1SLois Curfman McInnes } 463a40ed3dSBarry Smith PetscFunctionReturn(0); 47b49de8d1SLois Curfman McInnes } 48b49de8d1SLois Curfman McInnes 495615d1e5SSatish Balay #undef __FUNC__ 505615d1e5SSatish Balay #define __FUNC__ "MatGetValues_MPIDense" 518f6be9afSLois Curfman McInnes int MatGetValues_MPIDense(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v) 52b49de8d1SLois Curfman McInnes { 53b49de8d1SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 54b49de8d1SLois Curfman McInnes int ierr, i, j, rstart = mdn->rstart, rend = mdn->rend, row; 55b49de8d1SLois Curfman McInnes 563a40ed3dSBarry Smith PetscFunctionBegin; 57b49de8d1SLois Curfman McInnes for ( i=0; i<m; i++ ) { 58a8c6a408SBarry Smith if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 59a8c6a408SBarry Smith if (idxm[i] >= mdn->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 60b49de8d1SLois Curfman McInnes if (idxm[i] >= rstart && idxm[i] < rend) { 61b49de8d1SLois Curfman McInnes row = idxm[i] - rstart; 62b49de8d1SLois Curfman McInnes for ( j=0; j<n; j++ ) { 63a8c6a408SBarry Smith if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 64a8c6a408SBarry Smith if (idxn[j] >= mdn->N) { 65a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 66a8c6a408SBarry Smith } 67b49de8d1SLois Curfman McInnes ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j); CHKERRQ(ierr); 68b49de8d1SLois Curfman McInnes } 69a8c6a408SBarry Smith } else { 70a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Only local values currently supported"); 718965ea79SLois Curfman McInnes } 728965ea79SLois Curfman McInnes } 733a40ed3dSBarry Smith PetscFunctionReturn(0); 748965ea79SLois Curfman McInnes } 758965ea79SLois Curfman McInnes 765615d1e5SSatish Balay #undef __FUNC__ 775615d1e5SSatish Balay #define __FUNC__ "MatGetArray_MPIDense" 788f6be9afSLois Curfman McInnes int MatGetArray_MPIDense(Mat A,Scalar **array) 79ff14e315SSatish Balay { 80ff14e315SSatish Balay Mat_MPIDense *a = (Mat_MPIDense *) A->data; 81ff14e315SSatish Balay int ierr; 82ff14e315SSatish Balay 833a40ed3dSBarry Smith PetscFunctionBegin; 84ff14e315SSatish Balay ierr = MatGetArray(a->A,array); CHKERRQ(ierr); 853a40ed3dSBarry Smith PetscFunctionReturn(0); 86ff14e315SSatish Balay } 87ff14e315SSatish Balay 885615d1e5SSatish Balay #undef __FUNC__ 895615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_MPIDense" 908f6be9afSLois Curfman McInnes int MatRestoreArray_MPIDense(Mat A,Scalar **array) 91ff14e315SSatish Balay { 923a40ed3dSBarry Smith PetscFunctionBegin; 933a40ed3dSBarry Smith PetscFunctionReturn(0); 94ff14e315SSatish Balay } 95ff14e315SSatish Balay 965615d1e5SSatish Balay #undef __FUNC__ 975615d1e5SSatish Balay #define __FUNC__ "MatAssemblyBegin_MPIDense" 988f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 998965ea79SLois Curfman McInnes { 10039ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 1018965ea79SLois Curfman McInnes MPI_Comm comm = mat->comm; 10239ddd567SLois Curfman McInnes int size = mdn->size, *owners = mdn->rowners, rank = mdn->rank; 1038965ea79SLois Curfman McInnes int *nprocs,i,j,idx,*procs,nsends,nreceives,nmax,*work; 10439ddd567SLois Curfman McInnes int tag = mat->tag, *owner,*starts,count,ierr; 1058965ea79SLois Curfman McInnes InsertMode addv; 10639ddd567SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 1078965ea79SLois Curfman McInnes Scalar *rvalues,*svalues; 1088965ea79SLois Curfman McInnes 1093a40ed3dSBarry Smith PetscFunctionBegin; 1108965ea79SLois Curfman McInnes /* make sure all processors are either in INSERTMODE or ADDMODE */ 111ca161407SBarry Smith ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 1127056b6fcSBarry Smith if (addv == (ADD_VALUES|INSERT_VALUES)) { 113a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Cannot mix adds/inserts on different procs"); 1148965ea79SLois Curfman McInnes } 115e0fa3b82SLois Curfman McInnes mat->insertmode = addv; /* in case this processor had no cache */ 1168965ea79SLois Curfman McInnes 1178965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 1180452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 119cddf8d76SBarry Smith PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 1200452661fSBarry Smith owner = (int *) PetscMalloc( (mdn->stash.n+1)*sizeof(int) ); CHKPTRQ(owner); 12139ddd567SLois Curfman McInnes for ( i=0; i<mdn->stash.n; i++ ) { 12239ddd567SLois Curfman McInnes idx = mdn->stash.idx[i]; 1238965ea79SLois Curfman McInnes for ( j=0; j<size; j++ ) { 1248965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 1258965ea79SLois Curfman McInnes nprocs[j]++; procs[j] = 1; owner[i] = j; break; 1268965ea79SLois Curfman McInnes } 1278965ea79SLois Curfman McInnes } 1288965ea79SLois Curfman McInnes } 1298965ea79SLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 1308965ea79SLois Curfman McInnes 1318965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 1320452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 133ca161407SBarry Smith ierr = MPI_Allreduce(procs,work,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 1348965ea79SLois Curfman McInnes nreceives = work[rank]; 135a8c6a408SBarry Smith if (nreceives > size) SETERRQ(PETSC_ERR_PLIB,0,"Internal PETSc error"); 136ca161407SBarry Smith ierr = MPI_Allreduce(nprocs,work,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); 1378965ea79SLois Curfman McInnes nmax = work[rank]; 1380452661fSBarry Smith PetscFree(work); 1398965ea79SLois Curfman McInnes 1408965ea79SLois Curfman McInnes /* post receives: 1418965ea79SLois Curfman McInnes 1) each message will consist of ordered pairs 1428965ea79SLois Curfman McInnes (global index,value) we store the global index as a double 1438965ea79SLois Curfman McInnes to simplify the message passing. 1448965ea79SLois Curfman McInnes 2) since we don't know how long each individual message is we 1458965ea79SLois Curfman McInnes allocate the largest needed buffer for each receive. Potentially 1468965ea79SLois Curfman McInnes this is a lot of wasted space. 1478965ea79SLois Curfman McInnes 1488965ea79SLois Curfman McInnes This could be done better. 1498965ea79SLois Curfman McInnes */ 1503b2fbd54SBarry Smith rvalues = (Scalar *) PetscMalloc(3*(nreceives+1)*(nmax+1)*sizeof(Scalar));CHKPTRQ(rvalues); 1513b2fbd54SBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nreceives+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits); 1528965ea79SLois Curfman McInnes for ( i=0; i<nreceives; i++ ) { 153ca161407SBarry Smith ierr = MPI_Irecv(rvalues+3*nmax*i,3*nmax,MPIU_SCALAR,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 1548965ea79SLois Curfman McInnes } 1558965ea79SLois Curfman McInnes 1568965ea79SLois Curfman McInnes /* do sends: 1578965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 1588965ea79SLois Curfman McInnes the ith processor 1598965ea79SLois Curfman McInnes */ 1603b2fbd54SBarry Smith svalues = (Scalar *) PetscMalloc( 3*(mdn->stash.n+1)*sizeof(Scalar));CHKPTRQ(svalues); 1613b2fbd54SBarry Smith send_waits = (MPI_Request *) PetscMalloc((nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits); 1620452661fSBarry Smith starts = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(starts); 1638965ea79SLois Curfman McInnes starts[0] = 0; 1648965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 16539ddd567SLois Curfman McInnes for ( i=0; i<mdn->stash.n; i++ ) { 16639ddd567SLois Curfman McInnes svalues[3*starts[owner[i]]] = (Scalar) mdn->stash.idx[i]; 16739ddd567SLois Curfman McInnes svalues[3*starts[owner[i]]+1] = (Scalar) mdn->stash.idy[i]; 16839ddd567SLois Curfman McInnes svalues[3*(starts[owner[i]]++)+2] = mdn->stash.array[i]; 1698965ea79SLois Curfman McInnes } 1700452661fSBarry Smith PetscFree(owner); 1718965ea79SLois Curfman McInnes starts[0] = 0; 1728965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 1738965ea79SLois Curfman McInnes count = 0; 1748965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 1758965ea79SLois Curfman McInnes if (procs[i]) { 176ca161407SBarry Smith ierr = MPI_Isend(svalues+3*starts[i],3*nprocs[i],MPIU_SCALAR,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 1778965ea79SLois Curfman McInnes } 1788965ea79SLois Curfman McInnes } 1790452661fSBarry Smith PetscFree(starts); PetscFree(nprocs); 1808965ea79SLois Curfman McInnes 1818965ea79SLois Curfman McInnes /* Free cache space */ 182d2dc9b81SLois Curfman McInnes PLogInfo(mat,"MatAssemblyBegin_MPIDense:Number of off-processor values %d\n",mdn->stash.n); 18339ddd567SLois Curfman McInnes ierr = StashDestroy_Private(&mdn->stash); CHKERRQ(ierr); 1848965ea79SLois Curfman McInnes 18539ddd567SLois Curfman McInnes mdn->svalues = svalues; mdn->rvalues = rvalues; 18639ddd567SLois Curfman McInnes mdn->nsends = nsends; mdn->nrecvs = nreceives; 18739ddd567SLois Curfman McInnes mdn->send_waits = send_waits; mdn->recv_waits = recv_waits; 18839ddd567SLois Curfman McInnes mdn->rmax = nmax; 1898965ea79SLois Curfman McInnes 1903a40ed3dSBarry Smith PetscFunctionReturn(0); 1918965ea79SLois Curfman McInnes } 19239ddd567SLois Curfman McInnes extern int MatSetUpMultiply_MPIDense(Mat); 1938965ea79SLois Curfman McInnes 1945615d1e5SSatish Balay #undef __FUNC__ 1955615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_MPIDense" 1968f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 1978965ea79SLois Curfman McInnes { 19839ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 1998965ea79SLois Curfman McInnes MPI_Status *send_status,recv_status; 20039ddd567SLois Curfman McInnes int imdex, nrecvs=mdn->nrecvs, count=nrecvs, i, n, ierr, row, col; 2018965ea79SLois Curfman McInnes Scalar *values,val; 202e0fa3b82SLois Curfman McInnes InsertMode addv = mat->insertmode; 2038965ea79SLois Curfman McInnes 2043a40ed3dSBarry Smith PetscFunctionBegin; 2058965ea79SLois Curfman McInnes /* wait on receives */ 2068965ea79SLois Curfman McInnes while (count) { 207ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,mdn->recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 2088965ea79SLois Curfman McInnes /* unpack receives into our local space */ 20939ddd567SLois Curfman McInnes values = mdn->rvalues + 3*imdex*mdn->rmax; 210ca161407SBarry Smith ierr = MPI_Get_count(&recv_status,MPIU_SCALAR,&n);CHKERRQ(ierr); 2118965ea79SLois Curfman McInnes n = n/3; 2128965ea79SLois Curfman McInnes for ( i=0; i<n; i++ ) { 213227d817aSBarry Smith row = (int) PetscReal(values[3*i]) - mdn->rstart; 214227d817aSBarry Smith col = (int) PetscReal(values[3*i+1]); 2158965ea79SLois Curfman McInnes val = values[3*i+2]; 21639ddd567SLois Curfman McInnes if (col >= 0 && col < mdn->N) { 21739ddd567SLois Curfman McInnes MatSetValues(mdn->A,1,&row,1,&col,&val,addv); 2188965ea79SLois Curfman McInnes } 219a8c6a408SBarry Smith else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid column");} 2208965ea79SLois Curfman McInnes } 2218965ea79SLois Curfman McInnes count--; 2228965ea79SLois Curfman McInnes } 2230452661fSBarry Smith PetscFree(mdn->recv_waits); PetscFree(mdn->rvalues); 2248965ea79SLois Curfman McInnes 2258965ea79SLois Curfman McInnes /* wait on sends */ 22639ddd567SLois Curfman McInnes if (mdn->nsends) { 2277056b6fcSBarry Smith send_status = (MPI_Status *) PetscMalloc(mdn->nsends*sizeof(MPI_Status));CHKPTRQ(send_status); 228ca161407SBarry Smith ierr = MPI_Waitall(mdn->nsends,mdn->send_waits,send_status);CHKERRQ(ierr); 2290452661fSBarry Smith PetscFree(send_status); 2308965ea79SLois Curfman McInnes } 2310452661fSBarry Smith PetscFree(mdn->send_waits); PetscFree(mdn->svalues); 2328965ea79SLois Curfman McInnes 23339ddd567SLois Curfman McInnes ierr = MatAssemblyBegin(mdn->A,mode); CHKERRQ(ierr); 23439ddd567SLois Curfman McInnes ierr = MatAssemblyEnd(mdn->A,mode); CHKERRQ(ierr); 2358965ea79SLois Curfman McInnes 2366d4a8577SBarry Smith if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 23739ddd567SLois Curfman McInnes ierr = MatSetUpMultiply_MPIDense(mat); CHKERRQ(ierr); 2388965ea79SLois Curfman McInnes } 2393a40ed3dSBarry Smith PetscFunctionReturn(0); 2408965ea79SLois Curfman McInnes } 2418965ea79SLois Curfman McInnes 2425615d1e5SSatish Balay #undef __FUNC__ 2435615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_MPIDense" 2448f6be9afSLois Curfman McInnes int MatZeroEntries_MPIDense(Mat A) 2458965ea79SLois Curfman McInnes { 2463a40ed3dSBarry Smith int ierr; 24739ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense *) A->data; 2483a40ed3dSBarry Smith 2493a40ed3dSBarry Smith PetscFunctionBegin; 2503a40ed3dSBarry Smith ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 2513a40ed3dSBarry Smith PetscFunctionReturn(0); 2528965ea79SLois Curfman McInnes } 2538965ea79SLois Curfman McInnes 2545615d1e5SSatish Balay #undef __FUNC__ 2555615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_MPIDense" 2568f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIDense(Mat A,int *bs) 2574e220ebcSLois Curfman McInnes { 2583a40ed3dSBarry Smith PetscFunctionBegin; 2594e220ebcSLois Curfman McInnes *bs = 1; 2603a40ed3dSBarry Smith PetscFunctionReturn(0); 2614e220ebcSLois Curfman McInnes } 2624e220ebcSLois Curfman McInnes 2638965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the 2648965ea79SLois Curfman McInnes matrix is square and the column and row owerships are identical. 2658965ea79SLois Curfman McInnes This is a BUG. The only way to fix it seems to be to access 2663501a2bdSLois Curfman McInnes mdn->A and mdn->B directly and not through the MatZeroRows() 2678965ea79SLois Curfman McInnes routine. 2688965ea79SLois Curfman McInnes */ 2695615d1e5SSatish Balay #undef __FUNC__ 2705615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_MPIDense" 2718f6be9afSLois Curfman McInnes int MatZeroRows_MPIDense(Mat A,IS is,Scalar *diag) 2728965ea79SLois Curfman McInnes { 27339ddd567SLois Curfman McInnes Mat_MPIDense *l = (Mat_MPIDense *) A->data; 2748965ea79SLois Curfman McInnes int i,ierr,N, *rows,*owners = l->rowners,size = l->size; 2758965ea79SLois Curfman McInnes int *procs,*nprocs,j,found,idx,nsends,*work; 2768965ea79SLois Curfman McInnes int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 2778965ea79SLois Curfman McInnes int *rvalues,tag = A->tag,count,base,slen,n,*source; 2788965ea79SLois Curfman McInnes int *lens,imdex,*lrows,*values; 2798965ea79SLois Curfman McInnes MPI_Comm comm = A->comm; 2808965ea79SLois Curfman McInnes MPI_Request *send_waits,*recv_waits; 2818965ea79SLois Curfman McInnes MPI_Status recv_status,*send_status; 2828965ea79SLois Curfman McInnes IS istmp; 2838965ea79SLois Curfman McInnes 2843a40ed3dSBarry Smith PetscFunctionBegin; 28577c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 2868965ea79SLois Curfman McInnes ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 2878965ea79SLois Curfman McInnes 2888965ea79SLois Curfman McInnes /* first count number of contributors to each processor */ 2890452661fSBarry Smith nprocs = (int *) PetscMalloc( 2*size*sizeof(int) ); CHKPTRQ(nprocs); 290cddf8d76SBarry Smith PetscMemzero(nprocs,2*size*sizeof(int)); procs = nprocs + size; 2910452661fSBarry Smith owner = (int *) PetscMalloc((N+1)*sizeof(int)); CHKPTRQ(owner); /* see note*/ 2928965ea79SLois Curfman McInnes for ( i=0; i<N; i++ ) { 2938965ea79SLois Curfman McInnes idx = rows[i]; 2948965ea79SLois Curfman McInnes found = 0; 2958965ea79SLois Curfman McInnes for ( j=0; j<size; j++ ) { 2968965ea79SLois Curfman McInnes if (idx >= owners[j] && idx < owners[j+1]) { 2978965ea79SLois Curfman McInnes nprocs[j]++; procs[j] = 1; owner[i] = j; found = 1; break; 2988965ea79SLois Curfman McInnes } 2998965ea79SLois Curfman McInnes } 300a8c6a408SBarry Smith if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Index out of range"); 3018965ea79SLois Curfman McInnes } 3028965ea79SLois Curfman McInnes nsends = 0; for ( i=0; i<size; i++ ) { nsends += procs[i];} 3038965ea79SLois Curfman McInnes 3048965ea79SLois Curfman McInnes /* inform other processors of number of messages and max length*/ 3050452661fSBarry Smith work = (int *) PetscMalloc( size*sizeof(int) ); CHKPTRQ(work); 306ca161407SBarry Smith ierr = MPI_Allreduce( procs, work,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 3078965ea79SLois Curfman McInnes nrecvs = work[rank]; 308ca161407SBarry Smith ierr = MPI_Allreduce( nprocs, work,size,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); 3098965ea79SLois Curfman McInnes nmax = work[rank]; 3100452661fSBarry Smith PetscFree(work); 3118965ea79SLois Curfman McInnes 3128965ea79SLois Curfman McInnes /* post receives: */ 3133a40ed3dSBarry Smith rvalues = (int *) PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int));CHKPTRQ(rvalues); 3143a40ed3dSBarry Smith recv_waits = (MPI_Request *) PetscMalloc((nrecvs+1)*sizeof(MPI_Request));CHKPTRQ(recv_waits); 3158965ea79SLois Curfman McInnes for ( i=0; i<nrecvs; i++ ) { 316ca161407SBarry Smith ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 3178965ea79SLois Curfman McInnes } 3188965ea79SLois Curfman McInnes 3198965ea79SLois Curfman McInnes /* do sends: 3208965ea79SLois Curfman McInnes 1) starts[i] gives the starting index in svalues for stuff going to 3218965ea79SLois Curfman McInnes the ith processor 3228965ea79SLois Curfman McInnes */ 3230452661fSBarry Smith svalues = (int *) PetscMalloc( (N+1)*sizeof(int) ); CHKPTRQ(svalues); 3247056b6fcSBarry Smith send_waits = (MPI_Request *) PetscMalloc((nsends+1)*sizeof(MPI_Request));CHKPTRQ(send_waits); 3250452661fSBarry Smith starts = (int *) PetscMalloc( (size+1)*sizeof(int) ); CHKPTRQ(starts); 3268965ea79SLois Curfman McInnes starts[0] = 0; 3278965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3288965ea79SLois Curfman McInnes for ( i=0; i<N; i++ ) { 3298965ea79SLois Curfman McInnes svalues[starts[owner[i]]++] = rows[i]; 3308965ea79SLois Curfman McInnes } 3318965ea79SLois Curfman McInnes ISRestoreIndices(is,&rows); 3328965ea79SLois Curfman McInnes 3338965ea79SLois Curfman McInnes starts[0] = 0; 3348965ea79SLois Curfman McInnes for ( i=1; i<size+1; i++ ) { starts[i] = starts[i-1] + nprocs[i-1];} 3358965ea79SLois Curfman McInnes count = 0; 3368965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 3378965ea79SLois Curfman McInnes if (procs[i]) { 338ca161407SBarry Smith ierr = MPI_Isend(svalues+starts[i],nprocs[i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 3398965ea79SLois Curfman McInnes } 3408965ea79SLois Curfman McInnes } 3410452661fSBarry Smith PetscFree(starts); 3428965ea79SLois Curfman McInnes 3438965ea79SLois Curfman McInnes base = owners[rank]; 3448965ea79SLois Curfman McInnes 3458965ea79SLois Curfman McInnes /* wait on receives */ 3460452661fSBarry Smith lens = (int *) PetscMalloc( 2*(nrecvs+1)*sizeof(int) ); CHKPTRQ(lens); 3478965ea79SLois Curfman McInnes source = lens + nrecvs; 3488965ea79SLois Curfman McInnes count = nrecvs; slen = 0; 3498965ea79SLois Curfman McInnes while (count) { 350ca161407SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 3518965ea79SLois Curfman McInnes /* unpack receives into our local space */ 352ca161407SBarry Smith ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr); 3538965ea79SLois Curfman McInnes source[imdex] = recv_status.MPI_SOURCE; 3548965ea79SLois Curfman McInnes lens[imdex] = n; 3558965ea79SLois Curfman McInnes slen += n; 3568965ea79SLois Curfman McInnes count--; 3578965ea79SLois Curfman McInnes } 3580452661fSBarry Smith PetscFree(recv_waits); 3598965ea79SLois Curfman McInnes 3608965ea79SLois Curfman McInnes /* move the data into the send scatter */ 3610452661fSBarry Smith lrows = (int *) PetscMalloc( (slen+1)*sizeof(int) ); CHKPTRQ(lrows); 3628965ea79SLois Curfman McInnes count = 0; 3638965ea79SLois Curfman McInnes for ( i=0; i<nrecvs; i++ ) { 3648965ea79SLois Curfman McInnes values = rvalues + i*nmax; 3658965ea79SLois Curfman McInnes for ( j=0; j<lens[i]; j++ ) { 3668965ea79SLois Curfman McInnes lrows[count++] = values[j] - base; 3678965ea79SLois Curfman McInnes } 3688965ea79SLois Curfman McInnes } 3690452661fSBarry Smith PetscFree(rvalues); PetscFree(lens); 3700452661fSBarry Smith PetscFree(owner); PetscFree(nprocs); 3718965ea79SLois Curfman McInnes 3728965ea79SLois Curfman McInnes /* actually zap the local rows */ 373029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 3748965ea79SLois Curfman McInnes PLogObjectParent(A,istmp); 3750452661fSBarry Smith PetscFree(lrows); 3768965ea79SLois Curfman McInnes ierr = MatZeroRows(l->A,istmp,diag); CHKERRQ(ierr); 3778965ea79SLois Curfman McInnes ierr = ISDestroy(istmp); CHKERRQ(ierr); 3788965ea79SLois Curfman McInnes 3798965ea79SLois Curfman McInnes /* wait on sends */ 3808965ea79SLois Curfman McInnes if (nsends) { 3817056b6fcSBarry Smith send_status = (MPI_Status *) PetscMalloc(nsends*sizeof(MPI_Status));CHKPTRQ(send_status); 382ca161407SBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 3830452661fSBarry Smith PetscFree(send_status); 3848965ea79SLois Curfman McInnes } 3850452661fSBarry Smith PetscFree(send_waits); PetscFree(svalues); 3868965ea79SLois Curfman McInnes 3873a40ed3dSBarry Smith PetscFunctionReturn(0); 3888965ea79SLois Curfman McInnes } 3898965ea79SLois Curfman McInnes 3905615d1e5SSatish Balay #undef __FUNC__ 3915615d1e5SSatish Balay #define __FUNC__ "MatMult_MPIDense" 3928f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 3938965ea79SLois Curfman McInnes { 39439ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 3958965ea79SLois Curfman McInnes int ierr; 396c456f294SBarry Smith 3973a40ed3dSBarry Smith PetscFunctionBegin; 39843a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 39943a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 40044cd7ae7SLois Curfman McInnes ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy); CHKERRQ(ierr); 4013a40ed3dSBarry Smith PetscFunctionReturn(0); 4028965ea79SLois Curfman McInnes } 4038965ea79SLois Curfman McInnes 4045615d1e5SSatish Balay #undef __FUNC__ 4055615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_MPIDense" 4068f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 4078965ea79SLois Curfman McInnes { 40839ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 4098965ea79SLois Curfman McInnes int ierr; 410c456f294SBarry Smith 4113a40ed3dSBarry Smith PetscFunctionBegin; 41243a90d84SBarry Smith ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 41343a90d84SBarry Smith ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 41444cd7ae7SLois Curfman McInnes ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz); CHKERRQ(ierr); 4153a40ed3dSBarry Smith PetscFunctionReturn(0); 4168965ea79SLois Curfman McInnes } 4178965ea79SLois Curfman McInnes 4185615d1e5SSatish Balay #undef __FUNC__ 4195615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_MPIDense" 4208f6be9afSLois Curfman McInnes int MatMultTrans_MPIDense(Mat A,Vec xx,Vec yy) 421096963f5SLois Curfman McInnes { 422096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 423096963f5SLois Curfman McInnes int ierr; 4243501a2bdSLois Curfman McInnes Scalar zero = 0.0; 425096963f5SLois Curfman McInnes 4263a40ed3dSBarry Smith PetscFunctionBegin; 4273501a2bdSLois Curfman McInnes ierr = VecSet(&zero,yy); CHKERRQ(ierr); 42844cd7ae7SLois Curfman McInnes ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec); CHKERRQ(ierr); 429537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 430537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 4313a40ed3dSBarry Smith PetscFunctionReturn(0); 432096963f5SLois Curfman McInnes } 433096963f5SLois Curfman McInnes 4345615d1e5SSatish Balay #undef __FUNC__ 4355615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_MPIDense" 4368f6be9afSLois Curfman McInnes int MatMultTransAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 437096963f5SLois Curfman McInnes { 438096963f5SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 439096963f5SLois Curfman McInnes int ierr; 440096963f5SLois Curfman McInnes 4413a40ed3dSBarry Smith PetscFunctionBegin; 4423501a2bdSLois Curfman McInnes ierr = VecCopy(yy,zz); CHKERRQ(ierr); 44344cd7ae7SLois Curfman McInnes ierr = MatMultTrans_SeqDense(a->A,xx,a->lvec); CHKERRQ(ierr); 444537820f0SBarry Smith ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 445537820f0SBarry Smith ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx); CHKERRQ(ierr); 4463a40ed3dSBarry Smith PetscFunctionReturn(0); 447096963f5SLois Curfman McInnes } 448096963f5SLois Curfman McInnes 4495615d1e5SSatish Balay #undef __FUNC__ 4505615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_MPIDense" 4518f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v) 4528965ea79SLois Curfman McInnes { 45339ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 454096963f5SLois Curfman McInnes Mat_SeqDense *aloc = (Mat_SeqDense *) a->A->data; 45544cd7ae7SLois Curfman McInnes int ierr, len, i, n, m = a->m, radd; 45644cd7ae7SLois Curfman McInnes Scalar *x, zero = 0.0; 457ed3cc1f0SBarry Smith 4583a40ed3dSBarry Smith PetscFunctionBegin; 45944cd7ae7SLois Curfman McInnes VecSet(&zero,v); 460096963f5SLois Curfman McInnes ierr = VecGetArray(v,&x); CHKERRQ(ierr); 461096963f5SLois Curfman McInnes ierr = VecGetSize(v,&n); CHKERRQ(ierr); 462a8c6a408SBarry Smith if (n != a->M) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec"); 46344cd7ae7SLois Curfman McInnes len = PetscMin(aloc->m,aloc->n); 4647ddc982cSLois Curfman McInnes radd = a->rstart*m; 46544cd7ae7SLois Curfman McInnes for ( i=0; i<len; i++ ) { 466096963f5SLois Curfman McInnes x[i] = aloc->v[radd + i*m + i]; 467096963f5SLois Curfman McInnes } 4683a40ed3dSBarry Smith PetscFunctionReturn(0); 4698965ea79SLois Curfman McInnes } 4708965ea79SLois Curfman McInnes 4715615d1e5SSatish Balay #undef __FUNC__ 4725615d1e5SSatish Balay #define __FUNC__ "MatDestroy_MPIDense" 473e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat) 4748965ea79SLois Curfman McInnes { 4753501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 4768965ea79SLois Curfman McInnes int ierr; 477ed3cc1f0SBarry Smith 4783a40ed3dSBarry Smith PetscFunctionBegin; 4793a40ed3dSBarry Smith #if defined(USE_PETSC_LOG) 480e1311b90SBarry Smith PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mdn->M,mdn->N); 4818965ea79SLois Curfman McInnes #endif 4820452661fSBarry Smith PetscFree(mdn->rowners); 4833501a2bdSLois Curfman McInnes ierr = MatDestroy(mdn->A); CHKERRQ(ierr); 4843501a2bdSLois Curfman McInnes if (mdn->lvec) VecDestroy(mdn->lvec); 4853501a2bdSLois Curfman McInnes if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 486622d7880SLois Curfman McInnes if (mdn->factor) { 487622d7880SLois Curfman McInnes if (mdn->factor->temp) PetscFree(mdn->factor->temp); 488622d7880SLois Curfman McInnes if (mdn->factor->tag) PetscFree(mdn->factor->tag); 489622d7880SLois Curfman McInnes if (mdn->factor->pivots) PetscFree(mdn->factor->pivots); 490622d7880SLois Curfman McInnes PetscFree(mdn->factor); 491622d7880SLois Curfman McInnes } 4920452661fSBarry Smith PetscFree(mdn); 49390f02eecSBarry Smith if (mat->mapping) { 49490f02eecSBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr); 49590f02eecSBarry Smith } 4968965ea79SLois Curfman McInnes PLogObjectDestroy(mat); 4970452661fSBarry Smith PetscHeaderDestroy(mat); 4983a40ed3dSBarry Smith PetscFunctionReturn(0); 4998965ea79SLois Curfman McInnes } 50039ddd567SLois Curfman McInnes 5018965ea79SLois Curfman McInnes #include "pinclude/pviewer.h" 5028965ea79SLois Curfman McInnes 5035615d1e5SSatish Balay #undef __FUNC__ 5045615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense_Binary" 50539ddd567SLois Curfman McInnes static int MatView_MPIDense_Binary(Mat mat,Viewer viewer) 5068965ea79SLois Curfman McInnes { 50739ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 5088965ea79SLois Curfman McInnes int ierr; 5097056b6fcSBarry Smith 5103a40ed3dSBarry Smith PetscFunctionBegin; 51139ddd567SLois Curfman McInnes if (mdn->size == 1) { 51239ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer); CHKERRQ(ierr); 5138965ea79SLois Curfman McInnes } 514a8c6a408SBarry Smith else SETERRQ(PETSC_ERR_SUP,0,"Only uniprocessor output supported"); 5153a40ed3dSBarry Smith PetscFunctionReturn(0); 5168965ea79SLois Curfman McInnes } 5178965ea79SLois Curfman McInnes 5185615d1e5SSatish Balay #undef __FUNC__ 5195615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense_ASCII" 52039ddd567SLois Curfman McInnes static int MatView_MPIDense_ASCII(Mat mat,Viewer viewer) 5218965ea79SLois Curfman McInnes { 52239ddd567SLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data; 52339ddd567SLois Curfman McInnes int ierr, format; 5248965ea79SLois Curfman McInnes FILE *fd; 52519bcc07fSBarry Smith ViewerType vtype; 5268965ea79SLois Curfman McInnes 5273a40ed3dSBarry Smith PetscFunctionBegin; 5283a40ed3dSBarry Smith ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 52990ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 53090ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 531639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 5324e220ebcSLois Curfman McInnes int rank; 5334e220ebcSLois Curfman McInnes MatInfo info; 534096963f5SLois Curfman McInnes MPI_Comm_rank(mat->comm,&rank); 5354e220ebcSLois Curfman McInnes ierr = MatGetInfo(mat,MAT_LOCAL,&info); 53677c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 5374e220ebcSLois Curfman McInnes fprintf(fd," [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mdn->m, 5384e220ebcSLois Curfman McInnes (int)info.nz_used,(int)info.nz_allocated,(int)info.memory); 539096963f5SLois Curfman McInnes fflush(fd); 54077c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 5413501a2bdSLois Curfman McInnes ierr = VecScatterView(mdn->Mvctx,viewer); CHKERRQ(ierr); 5423a40ed3dSBarry Smith PetscFunctionReturn(0); 5433501a2bdSLois Curfman McInnes } 54402cad45dSBarry Smith else if (format == VIEWER_FORMAT_ASCII_INFO) { 5453a40ed3dSBarry Smith PetscFunctionReturn(0); 5468965ea79SLois Curfman McInnes } 54719bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER) { 54877c4ece6SBarry Smith PetscSequentialPhaseBegin(mat->comm,1); 54939ddd567SLois Curfman McInnes fprintf(fd,"[%d] rows %d starts %d ends %d cols %d\n", 55039ddd567SLois Curfman McInnes mdn->rank,mdn->m,mdn->rstart,mdn->rend,mdn->n); 55139ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer); CHKERRQ(ierr); 5528965ea79SLois Curfman McInnes fflush(fd); 55377c4ece6SBarry Smith PetscSequentialPhaseEnd(mat->comm,1); 5543a40ed3dSBarry Smith } else { 55539ddd567SLois Curfman McInnes int size = mdn->size, rank = mdn->rank; 5568965ea79SLois Curfman McInnes if (size == 1) { 55739ddd567SLois Curfman McInnes ierr = MatView(mdn->A,viewer); CHKERRQ(ierr); 5583a40ed3dSBarry Smith } else { 5598965ea79SLois Curfman McInnes /* assemble the entire matrix onto first processor. */ 5608965ea79SLois Curfman McInnes Mat A; 56139ddd567SLois Curfman McInnes int M = mdn->M, N = mdn->N,m,row,i, nz, *cols; 56239ddd567SLois Curfman McInnes Scalar *vals; 56339ddd567SLois Curfman McInnes Mat_SeqDense *Amdn = (Mat_SeqDense*) mdn->A->data; 5648965ea79SLois Curfman McInnes 5658965ea79SLois Curfman McInnes if (!rank) { 5660513a670SBarry Smith ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A); CHKERRQ(ierr); 5673a40ed3dSBarry Smith } else { 5680513a670SBarry Smith ierr = MatCreateMPIDense(mat->comm,0,N,M,N,PETSC_NULL,&A); CHKERRQ(ierr); 5698965ea79SLois Curfman McInnes } 5708965ea79SLois Curfman McInnes PLogObjectParent(mat,A); 5718965ea79SLois Curfman McInnes 57239ddd567SLois Curfman McInnes /* Copy the matrix ... This isn't the most efficient means, 57339ddd567SLois Curfman McInnes but it's quick for now */ 57439ddd567SLois Curfman McInnes row = mdn->rstart; m = Amdn->m; 5758965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 57639ddd567SLois Curfman McInnes ierr = MatGetRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr); 57739ddd567SLois Curfman McInnes ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES); CHKERRQ(ierr); 57839ddd567SLois Curfman McInnes ierr = MatRestoreRow(mat,row,&nz,&cols,&vals); CHKERRQ(ierr); 57939ddd567SLois Curfman McInnes row++; 5808965ea79SLois Curfman McInnes } 5818965ea79SLois Curfman McInnes 5826d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 5836d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 5848965ea79SLois Curfman McInnes if (!rank) { 58539ddd567SLois Curfman McInnes ierr = MatView(((Mat_MPIDense*)(A->data))->A,viewer); CHKERRQ(ierr); 5868965ea79SLois Curfman McInnes } 5878965ea79SLois Curfman McInnes ierr = MatDestroy(A); CHKERRQ(ierr); 5888965ea79SLois Curfman McInnes } 5898965ea79SLois Curfman McInnes } 5903a40ed3dSBarry Smith PetscFunctionReturn(0); 5918965ea79SLois Curfman McInnes } 5928965ea79SLois Curfman McInnes 5935615d1e5SSatish Balay #undef __FUNC__ 5945615d1e5SSatish Balay #define __FUNC__ "MatView_MPIDense" 595e1311b90SBarry Smith int MatView_MPIDense(Mat mat,Viewer viewer) 5968965ea79SLois Curfman McInnes { 59739ddd567SLois Curfman McInnes int ierr; 598bcd2baecSBarry Smith ViewerType vtype; 5998965ea79SLois Curfman McInnes 600bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 601bcd2baecSBarry Smith if (vtype == ASCII_FILE_VIEWER) { 60239ddd567SLois Curfman McInnes ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr); 6033a40ed3dSBarry Smith } else if (vtype == ASCII_FILES_VIEWER) { 60439ddd567SLois Curfman McInnes ierr = MatView_MPIDense_ASCII(mat,viewer); CHKERRQ(ierr); 6053a40ed3dSBarry Smith } else if (vtype == BINARY_FILE_VIEWER) { 6063a40ed3dSBarry Smith ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 6075cd90555SBarry Smith } else { 6085cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 6098965ea79SLois Curfman McInnes } 6103a40ed3dSBarry Smith PetscFunctionReturn(0); 6118965ea79SLois Curfman McInnes } 6128965ea79SLois Curfman McInnes 6135615d1e5SSatish Balay #undef __FUNC__ 6145615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_MPIDense" 6158f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 6168965ea79SLois Curfman McInnes { 6173501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 6183501a2bdSLois Curfman McInnes Mat mdn = mat->A; 6194e220ebcSLois Curfman McInnes int ierr; 6204e220ebcSLois Curfman McInnes double isend[5], irecv[5]; 6218965ea79SLois Curfman McInnes 6223a40ed3dSBarry Smith PetscFunctionBegin; 6234e220ebcSLois Curfman McInnes info->rows_global = (double)mat->M; 6244e220ebcSLois Curfman McInnes info->columns_global = (double)mat->N; 6254e220ebcSLois Curfman McInnes info->rows_local = (double)mat->m; 6264e220ebcSLois Curfman McInnes info->columns_local = (double)mat->N; 6274e220ebcSLois Curfman McInnes info->block_size = 1.0; 6284e220ebcSLois Curfman McInnes ierr = MatGetInfo(mdn,MAT_LOCAL,info); CHKERRQ(ierr); 6294e220ebcSLois Curfman McInnes isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 6304e220ebcSLois Curfman McInnes isend[3] = info->memory; isend[4] = info->mallocs; 6318965ea79SLois Curfman McInnes if (flag == MAT_LOCAL) { 6324e220ebcSLois Curfman McInnes info->nz_used = isend[0]; 6334e220ebcSLois Curfman McInnes info->nz_allocated = isend[1]; 6344e220ebcSLois Curfman McInnes info->nz_unneeded = isend[2]; 6354e220ebcSLois Curfman McInnes info->memory = isend[3]; 6364e220ebcSLois Curfman McInnes info->mallocs = isend[4]; 6378965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_MAX) { 638ca161407SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_INT,MPI_MAX,A->comm);CHKERRQ(ierr); 6394e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6404e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6414e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6424e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6434e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6448965ea79SLois Curfman McInnes } else if (flag == MAT_GLOBAL_SUM) { 645ca161407SBarry Smith ierr = MPI_Allreduce(isend,irecv,5,MPI_INT,MPI_SUM,A->comm);CHKERRQ(ierr); 6464e220ebcSLois Curfman McInnes info->nz_used = irecv[0]; 6474e220ebcSLois Curfman McInnes info->nz_allocated = irecv[1]; 6484e220ebcSLois Curfman McInnes info->nz_unneeded = irecv[2]; 6494e220ebcSLois Curfman McInnes info->memory = irecv[3]; 6504e220ebcSLois Curfman McInnes info->mallocs = irecv[4]; 6518965ea79SLois Curfman McInnes } 6524e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 6534e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 6544e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 6553a40ed3dSBarry Smith PetscFunctionReturn(0); 6568965ea79SLois Curfman McInnes } 6578965ea79SLois Curfman McInnes 6588c469469SLois Curfman McInnes /* extern int MatLUFactorSymbolic_MPIDense(Mat,IS,IS,double,Mat*); 6598aaee692SLois Curfman McInnes extern int MatLUFactorNumeric_MPIDense(Mat,Mat*); 6608aaee692SLois Curfman McInnes extern int MatLUFactor_MPIDense(Mat,IS,IS,double); 6618aaee692SLois Curfman McInnes extern int MatSolve_MPIDense(Mat,Vec,Vec); 6628c469469SLois Curfman McInnes extern int MatSolveAdd_MPIDense(Mat,Vec,Vec,Vec); 6638aaee692SLois Curfman McInnes extern int MatSolveTrans_MPIDense(Mat,Vec,Vec); 6648aaee692SLois Curfman McInnes extern int MatSolveTransAdd_MPIDense(Mat,Vec,Vec,Vec); */ 6658aaee692SLois Curfman McInnes 6665615d1e5SSatish Balay #undef __FUNC__ 6675615d1e5SSatish Balay #define __FUNC__ "MatSetOption_MPIDense" 6688f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op) 6698965ea79SLois Curfman McInnes { 67039ddd567SLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 6718965ea79SLois Curfman McInnes 6723a40ed3dSBarry Smith PetscFunctionBegin; 6736d4a8577SBarry Smith if (op == MAT_NO_NEW_NONZERO_LOCATIONS || 6746d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 675c2653b3dSLois Curfman McInnes op == MAT_NEW_NONZERO_LOCATION_ERROR || 67696854ed6SLois Curfman McInnes op == MAT_NEW_NONZERO_ALLOCATION_ERROR || 677219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_SORTED || 678219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED) { 679b1fbbac0SLois Curfman McInnes MatSetOption(a->A,op); 680b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROW_ORIENTED) { 681aeafbbfcSLois Curfman McInnes a->roworiented = 1; 6828965ea79SLois Curfman McInnes MatSetOption(a->A,op); 683b1fbbac0SLois Curfman McInnes } else if (op == MAT_ROWS_SORTED || 684219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 6856d4a8577SBarry Smith op == MAT_SYMMETRIC || 6866d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 687b51ba29fSSatish Balay op == MAT_YES_NEW_DIAGONALS || 688b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 689981c4779SBarry Smith PLogInfo(A,"MatSetOption_MPIDense:Option ignored\n"); 6903a40ed3dSBarry Smith } else if (op == MAT_COLUMN_ORIENTED) { 6913a40ed3dSBarry Smith a->roworiented = 0; MatSetOption(a->A,op); 6923a40ed3dSBarry Smith } else if (op == MAT_NO_NEW_DIAGONALS) { 6933a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 6943a40ed3dSBarry Smith } else { 6953a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 6963a40ed3dSBarry Smith } 6973a40ed3dSBarry Smith PetscFunctionReturn(0); 6988965ea79SLois Curfman McInnes } 6998965ea79SLois Curfman McInnes 7005615d1e5SSatish Balay #undef __FUNC__ 7015615d1e5SSatish Balay #define __FUNC__ "MatGetSize_MPIDense" 7028f6be9afSLois Curfman McInnes int MatGetSize_MPIDense(Mat A,int *m,int *n) 7038965ea79SLois Curfman McInnes { 7043501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 7053a40ed3dSBarry Smith 7063a40ed3dSBarry Smith PetscFunctionBegin; 7078965ea79SLois Curfman McInnes *m = mat->M; *n = mat->N; 7083a40ed3dSBarry Smith PetscFunctionReturn(0); 7098965ea79SLois Curfman McInnes } 7108965ea79SLois Curfman McInnes 7115615d1e5SSatish Balay #undef __FUNC__ 7125615d1e5SSatish Balay #define __FUNC__ "MatGetLocalSize_MPIDense" 7138f6be9afSLois Curfman McInnes int MatGetLocalSize_MPIDense(Mat A,int *m,int *n) 7148965ea79SLois Curfman McInnes { 7153501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 7163a40ed3dSBarry Smith 7173a40ed3dSBarry Smith PetscFunctionBegin; 7188965ea79SLois Curfman McInnes *m = mat->m; *n = mat->N; 7193a40ed3dSBarry Smith PetscFunctionReturn(0); 7208965ea79SLois Curfman McInnes } 7218965ea79SLois Curfman McInnes 7225615d1e5SSatish Balay #undef __FUNC__ 7235615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_MPIDense" 7248f6be9afSLois Curfman McInnes int MatGetOwnershipRange_MPIDense(Mat A,int *m,int *n) 7258965ea79SLois Curfman McInnes { 7263501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 7273a40ed3dSBarry Smith 7283a40ed3dSBarry Smith PetscFunctionBegin; 7298965ea79SLois Curfman McInnes *m = mat->rstart; *n = mat->rend; 7303a40ed3dSBarry Smith PetscFunctionReturn(0); 7318965ea79SLois Curfman McInnes } 7328965ea79SLois Curfman McInnes 7335615d1e5SSatish Balay #undef __FUNC__ 7345615d1e5SSatish Balay #define __FUNC__ "MatGetRow_MPIDense" 7358f6be9afSLois Curfman McInnes int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,Scalar **v) 7368965ea79SLois Curfman McInnes { 7373501a2bdSLois Curfman McInnes Mat_MPIDense *mat = (Mat_MPIDense *) A->data; 7383a40ed3dSBarry Smith int lrow, rstart = mat->rstart, rend = mat->rend,ierr; 7398965ea79SLois Curfman McInnes 7403a40ed3dSBarry Smith PetscFunctionBegin; 741a8c6a408SBarry Smith if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,0,"only local rows") 7428965ea79SLois Curfman McInnes lrow = row - rstart; 7433a40ed3dSBarry Smith ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr); 7443a40ed3dSBarry Smith PetscFunctionReturn(0); 7458965ea79SLois Curfman McInnes } 7468965ea79SLois Curfman McInnes 7475615d1e5SSatish Balay #undef __FUNC__ 7485615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_MPIDense" 7498f6be9afSLois Curfman McInnes int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,Scalar **v) 7508965ea79SLois Curfman McInnes { 7513a40ed3dSBarry Smith PetscFunctionBegin; 7520452661fSBarry Smith if (idx) PetscFree(*idx); 7530452661fSBarry Smith if (v) PetscFree(*v); 7543a40ed3dSBarry Smith PetscFunctionReturn(0); 7558965ea79SLois Curfman McInnes } 7568965ea79SLois Curfman McInnes 7575615d1e5SSatish Balay #undef __FUNC__ 7585615d1e5SSatish Balay #define __FUNC__ "MatNorm_MPIDense" 7598f6be9afSLois Curfman McInnes int MatNorm_MPIDense(Mat A,NormType type,double *norm) 760096963f5SLois Curfman McInnes { 7613501a2bdSLois Curfman McInnes Mat_MPIDense *mdn = (Mat_MPIDense *) A->data; 7623501a2bdSLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*) mdn->A->data; 7633501a2bdSLois Curfman McInnes int ierr, i, j; 7643501a2bdSLois Curfman McInnes double sum = 0.0; 7653501a2bdSLois Curfman McInnes Scalar *v = mat->v; 7663501a2bdSLois Curfman McInnes 7673a40ed3dSBarry Smith PetscFunctionBegin; 7683501a2bdSLois Curfman McInnes if (mdn->size == 1) { 7693501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,norm); CHKERRQ(ierr); 7703501a2bdSLois Curfman McInnes } else { 7713501a2bdSLois Curfman McInnes if (type == NORM_FROBENIUS) { 7723501a2bdSLois Curfman McInnes for (i=0; i<mat->n*mat->m; i++ ) { 7733a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 7743501a2bdSLois Curfman McInnes sum += real(conj(*v)*(*v)); v++; 7753501a2bdSLois Curfman McInnes #else 7763501a2bdSLois Curfman McInnes sum += (*v)*(*v); v++; 7773501a2bdSLois Curfman McInnes #endif 7783501a2bdSLois Curfman McInnes } 779ca161407SBarry Smith ierr = MPI_Allreduce(&sum,norm,1,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 7803501a2bdSLois Curfman McInnes *norm = sqrt(*norm); 7813501a2bdSLois Curfman McInnes PLogFlops(2*mat->n*mat->m); 7823a40ed3dSBarry Smith } else if (type == NORM_1) { 7833501a2bdSLois Curfman McInnes double *tmp, *tmp2; 7840452661fSBarry Smith tmp = (double *) PetscMalloc( 2*mdn->N*sizeof(double) ); CHKPTRQ(tmp); 7853501a2bdSLois Curfman McInnes tmp2 = tmp + mdn->N; 786cddf8d76SBarry Smith PetscMemzero(tmp,2*mdn->N*sizeof(double)); 787096963f5SLois Curfman McInnes *norm = 0.0; 7883501a2bdSLois Curfman McInnes v = mat->v; 7893501a2bdSLois Curfman McInnes for ( j=0; j<mat->n; j++ ) { 7903501a2bdSLois Curfman McInnes for ( i=0; i<mat->m; i++ ) { 79167e560aaSBarry Smith tmp[j] += PetscAbsScalar(*v); v++; 7923501a2bdSLois Curfman McInnes } 7933501a2bdSLois Curfman McInnes } 794ca161407SBarry Smith ierr = MPI_Allreduce(tmp,tmp2,mdn->N,MPI_DOUBLE,MPI_SUM,A->comm);CHKERRQ(ierr); 7953501a2bdSLois Curfman McInnes for ( j=0; j<mdn->N; j++ ) { 7963501a2bdSLois Curfman McInnes if (tmp2[j] > *norm) *norm = tmp2[j]; 7973501a2bdSLois Curfman McInnes } 7980452661fSBarry Smith PetscFree(tmp); 7993501a2bdSLois Curfman McInnes PLogFlops(mat->n*mat->m); 8003a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { /* max row norm */ 8013501a2bdSLois Curfman McInnes double ntemp; 8023501a2bdSLois Curfman McInnes ierr = MatNorm(mdn->A,type,&ntemp); CHKERRQ(ierr); 803ca161407SBarry Smith ierr = MPI_Allreduce(&ntemp,norm,1,MPI_DOUBLE,MPI_MAX,A->comm);CHKERRQ(ierr); 8043a40ed3dSBarry Smith } else { 805a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"No support for two norm"); 8063501a2bdSLois Curfman McInnes } 8073501a2bdSLois Curfman McInnes } 8083a40ed3dSBarry Smith PetscFunctionReturn(0); 8093501a2bdSLois Curfman McInnes } 8103501a2bdSLois Curfman McInnes 8115615d1e5SSatish Balay #undef __FUNC__ 8125615d1e5SSatish Balay #define __FUNC__ "MatTranspose_MPIDense" 8138f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout) 8143501a2bdSLois Curfman McInnes { 8153501a2bdSLois Curfman McInnes Mat_MPIDense *a = (Mat_MPIDense *) A->data; 8163501a2bdSLois Curfman McInnes Mat_SeqDense *Aloc = (Mat_SeqDense *) a->A->data; 8173501a2bdSLois Curfman McInnes Mat B; 8183501a2bdSLois Curfman McInnes int M = a->M, N = a->N, m, n, *rwork, rstart = a->rstart; 8193501a2bdSLois Curfman McInnes int j, i, ierr; 8203501a2bdSLois Curfman McInnes Scalar *v; 8213501a2bdSLois Curfman McInnes 8223a40ed3dSBarry Smith PetscFunctionBegin; 8237056b6fcSBarry Smith if (matout == PETSC_NULL && M != N) { 824a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Supports square matrix only in-place"); 8257056b6fcSBarry Smith } 8267056b6fcSBarry Smith ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr); 8273501a2bdSLois Curfman McInnes 8283501a2bdSLois Curfman McInnes m = Aloc->m; n = Aloc->n; v = Aloc->v; 8290452661fSBarry Smith rwork = (int *) PetscMalloc(n*sizeof(int)); CHKPTRQ(rwork); 8303501a2bdSLois Curfman McInnes for ( j=0; j<n; j++ ) { 8313501a2bdSLois Curfman McInnes for (i=0; i<m; i++) rwork[i] = rstart + i; 8323501a2bdSLois Curfman McInnes ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES); CHKERRQ(ierr); 8333501a2bdSLois Curfman McInnes v += m; 8343501a2bdSLois Curfman McInnes } 8350452661fSBarry Smith PetscFree(rwork); 8366d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 8376d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 8383638b69dSLois Curfman McInnes if (matout != PETSC_NULL) { 8393501a2bdSLois Curfman McInnes *matout = B; 8403501a2bdSLois Curfman McInnes } else { 841f830108cSBarry Smith PetscOps *Abops; 842f830108cSBarry Smith struct _MatOps *Aops; 843f830108cSBarry Smith 8443501a2bdSLois Curfman McInnes /* This isn't really an in-place transpose, but free data struct from a */ 8450452661fSBarry Smith PetscFree(a->rowners); 8463501a2bdSLois Curfman McInnes ierr = MatDestroy(a->A); CHKERRQ(ierr); 8473501a2bdSLois Curfman McInnes if (a->lvec) VecDestroy(a->lvec); 8483501a2bdSLois Curfman McInnes if (a->Mvctx) VecScatterDestroy(a->Mvctx); 8490452661fSBarry Smith PetscFree(a); 850f830108cSBarry Smith 851f830108cSBarry Smith /* 852f830108cSBarry Smith This is horrible, horrible code. We need to keep the 853f830108cSBarry Smith A pointers for the bops and ops but copy everything 854f830108cSBarry Smith else from C. 855f830108cSBarry Smith */ 856f830108cSBarry Smith Abops = A->bops; 857f830108cSBarry Smith Aops = A->ops; 858f09e8eb9SSatish Balay PetscMemcpy(A,B,sizeof(struct _p_Mat)); 859f830108cSBarry Smith A->bops = Abops; 860f830108cSBarry Smith A->ops = Aops; 861f830108cSBarry Smith 8620452661fSBarry Smith PetscHeaderDestroy(B); 8633501a2bdSLois Curfman McInnes } 8643a40ed3dSBarry Smith PetscFunctionReturn(0); 865096963f5SLois Curfman McInnes } 866096963f5SLois Curfman McInnes 86744cd7ae7SLois Curfman McInnes #include "pinclude/plapack.h" 8685615d1e5SSatish Balay #undef __FUNC__ 8695615d1e5SSatish Balay #define __FUNC__ "MatScale_MPIDense" 8708f6be9afSLois Curfman McInnes int MatScale_MPIDense(Scalar *alpha,Mat inA) 87144cd7ae7SLois Curfman McInnes { 87244cd7ae7SLois Curfman McInnes Mat_MPIDense *A = (Mat_MPIDense *) inA->data; 87344cd7ae7SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->A->data; 87444cd7ae7SLois Curfman McInnes int one = 1, nz; 87544cd7ae7SLois Curfman McInnes 8763a40ed3dSBarry Smith PetscFunctionBegin; 87744cd7ae7SLois Curfman McInnes nz = a->m*a->n; 87844cd7ae7SLois Curfman McInnes BLscal_( &nz, alpha, a->v, &one ); 87944cd7ae7SLois Curfman McInnes PLogFlops(nz); 8803a40ed3dSBarry Smith PetscFunctionReturn(0); 88144cd7ae7SLois Curfman McInnes } 88244cd7ae7SLois Curfman McInnes 8833d1612f7SBarry Smith static int MatConvertSameType_MPIDense(Mat,Mat *,int); 8848965ea79SLois Curfman McInnes 8858965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/ 88639ddd567SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_MPIDense, 88739ddd567SLois Curfman McInnes MatGetRow_MPIDense,MatRestoreRow_MPIDense, 88839ddd567SLois Curfman McInnes MatMult_MPIDense,MatMultAdd_MPIDense, 889096963f5SLois Curfman McInnes MatMultTrans_MPIDense,MatMultTransAdd_MPIDense, 890e7ca0642SLois Curfman McInnes 0,0, 89139ddd567SLois Curfman McInnes 0,0, 8928c469469SLois Curfman McInnes 0,0, 8938aaee692SLois Curfman McInnes 0,MatTranspose_MPIDense, 89439ddd567SLois Curfman McInnes MatGetInfo_MPIDense,0, 895096963f5SLois Curfman McInnes MatGetDiagonal_MPIDense,0,MatNorm_MPIDense, 89639ddd567SLois Curfman McInnes MatAssemblyBegin_MPIDense,MatAssemblyEnd_MPIDense, 8978965ea79SLois Curfman McInnes 0, 89839ddd567SLois Curfman McInnes MatSetOption_MPIDense,MatZeroEntries_MPIDense,MatZeroRows_MPIDense, 8993b2fbd54SBarry Smith 0,0, 9008aaee692SLois Curfman McInnes 0,0, 90139ddd567SLois Curfman McInnes MatGetSize_MPIDense,MatGetLocalSize_MPIDense, 90239ddd567SLois Curfman McInnes MatGetOwnershipRange_MPIDense, 903ff14e315SSatish Balay 0,0, MatGetArray_MPIDense, MatRestoreArray_MPIDense, 90494a9d846SBarry Smith MatConvertSameType_MPIDense, 905b49de8d1SLois Curfman McInnes 0,0,0,0,0, 9064e220ebcSLois Curfman McInnes 0,0,MatGetValues_MPIDense,0,0,MatScale_MPIDense, 9074e220ebcSLois Curfman McInnes 0,0,0,MatGetBlockSize_MPIDense}; 9088965ea79SLois Curfman McInnes 9095615d1e5SSatish Balay #undef __FUNC__ 9105615d1e5SSatish Balay #define __FUNC__ "MatCreateMPIDense" 9118965ea79SLois Curfman McInnes /*@C 91239ddd567SLois Curfman McInnes MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 9138965ea79SLois Curfman McInnes 914*db81eaa0SLois Curfman McInnes Collective on MPI_Comm 915*db81eaa0SLois Curfman McInnes 9168965ea79SLois Curfman McInnes Input Parameters: 917*db81eaa0SLois Curfman McInnes + comm - MPI communicator 9188965ea79SLois Curfman McInnes . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 919*db81eaa0SLois Curfman McInnes . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 9208965ea79SLois Curfman McInnes . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 921*db81eaa0SLois Curfman McInnes . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 922*db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 923dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 9248965ea79SLois Curfman McInnes 9258965ea79SLois Curfman McInnes Output Parameter: 926477f1c0bSLois Curfman McInnes . A - the matrix 9278965ea79SLois Curfman McInnes 928b259b22eSLois Curfman McInnes Notes: 92939ddd567SLois Curfman McInnes The dense format is fully compatible with standard Fortran 77 93039ddd567SLois Curfman McInnes storage by columns. 9318965ea79SLois Curfman McInnes 93218f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 93318f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 934b4fd4287SBarry Smith set data=PETSC_NULL. 93518f449edSLois Curfman McInnes 9368965ea79SLois Curfman McInnes The user MUST specify either the local or global matrix dimensions 9378965ea79SLois Curfman McInnes (possibly both). 9388965ea79SLois Curfman McInnes 9393501a2bdSLois Curfman McInnes Currently, the only parallel dense matrix decomposition is by rows, 9403501a2bdSLois Curfman McInnes so that n=N and each submatrix owns all of the global columns. 9413501a2bdSLois Curfman McInnes 94239ddd567SLois Curfman McInnes .keywords: matrix, dense, parallel 9438965ea79SLois Curfman McInnes 94439ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 9458965ea79SLois Curfman McInnes @*/ 946477f1c0bSLois Curfman McInnes int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,Scalar *data,Mat *A) 9478965ea79SLois Curfman McInnes { 9488965ea79SLois Curfman McInnes Mat mat; 94939ddd567SLois Curfman McInnes Mat_MPIDense *a; 95025cdf11fSBarry Smith int ierr, i,flg; 9518965ea79SLois Curfman McInnes 9523a40ed3dSBarry Smith PetscFunctionBegin; 953ed2daf61SLois Curfman McInnes /* Note: For now, when data is specified above, this assumes the user correctly 954ed2daf61SLois Curfman McInnes allocates the local dense storage space. We should add error checking. */ 95518f449edSLois Curfman McInnes 956477f1c0bSLois Curfman McInnes *A = 0; 957f830108cSBarry Smith PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,comm,MatDestroy,MatView); 9588965ea79SLois Curfman McInnes PLogObjectCreate(mat); 9590452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIDense)); CHKPTRQ(a); 960f830108cSBarry Smith PetscMemcpy(mat->ops,&MatOps,sizeof(struct _MatOps)); 961e1311b90SBarry Smith mat->ops->destroy = MatDestroy_MPIDense; 962e1311b90SBarry Smith mat->ops->view = MatView_MPIDense; 9638965ea79SLois Curfman McInnes mat->factor = 0; 96490f02eecSBarry Smith mat->mapping = 0; 9658965ea79SLois Curfman McInnes 966622d7880SLois Curfman McInnes a->factor = 0; 967e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 9688965ea79SLois Curfman McInnes MPI_Comm_rank(comm,&a->rank); 9698965ea79SLois Curfman McInnes MPI_Comm_size(comm,&a->size); 9708965ea79SLois Curfman McInnes 971ca161407SBarry Smith if (M == PETSC_DECIDE) {ierr = MPI_Allreduce(&m,&M,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr);} 9728965ea79SLois Curfman McInnes if (m == PETSC_DECIDE) {m = M/a->size + ((M % a->size) > a->rank);} 97339ddd567SLois Curfman McInnes 97439ddd567SLois Curfman McInnes /* each row stores all columns */ 97539ddd567SLois Curfman McInnes if (N == PETSC_DECIDE) N = n; 976d7e8b826SBarry Smith if (n == PETSC_DECIDE) {n = N/a->size + ((N % a->size) > a->rank);} 977a8c6a408SBarry Smith /* if (n != N) SETERRQ(PETSC_ERR_SUP,0,"For now, only n=N is supported"); */ 978aca0ad90SLois Curfman McInnes a->N = mat->N = N; 979aca0ad90SLois Curfman McInnes a->M = mat->M = M; 980aca0ad90SLois Curfman McInnes a->m = mat->m = m; 981aca0ad90SLois Curfman McInnes a->n = mat->n = n; 9828965ea79SLois Curfman McInnes 9838965ea79SLois Curfman McInnes /* build local table of row and column ownerships */ 984d7e8b826SBarry Smith a->rowners = (int *) PetscMalloc(2*(a->size+2)*sizeof(int)); CHKPTRQ(a->rowners); 985d7e8b826SBarry Smith a->cowners = a->rowners + a->size + 1; 986f09e8eb9SSatish Balay PLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 987ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,a->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 9888965ea79SLois Curfman McInnes a->rowners[0] = 0; 9898965ea79SLois Curfman McInnes for ( i=2; i<=a->size; i++ ) { 9908965ea79SLois Curfman McInnes a->rowners[i] += a->rowners[i-1]; 9918965ea79SLois Curfman McInnes } 9928965ea79SLois Curfman McInnes a->rstart = a->rowners[a->rank]; 9938965ea79SLois Curfman McInnes a->rend = a->rowners[a->rank+1]; 994ca161407SBarry Smith ierr = MPI_Allgather(&n,1,MPI_INT,a->cowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 995d7e8b826SBarry Smith a->cowners[0] = 0; 996d7e8b826SBarry Smith for ( i=2; i<=a->size; i++ ) { 997d7e8b826SBarry Smith a->cowners[i] += a->cowners[i-1]; 998d7e8b826SBarry Smith } 9998965ea79SLois Curfman McInnes 1000029af93fSBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,m,N,data,&a->A); CHKERRQ(ierr); 10018965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 10028965ea79SLois Curfman McInnes 10038965ea79SLois Curfman McInnes /* build cache for off array entries formed */ 10048965ea79SLois Curfman McInnes ierr = StashBuild_Private(&a->stash); CHKERRQ(ierr); 10058965ea79SLois Curfman McInnes 10068965ea79SLois Curfman McInnes /* stuff used for matrix vector multiply */ 10078965ea79SLois Curfman McInnes a->lvec = 0; 10088965ea79SLois Curfman McInnes a->Mvctx = 0; 100939b7565bSBarry Smith a->roworiented = 1; 10108965ea79SLois Curfman McInnes 1011477f1c0bSLois Curfman McInnes *A = mat; 101225cdf11fSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 101325cdf11fSBarry Smith if (flg) { 10148c469469SLois Curfman McInnes ierr = MatPrintHelp(mat); CHKERRQ(ierr); 10158c469469SLois Curfman McInnes } 10163a40ed3dSBarry Smith PetscFunctionReturn(0); 10178965ea79SLois Curfman McInnes } 10188965ea79SLois Curfman McInnes 10195615d1e5SSatish Balay #undef __FUNC__ 10205615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_MPIDense" 10213d1612f7SBarry Smith static int MatConvertSameType_MPIDense(Mat A,Mat *newmat,int cpvalues) 10228965ea79SLois Curfman McInnes { 10238965ea79SLois Curfman McInnes Mat mat; 10243501a2bdSLois Curfman McInnes Mat_MPIDense *a,*oldmat = (Mat_MPIDense *) A->data; 102539ddd567SLois Curfman McInnes int ierr; 10262ba99913SLois Curfman McInnes FactorCtx *factor; 10278965ea79SLois Curfman McInnes 10283a40ed3dSBarry Smith PetscFunctionBegin; 10298965ea79SLois Curfman McInnes *newmat = 0; 1030f830108cSBarry Smith PetscHeaderCreate(mat,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIDENSE,A->comm,MatDestroy,MatView); 10318965ea79SLois Curfman McInnes PLogObjectCreate(mat); 10320452661fSBarry Smith mat->data = (void *) (a = PetscNew(Mat_MPIDense)); CHKPTRQ(a); 1033f830108cSBarry Smith PetscMemcpy(mat->ops,&MatOps,sizeof(struct _MatOps)); 1034e1311b90SBarry Smith mat->ops->destroy = MatDestroy_MPIDense; 1035e1311b90SBarry Smith mat->ops->view = MatView_MPIDense; 10363501a2bdSLois Curfman McInnes mat->factor = A->factor; 1037c456f294SBarry Smith mat->assembled = PETSC_TRUE; 10388965ea79SLois Curfman McInnes 103944cd7ae7SLois Curfman McInnes a->m = mat->m = oldmat->m; 104044cd7ae7SLois Curfman McInnes a->n = mat->n = oldmat->n; 104144cd7ae7SLois Curfman McInnes a->M = mat->M = oldmat->M; 104244cd7ae7SLois Curfman McInnes a->N = mat->N = oldmat->N; 10432ba99913SLois Curfman McInnes if (oldmat->factor) { 10442ba99913SLois Curfman McInnes a->factor = (FactorCtx *) (factor = PetscNew(FactorCtx)); CHKPTRQ(factor); 10452ba99913SLois Curfman McInnes /* copy factor contents ... add this code! */ 10462ba99913SLois Curfman McInnes } else a->factor = 0; 10478965ea79SLois Curfman McInnes 10488965ea79SLois Curfman McInnes a->rstart = oldmat->rstart; 10498965ea79SLois Curfman McInnes a->rend = oldmat->rend; 10508965ea79SLois Curfman McInnes a->size = oldmat->size; 10518965ea79SLois Curfman McInnes a->rank = oldmat->rank; 1052e0fa3b82SLois Curfman McInnes mat->insertmode = NOT_SET_VALUES; 10538965ea79SLois Curfman McInnes 10540452661fSBarry Smith a->rowners = (int *) PetscMalloc((a->size+1)*sizeof(int)); CHKPTRQ(a->rowners); 1055f09e8eb9SSatish Balay PLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 10568965ea79SLois Curfman McInnes PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int)); 10578965ea79SLois Curfman McInnes ierr = StashInitialize_Private(&a->stash); CHKERRQ(ierr); 10588965ea79SLois Curfman McInnes 10598965ea79SLois Curfman McInnes ierr = VecDuplicate(oldmat->lvec,&a->lvec); CHKERRQ(ierr); 10608965ea79SLois Curfman McInnes PLogObjectParent(mat,a->lvec); 106155659b69SBarry Smith ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx); CHKERRQ(ierr); 10628965ea79SLois Curfman McInnes PLogObjectParent(mat,a->Mvctx); 10638965ea79SLois Curfman McInnes ierr = MatConvert(oldmat->A,MATSAME,&a->A); CHKERRQ(ierr); 10648965ea79SLois Curfman McInnes PLogObjectParent(mat,a->A); 10658965ea79SLois Curfman McInnes *newmat = mat; 10663a40ed3dSBarry Smith PetscFunctionReturn(0); 10678965ea79SLois Curfman McInnes } 10688965ea79SLois Curfman McInnes 106977c4ece6SBarry Smith #include "sys.h" 10708965ea79SLois Curfman McInnes 10715615d1e5SSatish Balay #undef __FUNC__ 10725615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense_DenseInFile" 107390ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M, int N, Mat *newmat) 107490ace30eSBarry Smith { 107540011551SBarry Smith int *rowners, i,size,rank,m,ierr,nz,j; 107690ace30eSBarry Smith Scalar *array,*vals,*vals_ptr; 107790ace30eSBarry Smith MPI_Status status; 107890ace30eSBarry Smith 10793a40ed3dSBarry Smith PetscFunctionBegin; 108090ace30eSBarry Smith MPI_Comm_rank(comm,&rank); 108190ace30eSBarry Smith MPI_Comm_size(comm,&size); 108290ace30eSBarry Smith 108390ace30eSBarry Smith /* determine ownership of all rows */ 108490ace30eSBarry Smith m = M/size + ((M % size) > rank); 108590ace30eSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners); 1086ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 108790ace30eSBarry Smith rowners[0] = 0; 108890ace30eSBarry Smith for ( i=2; i<=size; i++ ) { 108990ace30eSBarry Smith rowners[i] += rowners[i-1]; 109090ace30eSBarry Smith } 109190ace30eSBarry Smith 109290ace30eSBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 109390ace30eSBarry Smith ierr = MatGetArray(*newmat,&array); CHKERRQ(ierr); 109490ace30eSBarry Smith 109590ace30eSBarry Smith if (!rank) { 109690ace30eSBarry Smith vals = (Scalar *) PetscMalloc( m*N*sizeof(Scalar) ); CHKPTRQ(vals); 109790ace30eSBarry Smith 109890ace30eSBarry Smith /* read in my part of the matrix numerical values */ 10990752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR); CHKERRQ(ierr); 110090ace30eSBarry Smith 110190ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 110290ace30eSBarry Smith vals_ptr = vals; 110390ace30eSBarry Smith for ( i=0; i<m; i++ ) { 110490ace30eSBarry Smith for ( j=0; j<N; j++ ) { 110590ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 110690ace30eSBarry Smith } 110790ace30eSBarry Smith } 110890ace30eSBarry Smith 110990ace30eSBarry Smith /* read in other processors and ship out */ 111090ace30eSBarry Smith for ( i=1; i<size; i++ ) { 111190ace30eSBarry Smith nz = (rowners[i+1] - rowners[i])*N; 11120752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr); 1113ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 111490ace30eSBarry Smith } 11153a40ed3dSBarry Smith } else { 111690ace30eSBarry Smith /* receive numeric values */ 111790ace30eSBarry Smith vals = (Scalar*) PetscMalloc( m*N*sizeof(Scalar) ); CHKPTRQ(vals); 111890ace30eSBarry Smith 111990ace30eSBarry Smith /* receive message of values*/ 1120ca161407SBarry Smith ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 112190ace30eSBarry Smith 112290ace30eSBarry Smith /* insert into matrix-by row (this is why cannot directly read into array */ 112390ace30eSBarry Smith vals_ptr = vals; 112490ace30eSBarry Smith for ( i=0; i<m; i++ ) { 112590ace30eSBarry Smith for ( j=0; j<N; j++ ) { 112690ace30eSBarry Smith array[i + j*m] = *vals_ptr++; 112790ace30eSBarry Smith } 112890ace30eSBarry Smith } 112990ace30eSBarry Smith } 113090ace30eSBarry Smith PetscFree(rowners); 113190ace30eSBarry Smith PetscFree(vals); 11326d4a8577SBarry Smith ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 11336d4a8577SBarry Smith ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 11343a40ed3dSBarry Smith PetscFunctionReturn(0); 113590ace30eSBarry Smith } 113690ace30eSBarry Smith 113790ace30eSBarry Smith 11385615d1e5SSatish Balay #undef __FUNC__ 11395615d1e5SSatish Balay #define __FUNC__ "MatLoad_MPIDense" 114019bcc07fSBarry Smith int MatLoad_MPIDense(Viewer viewer,MatType type,Mat *newmat) 11418965ea79SLois Curfman McInnes { 11428965ea79SLois Curfman McInnes Mat A; 11438965ea79SLois Curfman McInnes Scalar *vals,*svals; 114419bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 11458965ea79SLois Curfman McInnes MPI_Status status; 11468965ea79SLois Curfman McInnes int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 11478965ea79SLois Curfman McInnes int *ourlens,*sndcounts = 0,*procsnz = 0, *offlens,jj,*mycols,*smycols; 114819bcc07fSBarry Smith int tag = ((PetscObject)viewer)->tag; 11493a40ed3dSBarry Smith int i, nz, ierr, j,rstart, rend, fd; 11508965ea79SLois Curfman McInnes 11513a40ed3dSBarry Smith PetscFunctionBegin; 11528965ea79SLois Curfman McInnes MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank); 11538965ea79SLois Curfman McInnes if (!rank) { 115490ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 11550752156aSBarry Smith ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT); CHKERRQ(ierr); 1156a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object"); 11578965ea79SLois Curfman McInnes } 11588965ea79SLois Curfman McInnes 1159ca161407SBarry Smith ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr); 116090ace30eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 116190ace30eSBarry Smith 116290ace30eSBarry Smith /* 116390ace30eSBarry Smith Handle case where matrix is stored on disk as a dense matrix 116490ace30eSBarry Smith */ 116590ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { 11663a40ed3dSBarry Smith ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr); 11673a40ed3dSBarry Smith PetscFunctionReturn(0); 116890ace30eSBarry Smith } 116990ace30eSBarry Smith 11708965ea79SLois Curfman McInnes /* determine ownership of all rows */ 11718965ea79SLois Curfman McInnes m = M/size + ((M % size) > rank); 11720452661fSBarry Smith rowners = (int *) PetscMalloc((size+2)*sizeof(int)); CHKPTRQ(rowners); 1173ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 11748965ea79SLois Curfman McInnes rowners[0] = 0; 11758965ea79SLois Curfman McInnes for ( i=2; i<=size; i++ ) { 11768965ea79SLois Curfman McInnes rowners[i] += rowners[i-1]; 11778965ea79SLois Curfman McInnes } 11788965ea79SLois Curfman McInnes rstart = rowners[rank]; 11798965ea79SLois Curfman McInnes rend = rowners[rank+1]; 11808965ea79SLois Curfman McInnes 11818965ea79SLois Curfman McInnes /* distribute row lengths to all processors */ 11820452661fSBarry Smith ourlens = (int*) PetscMalloc( 2*(rend-rstart)*sizeof(int) ); CHKPTRQ(ourlens); 11838965ea79SLois Curfman McInnes offlens = ourlens + (rend-rstart); 11848965ea79SLois Curfman McInnes if (!rank) { 11850452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 11860752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr); 11870452661fSBarry Smith sndcounts = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(sndcounts); 11888965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) sndcounts[i] = rowners[i+1] - rowners[i]; 1189ca161407SBarry Smith ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 11900452661fSBarry Smith PetscFree(sndcounts); 1191ca161407SBarry Smith } else { 1192ca161407SBarry Smith ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT, 0,comm);CHKERRQ(ierr); 11938965ea79SLois Curfman McInnes } 11948965ea79SLois Curfman McInnes 11958965ea79SLois Curfman McInnes if (!rank) { 11968965ea79SLois Curfman McInnes /* calculate the number of nonzeros on each processor */ 11970452661fSBarry Smith procsnz = (int*) PetscMalloc( size*sizeof(int) ); CHKPTRQ(procsnz); 1198cddf8d76SBarry Smith PetscMemzero(procsnz,size*sizeof(int)); 11998965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 12008965ea79SLois Curfman McInnes for ( j=rowners[i]; j< rowners[i+1]; j++ ) { 12018965ea79SLois Curfman McInnes procsnz[i] += rowlengths[j]; 12028965ea79SLois Curfman McInnes } 12038965ea79SLois Curfman McInnes } 12040452661fSBarry Smith PetscFree(rowlengths); 12058965ea79SLois Curfman McInnes 12068965ea79SLois Curfman McInnes /* determine max buffer needed and allocate it */ 12078965ea79SLois Curfman McInnes maxnz = 0; 12088965ea79SLois Curfman McInnes for ( i=0; i<size; i++ ) { 12090452661fSBarry Smith maxnz = PetscMax(maxnz,procsnz[i]); 12108965ea79SLois Curfman McInnes } 12110452661fSBarry Smith cols = (int *) PetscMalloc( maxnz*sizeof(int) ); CHKPTRQ(cols); 12128965ea79SLois Curfman McInnes 12138965ea79SLois Curfman McInnes /* read in my part of the matrix column indices */ 12148965ea79SLois Curfman McInnes nz = procsnz[0]; 12150452661fSBarry Smith mycols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 12160752156aSBarry Smith ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT); CHKERRQ(ierr); 12178965ea79SLois Curfman McInnes 12188965ea79SLois Curfman McInnes /* read in every one elses and ship off */ 12198965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { 12208965ea79SLois Curfman McInnes nz = procsnz[i]; 12210752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT); CHKERRQ(ierr); 1222ca161407SBarry Smith ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr); 12238965ea79SLois Curfman McInnes } 12240452661fSBarry Smith PetscFree(cols); 12253a40ed3dSBarry Smith } else { 12268965ea79SLois Curfman McInnes /* determine buffer space needed for message */ 12278965ea79SLois Curfman McInnes nz = 0; 12288965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12298965ea79SLois Curfman McInnes nz += ourlens[i]; 12308965ea79SLois Curfman McInnes } 12310452661fSBarry Smith mycols = (int*) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(mycols); 12328965ea79SLois Curfman McInnes 12338965ea79SLois Curfman McInnes /* receive message of column indices*/ 1234ca161407SBarry Smith ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr); 1235ca161407SBarry Smith ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr); 1236a8c6a408SBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 12378965ea79SLois Curfman McInnes } 12388965ea79SLois Curfman McInnes 12398965ea79SLois Curfman McInnes /* loop over local rows, determining number of off diagonal entries */ 1240cddf8d76SBarry Smith PetscMemzero(offlens,m*sizeof(int)); 12418965ea79SLois Curfman McInnes jj = 0; 12428965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12438965ea79SLois Curfman McInnes for ( j=0; j<ourlens[i]; j++ ) { 12448965ea79SLois Curfman McInnes if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 12458965ea79SLois Curfman McInnes jj++; 12468965ea79SLois Curfman McInnes } 12478965ea79SLois Curfman McInnes } 12488965ea79SLois Curfman McInnes 12498965ea79SLois Curfman McInnes /* create our matrix */ 12508965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12518965ea79SLois Curfman McInnes ourlens[i] -= offlens[i]; 12528965ea79SLois Curfman McInnes } 1253b4fd4287SBarry Smith ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr); 12548965ea79SLois Curfman McInnes A = *newmat; 12558965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12568965ea79SLois Curfman McInnes ourlens[i] += offlens[i]; 12578965ea79SLois Curfman McInnes } 12588965ea79SLois Curfman McInnes 12598965ea79SLois Curfman McInnes if (!rank) { 12600452661fSBarry Smith vals = (Scalar *) PetscMalloc( maxnz*sizeof(Scalar) ); CHKPTRQ(vals); 12618965ea79SLois Curfman McInnes 12628965ea79SLois Curfman McInnes /* read in my part of the matrix numerical values */ 12638965ea79SLois Curfman McInnes nz = procsnz[0]; 12640752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr); 12658965ea79SLois Curfman McInnes 12668965ea79SLois Curfman McInnes /* insert into matrix */ 12678965ea79SLois Curfman McInnes jj = rstart; 12688965ea79SLois Curfman McInnes smycols = mycols; 12698965ea79SLois Curfman McInnes svals = vals; 12708965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12718965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 12728965ea79SLois Curfman McInnes smycols += ourlens[i]; 12738965ea79SLois Curfman McInnes svals += ourlens[i]; 12748965ea79SLois Curfman McInnes jj++; 12758965ea79SLois Curfman McInnes } 12768965ea79SLois Curfman McInnes 12778965ea79SLois Curfman McInnes /* read in other processors and ship out */ 12788965ea79SLois Curfman McInnes for ( i=1; i<size; i++ ) { 12798965ea79SLois Curfman McInnes nz = procsnz[i]; 12800752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr); 1281ca161407SBarry Smith ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 12828965ea79SLois Curfman McInnes } 12830452661fSBarry Smith PetscFree(procsnz); 12843a40ed3dSBarry Smith } else { 12858965ea79SLois Curfman McInnes /* receive numeric values */ 12860452661fSBarry Smith vals = (Scalar*) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals); 12878965ea79SLois Curfman McInnes 12888965ea79SLois Curfman McInnes /* receive message of values*/ 1289ca161407SBarry Smith ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1290ca161407SBarry Smith ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 1291a8c6a408SBarry Smith if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"something is wrong with file"); 12928965ea79SLois Curfman McInnes 12938965ea79SLois Curfman McInnes /* insert into matrix */ 12948965ea79SLois Curfman McInnes jj = rstart; 12958965ea79SLois Curfman McInnes smycols = mycols; 12968965ea79SLois Curfman McInnes svals = vals; 12978965ea79SLois Curfman McInnes for ( i=0; i<m; i++ ) { 12988965ea79SLois Curfman McInnes ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 12998965ea79SLois Curfman McInnes smycols += ourlens[i]; 13008965ea79SLois Curfman McInnes svals += ourlens[i]; 13018965ea79SLois Curfman McInnes jj++; 13028965ea79SLois Curfman McInnes } 13038965ea79SLois Curfman McInnes } 13040452661fSBarry Smith PetscFree(ourlens); PetscFree(vals); PetscFree(mycols); PetscFree(rowners); 13058965ea79SLois Curfman McInnes 13066d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 13076d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 13083a40ed3dSBarry Smith PetscFunctionReturn(0); 13098965ea79SLois Curfman McInnes } 131090ace30eSBarry Smith 131190ace30eSBarry Smith 131290ace30eSBarry Smith 131390ace30eSBarry Smith 131490ace30eSBarry Smith 1315