xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision 5c5985e7f6792cd82332f5d47aeab782735e3201)
173f4d377SMatthew Knepley /*$Id: mpidense.c,v 1.159 2001/08/10 03:30:41 bsmith Exp $*/
28965ea79SLois Curfman McInnes 
3ed3cc1f0SBarry Smith /*
4ed3cc1f0SBarry Smith    Basic functions for basic parallel dense matrices.
5ed3cc1f0SBarry Smith */
6ed3cc1f0SBarry Smith 
770f55243SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h"
88965ea79SLois Curfman McInnes 
90de54da6SSatish Balay EXTERN_C_BEGIN
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense"
120de54da6SSatish Balay int MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B)
130de54da6SSatish Balay {
140de54da6SSatish Balay   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
15cfce73b9SSatish Balay   int          m = A->m,rstart = mdn->rstart,ierr;
1687828ca2SBarry Smith   PetscScalar  *array;
170de54da6SSatish Balay   MPI_Comm     comm;
180de54da6SSatish Balay 
190de54da6SSatish Balay   PetscFunctionBegin;
20273d9f13SBarry Smith   if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported.");
210de54da6SSatish Balay 
220de54da6SSatish Balay   /* The reuse aspect is not implemented efficiently */
230de54da6SSatish Balay   if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);}
240de54da6SSatish Balay 
250de54da6SSatish Balay   ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr);
260de54da6SSatish Balay   ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr);
27*5c5985e7SKris Buschelman   ierr = MatCreate(comm,m,m,m,m,B);CHKERRQ(ierr);
28*5c5985e7SKris Buschelman   ierr = MatSetType(*B,mdn->A->type_name);CHKERRQ(ierr);
29*5c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr);
300de54da6SSatish Balay   ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr);
310de54da6SSatish Balay   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
320de54da6SSatish Balay   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
330de54da6SSatish Balay 
340de54da6SSatish Balay   *iscopy = PETSC_TRUE;
350de54da6SSatish Balay   PetscFunctionReturn(0);
360de54da6SSatish Balay }
370de54da6SSatish Balay EXTERN_C_END
380de54da6SSatish Balay 
394a2ae208SSatish Balay #undef __FUNCT__
404a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense"
41f15d580aSBarry Smith int MatSetValues_MPIDense(Mat mat,int m,const int idxm[],int n,const int idxn[],const PetscScalar v[],InsertMode addv)
428965ea79SLois Curfman McInnes {
4339b7565bSBarry Smith   Mat_MPIDense *A = (Mat_MPIDense*)mat->data;
4439b7565bSBarry Smith   int          ierr,i,j,rstart = A->rstart,rend = A->rend,row;
45273d9f13SBarry Smith   PetscTruth   roworiented = A->roworiented;
468965ea79SLois Curfman McInnes 
473a40ed3dSBarry Smith   PetscFunctionBegin;
488965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
495ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
50273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
518965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
528965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
5339b7565bSBarry Smith       if (roworiented) {
5439b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
553a40ed3dSBarry Smith       } else {
568965ea79SLois Curfman McInnes         for (j=0; j<n; j++) {
575ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
58273d9f13SBarry Smith           if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
5939b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
6039b7565bSBarry Smith         }
618965ea79SLois Curfman McInnes       }
623a40ed3dSBarry Smith     } else {
633782ba37SSatish Balay       if (!A->donotstash) {
6439b7565bSBarry Smith         if (roworiented) {
658798bf22SSatish Balay           ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr);
66d36fbae8SSatish Balay         } else {
678798bf22SSatish Balay           ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr);
6839b7565bSBarry Smith         }
69b49de8d1SLois Curfman McInnes       }
70b49de8d1SLois Curfman McInnes     }
713782ba37SSatish Balay   }
723a40ed3dSBarry Smith   PetscFunctionReturn(0);
73b49de8d1SLois Curfman McInnes }
74b49de8d1SLois Curfman McInnes 
754a2ae208SSatish Balay #undef __FUNCT__
764a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense"
77f15d580aSBarry Smith int MatGetValues_MPIDense(Mat mat,int m,const int idxm[],int n,const int idxn[],PetscScalar v[])
78b49de8d1SLois Curfman McInnes {
79b49de8d1SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
80b49de8d1SLois Curfman McInnes   int          ierr,i,j,rstart = mdn->rstart,rend = mdn->rend,row;
81b49de8d1SLois Curfman McInnes 
823a40ed3dSBarry Smith   PetscFunctionBegin;
83b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
8429bbc08cSBarry Smith     if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
85273d9f13SBarry Smith     if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
86b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
87b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
88b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
8929bbc08cSBarry Smith         if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
90273d9f13SBarry Smith         if (idxn[j] >= mat->N) {
9129bbc08cSBarry Smith           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
92a8c6a408SBarry Smith         }
93b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
94b49de8d1SLois Curfman McInnes       }
95a8c6a408SBarry Smith     } else {
9629bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"Only local values currently supported");
978965ea79SLois Curfman McInnes     }
988965ea79SLois Curfman McInnes   }
993a40ed3dSBarry Smith   PetscFunctionReturn(0);
1008965ea79SLois Curfman McInnes }
1018965ea79SLois Curfman McInnes 
1024a2ae208SSatish Balay #undef __FUNCT__
1034a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense"
1044e7234bfSBarry Smith int MatGetArray_MPIDense(Mat A,PetscScalar *array[])
105ff14e315SSatish Balay {
106ff14e315SSatish Balay   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
107ff14e315SSatish Balay   int          ierr;
108ff14e315SSatish Balay 
1093a40ed3dSBarry Smith   PetscFunctionBegin;
110ff14e315SSatish Balay   ierr = MatGetArray(a->A,array);CHKERRQ(ierr);
1113a40ed3dSBarry Smith   PetscFunctionReturn(0);
112ff14e315SSatish Balay }
113ff14e315SSatish Balay 
1144a2ae208SSatish Balay #undef __FUNCT__
1154a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense"
116ca3fa75bSLois Curfman McInnes static int MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
117ca3fa75bSLois Curfman McInnes {
118ca3fa75bSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd;
119ca3fa75bSLois Curfman McInnes   Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data;
120cfce73b9SSatish Balay   int          i,j,ierr,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols;
12187828ca2SBarry Smith   PetscScalar  *av,*bv,*v = lmat->v;
122ca3fa75bSLois Curfman McInnes   Mat          newmat;
123ca3fa75bSLois Curfman McInnes 
124ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
125ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
126ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
127b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
128b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
129ca3fa75bSLois Curfman McInnes 
130ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
1317eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
1327eba5e9cSLois Curfman McInnes      original matrix! */
133ca3fa75bSLois Curfman McInnes 
134ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
1357eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
136ca3fa75bSLois Curfman McInnes 
137ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
138ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
13929bbc08cSBarry Smith     /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
1407eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
141ca3fa75bSLois Curfman McInnes     newmat = *B;
142ca3fa75bSLois Curfman McInnes   } else {
143ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
14432828cfdSBarry Smith     ierr = MatCreateMPIDense(A->comm,nrows,cs,PETSC_DECIDE,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr);
145ca3fa75bSLois Curfman McInnes   }
146ca3fa75bSLois Curfman McInnes 
147ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
148ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
149ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense *)newmatd->A->data)->v;
150ca3fa75bSLois Curfman McInnes 
151ca3fa75bSLois Curfman McInnes   for (i=0; i<ncols; i++) {
152ca3fa75bSLois Curfman McInnes     av = v + nlrows*icol[i];
153ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
1547eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
155ca3fa75bSLois Curfman McInnes     }
156ca3fa75bSLois Curfman McInnes   }
157ca3fa75bSLois Curfman McInnes 
158ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
159ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
160ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
161ca3fa75bSLois Curfman McInnes 
162ca3fa75bSLois Curfman McInnes   /* Free work space */
163ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
164ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
165ca3fa75bSLois Curfman McInnes   *B = newmat;
166ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
167ca3fa75bSLois Curfman McInnes }
168ca3fa75bSLois Curfman McInnes 
1694a2ae208SSatish Balay #undef __FUNCT__
1704a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense"
1714e7234bfSBarry Smith int MatRestoreArray_MPIDense(Mat A,PetscScalar *array[])
172ff14e315SSatish Balay {
1733a40ed3dSBarry Smith   PetscFunctionBegin;
1743a40ed3dSBarry Smith   PetscFunctionReturn(0);
175ff14e315SSatish Balay }
176ff14e315SSatish Balay 
1774a2ae208SSatish Balay #undef __FUNCT__
1784a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense"
1798f6be9afSLois Curfman McInnes int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
1808965ea79SLois Curfman McInnes {
18139ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
1828965ea79SLois Curfman McInnes   MPI_Comm     comm = mat->comm;
183d36fbae8SSatish Balay   int          ierr,nstash,reallocs;
1848965ea79SLois Curfman McInnes   InsertMode   addv;
1858965ea79SLois Curfman McInnes 
1863a40ed3dSBarry Smith   PetscFunctionBegin;
1878965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
188ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
1897056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
19029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
1918965ea79SLois Curfman McInnes   }
192e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
1938965ea79SLois Curfman McInnes 
1948798bf22SSatish Balay   ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr);
1958798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
196b0a32e0cSBarry Smith   PetscLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %d entries, uses %d mallocs.\n",nstash,reallocs);
1973a40ed3dSBarry Smith   PetscFunctionReturn(0);
1988965ea79SLois Curfman McInnes }
1998965ea79SLois Curfman McInnes 
2004a2ae208SSatish Balay #undef __FUNCT__
2014a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense"
2028f6be9afSLois Curfman McInnes int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2038965ea79SLois Curfman McInnes {
20439ddd567SLois Curfman McInnes   Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data;
2057ef1d9bdSSatish Balay   int          i,n,ierr,*row,*col,flg,j,rstart,ncols;
20687828ca2SBarry Smith   PetscScalar  *val;
207e0fa3b82SLois Curfman McInnes   InsertMode   addv=mat->insertmode;
2088965ea79SLois Curfman McInnes 
2093a40ed3dSBarry Smith   PetscFunctionBegin;
2108965ea79SLois Curfman McInnes   /*  wait on receives */
2117ef1d9bdSSatish Balay   while (1) {
2128798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2137ef1d9bdSSatish Balay     if (!flg) break;
2148965ea79SLois Curfman McInnes 
2157ef1d9bdSSatish Balay     for (i=0; i<n;) {
2167ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2177ef1d9bdSSatish Balay       for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
2187ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2197ef1d9bdSSatish Balay       else       ncols = n-i;
2207ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
2217ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
2227ef1d9bdSSatish Balay       i = j;
2238965ea79SLois Curfman McInnes     }
2247ef1d9bdSSatish Balay   }
2258798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
2268965ea79SLois Curfman McInnes 
22739ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
22839ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
2298965ea79SLois Curfman McInnes 
2306d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
23139ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
2328965ea79SLois Curfman McInnes   }
2333a40ed3dSBarry Smith   PetscFunctionReturn(0);
2348965ea79SLois Curfman McInnes }
2358965ea79SLois Curfman McInnes 
2364a2ae208SSatish Balay #undef __FUNCT__
2374a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense"
2388f6be9afSLois Curfman McInnes int MatZeroEntries_MPIDense(Mat A)
2398965ea79SLois Curfman McInnes {
2403a40ed3dSBarry Smith   int          ierr;
24139ddd567SLois Curfman McInnes   Mat_MPIDense *l = (Mat_MPIDense*)A->data;
2423a40ed3dSBarry Smith 
2433a40ed3dSBarry Smith   PetscFunctionBegin;
2443a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
2453a40ed3dSBarry Smith   PetscFunctionReturn(0);
2468965ea79SLois Curfman McInnes }
2478965ea79SLois Curfman McInnes 
2484a2ae208SSatish Balay #undef __FUNCT__
2494a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_MPIDense"
2508f6be9afSLois Curfman McInnes int MatGetBlockSize_MPIDense(Mat A,int *bs)
2514e220ebcSLois Curfman McInnes {
2523a40ed3dSBarry Smith   PetscFunctionBegin;
2534e220ebcSLois Curfman McInnes   *bs = 1;
2543a40ed3dSBarry Smith   PetscFunctionReturn(0);
2554e220ebcSLois Curfman McInnes }
2564e220ebcSLois Curfman McInnes 
2578965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
2588965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
2598965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
2603501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
2618965ea79SLois Curfman McInnes    routine.
2628965ea79SLois Curfman McInnes */
2634a2ae208SSatish Balay #undef __FUNCT__
2644a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense"
265268466fbSBarry Smith int MatZeroRows_MPIDense(Mat A,IS is,const PetscScalar *diag)
2668965ea79SLois Curfman McInnes {
26739ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
2688965ea79SLois Curfman McInnes   int            i,ierr,N,*rows,*owners = l->rowners,size = l->size;
269c1dc657dSBarry Smith   int            *nprocs,j,idx,nsends;
2708965ea79SLois Curfman McInnes   int            nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank;
2718965ea79SLois Curfman McInnes   int            *rvalues,tag = A->tag,count,base,slen,n,*source;
2728965ea79SLois Curfman McInnes   int            *lens,imdex,*lrows,*values;
2738965ea79SLois Curfman McInnes   MPI_Comm       comm = A->comm;
2748965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
2758965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
2768965ea79SLois Curfman McInnes   IS             istmp;
27735d8aa7fSBarry Smith   PetscTruth     found;
2788965ea79SLois Curfman McInnes 
2793a40ed3dSBarry Smith   PetscFunctionBegin;
280b9b97703SBarry Smith   ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr);
2818965ea79SLois Curfman McInnes   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
2828965ea79SLois Curfman McInnes 
2838965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
284b0a32e0cSBarry Smith   ierr  = PetscMalloc(2*size*sizeof(int),&nprocs);CHKERRQ(ierr);
285549d3d68SSatish Balay   ierr  = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr);
286b0a32e0cSBarry Smith   ierr  = PetscMalloc((N+1)*sizeof(int),&owner);CHKERRQ(ierr); /* see note*/
2878965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
2888965ea79SLois Curfman McInnes     idx = rows[i];
28935d8aa7fSBarry Smith     found = PETSC_FALSE;
2908965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
2918965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
292c1dc657dSBarry Smith         nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
2938965ea79SLois Curfman McInnes       }
2948965ea79SLois Curfman McInnes     }
29529bbc08cSBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
2968965ea79SLois Curfman McInnes   }
297c1dc657dSBarry Smith   nsends = 0;  for (i=0; i<size; i++) { nsends += nprocs[2*i+1];}
2988965ea79SLois Curfman McInnes 
2998965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
300c1dc657dSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
3018965ea79SLois Curfman McInnes 
3028965ea79SLois Curfman McInnes   /* post receives:   */
303b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int),&rvalues);CHKERRQ(ierr);
304b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr);
3058965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
306ca161407SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3078965ea79SLois Curfman McInnes   }
3088965ea79SLois Curfman McInnes 
3098965ea79SLois Curfman McInnes   /* do sends:
3108965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3118965ea79SLois Curfman McInnes          the ith processor
3128965ea79SLois Curfman McInnes   */
313b0a32e0cSBarry Smith   ierr = PetscMalloc((N+1)*sizeof(int),&svalues);CHKERRQ(ierr);
314b0a32e0cSBarry Smith   ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr);
315b0a32e0cSBarry Smith   ierr = PetscMalloc((size+1)*sizeof(int),&starts);CHKERRQ(ierr);
3168965ea79SLois Curfman McInnes   starts[0]  = 0;
317c1dc657dSBarry Smith   for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3188965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3198965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3208965ea79SLois Curfman McInnes   }
3218965ea79SLois Curfman McInnes   ISRestoreIndices(is,&rows);
3228965ea79SLois Curfman McInnes 
3238965ea79SLois Curfman McInnes   starts[0] = 0;
324c1dc657dSBarry Smith   for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3258965ea79SLois Curfman McInnes   count = 0;
3268965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
327c1dc657dSBarry Smith     if (nprocs[2*i+1]) {
328c1dc657dSBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3298965ea79SLois Curfman McInnes     }
3308965ea79SLois Curfman McInnes   }
331606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
3328965ea79SLois Curfman McInnes 
3338965ea79SLois Curfman McInnes   base = owners[rank];
3348965ea79SLois Curfman McInnes 
3358965ea79SLois Curfman McInnes   /*  wait on receives */
336b0a32e0cSBarry Smith   ierr   = PetscMalloc(2*(nrecvs+1)*sizeof(int),&lens);CHKERRQ(ierr);
3378965ea79SLois Curfman McInnes   source = lens + nrecvs;
3388965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
3398965ea79SLois Curfman McInnes   while (count) {
340ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
3418965ea79SLois Curfman McInnes     /* unpack receives into our local space */
342ca161407SBarry Smith     ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr);
3438965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
3448965ea79SLois Curfman McInnes     lens[imdex]  = n;
3458965ea79SLois Curfman McInnes     slen += n;
3468965ea79SLois Curfman McInnes     count--;
3478965ea79SLois Curfman McInnes   }
348606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
3498965ea79SLois Curfman McInnes 
3508965ea79SLois Curfman McInnes   /* move the data into the send scatter */
351b0a32e0cSBarry Smith   ierr = PetscMalloc((slen+1)*sizeof(int),&lrows);CHKERRQ(ierr);
3528965ea79SLois Curfman McInnes   count = 0;
3538965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
3548965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
3558965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
3568965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
3578965ea79SLois Curfman McInnes     }
3588965ea79SLois Curfman McInnes   }
359606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
360606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
361606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
362606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
3638965ea79SLois Curfman McInnes 
3648965ea79SLois Curfman McInnes   /* actually zap the local rows */
365029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr);
366b0a32e0cSBarry Smith   PetscLogObjectParent(A,istmp);
367606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
3688965ea79SLois Curfman McInnes   ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr);
3698965ea79SLois Curfman McInnes   ierr = ISDestroy(istmp);CHKERRQ(ierr);
3708965ea79SLois Curfman McInnes 
3718965ea79SLois Curfman McInnes   /* wait on sends */
3728965ea79SLois Curfman McInnes   if (nsends) {
373b0a32e0cSBarry Smith     ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
374ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
375606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
3768965ea79SLois Curfman McInnes   }
377606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
378606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
3798965ea79SLois Curfman McInnes 
3803a40ed3dSBarry Smith   PetscFunctionReturn(0);
3818965ea79SLois Curfman McInnes }
3828965ea79SLois Curfman McInnes 
3834a2ae208SSatish Balay #undef __FUNCT__
3844a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense"
3858f6be9afSLois Curfman McInnes int MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
3868965ea79SLois Curfman McInnes {
38739ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
3888965ea79SLois Curfman McInnes   int          ierr;
389c456f294SBarry Smith 
3903a40ed3dSBarry Smith   PetscFunctionBegin;
39143a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
39243a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
39344cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
3943a40ed3dSBarry Smith   PetscFunctionReturn(0);
3958965ea79SLois Curfman McInnes }
3968965ea79SLois Curfman McInnes 
3974a2ae208SSatish Balay #undef __FUNCT__
3984a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense"
3998f6be9afSLois Curfman McInnes int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4008965ea79SLois Curfman McInnes {
40139ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
4028965ea79SLois Curfman McInnes   int          ierr;
403c456f294SBarry Smith 
4043a40ed3dSBarry Smith   PetscFunctionBegin;
40543a90d84SBarry Smith   ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
40643a90d84SBarry Smith   ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
40744cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4083a40ed3dSBarry Smith   PetscFunctionReturn(0);
4098965ea79SLois Curfman McInnes }
4108965ea79SLois Curfman McInnes 
4114a2ae208SSatish Balay #undef __FUNCT__
4124a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense"
4137c922b88SBarry Smith int MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
414096963f5SLois Curfman McInnes {
415096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
416096963f5SLois Curfman McInnes   int          ierr;
41787828ca2SBarry Smith   PetscScalar  zero = 0.0;
418096963f5SLois Curfman McInnes 
4193a40ed3dSBarry Smith   PetscFunctionBegin;
4203501a2bdSLois Curfman McInnes   ierr = VecSet(&zero,yy);CHKERRQ(ierr);
4217c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
422537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
423537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4243a40ed3dSBarry Smith   PetscFunctionReturn(0);
425096963f5SLois Curfman McInnes }
426096963f5SLois Curfman McInnes 
4274a2ae208SSatish Balay #undef __FUNCT__
4284a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense"
4297c922b88SBarry Smith int MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
430096963f5SLois Curfman McInnes {
431096963f5SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
432096963f5SLois Curfman McInnes   int          ierr;
433096963f5SLois Curfman McInnes 
4343a40ed3dSBarry Smith   PetscFunctionBegin;
4353501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
4367c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
437537820f0SBarry Smith   ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
438537820f0SBarry Smith   ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr);
4393a40ed3dSBarry Smith   PetscFunctionReturn(0);
440096963f5SLois Curfman McInnes }
441096963f5SLois Curfman McInnes 
4424a2ae208SSatish Balay #undef __FUNCT__
4434a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense"
4448f6be9afSLois Curfman McInnes int MatGetDiagonal_MPIDense(Mat A,Vec v)
4458965ea79SLois Curfman McInnes {
44639ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
447096963f5SLois Curfman McInnes   Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data;
448273d9f13SBarry Smith   int          ierr,len,i,n,m = A->m,radd;
44987828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
450ed3cc1f0SBarry Smith 
4513a40ed3dSBarry Smith   PetscFunctionBegin;
452273d9f13SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
453b1d4fb26SBarry Smith   ierr = VecGetArrayFast(v,&x);CHKERRQ(ierr);
454096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
455273d9f13SBarry Smith   if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
456273d9f13SBarry Smith   len  = PetscMin(a->A->m,a->A->n);
4577ddc982cSLois Curfman McInnes   radd = a->rstart*m;
45844cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
459096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
460096963f5SLois Curfman McInnes   }
461b1d4fb26SBarry Smith   ierr = VecRestoreArrayFast(v,&x);CHKERRQ(ierr);
4623a40ed3dSBarry Smith   PetscFunctionReturn(0);
4638965ea79SLois Curfman McInnes }
4648965ea79SLois Curfman McInnes 
4654a2ae208SSatish Balay #undef __FUNCT__
4664a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense"
467e1311b90SBarry Smith int MatDestroy_MPIDense(Mat mat)
4688965ea79SLois Curfman McInnes {
4693501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
4708965ea79SLois Curfman McInnes   int          ierr;
471ed3cc1f0SBarry Smith 
4723a40ed3dSBarry Smith   PetscFunctionBegin;
47394d884c6SBarry Smith 
474aa482453SBarry Smith #if defined(PETSC_USE_LOG)
475b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mat->M,mat->N);
4768965ea79SLois Curfman McInnes #endif
4778798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
478606d414cSSatish Balay   ierr = PetscFree(mdn->rowners);CHKERRQ(ierr);
4793501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A);CHKERRQ(ierr);
4803501a2bdSLois Curfman McInnes   if (mdn->lvec)   VecDestroy(mdn->lvec);
4813501a2bdSLois Curfman McInnes   if (mdn->Mvctx)  VecScatterDestroy(mdn->Mvctx);
482622d7880SLois Curfman McInnes   if (mdn->factor) {
483606d414cSSatish Balay     if (mdn->factor->temp)   {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);}
484606d414cSSatish Balay     if (mdn->factor->tag)    {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);}
485606d414cSSatish Balay     if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);}
486606d414cSSatish Balay     ierr = PetscFree(mdn->factor);CHKERRQ(ierr);
487622d7880SLois Curfman McInnes   }
488606d414cSSatish Balay   ierr = PetscFree(mdn);CHKERRQ(ierr);
4893a40ed3dSBarry Smith   PetscFunctionReturn(0);
4908965ea79SLois Curfman McInnes }
49139ddd567SLois Curfman McInnes 
4924a2ae208SSatish Balay #undef __FUNCT__
4934a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary"
494b0a32e0cSBarry Smith static int MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
4958965ea79SLois Curfman McInnes {
49639ddd567SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
4978965ea79SLois Curfman McInnes   int          ierr;
4987056b6fcSBarry Smith 
4993a40ed3dSBarry Smith   PetscFunctionBegin;
50039ddd567SLois Curfman McInnes   if (mdn->size == 1) {
50139ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5028965ea79SLois Curfman McInnes   }
50329bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported");
5043a40ed3dSBarry Smith   PetscFunctionReturn(0);
5058965ea79SLois Curfman McInnes }
5068965ea79SLois Curfman McInnes 
5074a2ae208SSatish Balay #undef __FUNCT__
5084a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket"
509b0a32e0cSBarry Smith static int MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
5108965ea79SLois Curfman McInnes {
51139ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
512fb9695e5SSatish Balay   int               ierr,size = mdn->size,rank = mdn->rank;
513b0a32e0cSBarry Smith   PetscViewerType   vtype;
514f1af5d2fSBarry Smith   PetscTruth        isascii,isdraw;
515b0a32e0cSBarry Smith   PetscViewer       sviewer;
516f3ef73ceSBarry Smith   PetscViewerFormat format;
5178965ea79SLois Curfman McInnes 
5183a40ed3dSBarry Smith   PetscFunctionBegin;
519b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
520fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
521f1af5d2fSBarry Smith   if (isascii) {
522b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
523b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
524456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
5254e220ebcSLois Curfman McInnes       MatInfo info;
526888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
527b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mat->m,
5286831982aSBarry Smith                    (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);CHKERRQ(ierr);
529b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5303501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
5313a40ed3dSBarry Smith       PetscFunctionReturn(0);
532fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
5333a40ed3dSBarry Smith       PetscFunctionReturn(0);
5348965ea79SLois Curfman McInnes     }
535f1af5d2fSBarry Smith   } else if (isdraw) {
536b0a32e0cSBarry Smith     PetscDraw       draw;
537f1af5d2fSBarry Smith     PetscTruth isnull;
538f1af5d2fSBarry Smith 
539b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
540b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
541f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
542f1af5d2fSBarry Smith   }
54377ed5343SBarry Smith 
5448965ea79SLois Curfman McInnes   if (size == 1) {
54539ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
5463a40ed3dSBarry Smith   } else {
5478965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
5488965ea79SLois Curfman McInnes     Mat          A;
549273d9f13SBarry Smith     int          M = mat->M,N = mat->N,m,row,i,nz,*cols;
55087828ca2SBarry Smith     PetscScalar  *vals;
5518965ea79SLois Curfman McInnes 
5528965ea79SLois Curfman McInnes     if (!rank) {
553f1af5d2fSBarry Smith       ierr = MatCreateMPIDense(mat->comm,M,N,M,N,PETSC_NULL,&A);CHKERRQ(ierr);
5543a40ed3dSBarry Smith     } else {
555f1af5d2fSBarry Smith       ierr = MatCreateMPIDense(mat->comm,0,0,M,N,PETSC_NULL,&A);CHKERRQ(ierr);
5568965ea79SLois Curfman McInnes     }
557b0a32e0cSBarry Smith     PetscLogObjectParent(mat,A);
5588965ea79SLois Curfman McInnes 
55939ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
56039ddd567SLois Curfman McInnes        but it's quick for now */
561273d9f13SBarry Smith     row = mdn->rstart; m = mdn->A->m;
5628965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
56339ddd567SLois Curfman McInnes       ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
56439ddd567SLois Curfman McInnes       ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
56539ddd567SLois Curfman McInnes       ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
56639ddd567SLois Curfman McInnes       row++;
5678965ea79SLois Curfman McInnes     }
5688965ea79SLois Curfman McInnes 
5696d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5706d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
571b0a32e0cSBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
572b9b97703SBarry Smith     if (!rank) {
5736831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
5748965ea79SLois Curfman McInnes     }
575b0a32e0cSBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
576b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5778965ea79SLois Curfman McInnes     ierr = MatDestroy(A);CHKERRQ(ierr);
5788965ea79SLois Curfman McInnes   }
5793a40ed3dSBarry Smith   PetscFunctionReturn(0);
5808965ea79SLois Curfman McInnes }
5818965ea79SLois Curfman McInnes 
5824a2ae208SSatish Balay #undef __FUNCT__
5834a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense"
584b0a32e0cSBarry Smith int MatView_MPIDense(Mat mat,PetscViewer viewer)
5858965ea79SLois Curfman McInnes {
58639ddd567SLois Curfman McInnes   int        ierr;
587f1af5d2fSBarry Smith   PetscTruth isascii,isbinary,isdraw,issocket;
5888965ea79SLois Curfman McInnes 
589433994e6SBarry Smith   PetscFunctionBegin;
5900f5bd95cSBarry Smith 
591b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
592fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
593b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
594fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
5950f5bd95cSBarry Smith 
596f1af5d2fSBarry Smith   if (isascii || issocket || isdraw) {
597f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
5980f5bd95cSBarry Smith   } else if (isbinary) {
5993a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
6005cd90555SBarry Smith   } else {
60129bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
6028965ea79SLois Curfman McInnes   }
6033a40ed3dSBarry Smith   PetscFunctionReturn(0);
6048965ea79SLois Curfman McInnes }
6058965ea79SLois Curfman McInnes 
6064a2ae208SSatish Balay #undef __FUNCT__
6074a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense"
6088f6be9afSLois Curfman McInnes int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
6098965ea79SLois Curfman McInnes {
6103501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
6113501a2bdSLois Curfman McInnes   Mat          mdn = mat->A;
6124e220ebcSLois Curfman McInnes   int          ierr;
613329f5518SBarry Smith   PetscReal    isend[5],irecv[5];
6148965ea79SLois Curfman McInnes 
6153a40ed3dSBarry Smith   PetscFunctionBegin;
616273d9f13SBarry Smith   info->rows_global    = (double)A->M;
617273d9f13SBarry Smith   info->columns_global = (double)A->N;
618273d9f13SBarry Smith   info->rows_local     = (double)A->m;
619273d9f13SBarry Smith   info->columns_local  = (double)A->N;
6204e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
6214e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
6224e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
6234e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
6248965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
6254e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
6264e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
6274e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
6284e220ebcSLois Curfman McInnes     info->memory       = isend[3];
6294e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
6308965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
631d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
6324e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
6334e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
6344e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
6354e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
6364e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
6378965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
638d7d1e502SBarry Smith     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,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   }
6454e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
6464e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
6474e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
6483a40ed3dSBarry Smith   PetscFunctionReturn(0);
6498965ea79SLois Curfman McInnes }
6508965ea79SLois Curfman McInnes 
6514a2ae208SSatish Balay #undef __FUNCT__
6524a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense"
6538f6be9afSLois Curfman McInnes int MatSetOption_MPIDense(Mat A,MatOption op)
6548965ea79SLois Curfman McInnes {
65539ddd567SLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
656273d9f13SBarry Smith   int          ierr;
6578965ea79SLois Curfman McInnes 
6583a40ed3dSBarry Smith   PetscFunctionBegin;
65912c028f9SKris Buschelman   switch (op) {
66012c028f9SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
66112c028f9SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
66212c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
66312c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
66412c028f9SKris Buschelman   case MAT_COLUMNS_SORTED:
66512c028f9SKris Buschelman   case MAT_COLUMNS_UNSORTED:
666273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
66712c028f9SKris Buschelman     break;
66812c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
669273d9f13SBarry Smith     a->roworiented = PETSC_TRUE;
670273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
67112c028f9SKris Buschelman     break;
67212c028f9SKris Buschelman   case MAT_ROWS_SORTED:
67312c028f9SKris Buschelman   case MAT_ROWS_UNSORTED:
67412c028f9SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
67512c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
676b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_MPIDense:Option ignored\n");
67712c028f9SKris Buschelman     break;
67812c028f9SKris Buschelman   case MAT_COLUMN_ORIENTED:
679273d9f13SBarry Smith     a->roworiented = PETSC_FALSE;
680273d9f13SBarry Smith     ierr = MatSetOption(a->A,op);CHKERRQ(ierr);
68112c028f9SKris Buschelman     break;
68212c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
683273d9f13SBarry Smith     a->donotstash = PETSC_TRUE;
68412c028f9SKris Buschelman     break;
68512c028f9SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
68629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
68777e54ba9SKris Buschelman   case MAT_SYMMETRIC:
68877e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
6899a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
6909a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
6919a4540c5SBarry Smith   case MAT_HERMITIAN:
6929a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
6939a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
6949a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
69577e54ba9SKris Buschelman     break;
69612c028f9SKris Buschelman   default:
69729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
6983a40ed3dSBarry Smith   }
6993a40ed3dSBarry Smith   PetscFunctionReturn(0);
7008965ea79SLois Curfman McInnes }
7018965ea79SLois Curfman McInnes 
7024a2ae208SSatish Balay #undef __FUNCT__
7034a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIDense"
70487828ca2SBarry Smith int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,PetscScalar **v)
7058965ea79SLois Curfman McInnes {
7063501a2bdSLois Curfman McInnes   Mat_MPIDense *mat = (Mat_MPIDense*)A->data;
7073a40ed3dSBarry Smith   int          lrow,rstart = mat->rstart,rend = mat->rend,ierr;
7088965ea79SLois Curfman McInnes 
7093a40ed3dSBarry Smith   PetscFunctionBegin;
71029bbc08cSBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows")
7118965ea79SLois Curfman McInnes   lrow = row - rstart;
7123a40ed3dSBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr);
7133a40ed3dSBarry Smith   PetscFunctionReturn(0);
7148965ea79SLois Curfman McInnes }
7158965ea79SLois Curfman McInnes 
7164a2ae208SSatish Balay #undef __FUNCT__
7174a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIDense"
71887828ca2SBarry Smith int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,PetscScalar **v)
7198965ea79SLois Curfman McInnes {
720606d414cSSatish Balay   int ierr;
721606d414cSSatish Balay 
7223a40ed3dSBarry Smith   PetscFunctionBegin;
723606d414cSSatish Balay   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
724606d414cSSatish Balay   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
7253a40ed3dSBarry Smith   PetscFunctionReturn(0);
7268965ea79SLois Curfman McInnes }
7278965ea79SLois Curfman McInnes 
7284a2ae208SSatish Balay #undef __FUNCT__
7294a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense"
7305b2fa520SLois Curfman McInnes int MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
7315b2fa520SLois Curfman McInnes {
7325b2fa520SLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7335b2fa520SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
73487828ca2SBarry Smith   PetscScalar  *l,*r,x,*v;
735273d9f13SBarry Smith   int          ierr,i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n;
7365b2fa520SLois Curfman McInnes 
7375b2fa520SLois Curfman McInnes   PetscFunctionBegin;
73872d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
7395b2fa520SLois Curfman McInnes   if (ll) {
74072d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
74129bbc08cSBarry Smith     if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size");
742b1d4fb26SBarry Smith     ierr = VecGetArrayFast(ll,&l);CHKERRQ(ierr);
7435b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
7445b2fa520SLois Curfman McInnes       x = l[i];
7455b2fa520SLois Curfman McInnes       v = mat->v + i;
7465b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
7475b2fa520SLois Curfman McInnes     }
748b1d4fb26SBarry Smith     ierr = VecRestoreArrayFast(ll,&l);CHKERRQ(ierr);
749b0a32e0cSBarry Smith     PetscLogFlops(n*m);
7505b2fa520SLois Curfman McInnes   }
7515b2fa520SLois Curfman McInnes   if (rr) {
75272d926a5SLois Curfman McInnes     ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr);
75329bbc08cSBarry Smith     if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size");
7545b2fa520SLois Curfman McInnes     ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
7555b2fa520SLois Curfman McInnes     ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr);
756b1d4fb26SBarry Smith     ierr = VecGetArrayFast(mdn->lvec,&r);CHKERRQ(ierr);
7575b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
7585b2fa520SLois Curfman McInnes       x = r[i];
7595b2fa520SLois Curfman McInnes       v = mat->v + i*m;
7605b2fa520SLois Curfman McInnes       for (j=0; j<m; j++) { (*v++) *= x;}
7615b2fa520SLois Curfman McInnes     }
762b1d4fb26SBarry Smith     ierr = VecRestoreArrayFast(mdn->lvec,&r);CHKERRQ(ierr);
763b0a32e0cSBarry Smith     PetscLogFlops(n*m);
7645b2fa520SLois Curfman McInnes   }
7655b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
7665b2fa520SLois Curfman McInnes }
7675b2fa520SLois Curfman McInnes 
7684a2ae208SSatish Balay #undef __FUNCT__
7694a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense"
770064f8208SBarry Smith int MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
771096963f5SLois Curfman McInnes {
7723501a2bdSLois Curfman McInnes   Mat_MPIDense *mdn = (Mat_MPIDense*)A->data;
7733501a2bdSLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data;
7743501a2bdSLois Curfman McInnes   int          ierr,i,j;
775329f5518SBarry Smith   PetscReal    sum = 0.0;
77687828ca2SBarry Smith   PetscScalar  *v = mat->v;
7773501a2bdSLois Curfman McInnes 
7783a40ed3dSBarry Smith   PetscFunctionBegin;
7793501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
780064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
7813501a2bdSLois Curfman McInnes   } else {
7823501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
783273d9f13SBarry Smith       for (i=0; i<mdn->A->n*mdn->A->m; i++) {
784aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
785329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
7863501a2bdSLois Curfman McInnes #else
7873501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
7883501a2bdSLois Curfman McInnes #endif
7893501a2bdSLois Curfman McInnes       }
790064f8208SBarry Smith       ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
791064f8208SBarry Smith       *nrm = sqrt(*nrm);
792b0a32e0cSBarry Smith       PetscLogFlops(2*mdn->A->n*mdn->A->m);
7933a40ed3dSBarry Smith     } else if (type == NORM_1) {
794329f5518SBarry Smith       PetscReal *tmp,*tmp2;
795b0a32e0cSBarry Smith       ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
796273d9f13SBarry Smith       tmp2 = tmp + A->N;
797273d9f13SBarry Smith       ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr);
798064f8208SBarry Smith       *nrm = 0.0;
7993501a2bdSLois Curfman McInnes       v = mat->v;
800273d9f13SBarry Smith       for (j=0; j<mdn->A->n; j++) {
801273d9f13SBarry Smith         for (i=0; i<mdn->A->m; i++) {
80267e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
8033501a2bdSLois Curfman McInnes         }
8043501a2bdSLois Curfman McInnes       }
805d7d1e502SBarry Smith       ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr);
806273d9f13SBarry Smith       for (j=0; j<A->N; j++) {
807064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
8083501a2bdSLois Curfman McInnes       }
809606d414cSSatish Balay       ierr = PetscFree(tmp);CHKERRQ(ierr);
810b0a32e0cSBarry Smith       PetscLogFlops(A->n*A->m);
8113a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
812329f5518SBarry Smith       PetscReal ntemp;
8133501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
814064f8208SBarry Smith       ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr);
8153a40ed3dSBarry Smith     } else {
81629bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"No support for two norm");
8173501a2bdSLois Curfman McInnes     }
8183501a2bdSLois Curfman McInnes   }
8193a40ed3dSBarry Smith   PetscFunctionReturn(0);
8203501a2bdSLois Curfman McInnes }
8213501a2bdSLois Curfman McInnes 
8224a2ae208SSatish Balay #undef __FUNCT__
8234a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense"
8248f6be9afSLois Curfman McInnes int MatTranspose_MPIDense(Mat A,Mat *matout)
8253501a2bdSLois Curfman McInnes {
8263501a2bdSLois Curfman McInnes   Mat_MPIDense *a = (Mat_MPIDense*)A->data;
8273501a2bdSLois Curfman McInnes   Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data;
8283501a2bdSLois Curfman McInnes   Mat          B;
829273d9f13SBarry Smith   int          M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart;
8303501a2bdSLois Curfman McInnes   int          j,i,ierr;
83187828ca2SBarry Smith   PetscScalar  *v;
8323501a2bdSLois Curfman McInnes 
8333a40ed3dSBarry Smith   PetscFunctionBegin;
8347c922b88SBarry Smith   if (!matout && M != N) {
83529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place");
8367056b6fcSBarry Smith   }
8377056b6fcSBarry Smith   ierr = MatCreateMPIDense(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,PETSC_NULL,&B);CHKERRQ(ierr);
8383501a2bdSLois Curfman McInnes 
839273d9f13SBarry Smith   m = a->A->m; n = a->A->n; v = Aloc->v;
840b0a32e0cSBarry Smith   ierr = PetscMalloc(n*sizeof(int),&rwork);CHKERRQ(ierr);
8413501a2bdSLois Curfman McInnes   for (j=0; j<n; j++) {
8423501a2bdSLois Curfman McInnes     for (i=0; i<m; i++) rwork[i] = rstart + i;
8433501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
8443501a2bdSLois Curfman McInnes     v   += m;
8453501a2bdSLois Curfman McInnes   }
846606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
8476d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8486d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8497c922b88SBarry Smith   if (matout) {
8503501a2bdSLois Curfman McInnes     *matout = B;
8513501a2bdSLois Curfman McInnes   } else {
852273d9f13SBarry Smith     ierr = MatHeaderCopy(A,B);CHKERRQ(ierr);
8533501a2bdSLois Curfman McInnes   }
8543a40ed3dSBarry Smith   PetscFunctionReturn(0);
855096963f5SLois Curfman McInnes }
856096963f5SLois Curfman McInnes 
857d9eff348SSatish Balay #include "petscblaslapack.h"
8584a2ae208SSatish Balay #undef __FUNCT__
8594a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense"
860268466fbSBarry Smith int MatScale_MPIDense(const PetscScalar *alpha,Mat inA)
86144cd7ae7SLois Curfman McInnes {
86244cd7ae7SLois Curfman McInnes   Mat_MPIDense *A = (Mat_MPIDense*)inA->data;
86344cd7ae7SLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->A->data;
86444cd7ae7SLois Curfman McInnes   int          one = 1,nz;
86544cd7ae7SLois Curfman McInnes 
8663a40ed3dSBarry Smith   PetscFunctionBegin;
867273d9f13SBarry Smith   nz = inA->m*inA->N;
868268466fbSBarry Smith   BLscal_(&nz,(PetscScalar*)alpha,a->v,&one);
869b0a32e0cSBarry Smith   PetscLogFlops(nz);
8703a40ed3dSBarry Smith   PetscFunctionReturn(0);
87144cd7ae7SLois Curfman McInnes }
87244cd7ae7SLois Curfman McInnes 
8735609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
8748965ea79SLois Curfman McInnes 
8754a2ae208SSatish Balay #undef __FUNCT__
8764a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense"
877273d9f13SBarry Smith int MatSetUpPreallocation_MPIDense(Mat A)
878273d9f13SBarry Smith {
879273d9f13SBarry Smith   int        ierr;
880273d9f13SBarry Smith 
881273d9f13SBarry Smith   PetscFunctionBegin;
882273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
883273d9f13SBarry Smith   PetscFunctionReturn(0);
884273d9f13SBarry Smith }
885273d9f13SBarry Smith 
8868965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
88709dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
88809dc0095SBarry Smith        MatGetRow_MPIDense,
88909dc0095SBarry Smith        MatRestoreRow_MPIDense,
89009dc0095SBarry Smith        MatMult_MPIDense,
89197304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense,
8927c922b88SBarry Smith        MatMultTranspose_MPIDense,
8937c922b88SBarry Smith        MatMultTransposeAdd_MPIDense,
8948965ea79SLois Curfman McInnes        0,
89509dc0095SBarry Smith        0,
89609dc0095SBarry Smith        0,
89797304618SKris Buschelman /*10*/ 0,
89809dc0095SBarry Smith        0,
89909dc0095SBarry Smith        0,
90009dc0095SBarry Smith        0,
90109dc0095SBarry Smith        MatTranspose_MPIDense,
90297304618SKris Buschelman /*15*/ MatGetInfo_MPIDense,
90397304618SKris Buschelman        0,
90409dc0095SBarry Smith        MatGetDiagonal_MPIDense,
9055b2fa520SLois Curfman McInnes        MatDiagonalScale_MPIDense,
90609dc0095SBarry Smith        MatNorm_MPIDense,
90797304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense,
90809dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
90909dc0095SBarry Smith        0,
91009dc0095SBarry Smith        MatSetOption_MPIDense,
91109dc0095SBarry Smith        MatZeroEntries_MPIDense,
91297304618SKris Buschelman /*25*/ MatZeroRows_MPIDense,
91309dc0095SBarry Smith        0,
91409dc0095SBarry Smith        0,
91509dc0095SBarry Smith        0,
91609dc0095SBarry Smith        0,
91797304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense,
918273d9f13SBarry Smith        0,
91909dc0095SBarry Smith        0,
92009dc0095SBarry Smith        MatGetArray_MPIDense,
92109dc0095SBarry Smith        MatRestoreArray_MPIDense,
92297304618SKris Buschelman /*35*/ MatDuplicate_MPIDense,
92309dc0095SBarry Smith        0,
92409dc0095SBarry Smith        0,
92509dc0095SBarry Smith        0,
92609dc0095SBarry Smith        0,
92797304618SKris Buschelman /*40*/ 0,
9282ce60cd0SSatish Balay        MatGetSubMatrices_MPIDense,
92909dc0095SBarry Smith        0,
93009dc0095SBarry Smith        MatGetValues_MPIDense,
93109dc0095SBarry Smith        0,
93297304618SKris Buschelman /*45*/ 0,
93309dc0095SBarry Smith        MatScale_MPIDense,
93409dc0095SBarry Smith        0,
93509dc0095SBarry Smith        0,
93609dc0095SBarry Smith        0,
93797304618SKris Buschelman /*50*/ MatGetBlockSize_MPIDense,
93809dc0095SBarry Smith        0,
93909dc0095SBarry Smith        0,
94009dc0095SBarry Smith        0,
94109dc0095SBarry Smith        0,
94297304618SKris Buschelman /*55*/ 0,
94309dc0095SBarry Smith        0,
94409dc0095SBarry Smith        0,
94509dc0095SBarry Smith        0,
94609dc0095SBarry Smith        0,
94797304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense,
948b9b97703SBarry Smith        MatDestroy_MPIDense,
949b9b97703SBarry Smith        MatView_MPIDense,
95097304618SKris Buschelman        MatGetPetscMaps_Petsc,
95197304618SKris Buschelman        0,
95297304618SKris Buschelman /*65*/ 0,
95397304618SKris Buschelman        0,
95497304618SKris Buschelman        0,
95597304618SKris Buschelman        0,
95697304618SKris Buschelman        0,
95797304618SKris Buschelman /*70*/ 0,
95897304618SKris Buschelman        0,
95997304618SKris Buschelman        0,
96097304618SKris Buschelman        0,
96197304618SKris Buschelman        0,
96297304618SKris Buschelman /*75*/ 0,
96397304618SKris Buschelman        0,
96497304618SKris Buschelman        0,
96597304618SKris Buschelman        0,
96697304618SKris Buschelman        0,
96797304618SKris Buschelman /*80*/ 0,
96897304618SKris Buschelman        0,
96997304618SKris Buschelman        0,
97097304618SKris Buschelman        0,
97197304618SKris Buschelman /*85*/ MatLoad_MPIDense};
9728965ea79SLois Curfman McInnes 
973273d9f13SBarry Smith EXTERN_C_BEGIN
9744a2ae208SSatish Balay #undef __FUNCT__
975a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense"
976a23d5eceSKris Buschelman int MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
977a23d5eceSKris Buschelman {
978a23d5eceSKris Buschelman   Mat_MPIDense *a;
979a23d5eceSKris Buschelman   int          ierr;
980a23d5eceSKris Buschelman 
981a23d5eceSKris Buschelman   PetscFunctionBegin;
982a23d5eceSKris Buschelman   mat->preallocated = PETSC_TRUE;
983a23d5eceSKris Buschelman   /* Note:  For now, when data is specified above, this assumes the user correctly
984a23d5eceSKris Buschelman    allocates the local dense storage space.  We should add error checking. */
985a23d5eceSKris Buschelman 
986a23d5eceSKris Buschelman   a    = (Mat_MPIDense*)mat->data;
987*5c5985e7SKris Buschelman   ierr = MatCreate(PETSC_COMM_SELF,mat->m,mat->N,mat->m,mat->N,&a->A);CHKERRQ(ierr);
988*5c5985e7SKris Buschelman   ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr);
989*5c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr);
990a23d5eceSKris Buschelman   PetscLogObjectParent(mat,a->A);
991a23d5eceSKris Buschelman   PetscFunctionReturn(0);
992a23d5eceSKris Buschelman }
993a23d5eceSKris Buschelman EXTERN_C_END
994a23d5eceSKris Buschelman 
9950bad9183SKris Buschelman /*MC
996fafad747SKris Buschelman    MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices.
9970bad9183SKris Buschelman 
9980bad9183SKris Buschelman    Options Database Keys:
9990bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions()
10000bad9183SKris Buschelman 
10010bad9183SKris Buschelman   Level: beginner
10020bad9183SKris Buschelman 
10030bad9183SKris Buschelman .seealso: MatCreateMPIDense
10040bad9183SKris Buschelman M*/
10050bad9183SKris Buschelman 
1006a23d5eceSKris Buschelman EXTERN_C_BEGIN
1007a23d5eceSKris Buschelman #undef __FUNCT__
10084a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense"
1009273d9f13SBarry Smith int MatCreate_MPIDense(Mat mat)
1010273d9f13SBarry Smith {
1011273d9f13SBarry Smith   Mat_MPIDense *a;
1012273d9f13SBarry Smith   int          ierr,i;
1013273d9f13SBarry Smith 
1014273d9f13SBarry Smith   PetscFunctionBegin;
1015b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
1016b0a32e0cSBarry Smith   mat->data         = (void*)a;
1017273d9f13SBarry Smith   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1018273d9f13SBarry Smith   mat->factor       = 0;
1019273d9f13SBarry Smith   mat->mapping      = 0;
1020273d9f13SBarry Smith 
1021273d9f13SBarry Smith   a->factor       = 0;
1022273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
1023273d9f13SBarry Smith   ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr);
1024273d9f13SBarry Smith   ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr);
1025273d9f13SBarry Smith 
1026273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr);
1027273d9f13SBarry Smith   ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr);
1028273d9f13SBarry Smith   a->nvec = mat->n;
1029273d9f13SBarry Smith 
1030273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
1031273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
10328a124369SBarry Smith   ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr);
10338a124369SBarry Smith   ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr);
1034273d9f13SBarry Smith 
1035273d9f13SBarry Smith   /* build local table of row and column ownerships */
1036b0a32e0cSBarry Smith   ierr       = PetscMalloc(2*(a->size+2)*sizeof(int),&a->rowners);CHKERRQ(ierr);
1037273d9f13SBarry Smith   a->cowners = a->rowners + a->size + 1;
1038b0a32e0cSBarry Smith   PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1039273d9f13SBarry Smith   ierr = MPI_Allgather(&mat->m,1,MPI_INT,a->rowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr);
1040273d9f13SBarry Smith   a->rowners[0] = 0;
1041273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1042273d9f13SBarry Smith     a->rowners[i] += a->rowners[i-1];
1043273d9f13SBarry Smith   }
1044273d9f13SBarry Smith   a->rstart = a->rowners[a->rank];
1045273d9f13SBarry Smith   a->rend   = a->rowners[a->rank+1];
1046273d9f13SBarry Smith   ierr      = MPI_Allgather(&mat->n,1,MPI_INT,a->cowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr);
1047273d9f13SBarry Smith   a->cowners[0] = 0;
1048273d9f13SBarry Smith   for (i=2; i<=a->size; i++) {
1049273d9f13SBarry Smith     a->cowners[i] += a->cowners[i-1];
1050273d9f13SBarry Smith   }
1051273d9f13SBarry Smith 
1052273d9f13SBarry Smith   /* build cache for off array entries formed */
1053273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
1054273d9f13SBarry Smith   ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr);
1055273d9f13SBarry Smith 
1056273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1057273d9f13SBarry Smith   a->lvec        = 0;
1058273d9f13SBarry Smith   a->Mvctx       = 0;
1059273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1060273d9f13SBarry Smith 
1061273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1062273d9f13SBarry Smith                                      "MatGetDiagonalBlock_MPIDense",
1063273d9f13SBarry Smith                                      MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1064a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C",
1065a23d5eceSKris Buschelman                                      "MatMPIDenseSetPreallocation_MPIDense",
1066a23d5eceSKris Buschelman                                      MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr);
1067273d9f13SBarry Smith   PetscFunctionReturn(0);
1068273d9f13SBarry Smith }
1069273d9f13SBarry Smith EXTERN_C_END
1070273d9f13SBarry Smith 
1071209238afSKris Buschelman /*MC
1072002d173eSKris Buschelman    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1073209238afSKris Buschelman 
1074209238afSKris Buschelman    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1075209238afSKris Buschelman    and MATMPIDENSE otherwise.
1076209238afSKris Buschelman 
1077209238afSKris Buschelman    Options Database Keys:
1078209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1079209238afSKris Buschelman 
1080209238afSKris Buschelman   Level: beginner
1081209238afSKris Buschelman 
1082209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1083209238afSKris Buschelman M*/
1084209238afSKris Buschelman 
1085209238afSKris Buschelman EXTERN_C_BEGIN
1086209238afSKris Buschelman #undef __FUNCT__
1087209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense"
1088209238afSKris Buschelman int MatCreate_Dense(Mat A) {
1089209238afSKris Buschelman   int ierr,size;
1090209238afSKris Buschelman 
1091209238afSKris Buschelman   PetscFunctionBegin;
1092209238afSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr);
1093209238afSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
1094209238afSKris Buschelman   if (size == 1) {
1095209238afSKris Buschelman     ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr);
1096209238afSKris Buschelman   } else {
1097209238afSKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
1098209238afSKris Buschelman   }
1099209238afSKris Buschelman   PetscFunctionReturn(0);
1100209238afSKris Buschelman }
1101209238afSKris Buschelman EXTERN_C_END
1102209238afSKris Buschelman 
11034a2ae208SSatish Balay #undef __FUNCT__
11044a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation"
1105273d9f13SBarry Smith /*@C
1106273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1107273d9f13SBarry Smith 
1108273d9f13SBarry Smith    Not collective
1109273d9f13SBarry Smith 
1110273d9f13SBarry Smith    Input Parameters:
1111273d9f13SBarry Smith .  A - the matrix
1112273d9f13SBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1113273d9f13SBarry Smith    to control all matrix memory allocation.
1114273d9f13SBarry Smith 
1115273d9f13SBarry Smith    Notes:
1116273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1117273d9f13SBarry Smith    storage by columns.
1118273d9f13SBarry Smith 
1119273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1120273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1121273d9f13SBarry Smith    set data=PETSC_NULL.
1122273d9f13SBarry Smith 
1123273d9f13SBarry Smith    Level: intermediate
1124273d9f13SBarry Smith 
1125273d9f13SBarry Smith .keywords: matrix,dense, parallel
1126273d9f13SBarry Smith 
1127273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1128273d9f13SBarry Smith @*/
112987828ca2SBarry Smith int MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data)
1130273d9f13SBarry Smith {
1131a23d5eceSKris Buschelman   int ierr,(*f)(Mat,PetscScalar *);
1132273d9f13SBarry Smith 
1133273d9f13SBarry Smith   PetscFunctionBegin;
1134565567f0SKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1135a23d5eceSKris Buschelman   if (f) {
1136a23d5eceSKris Buschelman     ierr = (*f)(mat,data);CHKERRQ(ierr);
1137a23d5eceSKris Buschelman   }
1138273d9f13SBarry Smith   PetscFunctionReturn(0);
1139273d9f13SBarry Smith }
1140273d9f13SBarry Smith 
11414a2ae208SSatish Balay #undef __FUNCT__
11424a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense"
11438965ea79SLois Curfman McInnes /*@C
114439ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
11458965ea79SLois Curfman McInnes 
1146db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1147db81eaa0SLois Curfman McInnes 
11488965ea79SLois Curfman McInnes    Input Parameters:
1149db81eaa0SLois Curfman McInnes +  comm - MPI communicator
11508965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1151db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
11528965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1153db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
11547f5ff6fdSBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1155dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
11568965ea79SLois Curfman McInnes 
11578965ea79SLois Curfman McInnes    Output Parameter:
1158477f1c0bSLois Curfman McInnes .  A - the matrix
11598965ea79SLois Curfman McInnes 
1160b259b22eSLois Curfman McInnes    Notes:
116139ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
116239ddd567SLois Curfman McInnes    storage by columns.
11638965ea79SLois Curfman McInnes 
116418f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
116518f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
11667f5ff6fdSBarry Smith    set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users).
116718f449edSLois Curfman McInnes 
11688965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
11698965ea79SLois Curfman McInnes    (possibly both).
11708965ea79SLois Curfman McInnes 
1171027ccd11SLois Curfman McInnes    Level: intermediate
1172027ccd11SLois Curfman McInnes 
117339ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
11748965ea79SLois Curfman McInnes 
117539ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
11768965ea79SLois Curfman McInnes @*/
117787828ca2SBarry Smith int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,PetscScalar *data,Mat *A)
11788965ea79SLois Curfman McInnes {
1179273d9f13SBarry Smith   int ierr,size;
11808965ea79SLois Curfman McInnes 
11813a40ed3dSBarry Smith   PetscFunctionBegin;
1182273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,M,N,A);CHKERRQ(ierr);
1183273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1184273d9f13SBarry Smith   if (size > 1) {
1185273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1186273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1187273d9f13SBarry Smith   } else {
1188273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1189273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
11908c469469SLois Curfman McInnes   }
11913a40ed3dSBarry Smith   PetscFunctionReturn(0);
11928965ea79SLois Curfman McInnes }
11938965ea79SLois Curfman McInnes 
11944a2ae208SSatish Balay #undef __FUNCT__
11954a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense"
11965609ef8eSBarry Smith static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
11978965ea79SLois Curfman McInnes {
11988965ea79SLois Curfman McInnes   Mat          mat;
11993501a2bdSLois Curfman McInnes   Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data;
120039ddd567SLois Curfman McInnes   int          ierr;
12018965ea79SLois Curfman McInnes 
12023a40ed3dSBarry Smith   PetscFunctionBegin;
12038965ea79SLois Curfman McInnes   *newmat       = 0;
1204273d9f13SBarry Smith   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&mat);CHKERRQ(ierr);
1205273d9f13SBarry Smith   ierr = MatSetType(mat,MATMPIDENSE);CHKERRQ(ierr);
1206b0a32e0cSBarry Smith   ierr              = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr);
1207b0a32e0cSBarry Smith   mat->data         = (void*)a;
1208549d3d68SSatish Balay   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
12093501a2bdSLois Curfman McInnes   mat->factor       = A->factor;
1210c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1211273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
12128965ea79SLois Curfman McInnes 
12138965ea79SLois Curfman McInnes   a->rstart       = oldmat->rstart;
12148965ea79SLois Curfman McInnes   a->rend         = oldmat->rend;
12158965ea79SLois Curfman McInnes   a->size         = oldmat->size;
12168965ea79SLois Curfman McInnes   a->rank         = oldmat->rank;
1217e0fa3b82SLois Curfman McInnes   mat->insertmode = NOT_SET_VALUES;
1218b9b97703SBarry Smith   a->nvec         = oldmat->nvec;
12193782ba37SSatish Balay   a->donotstash   = oldmat->donotstash;
1220b0a32e0cSBarry Smith   ierr            = PetscMalloc((a->size+1)*sizeof(int),&a->rowners);CHKERRQ(ierr);
1221b0a32e0cSBarry Smith   PetscLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense));
1222549d3d68SSatish Balay   ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr);
12238798bf22SSatish Balay   ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr);
12248965ea79SLois Curfman McInnes 
1225329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
12265609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
1227b0a32e0cSBarry Smith   PetscLogObjectParent(mat,a->A);
12288965ea79SLois Curfman McInnes   *newmat = mat;
12293a40ed3dSBarry Smith   PetscFunctionReturn(0);
12308965ea79SLois Curfman McInnes }
12318965ea79SLois Curfman McInnes 
1232e090d566SSatish Balay #include "petscsys.h"
12338965ea79SLois Curfman McInnes 
12344a2ae208SSatish Balay #undef __FUNCT__
12354a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile"
123690ace30eSBarry Smith int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M,int N,Mat *newmat)
123790ace30eSBarry Smith {
123840011551SBarry Smith   int          *rowners,i,size,rank,m,ierr,nz,j;
123987828ca2SBarry Smith   PetscScalar  *array,*vals,*vals_ptr;
124090ace30eSBarry Smith   MPI_Status   status;
124190ace30eSBarry Smith 
12423a40ed3dSBarry Smith   PetscFunctionBegin;
1243d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1244d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
124590ace30eSBarry Smith 
124690ace30eSBarry Smith   /* determine ownership of all rows */
124790ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
1248b0a32e0cSBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr);
1249ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
125090ace30eSBarry Smith   rowners[0] = 0;
125190ace30eSBarry Smith   for (i=2; i<=size; i++) {
125290ace30eSBarry Smith     rowners[i] += rowners[i-1];
125390ace30eSBarry Smith   }
125490ace30eSBarry Smith 
125590ace30eSBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
125690ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr);
125790ace30eSBarry Smith 
125890ace30eSBarry Smith   if (!rank) {
125987828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
126090ace30eSBarry Smith 
126190ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
12620752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
126390ace30eSBarry Smith 
126490ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
126590ace30eSBarry Smith     vals_ptr = vals;
126690ace30eSBarry Smith     for (i=0; i<m; i++) {
126790ace30eSBarry Smith       for (j=0; j<N; j++) {
126890ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
126990ace30eSBarry Smith       }
127090ace30eSBarry Smith     }
127190ace30eSBarry Smith 
127290ace30eSBarry Smith     /* read in other processors and ship out */
127390ace30eSBarry Smith     for (i=1; i<size; i++) {
127490ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
12750752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1276ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr);
127790ace30eSBarry Smith     }
12783a40ed3dSBarry Smith   } else {
127990ace30eSBarry Smith     /* receive numeric values */
128087828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
128190ace30eSBarry Smith 
128290ace30eSBarry Smith     /* receive message of values*/
1283ca161407SBarry Smith     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr);
128490ace30eSBarry Smith 
128590ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
128690ace30eSBarry Smith     vals_ptr = vals;
128790ace30eSBarry Smith     for (i=0; i<m; i++) {
128890ace30eSBarry Smith       for (j=0; j<N; j++) {
128990ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
129090ace30eSBarry Smith       }
129190ace30eSBarry Smith     }
129290ace30eSBarry Smith   }
1293606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
1294606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
12956d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12966d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12973a40ed3dSBarry Smith   PetscFunctionReturn(0);
129890ace30eSBarry Smith }
129990ace30eSBarry Smith 
13004a2ae208SSatish Balay #undef __FUNCT__
13014a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense"
13028e9aea5cSBarry Smith int MatLoad_MPIDense(PetscViewer viewer,const MatType type,Mat *newmat)
13038965ea79SLois Curfman McInnes {
13048965ea79SLois Curfman McInnes   Mat          A;
130587828ca2SBarry Smith   PetscScalar  *vals,*svals;
130619bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
13078965ea79SLois Curfman McInnes   MPI_Status   status;
13088965ea79SLois Curfman McInnes   int          header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols;
13098965ea79SLois Curfman McInnes   int          *ourlens,*sndcounts = 0,*procsnz = 0,*offlens,jj,*mycols,*smycols;
131019bcc07fSBarry Smith   int          tag = ((PetscObject)viewer)->tag;
13113a40ed3dSBarry Smith   int          i,nz,ierr,j,rstart,rend,fd;
13128965ea79SLois Curfman McInnes 
13133a40ed3dSBarry Smith   PetscFunctionBegin;
1314d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1315d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13168965ea79SLois Curfman McInnes   if (!rank) {
1317b0a32e0cSBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
13180752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
1319552e946dSBarry Smith     if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
13208965ea79SLois Curfman McInnes   }
13218965ea79SLois Curfman McInnes 
1322ca161407SBarry Smith   ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr);
132390ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
132490ace30eSBarry Smith 
132590ace30eSBarry Smith   /*
132690ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
132790ace30eSBarry Smith   */
132890ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
13293a40ed3dSBarry Smith     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,newmat);CHKERRQ(ierr);
13303a40ed3dSBarry Smith     PetscFunctionReturn(0);
133190ace30eSBarry Smith   }
133290ace30eSBarry Smith 
13338965ea79SLois Curfman McInnes   /* determine ownership of all rows */
13348965ea79SLois Curfman McInnes   m          = M/size + ((M % size) > rank);
1335b0a32e0cSBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr);
1336ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
13378965ea79SLois Curfman McInnes   rowners[0] = 0;
13388965ea79SLois Curfman McInnes   for (i=2; i<=size; i++) {
13398965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
13408965ea79SLois Curfman McInnes   }
13418965ea79SLois Curfman McInnes   rstart = rowners[rank];
13428965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
13438965ea79SLois Curfman McInnes 
13448965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
1345b0a32e0cSBarry Smith   ierr    = PetscMalloc(2*(rend-rstart)*sizeof(int),&ourlens);CHKERRQ(ierr);
13468965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
13478965ea79SLois Curfman McInnes   if (!rank) {
1348b0a32e0cSBarry Smith     ierr = PetscMalloc(M*sizeof(int),&rowlengths);CHKERRQ(ierr);
13490752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
1350b0a32e0cSBarry Smith     ierr = PetscMalloc(size*sizeof(int),&sndcounts);CHKERRQ(ierr);
13518965ea79SLois Curfman McInnes     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
1352ca161407SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
1353606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1354ca161407SBarry Smith   } else {
1355ca161407SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr);
13568965ea79SLois Curfman McInnes   }
13578965ea79SLois Curfman McInnes 
13588965ea79SLois Curfman McInnes   if (!rank) {
13598965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
1360b0a32e0cSBarry Smith     ierr = PetscMalloc(size*sizeof(int),&procsnz);CHKERRQ(ierr);
1361549d3d68SSatish Balay     ierr = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr);
13628965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
13638965ea79SLois Curfman McInnes       for (j=rowners[i]; j< rowners[i+1]; j++) {
13648965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
13658965ea79SLois Curfman McInnes       }
13668965ea79SLois Curfman McInnes     }
1367606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
13688965ea79SLois Curfman McInnes 
13698965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
13708965ea79SLois Curfman McInnes     maxnz = 0;
13718965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
13720452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
13738965ea79SLois Curfman McInnes     }
1374b0a32e0cSBarry Smith     ierr = PetscMalloc(maxnz*sizeof(int),&cols);CHKERRQ(ierr);
13758965ea79SLois Curfman McInnes 
13768965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
13778965ea79SLois Curfman McInnes     nz = procsnz[0];
1378b0a32e0cSBarry Smith     ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr);
13790752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
13808965ea79SLois Curfman McInnes 
13818965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
13828965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
13838965ea79SLois Curfman McInnes       nz   = procsnz[i];
13840752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
1385ca161407SBarry Smith       ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr);
13868965ea79SLois Curfman McInnes     }
1387606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
13883a40ed3dSBarry Smith   } else {
13898965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
13908965ea79SLois Curfman McInnes     nz = 0;
13918965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
13928965ea79SLois Curfman McInnes       nz += ourlens[i];
13938965ea79SLois Curfman McInnes     }
1394b0a32e0cSBarry Smith     ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr);
13958965ea79SLois Curfman McInnes 
13968965ea79SLois Curfman McInnes     /* receive message of column indices*/
1397ca161407SBarry Smith     ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr);
1398ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr);
139929bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
14008965ea79SLois Curfman McInnes   }
14018965ea79SLois Curfman McInnes 
14028965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
1403549d3d68SSatish Balay   ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr);
14048965ea79SLois Curfman McInnes   jj = 0;
14058965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14068965ea79SLois Curfman McInnes     for (j=0; j<ourlens[i]; j++) {
14078965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
14088965ea79SLois Curfman McInnes       jj++;
14098965ea79SLois Curfman McInnes     }
14108965ea79SLois Curfman McInnes   }
14118965ea79SLois Curfman McInnes 
14128965ea79SLois Curfman McInnes   /* create our matrix */
14138965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14148965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
14158965ea79SLois Curfman McInnes   }
1416b4fd4287SBarry Smith   ierr = MatCreateMPIDense(comm,m,PETSC_DECIDE,M,N,PETSC_NULL,newmat);CHKERRQ(ierr);
14178965ea79SLois Curfman McInnes   A = *newmat;
14188965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
14198965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
14208965ea79SLois Curfman McInnes   }
14218965ea79SLois Curfman McInnes 
14228965ea79SLois Curfman McInnes   if (!rank) {
142387828ca2SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
14248965ea79SLois Curfman McInnes 
14258965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
14268965ea79SLois Curfman McInnes     nz = procsnz[0];
14270752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
14288965ea79SLois Curfman McInnes 
14298965ea79SLois Curfman McInnes     /* insert into matrix */
14308965ea79SLois Curfman McInnes     jj      = rstart;
14318965ea79SLois Curfman McInnes     smycols = mycols;
14328965ea79SLois Curfman McInnes     svals   = vals;
14338965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
14348965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
14358965ea79SLois Curfman McInnes       smycols += ourlens[i];
14368965ea79SLois Curfman McInnes       svals   += ourlens[i];
14378965ea79SLois Curfman McInnes       jj++;
14388965ea79SLois Curfman McInnes     }
14398965ea79SLois Curfman McInnes 
14408965ea79SLois Curfman McInnes     /* read in other processors and ship out */
14418965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
14428965ea79SLois Curfman McInnes       nz   = procsnz[i];
14430752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
1444ca161407SBarry Smith       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr);
14458965ea79SLois Curfman McInnes     }
1446606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
14473a40ed3dSBarry Smith   } else {
14488965ea79SLois Curfman McInnes     /* receive numeric values */
144987828ca2SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
14508965ea79SLois Curfman McInnes 
14518965ea79SLois Curfman McInnes     /* receive message of values*/
1452ca161407SBarry Smith     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr);
1453ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
145429bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
14558965ea79SLois Curfman McInnes 
14568965ea79SLois Curfman McInnes     /* insert into matrix */
14578965ea79SLois Curfman McInnes     jj      = rstart;
14588965ea79SLois Curfman McInnes     smycols = mycols;
14598965ea79SLois Curfman McInnes     svals   = vals;
14608965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
14618965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
14628965ea79SLois Curfman McInnes       smycols += ourlens[i];
14638965ea79SLois Curfman McInnes       svals   += ourlens[i];
14648965ea79SLois Curfman McInnes       jj++;
14658965ea79SLois Curfman McInnes     }
14668965ea79SLois Curfman McInnes   }
1467606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1468606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1469606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1470606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
14718965ea79SLois Curfman McInnes 
14726d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14736d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14743a40ed3dSBarry Smith   PetscFunctionReturn(0);
14758965ea79SLois Curfman McInnes }
147690ace30eSBarry Smith 
147790ace30eSBarry Smith 
147890ace30eSBarry Smith 
147990ace30eSBarry Smith 
148090ace30eSBarry Smith 
1481