xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision aeccfd6fd46a34231537ae422040d581ee0e5930)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
3ed3cc1f0SBarry Smith /*
4ed3cc1f0SBarry Smith    Basic functions for basic parallel dense matrices.
5ed3cc1f0SBarry Smith */
6ed3cc1f0SBarry Smith 
704fea9ffSBarry Smith 
889ac3759SBarry Smith #include "src/mat/impls/dense/mpi/mpidense.h"    /*I   "petscmat.h"  I*/
904fea9ffSBarry Smith PetscInt Plapack_nprows,Plapack_npcols;
1004fea9ffSBarry Smith MPI_Comm Plapack_comm_2d;
118965ea79SLois Curfman McInnes 
12ba8c8a56SBarry Smith #undef __FUNCT__
13ab92ecdeSBarry Smith #define __FUNCT__ "MatDenseGetLocalMatrix"
14ab92ecdeSBarry Smith /*@
15ab92ecdeSBarry Smith 
16ab92ecdeSBarry Smith       MatDenseGetLocalMatrix - For a MATMPIDENSE or MATSEQDENSE matrix returns the sequential
17ab92ecdeSBarry Smith               matrix that represents the operator. For sequential matrices it returns itself.
18ab92ecdeSBarry Smith 
19ab92ecdeSBarry Smith     Input Parameter:
20ab92ecdeSBarry Smith .      A - the Seq or MPI dense matrix
21ab92ecdeSBarry Smith 
22ab92ecdeSBarry Smith     Output Parameter:
23ab92ecdeSBarry Smith .      B - the inner matrix
24ab92ecdeSBarry Smith 
258e6c10adSSatish Balay     Level: intermediate
268e6c10adSSatish Balay 
27ab92ecdeSBarry Smith @*/
28ab92ecdeSBarry Smith PetscErrorCode MatDenseGetLocalMatrix(Mat A,Mat *B)
29ab92ecdeSBarry Smith {
30ab92ecdeSBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
31ab92ecdeSBarry Smith   PetscErrorCode ierr;
32ab92ecdeSBarry Smith   PetscTruth     flg;
33ab92ecdeSBarry Smith 
34ab92ecdeSBarry Smith   PetscFunctionBegin;
35ab92ecdeSBarry Smith   ierr = PetscTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr);
36ab92ecdeSBarry Smith   if (flg) {
37ab92ecdeSBarry Smith     *B = mat->A;
38ab92ecdeSBarry Smith   } else {
39ab92ecdeSBarry Smith     *B = A;
40ab92ecdeSBarry Smith   }
41ab92ecdeSBarry Smith   PetscFunctionReturn(0);
42ab92ecdeSBarry Smith }
43ab92ecdeSBarry Smith 
44ab92ecdeSBarry Smith #undef __FUNCT__
45ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense"
46ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
47ba8c8a56SBarry Smith {
48ba8c8a56SBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
49ba8c8a56SBarry Smith   PetscErrorCode ierr;
50899cda47SBarry Smith   PetscInt       lrow,rstart = A->rmap.rstart,rend = A->rmap.rend;
51ba8c8a56SBarry Smith 
52ba8c8a56SBarry Smith   PetscFunctionBegin;
53ba8c8a56SBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows")
54ba8c8a56SBarry Smith   lrow = row - rstart;
55ba8c8a56SBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt **)idx,(const PetscScalar **)v);CHKERRQ(ierr);
56ba8c8a56SBarry Smith   PetscFunctionReturn(0);
57ba8c8a56SBarry Smith }
58ba8c8a56SBarry Smith 
59ba8c8a56SBarry Smith #undef __FUNCT__
60ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense"
61ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
62ba8c8a56SBarry Smith {
63ba8c8a56SBarry Smith   PetscErrorCode ierr;
64ba8c8a56SBarry Smith 
65ba8c8a56SBarry Smith   PetscFunctionBegin;
66ba8c8a56SBarry Smith   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
67ba8c8a56SBarry Smith   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
68ba8c8a56SBarry Smith   PetscFunctionReturn(0);
69ba8c8a56SBarry Smith }
70ba8c8a56SBarry Smith 
710de54da6SSatish Balay EXTERN_C_BEGIN
724a2ae208SSatish Balay #undef __FUNCT__
734a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense"
74be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B)
750de54da6SSatish Balay {
760de54da6SSatish Balay   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
776849ba73SBarry Smith   PetscErrorCode ierr;
78899cda47SBarry Smith   PetscInt       m = A->rmap.n,rstart = A->rmap.rstart;
7987828ca2SBarry Smith   PetscScalar    *array;
800de54da6SSatish Balay   MPI_Comm       comm;
810de54da6SSatish Balay 
820de54da6SSatish Balay   PetscFunctionBegin;
83899cda47SBarry Smith   if (A->rmap.N != A->cmap.N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported.");
840de54da6SSatish Balay 
850de54da6SSatish Balay   /* The reuse aspect is not implemented efficiently */
860de54da6SSatish Balay   if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);}
870de54da6SSatish Balay 
880de54da6SSatish Balay   ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr);
890de54da6SSatish Balay   ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr);
90f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,B);CHKERRQ(ierr);
91f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*B,m,m,m,m);CHKERRQ(ierr);
927adad957SLisandro Dalcin   ierr = MatSetType(*B,((PetscObject)mdn->A)->type_name);CHKERRQ(ierr);
935c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr);
940de54da6SSatish Balay   ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr);
950de54da6SSatish Balay   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
960de54da6SSatish Balay   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
970de54da6SSatish Balay 
980de54da6SSatish Balay   *iscopy = PETSC_TRUE;
990de54da6SSatish Balay   PetscFunctionReturn(0);
1000de54da6SSatish Balay }
1010de54da6SSatish Balay EXTERN_C_END
1020de54da6SSatish Balay 
1034a2ae208SSatish Balay #undef __FUNCT__
1044a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense"
10513f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
1068965ea79SLois Curfman McInnes {
10739b7565bSBarry Smith   Mat_MPIDense   *A = (Mat_MPIDense*)mat->data;
108dfbe8321SBarry Smith   PetscErrorCode ierr;
109899cda47SBarry Smith   PetscInt       i,j,rstart = mat->rmap.rstart,rend = mat->rmap.rend,row;
110273d9f13SBarry Smith   PetscTruth     roworiented = A->roworiented;
1118965ea79SLois Curfman McInnes 
1123a40ed3dSBarry Smith   PetscFunctionBegin;
1138965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
1145ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
115899cda47SBarry Smith     if (idxm[i] >= mat->rmap.N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
1168965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
1178965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
11839b7565bSBarry Smith       if (roworiented) {
11939b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
1203a40ed3dSBarry Smith       } else {
1218965ea79SLois Curfman McInnes         for (j=0; j<n; j++) {
1225ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
123899cda47SBarry Smith           if (idxn[j] >= mat->cmap.N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
12439b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
12539b7565bSBarry Smith         }
1268965ea79SLois Curfman McInnes       }
1273a40ed3dSBarry Smith     } else {
1283782ba37SSatish Balay       if (!A->donotstash) {
12939b7565bSBarry Smith         if (roworiented) {
1308798bf22SSatish Balay           ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr);
131d36fbae8SSatish Balay         } else {
1328798bf22SSatish Balay           ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr);
13339b7565bSBarry Smith         }
134b49de8d1SLois Curfman McInnes       }
135b49de8d1SLois Curfman McInnes     }
1363782ba37SSatish Balay   }
1373a40ed3dSBarry Smith   PetscFunctionReturn(0);
138b49de8d1SLois Curfman McInnes }
139b49de8d1SLois Curfman McInnes 
1404a2ae208SSatish Balay #undef __FUNCT__
1414a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense"
14213f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
143b49de8d1SLois Curfman McInnes {
144b49de8d1SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
145dfbe8321SBarry Smith   PetscErrorCode ierr;
146899cda47SBarry Smith   PetscInt       i,j,rstart = mat->rmap.rstart,rend = mat->rmap.rend,row;
147b49de8d1SLois Curfman McInnes 
1483a40ed3dSBarry Smith   PetscFunctionBegin;
149b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
15097e567efSBarry Smith     if (idxm[i] < 0) continue; /* SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); */
151899cda47SBarry Smith     if (idxm[i] >= mat->rmap.N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
152b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
153b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
154b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
15597e567efSBarry Smith         if (idxn[j] < 0) continue; /* SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); */
156899cda47SBarry Smith         if (idxn[j] >= mat->cmap.N) {
15729bbc08cSBarry Smith           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
158a8c6a408SBarry Smith         }
159b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
160b49de8d1SLois Curfman McInnes       }
161a8c6a408SBarry Smith     } else {
16229bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"Only local values currently supported");
1638965ea79SLois Curfman McInnes     }
1648965ea79SLois Curfman McInnes   }
1653a40ed3dSBarry Smith   PetscFunctionReturn(0);
1668965ea79SLois Curfman McInnes }
1678965ea79SLois Curfman McInnes 
1684a2ae208SSatish Balay #undef __FUNCT__
1694a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense"
170dfbe8321SBarry Smith PetscErrorCode MatGetArray_MPIDense(Mat A,PetscScalar *array[])
171ff14e315SSatish Balay {
172ff14e315SSatish Balay   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
173dfbe8321SBarry Smith   PetscErrorCode ierr;
174ff14e315SSatish Balay 
1753a40ed3dSBarry Smith   PetscFunctionBegin;
176ff14e315SSatish Balay   ierr = MatGetArray(a->A,array);CHKERRQ(ierr);
1773a40ed3dSBarry Smith   PetscFunctionReturn(0);
178ff14e315SSatish Balay }
179ff14e315SSatish Balay 
1804a2ae208SSatish Balay #undef __FUNCT__
1814a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense"
18213f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
183ca3fa75bSLois Curfman McInnes {
184ca3fa75bSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data,*newmatd;
185ca3fa75bSLois Curfman McInnes   Mat_SeqDense   *lmat = (Mat_SeqDense*)mat->A->data;
1866849ba73SBarry Smith   PetscErrorCode ierr;
18713f74950SBarry Smith   PetscInt       i,j,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols;
18887828ca2SBarry Smith   PetscScalar    *av,*bv,*v = lmat->v;
189ca3fa75bSLois Curfman McInnes   Mat            newmat;
190ca3fa75bSLois Curfman McInnes 
191ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
192ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
193ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
194b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
195b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
196ca3fa75bSLois Curfman McInnes 
197ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
1987eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
1997eba5e9cSLois Curfman McInnes      original matrix! */
200ca3fa75bSLois Curfman McInnes 
201ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
2027eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
203ca3fa75bSLois Curfman McInnes 
204ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
205ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
20629bbc08cSBarry Smith     /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
2077eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
208ca3fa75bSLois Curfman McInnes     newmat = *B;
209ca3fa75bSLois Curfman McInnes   } else {
210ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
2117adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr);
212f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,cs,PETSC_DECIDE,ncols);CHKERRQ(ierr);
2137adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
214878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
215ca3fa75bSLois Curfman McInnes   }
216ca3fa75bSLois Curfman McInnes 
217ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
218ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
219ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense *)newmatd->A->data)->v;
220ca3fa75bSLois Curfman McInnes 
221ca3fa75bSLois Curfman McInnes   for (i=0; i<ncols; i++) {
222ca3fa75bSLois Curfman McInnes     av = v + nlrows*icol[i];
223ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
2247eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
225ca3fa75bSLois Curfman McInnes     }
226ca3fa75bSLois Curfman McInnes   }
227ca3fa75bSLois Curfman McInnes 
228ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
229ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
230ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
231ca3fa75bSLois Curfman McInnes 
232ca3fa75bSLois Curfman McInnes   /* Free work space */
233ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
234ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
235ca3fa75bSLois Curfman McInnes   *B = newmat;
236ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
237ca3fa75bSLois Curfman McInnes }
238ca3fa75bSLois Curfman McInnes 
2394a2ae208SSatish Balay #undef __FUNCT__
2404a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense"
241dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_MPIDense(Mat A,PetscScalar *array[])
242ff14e315SSatish Balay {
2433a40ed3dSBarry Smith   PetscFunctionBegin;
2443a40ed3dSBarry Smith   PetscFunctionReturn(0);
245ff14e315SSatish Balay }
246ff14e315SSatish Balay 
2474a2ae208SSatish Balay #undef __FUNCT__
2484a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense"
249dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
2508965ea79SLois Curfman McInnes {
25139ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
2527adad957SLisandro Dalcin   MPI_Comm       comm = ((PetscObject)mat)->comm;
253dfbe8321SBarry Smith   PetscErrorCode ierr;
25413f74950SBarry Smith   PetscInt       nstash,reallocs;
2558965ea79SLois Curfman McInnes   InsertMode     addv;
2568965ea79SLois Curfman McInnes 
2573a40ed3dSBarry Smith   PetscFunctionBegin;
2588965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
259ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
2607056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
26129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
2628965ea79SLois Curfman McInnes   }
263e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
2648965ea79SLois Curfman McInnes 
2651e2582c4SBarry Smith   ierr = MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap.range);CHKERRQ(ierr);
2668798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
267ae15b995SBarry Smith   ierr = PetscInfo2(mdn->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);CHKERRQ(ierr);
2683a40ed3dSBarry Smith   PetscFunctionReturn(0);
2698965ea79SLois Curfman McInnes }
2708965ea79SLois Curfman McInnes 
2714a2ae208SSatish Balay #undef __FUNCT__
2724a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense"
273dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2748965ea79SLois Curfman McInnes {
27539ddd567SLois Curfman McInnes   Mat_MPIDense    *mdn=(Mat_MPIDense*)mat->data;
2766849ba73SBarry Smith   PetscErrorCode  ierr;
27713f74950SBarry Smith   PetscInt        i,*row,*col,flg,j,rstart,ncols;
27813f74950SBarry Smith   PetscMPIInt     n;
27987828ca2SBarry Smith   PetscScalar     *val;
280e0fa3b82SLois Curfman McInnes   InsertMode      addv=mat->insertmode;
2818965ea79SLois Curfman McInnes 
2823a40ed3dSBarry Smith   PetscFunctionBegin;
2838965ea79SLois Curfman McInnes   /*  wait on receives */
2847ef1d9bdSSatish Balay   while (1) {
2858798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2867ef1d9bdSSatish Balay     if (!flg) break;
2878965ea79SLois Curfman McInnes 
2887ef1d9bdSSatish Balay     for (i=0; i<n;) {
2897ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2907ef1d9bdSSatish Balay       for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
2917ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2927ef1d9bdSSatish Balay       else       ncols = n-i;
2937ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
2947ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
2957ef1d9bdSSatish Balay       i = j;
2968965ea79SLois Curfman McInnes     }
2977ef1d9bdSSatish Balay   }
2988798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
2998965ea79SLois Curfman McInnes 
30039ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
30139ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
3028965ea79SLois Curfman McInnes 
3036d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
30439ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
3058965ea79SLois Curfman McInnes   }
3063a40ed3dSBarry Smith   PetscFunctionReturn(0);
3078965ea79SLois Curfman McInnes }
3088965ea79SLois Curfman McInnes 
3094a2ae208SSatish Balay #undef __FUNCT__
3104a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense"
311dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A)
3128965ea79SLois Curfman McInnes {
313dfbe8321SBarry Smith   PetscErrorCode ierr;
31439ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
3153a40ed3dSBarry Smith 
3163a40ed3dSBarry Smith   PetscFunctionBegin;
3173a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
3183a40ed3dSBarry Smith   PetscFunctionReturn(0);
3198965ea79SLois Curfman McInnes }
3208965ea79SLois Curfman McInnes 
3218965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
3228965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
3238965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
3243501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
3258965ea79SLois Curfman McInnes    routine.
3268965ea79SLois Curfman McInnes */
3274a2ae208SSatish Balay #undef __FUNCT__
3284a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense"
329f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
3308965ea79SLois Curfman McInnes {
33139ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
3326849ba73SBarry Smith   PetscErrorCode ierr;
333899cda47SBarry Smith   PetscInt       i,*owners = A->rmap.range;
33413f74950SBarry Smith   PetscInt       *nprocs,j,idx,nsends;
33513f74950SBarry Smith   PetscInt       nmax,*svalues,*starts,*owner,nrecvs;
3367adad957SLisandro Dalcin   PetscInt       *rvalues,tag = ((PetscObject)A)->tag,count,base,slen,*source;
33713f74950SBarry Smith   PetscInt       *lens,*lrows,*values;
33813f74950SBarry Smith   PetscMPIInt    n,imdex,rank = l->rank,size = l->size;
3397adad957SLisandro Dalcin   MPI_Comm       comm = ((PetscObject)A)->comm;
3408965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
3418965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
34235d8aa7fSBarry Smith   PetscTruth     found;
3438965ea79SLois Curfman McInnes 
3443a40ed3dSBarry Smith   PetscFunctionBegin;
3458965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
34613f74950SBarry Smith   ierr  = PetscMalloc(2*size*sizeof(PetscInt),&nprocs);CHKERRQ(ierr);
34713f74950SBarry Smith   ierr  = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr);
34813f74950SBarry Smith   ierr  = PetscMalloc((N+1)*sizeof(PetscInt),&owner);CHKERRQ(ierr); /* see note*/
3498965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3508965ea79SLois Curfman McInnes     idx = rows[i];
35135d8aa7fSBarry Smith     found = PETSC_FALSE;
3528965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
3538965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
354c1dc657dSBarry Smith         nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
3558965ea79SLois Curfman McInnes       }
3568965ea79SLois Curfman McInnes     }
35729bbc08cSBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
3588965ea79SLois Curfman McInnes   }
359c1dc657dSBarry Smith   nsends = 0;  for (i=0; i<size; i++) { nsends += nprocs[2*i+1];}
3608965ea79SLois Curfman McInnes 
3618965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
362c1dc657dSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
3638965ea79SLois Curfman McInnes 
3648965ea79SLois Curfman McInnes   /* post receives:   */
36513f74950SBarry Smith   ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);CHKERRQ(ierr);
366b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr);
3678965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
36813f74950SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3698965ea79SLois Curfman McInnes   }
3708965ea79SLois Curfman McInnes 
3718965ea79SLois Curfman McInnes   /* do sends:
3728965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3738965ea79SLois Curfman McInnes          the ith processor
3748965ea79SLois Curfman McInnes   */
37513f74950SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&svalues);CHKERRQ(ierr);
376b0a32e0cSBarry Smith   ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr);
37713f74950SBarry Smith   ierr = PetscMalloc((size+1)*sizeof(PetscInt),&starts);CHKERRQ(ierr);
3788965ea79SLois Curfman McInnes   starts[0]  = 0;
379c1dc657dSBarry Smith   for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3808965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3818965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3828965ea79SLois Curfman McInnes   }
3838965ea79SLois Curfman McInnes 
3848965ea79SLois Curfman McInnes   starts[0] = 0;
385c1dc657dSBarry Smith   for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3868965ea79SLois Curfman McInnes   count = 0;
3878965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
388c1dc657dSBarry Smith     if (nprocs[2*i+1]) {
38913f74950SBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3908965ea79SLois Curfman McInnes     }
3918965ea79SLois Curfman McInnes   }
392606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
3938965ea79SLois Curfman McInnes 
3948965ea79SLois Curfman McInnes   base = owners[rank];
3958965ea79SLois Curfman McInnes 
3968965ea79SLois Curfman McInnes   /*  wait on receives */
39713f74950SBarry Smith   ierr   = PetscMalloc(2*(nrecvs+1)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
3988965ea79SLois Curfman McInnes   source = lens + nrecvs;
3998965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
4008965ea79SLois Curfman McInnes   while (count) {
401ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
4028965ea79SLois Curfman McInnes     /* unpack receives into our local space */
40313f74950SBarry Smith     ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr);
4048965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
4058965ea79SLois Curfman McInnes     lens[imdex]    = n;
4068965ea79SLois Curfman McInnes     slen += n;
4078965ea79SLois Curfman McInnes     count--;
4088965ea79SLois Curfman McInnes   }
409606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
4108965ea79SLois Curfman McInnes 
4118965ea79SLois Curfman McInnes   /* move the data into the send scatter */
41213f74950SBarry Smith   ierr = PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);CHKERRQ(ierr);
4138965ea79SLois Curfman McInnes   count = 0;
4148965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
4158965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
4168965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
4178965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
4188965ea79SLois Curfman McInnes     }
4198965ea79SLois Curfman McInnes   }
420606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
421606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
422606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
423606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
4248965ea79SLois Curfman McInnes 
4258965ea79SLois Curfman McInnes   /* actually zap the local rows */
426f4df32b1SMatthew Knepley   ierr = MatZeroRows(l->A,slen,lrows,diag);CHKERRQ(ierr);
427606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
4288965ea79SLois Curfman McInnes 
4298965ea79SLois Curfman McInnes   /* wait on sends */
4308965ea79SLois Curfman McInnes   if (nsends) {
431b0a32e0cSBarry Smith     ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
432ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
433606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
4348965ea79SLois Curfman McInnes   }
435606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
436606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
4378965ea79SLois Curfman McInnes 
4383a40ed3dSBarry Smith   PetscFunctionReturn(0);
4398965ea79SLois Curfman McInnes }
4408965ea79SLois Curfman McInnes 
4414a2ae208SSatish Balay #undef __FUNCT__
4424a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense"
443dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
4448965ea79SLois Curfman McInnes {
44539ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
446dfbe8321SBarry Smith   PetscErrorCode ierr;
447c456f294SBarry Smith 
4483a40ed3dSBarry Smith   PetscFunctionBegin;
449ca9f406cSSatish Balay   ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
450ca9f406cSSatish Balay   ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
45144cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4523a40ed3dSBarry Smith   PetscFunctionReturn(0);
4538965ea79SLois Curfman McInnes }
4548965ea79SLois Curfman McInnes 
4554a2ae208SSatish Balay #undef __FUNCT__
4564a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense"
457dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4588965ea79SLois Curfman McInnes {
45939ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
460dfbe8321SBarry Smith   PetscErrorCode ierr;
461c456f294SBarry Smith 
4623a40ed3dSBarry Smith   PetscFunctionBegin;
463ca9f406cSSatish Balay   ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
464ca9f406cSSatish Balay   ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
46544cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4663a40ed3dSBarry Smith   PetscFunctionReturn(0);
4678965ea79SLois Curfman McInnes }
4688965ea79SLois Curfman McInnes 
4694a2ae208SSatish Balay #undef __FUNCT__
4704a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense"
471dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
472096963f5SLois Curfman McInnes {
473096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
474dfbe8321SBarry Smith   PetscErrorCode ierr;
47587828ca2SBarry Smith   PetscScalar    zero = 0.0;
476096963f5SLois Curfman McInnes 
4773a40ed3dSBarry Smith   PetscFunctionBegin;
4782dcb1b2aSMatthew Knepley   ierr = VecSet(yy,zero);CHKERRQ(ierr);
4797c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
480ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
481ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4823a40ed3dSBarry Smith   PetscFunctionReturn(0);
483096963f5SLois Curfman McInnes }
484096963f5SLois Curfman McInnes 
4854a2ae208SSatish Balay #undef __FUNCT__
4864a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense"
487dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
488096963f5SLois Curfman McInnes {
489096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
490dfbe8321SBarry Smith   PetscErrorCode ierr;
491096963f5SLois Curfman McInnes 
4923a40ed3dSBarry Smith   PetscFunctionBegin;
4933501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
4947c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
495ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
496ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4973a40ed3dSBarry Smith   PetscFunctionReturn(0);
498096963f5SLois Curfman McInnes }
499096963f5SLois Curfman McInnes 
5004a2ae208SSatish Balay #undef __FUNCT__
5014a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense"
502dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v)
5038965ea79SLois Curfman McInnes {
50439ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
505096963f5SLois Curfman McInnes   Mat_SeqDense   *aloc = (Mat_SeqDense*)a->A->data;
506dfbe8321SBarry Smith   PetscErrorCode ierr;
507899cda47SBarry Smith   PetscInt       len,i,n,m = A->rmap.n,radd;
50887828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
509ed3cc1f0SBarry Smith 
5103a40ed3dSBarry Smith   PetscFunctionBegin;
5112dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
5121ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
513096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
514899cda47SBarry Smith   if (n != A->rmap.N) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
515899cda47SBarry Smith   len  = PetscMin(a->A->rmap.n,a->A->cmap.n);
516899cda47SBarry Smith   radd = A->rmap.rstart*m;
51744cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
518096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
519096963f5SLois Curfman McInnes   }
5201ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
5213a40ed3dSBarry Smith   PetscFunctionReturn(0);
5228965ea79SLois Curfman McInnes }
5238965ea79SLois Curfman McInnes 
5244a2ae208SSatish Balay #undef __FUNCT__
5254a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense"
526dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat)
5278965ea79SLois Curfman McInnes {
5283501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
529dfbe8321SBarry Smith   PetscErrorCode ierr;
53001b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
53101b82886SBarry Smith   Mat_Plapack   *lu=(Mat_Plapack*)(mat->spptr);
53201b82886SBarry Smith #endif
533ed3cc1f0SBarry Smith 
5343a40ed3dSBarry Smith   PetscFunctionBegin;
53594d884c6SBarry Smith 
536aa482453SBarry Smith #if defined(PETSC_USE_LOG)
537899cda47SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->rmap.N,mat->cmap.N);
5388965ea79SLois Curfman McInnes #endif
5398798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
5403501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A);CHKERRQ(ierr);
54105b42c5fSBarry Smith   if (mdn->lvec)   {ierr = VecDestroy(mdn->lvec);CHKERRQ(ierr);}
54205b42c5fSBarry Smith   if (mdn->Mvctx)  {ierr = VecScatterDestroy(mdn->Mvctx);CHKERRQ(ierr);}
54301b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
54401b82886SBarry Smith   if (lu->CleanUpPlapack) {
54562b4c0b3SBarry Smith     ierr = PLA_Obj_free(&lu->A);CHKERRQ(ierr);
54662b4c0b3SBarry Smith     ierr = PLA_Obj_free (&lu->pivots);CHKERRQ(ierr);
54762b4c0b3SBarry Smith     ierr = PLA_Temp_free(&lu->templ);CHKERRQ(ierr);
54809d27a7eSBarry Smith 
54901b82886SBarry Smith 
55001b82886SBarry Smith     ierr = ISDestroy(lu->is_pla);CHKERRQ(ierr);
55101b82886SBarry Smith     ierr = ISDestroy(lu->is_petsc);CHKERRQ(ierr);
55201b82886SBarry Smith     ierr = VecScatterDestroy(lu->ctx);CHKERRQ(ierr);
553622d7880SLois Curfman McInnes   }
55401b82886SBarry Smith #endif
55501b82886SBarry Smith 
556606d414cSSatish Balay   ierr = PetscFree(mdn);CHKERRQ(ierr);
557dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
558901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);CHKERRQ(ierr);
559901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
5604ae313f4SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C","",PETSC_NULL);CHKERRQ(ierr);
5614ae313f4SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C","",PETSC_NULL);CHKERRQ(ierr);
5624ae313f4SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C","",PETSC_NULL);CHKERRQ(ierr);
5633a40ed3dSBarry Smith   PetscFunctionReturn(0);
5648965ea79SLois Curfman McInnes }
56539ddd567SLois Curfman McInnes 
5664a2ae208SSatish Balay #undef __FUNCT__
5674a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary"
5686849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
5698965ea79SLois Curfman McInnes {
57039ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
571dfbe8321SBarry Smith   PetscErrorCode    ierr;
572aa05aa95SBarry Smith   PetscViewerFormat format;
573aa05aa95SBarry Smith   int               fd;
574aa05aa95SBarry Smith   PetscInt          header[4],mmax,N = mat->cmap.N,i,j,m,k;
575aa05aa95SBarry Smith   PetscMPIInt       rank,tag  = ((PetscObject)viewer)->tag,size;
576578230a0SSatish Balay   PetscScalar       *work,*v,*vv;
577aa05aa95SBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)mdn->A->data;
578aa05aa95SBarry Smith   MPI_Status        status;
5797056b6fcSBarry Smith 
5803a40ed3dSBarry Smith   PetscFunctionBegin;
58139ddd567SLois Curfman McInnes   if (mdn->size == 1) {
58239ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
583aa05aa95SBarry Smith   } else {
584aa05aa95SBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
5857adad957SLisandro Dalcin     ierr = MPI_Comm_rank(((PetscObject)mat)->comm,&rank);CHKERRQ(ierr);
5867adad957SLisandro Dalcin     ierr = MPI_Comm_size(((PetscObject)mat)->comm,&size);CHKERRQ(ierr);
587aa05aa95SBarry Smith 
588aa05aa95SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
589aa05aa95SBarry Smith     if (format == PETSC_VIEWER_BINARY_NATIVE) {
590aa05aa95SBarry Smith 
591aa05aa95SBarry Smith       if (!rank) {
592aa05aa95SBarry Smith 	/* store the matrix as a dense matrix */
593aa05aa95SBarry Smith 	header[0] = MAT_FILE_COOKIE;
594aa05aa95SBarry Smith 	header[1] = mat->rmap.N;
595aa05aa95SBarry Smith 	header[2] = N;
596aa05aa95SBarry Smith 	header[3] = MATRIX_BINARY_FORMAT_DENSE;
597aa05aa95SBarry Smith 	ierr = PetscBinaryWrite(fd,header,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
598aa05aa95SBarry Smith 
599aa05aa95SBarry Smith 	/* get largest work array needed for transposing array */
600aa05aa95SBarry Smith         mmax = mat->rmap.n;
601aa05aa95SBarry Smith         for (i=1; i<size; i++) {
602aa05aa95SBarry Smith           mmax = PetscMax(mmax,mat->rmap.range[i+1] - mat->rmap.range[i]);
6038965ea79SLois Curfman McInnes         }
604aa05aa95SBarry Smith 	ierr = PetscMalloc(mmax*N*sizeof(PetscScalar),&work);CHKERRQ(ierr);
605aa05aa95SBarry Smith 
606aa05aa95SBarry Smith 	/* write out local array, by rows */
607aa05aa95SBarry Smith         m    = mat->rmap.n;
608aa05aa95SBarry Smith 	v    = a->v;
609aa05aa95SBarry Smith         for (j=0; j<N; j++) {
610aa05aa95SBarry Smith           for (i=0; i<m; i++) {
611578230a0SSatish Balay 	    work[j + i*N] = *v++;
612aa05aa95SBarry Smith 	  }
613aa05aa95SBarry Smith 	}
614aa05aa95SBarry Smith 	ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
615aa05aa95SBarry Smith         /* get largest work array to receive messages from other processes, excludes process zero */
616aa05aa95SBarry Smith         mmax = 0;
617aa05aa95SBarry Smith         for (i=1; i<size; i++) {
618aa05aa95SBarry Smith           mmax = PetscMax(mmax,mat->rmap.range[i+1] - mat->rmap.range[i]);
619aa05aa95SBarry Smith         }
620578230a0SSatish Balay 	ierr = PetscMalloc(mmax*N*sizeof(PetscScalar),&vv);CHKERRQ(ierr);
621578230a0SSatish Balay         v = vv;
622aa05aa95SBarry Smith         for (k=1; k<size; k++) {
623aa05aa95SBarry Smith           m    = mat->rmap.range[k+1] - mat->rmap.range[k];
6247adad957SLisandro Dalcin           ierr = MPI_Recv(v,m*N,MPIU_SCALAR,k,tag,((PetscObject)mat)->comm,&status);CHKERRQ(ierr);
625aa05aa95SBarry Smith 
626aa05aa95SBarry Smith           for (j=0; j<N; j++) {
627aa05aa95SBarry Smith             for (i=0; i<m; i++) {
628578230a0SSatish Balay               work[j + i*N] = *v++;
629aa05aa95SBarry Smith 	    }
630aa05aa95SBarry Smith 	  }
631aa05aa95SBarry Smith 	  ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
632aa05aa95SBarry Smith         }
633aa05aa95SBarry Smith         ierr = PetscFree(work);CHKERRQ(ierr);
634578230a0SSatish Balay         ierr = PetscFree(vv);CHKERRQ(ierr);
635aa05aa95SBarry Smith       } else {
6367adad957SLisandro Dalcin         ierr = MPI_Send(a->v,mat->rmap.n*mat->cmap.N,MPIU_SCALAR,0,tag,((PetscObject)mat)->comm);CHKERRQ(ierr);
637aa05aa95SBarry Smith       }
638aa05aa95SBarry Smith     }
639aa05aa95SBarry Smith   }
6403a40ed3dSBarry Smith   PetscFunctionReturn(0);
6418965ea79SLois Curfman McInnes }
6428965ea79SLois Curfman McInnes 
6434a2ae208SSatish Balay #undef __FUNCT__
6444a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket"
6456849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
6468965ea79SLois Curfman McInnes {
64739ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
648dfbe8321SBarry Smith   PetscErrorCode    ierr;
64932dcc486SBarry Smith   PetscMPIInt       size = mdn->size,rank = mdn->rank;
650b0a32e0cSBarry Smith   PetscViewerType   vtype;
65132077d6dSBarry Smith   PetscTruth        iascii,isdraw;
652b0a32e0cSBarry Smith   PetscViewer       sviewer;
653f3ef73ceSBarry Smith   PetscViewerFormat format;
65401b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
65501b82886SBarry Smith   Mat_Plapack       *lu=(Mat_Plapack*)(mat->spptr);
65601b82886SBarry Smith #endif
6578965ea79SLois Curfman McInnes 
6583a40ed3dSBarry Smith   PetscFunctionBegin;
65932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
660fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
66132077d6dSBarry Smith   if (iascii) {
662b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
663b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
664456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
6654e220ebcSLois Curfman McInnes       MatInfo info;
666888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
667899cda47SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->rmap.n,
66877431f27SBarry Smith                    (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr);
669b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6703501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
67101b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
67201b82886SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"PLAPACK run parameters:\n");CHKERRQ(ierr);
67304fea9ffSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Processor mesh: nprows %d, npcols %d\n",Plapack_nprows, Plapack_npcols);CHKERRQ(ierr);
67401b82886SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Distr. block size nb: %d \n",lu->nb);CHKERRQ(ierr);
67501b82886SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Error checking: %d\n",lu->ierror);CHKERRQ(ierr);
67601b82886SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Algorithmic block size: %d\n",lu->nb_alg);CHKERRQ(ierr);
67701b82886SBarry Smith #endif
6783a40ed3dSBarry Smith       PetscFunctionReturn(0);
679fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
6803a40ed3dSBarry Smith       PetscFunctionReturn(0);
6818965ea79SLois Curfman McInnes     }
682f1af5d2fSBarry Smith   } else if (isdraw) {
683b0a32e0cSBarry Smith     PetscDraw  draw;
684f1af5d2fSBarry Smith     PetscTruth isnull;
685f1af5d2fSBarry Smith 
686b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
687b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
688f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
689f1af5d2fSBarry Smith   }
69077ed5343SBarry Smith 
6918965ea79SLois Curfman McInnes   if (size == 1) {
69239ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
6933a40ed3dSBarry Smith   } else {
6948965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
6958965ea79SLois Curfman McInnes     Mat         A;
696899cda47SBarry Smith     PetscInt    M = mat->rmap.N,N = mat->cmap.N,m,row,i,nz;
697ba8c8a56SBarry Smith     PetscInt    *cols;
698ba8c8a56SBarry Smith     PetscScalar *vals;
6998965ea79SLois Curfman McInnes 
7007adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)mat)->comm,&A);CHKERRQ(ierr);
7018965ea79SLois Curfman McInnes     if (!rank) {
702f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr);
7033a40ed3dSBarry Smith     } else {
704f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr);
7058965ea79SLois Curfman McInnes     }
7067adad957SLisandro Dalcin     /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */
707878740d9SKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
708878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL);
70952e6d16bSBarry Smith     ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr);
7108965ea79SLois Curfman McInnes 
71139ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
71239ddd567SLois Curfman McInnes        but it's quick for now */
71351022da4SBarry Smith     A->insertmode = INSERT_VALUES;
714899cda47SBarry Smith     row = mat->rmap.rstart; m = mdn->A->rmap.n;
7158965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
716ba8c8a56SBarry Smith       ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
717ba8c8a56SBarry Smith       ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
718ba8c8a56SBarry Smith       ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
71939ddd567SLois Curfman McInnes       row++;
7208965ea79SLois Curfman McInnes     }
7218965ea79SLois Curfman McInnes 
7226d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7236d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
724b0a32e0cSBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
725b9b97703SBarry Smith     if (!rank) {
7266831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
7278965ea79SLois Curfman McInnes     }
728b0a32e0cSBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
729b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
7308965ea79SLois Curfman McInnes     ierr = MatDestroy(A);CHKERRQ(ierr);
7318965ea79SLois Curfman McInnes   }
7323a40ed3dSBarry Smith   PetscFunctionReturn(0);
7338965ea79SLois Curfman McInnes }
7348965ea79SLois Curfman McInnes 
7354a2ae208SSatish Balay #undef __FUNCT__
7364a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense"
737dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
7388965ea79SLois Curfman McInnes {
739dfbe8321SBarry Smith   PetscErrorCode ierr;
74032077d6dSBarry Smith   PetscTruth     iascii,isbinary,isdraw,issocket;
7418965ea79SLois Curfman McInnes 
742433994e6SBarry Smith   PetscFunctionBegin;
7430f5bd95cSBarry Smith 
74432077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
745fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
746b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
747fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
7480f5bd95cSBarry Smith 
74932077d6dSBarry Smith   if (iascii || issocket || isdraw) {
750f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
7510f5bd95cSBarry Smith   } else if (isbinary) {
7523a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
7535cd90555SBarry Smith   } else {
754958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
7558965ea79SLois Curfman McInnes   }
7563a40ed3dSBarry Smith   PetscFunctionReturn(0);
7578965ea79SLois Curfman McInnes }
7588965ea79SLois Curfman McInnes 
7594a2ae208SSatish Balay #undef __FUNCT__
7604a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense"
761dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
7628965ea79SLois Curfman McInnes {
7633501a2bdSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
7643501a2bdSLois Curfman McInnes   Mat            mdn = mat->A;
765dfbe8321SBarry Smith   PetscErrorCode ierr;
766329f5518SBarry Smith   PetscReal      isend[5],irecv[5];
7678965ea79SLois Curfman McInnes 
7683a40ed3dSBarry Smith   PetscFunctionBegin;
769899cda47SBarry Smith   info->rows_global    = (double)A->rmap.N;
770899cda47SBarry Smith   info->columns_global = (double)A->cmap.N;
771899cda47SBarry Smith   info->rows_local     = (double)A->rmap.n;
772899cda47SBarry Smith   info->columns_local  = (double)A->cmap.N;
7734e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
7744e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
7754e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
7764e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
7778965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
7784e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
7794e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
7804e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
7814e220ebcSLois Curfman McInnes     info->memory       = isend[3];
7824e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
7838965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
7847adad957SLisandro Dalcin     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,((PetscObject)A)->comm);CHKERRQ(ierr);
7854e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7864e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7874e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7884e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7894e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7908965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
7917adad957SLisandro Dalcin     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
7924e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7934e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7944e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7954e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7964e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7978965ea79SLois Curfman McInnes   }
7984e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
7994e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
8004e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
8013a40ed3dSBarry Smith   PetscFunctionReturn(0);
8028965ea79SLois Curfman McInnes }
8038965ea79SLois Curfman McInnes 
8044a2ae208SSatish Balay #undef __FUNCT__
8054a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense"
8064e0d8c25SBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op,PetscTruth flg)
8078965ea79SLois Curfman McInnes {
80839ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
809dfbe8321SBarry Smith   PetscErrorCode ierr;
8108965ea79SLois Curfman McInnes 
8113a40ed3dSBarry Smith   PetscFunctionBegin;
81212c028f9SKris Buschelman   switch (op) {
813512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
81412c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
81512c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
8164e0d8c25SBarry Smith     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
81712c028f9SKris Buschelman     break;
81812c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
8194e0d8c25SBarry Smith     a->roworiented = flg;
8204e0d8c25SBarry Smith     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
82112c028f9SKris Buschelman     break;
8224e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
82312c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
824290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
82512c028f9SKris Buschelman     break;
82612c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
8274e0d8c25SBarry Smith     a->donotstash = flg;
82812c028f9SKris Buschelman     break;
82977e54ba9SKris Buschelman   case MAT_SYMMETRIC:
83077e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
8319a4540c5SBarry Smith   case MAT_HERMITIAN:
8329a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
833600fe468SBarry Smith   case MAT_IGNORE_LOWER_TRIANGULAR:
834290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
83577e54ba9SKris Buschelman     break;
83612c028f9SKris Buschelman   default:
837600fe468SBarry Smith     SETERRQ1(PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
8383a40ed3dSBarry Smith   }
8393a40ed3dSBarry Smith   PetscFunctionReturn(0);
8408965ea79SLois Curfman McInnes }
8418965ea79SLois Curfman McInnes 
8428965ea79SLois Curfman McInnes 
8434a2ae208SSatish Balay #undef __FUNCT__
8444a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense"
845dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
8465b2fa520SLois Curfman McInnes {
8475b2fa520SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8485b2fa520SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
84987828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
850dfbe8321SBarry Smith   PetscErrorCode ierr;
851899cda47SBarry Smith   PetscInt       i,j,s2a,s3a,s2,s3,m=mdn->A->rmap.n,n=mdn->A->cmap.n;
8525b2fa520SLois Curfman McInnes 
8535b2fa520SLois Curfman McInnes   PetscFunctionBegin;
85472d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
8555b2fa520SLois Curfman McInnes   if (ll) {
85672d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
857175be7b4SMatthew Knepley     if (s2a != s2) SETERRQ2(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2);
8581ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
8595b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
8605b2fa520SLois Curfman McInnes       x = l[i];
8615b2fa520SLois Curfman McInnes       v = mat->v + i;
8625b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
8635b2fa520SLois Curfman McInnes     }
8641ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
865efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8665b2fa520SLois Curfman McInnes   }
8675b2fa520SLois Curfman McInnes   if (rr) {
868175be7b4SMatthew Knepley     ierr = VecGetLocalSize(rr,&s3a);CHKERRQ(ierr);
869175be7b4SMatthew Knepley     if (s3a != s3) SETERRQ2(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3);
870ca9f406cSSatish Balay     ierr = VecScatterBegin(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
871ca9f406cSSatish Balay     ierr = VecScatterEnd(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8721ebc52fbSHong Zhang     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
8735b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
8745b2fa520SLois Curfman McInnes       x = r[i];
8755b2fa520SLois Curfman McInnes       v = mat->v + i*m;
8765b2fa520SLois Curfman McInnes       for (j=0; j<m; j++) { (*v++) *= x;}
8775b2fa520SLois Curfman McInnes     }
8781ebc52fbSHong Zhang     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
879efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8805b2fa520SLois Curfman McInnes   }
8815b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
8825b2fa520SLois Curfman McInnes }
8835b2fa520SLois Curfman McInnes 
8844a2ae208SSatish Balay #undef __FUNCT__
8854a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense"
886dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
887096963f5SLois Curfman McInnes {
8883501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8893501a2bdSLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
890dfbe8321SBarry Smith   PetscErrorCode ierr;
89113f74950SBarry Smith   PetscInt       i,j;
892329f5518SBarry Smith   PetscReal      sum = 0.0;
89387828ca2SBarry Smith   PetscScalar    *v = mat->v;
8943501a2bdSLois Curfman McInnes 
8953a40ed3dSBarry Smith   PetscFunctionBegin;
8963501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
897064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
8983501a2bdSLois Curfman McInnes   } else {
8993501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
900899cda47SBarry Smith       for (i=0; i<mdn->A->cmap.n*mdn->A->rmap.n; i++) {
901aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
902329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
9033501a2bdSLois Curfman McInnes #else
9043501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
9053501a2bdSLois Curfman McInnes #endif
9063501a2bdSLois Curfman McInnes       }
9077adad957SLisandro Dalcin       ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
908064f8208SBarry Smith       *nrm = sqrt(*nrm);
909899cda47SBarry Smith       ierr = PetscLogFlops(2*mdn->A->cmap.n*mdn->A->rmap.n);CHKERRQ(ierr);
9103a40ed3dSBarry Smith     } else if (type == NORM_1) {
911329f5518SBarry Smith       PetscReal *tmp,*tmp2;
912899cda47SBarry Smith       ierr = PetscMalloc(2*A->cmap.N*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
913899cda47SBarry Smith       tmp2 = tmp + A->cmap.N;
914899cda47SBarry Smith       ierr = PetscMemzero(tmp,2*A->cmap.N*sizeof(PetscReal));CHKERRQ(ierr);
915064f8208SBarry Smith       *nrm = 0.0;
9163501a2bdSLois Curfman McInnes       v = mat->v;
917899cda47SBarry Smith       for (j=0; j<mdn->A->cmap.n; j++) {
918899cda47SBarry Smith         for (i=0; i<mdn->A->rmap.n; i++) {
91967e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
9203501a2bdSLois Curfman McInnes         }
9213501a2bdSLois Curfman McInnes       }
9227adad957SLisandro Dalcin       ierr = MPI_Allreduce(tmp,tmp2,A->cmap.N,MPIU_REAL,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
923899cda47SBarry Smith       for (j=0; j<A->cmap.N; j++) {
924064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
9253501a2bdSLois Curfman McInnes       }
926606d414cSSatish Balay       ierr = PetscFree(tmp);CHKERRQ(ierr);
927899cda47SBarry Smith       ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr);
9283a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
929329f5518SBarry Smith       PetscReal ntemp;
9303501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
9317adad957SLisandro Dalcin       ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,((PetscObject)A)->comm);CHKERRQ(ierr);
9323a40ed3dSBarry Smith     } else {
93329bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"No support for two norm");
9343501a2bdSLois Curfman McInnes     }
9353501a2bdSLois Curfman McInnes   }
9363a40ed3dSBarry Smith   PetscFunctionReturn(0);
9373501a2bdSLois Curfman McInnes }
9383501a2bdSLois Curfman McInnes 
9394a2ae208SSatish Balay #undef __FUNCT__
9404a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense"
941fc4dec0aSBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,MatReuse reuse,Mat *matout)
9423501a2bdSLois Curfman McInnes {
9433501a2bdSLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
9443501a2bdSLois Curfman McInnes   Mat_SeqDense   *Aloc = (Mat_SeqDense*)a->A->data;
9453501a2bdSLois Curfman McInnes   Mat            B;
946899cda47SBarry Smith   PetscInt       M = A->rmap.N,N = A->cmap.N,m,n,*rwork,rstart = A->rmap.rstart;
9476849ba73SBarry Smith   PetscErrorCode ierr;
94813f74950SBarry Smith   PetscInt       j,i;
94987828ca2SBarry Smith   PetscScalar    *v;
9503501a2bdSLois Curfman McInnes 
9513a40ed3dSBarry Smith   PetscFunctionBegin;
952fc4dec0aSBarry Smith   if (A == *matout && M != N) {
95329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place");
9547056b6fcSBarry Smith   }
955fc4dec0aSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || A == *matout) {
9567adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
957f69a0ea3SMatthew Knepley     ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,M);CHKERRQ(ierr);
9587adad957SLisandro Dalcin     ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
959878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr);
960fc4dec0aSBarry Smith   } else {
961fc4dec0aSBarry Smith     B = *matout;
962fc4dec0aSBarry Smith   }
9633501a2bdSLois Curfman McInnes 
964899cda47SBarry Smith   m = a->A->rmap.n; n = a->A->cmap.n; v = Aloc->v;
9651acff37aSSatish Balay   ierr = PetscMalloc(m*sizeof(PetscInt),&rwork);CHKERRQ(ierr);
9663501a2bdSLois Curfman McInnes   for (i=0; i<m; i++) rwork[i] = rstart + i;
9671acff37aSSatish Balay   for (j=0; j<n; j++) {
9683501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
9693501a2bdSLois Curfman McInnes     v   += m;
9703501a2bdSLois Curfman McInnes   }
971606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
9726d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9736d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
974fc4dec0aSBarry Smith   if (*matout != A) {
9753501a2bdSLois Curfman McInnes     *matout = B;
9763501a2bdSLois Curfman McInnes   } else {
977273d9f13SBarry Smith     ierr = MatHeaderCopy(A,B);CHKERRQ(ierr);
9783501a2bdSLois Curfman McInnes   }
9793a40ed3dSBarry Smith   PetscFunctionReturn(0);
980096963f5SLois Curfman McInnes }
981096963f5SLois Curfman McInnes 
982d9eff348SSatish Balay #include "petscblaslapack.h"
9834a2ae208SSatish Balay #undef __FUNCT__
9844a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense"
985f4df32b1SMatthew Knepley PetscErrorCode MatScale_MPIDense(Mat inA,PetscScalar alpha)
98644cd7ae7SLois Curfman McInnes {
98744cd7ae7SLois Curfman McInnes   Mat_MPIDense   *A = (Mat_MPIDense*)inA->data;
98844cd7ae7SLois Curfman McInnes   Mat_SeqDense   *a = (Mat_SeqDense*)A->A->data;
989f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
990efee365bSSatish Balay   PetscErrorCode ierr;
9910805154bSBarry Smith   PetscBLASInt   one = 1,nz = PetscBLASIntCast(inA->rmap.n*inA->cmap.N);
99244cd7ae7SLois Curfman McInnes 
9933a40ed3dSBarry Smith   PetscFunctionBegin;
994f4df32b1SMatthew Knepley   BLASscal_(&nz,&oalpha,a->v,&one);
995efee365bSSatish Balay   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
9963a40ed3dSBarry Smith   PetscFunctionReturn(0);
99744cd7ae7SLois Curfman McInnes }
99844cd7ae7SLois Curfman McInnes 
9996849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
10008965ea79SLois Curfman McInnes 
10014a2ae208SSatish Balay #undef __FUNCT__
10024a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense"
1003dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_MPIDense(Mat A)
1004273d9f13SBarry Smith {
1005dfbe8321SBarry Smith   PetscErrorCode ierr;
1006273d9f13SBarry Smith 
1007273d9f13SBarry Smith   PetscFunctionBegin;
1008273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
1009273d9f13SBarry Smith   PetscFunctionReturn(0);
1010273d9f13SBarry Smith }
1011273d9f13SBarry Smith 
10124ae313f4SHong Zhang #undef __FUNCT__
10134ae313f4SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIDense_MPIDense"
10144ae313f4SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
10154ae313f4SHong Zhang {
10164ae313f4SHong Zhang   PetscErrorCode ierr;
10174ae313f4SHong Zhang   PetscInt       m=A->rmap.n,n=B->cmap.n;
10184ae313f4SHong Zhang   Mat            Cmat;
10194ae313f4SHong Zhang 
10204ae313f4SHong Zhang   PetscFunctionBegin;
10214ae313f4SHong Zhang   if (A->cmap.n != B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"A->cmap.n %d != B->rmap.n %d\n",A->cmap.n,B->rmap.n);
10227adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)B)->comm,&Cmat);CHKERRQ(ierr);
10234ae313f4SHong Zhang   ierr = MatSetSizes(Cmat,m,n,A->rmap.N,B->cmap.N);CHKERRQ(ierr);
10244ae313f4SHong Zhang   ierr = MatSetType(Cmat,MATMPIDENSE);CHKERRQ(ierr);
10254ae313f4SHong Zhang   ierr = MatMPIDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
10262a6744ebSBarry Smith   ierr = MatAssemblyBegin(Cmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10270a71e23eSMatthew Knepley   ierr = MatAssemblyEnd(Cmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10284ae313f4SHong Zhang   *C = Cmat;
10294ae313f4SHong Zhang   PetscFunctionReturn(0);
10304ae313f4SHong Zhang }
10314ae313f4SHong Zhang 
103201b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
103301b82886SBarry Smith #undef __FUNCT__
103401b82886SBarry Smith #define __FUNCT__ "MatSolve_MPIDense"
103501b82886SBarry Smith PetscErrorCode MatSolve_MPIDense(Mat A,Vec b,Vec x)
103601b82886SBarry Smith {
103704fea9ffSBarry Smith   MPI_Comm       comm = ((PetscObject)A)->comm,dummy_comm;
103801b82886SBarry Smith   Mat_Plapack    *lu = (Mat_Plapack*)A->spptr;
103901b82886SBarry Smith   PetscErrorCode ierr;
104001b82886SBarry Smith   PetscInt       M=A->rmap.N,m=A->rmap.n,rstart,i,j,*idx_pla,*idx_petsc,loc_m,loc_stride;
104101b82886SBarry Smith   PetscScalar    *array;
104201b82886SBarry Smith   PetscReal      one = 1.0;
1043*aeccfd6fSBarry Smith   PetscMPIInt    size,r_rank,r_nproc,c_rank,c_nproc;;
104401b82886SBarry Smith   PLA_Obj        v_pla = NULL;
104501b82886SBarry Smith   PetscScalar    *loc_buf;
104601b82886SBarry Smith   Vec            loc_x;
104701b82886SBarry Smith 
104801b82886SBarry Smith   PetscFunctionBegin;
104901b82886SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
105001b82886SBarry Smith 
105101b82886SBarry Smith   /* Create PLAPACK vector objects, then copy b into PLAPACK b */
105262b4c0b3SBarry Smith   ierr = PLA_Mvector_create(lu->datatype,M,1,lu->templ,PLA_ALIGN_FIRST,&v_pla);CHKERRQ(ierr);
105362b4c0b3SBarry Smith   ierr = PLA_Obj_set_to_zero(v_pla);CHKERRQ(ierr);
105401b82886SBarry Smith 
105501b82886SBarry Smith   /* Copy b into rhs_pla */
105662b4c0b3SBarry Smith   ierr = PLA_API_begin();CHKERRQ(ierr);
105762b4c0b3SBarry Smith   ierr = PLA_Obj_API_open(v_pla);CHKERRQ(ierr);
105801b82886SBarry Smith   ierr = VecGetArray(b,&array);CHKERRQ(ierr);
105962b4c0b3SBarry Smith   ierr = PLA_API_axpy_vector_to_global(m,&one,(void *)array,1,v_pla,lu->rstart);CHKERRQ(ierr);
106001b82886SBarry Smith   ierr = VecRestoreArray(b,&array);CHKERRQ(ierr);
106162b4c0b3SBarry Smith   ierr = PLA_Obj_API_close(v_pla);CHKERRQ(ierr);
106262b4c0b3SBarry Smith   ierr = PLA_API_end();CHKERRQ(ierr);
106301b82886SBarry Smith 
106401b82886SBarry Smith   if (A->factor == FACTOR_LU){
106501b82886SBarry Smith     /* Apply the permutations to the right hand sides */
106662b4c0b3SBarry Smith     ierr = PLA_Apply_pivots_to_rows (v_pla,lu->pivots);CHKERRQ(ierr);
106701b82886SBarry Smith 
106801b82886SBarry Smith     /* Solve L y = b, overwriting b with y */
106962b4c0b3SBarry Smith     ierr = PLA_Trsv( PLA_LOWER_TRIANGULAR,PLA_NO_TRANSPOSE,PLA_UNIT_DIAG,lu->A,v_pla );CHKERRQ(ierr);
107001b82886SBarry Smith 
107101b82886SBarry Smith     /* Solve U x = y (=b), overwriting b with x */
107262b4c0b3SBarry Smith     ierr = PLA_Trsv( PLA_UPPER_TRIANGULAR,PLA_NO_TRANSPOSE,PLA_NONUNIT_DIAG,lu->A,v_pla );CHKERRQ(ierr);
107301b82886SBarry Smith   } else { /* FACTOR_CHOLESKY */
107462b4c0b3SBarry Smith     ierr = PLA_Trsv( PLA_LOWER_TRIANGULAR,PLA_NO_TRANSPOSE,PLA_NONUNIT_DIAG,lu->A,v_pla);CHKERRQ(ierr);
107562b4c0b3SBarry Smith     ierr = PLA_Trsv( PLA_LOWER_TRIANGULAR,(lu->datatype == MPI_DOUBLE ? PLA_TRANSPOSE : PLA_CONJUGATE_TRANSPOSE),PLA_NONUNIT_DIAG,lu->A,v_pla);CHKERRQ(ierr);
107601b82886SBarry Smith   }
107701b82886SBarry Smith 
107801b82886SBarry Smith   /* Copy PLAPACK x into Petsc vector x  */
107962b4c0b3SBarry Smith   ierr = PLA_Obj_local_length(v_pla, &loc_m);CHKERRQ(ierr);
108062b4c0b3SBarry Smith   ierr = PLA_Obj_local_buffer(v_pla, (void**)&loc_buf);CHKERRQ(ierr);
108162b4c0b3SBarry Smith   ierr = PLA_Obj_local_stride(v_pla, &loc_stride);CHKERRQ(ierr);
108201b82886SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,loc_m*loc_stride,loc_buf,&loc_x);CHKERRQ(ierr);
108301b82886SBarry Smith   if (!lu->pla_solved){
108401b82886SBarry Smith 
108504fea9ffSBarry Smith     ierr = PLA_Temp_comm_row_info(lu->templ,&dummy_comm,&r_rank,&r_nproc);CHKERRQ(ierr);
108604fea9ffSBarry Smith     ierr = PLA_Temp_comm_col_info(lu->templ,&dummy_comm,&c_rank,&c_nproc);CHKERRQ(ierr);
108701b82886SBarry Smith     /* printf(" [%d] rank: %d %d, nproc: %d %d\n",rank,r_rank,c_rank,r_nproc,c_nproc); */
108801b82886SBarry Smith 
108901b82886SBarry Smith     /* Create IS and cts for VecScatterring */
109062b4c0b3SBarry Smith     ierr = PLA_Obj_local_length(v_pla, &loc_m);CHKERRQ(ierr);
109162b4c0b3SBarry Smith     ierr = PLA_Obj_local_stride(v_pla, &loc_stride);CHKERRQ(ierr);
109201b82886SBarry Smith     ierr = PetscMalloc((2*loc_m+1)*sizeof(PetscInt),&idx_pla);CHKERRQ(ierr);
109301b82886SBarry Smith     idx_petsc = idx_pla + loc_m;
109401b82886SBarry Smith 
109501b82886SBarry Smith     rstart = (r_rank*c_nproc+c_rank)*lu->nb;
109601b82886SBarry Smith     for (i=0; i<loc_m; i+=lu->nb){
109701b82886SBarry Smith       j = 0;
109801b82886SBarry Smith       while (j < lu->nb && i+j < loc_m){
109901b82886SBarry Smith         idx_petsc[i+j] = rstart + j; j++;
110001b82886SBarry Smith       }
110101b82886SBarry Smith       rstart += size*lu->nb;
110201b82886SBarry Smith     }
110301b82886SBarry Smith 
110401b82886SBarry Smith     for (i=0; i<loc_m; i++) idx_pla[i] = i*loc_stride;
110501b82886SBarry Smith 
110601b82886SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,loc_m,idx_pla,&lu->is_pla);CHKERRQ(ierr);
110701b82886SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,loc_m,idx_petsc,&lu->is_petsc);CHKERRQ(ierr);
110801b82886SBarry Smith     ierr = PetscFree(idx_pla);CHKERRQ(ierr);
110901b82886SBarry Smith     ierr = VecScatterCreate(loc_x,lu->is_pla,x,lu->is_petsc,&lu->ctx);CHKERRQ(ierr);
111001b82886SBarry Smith   }
111101b82886SBarry Smith   ierr = VecScatterBegin(lu->ctx,loc_x,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
111201b82886SBarry Smith   ierr = VecScatterEnd(lu->ctx,loc_x,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
111301b82886SBarry Smith 
111401b82886SBarry Smith   /* Free data */
111501b82886SBarry Smith   ierr = VecDestroy(loc_x);CHKERRQ(ierr);
111662b4c0b3SBarry Smith   ierr = PLA_Obj_free(&v_pla);CHKERRQ(ierr);
111701b82886SBarry Smith 
111801b82886SBarry Smith   lu->pla_solved = PETSC_TRUE;
111901b82886SBarry Smith   PetscFunctionReturn(0);
112001b82886SBarry Smith }
112101b82886SBarry Smith 
112201b82886SBarry Smith #undef __FUNCT__
112301b82886SBarry Smith #define __FUNCT__ "MatLUFactorNumeric_MPIDense"
112401b82886SBarry Smith PetscErrorCode MatLUFactorNumeric_MPIDense(Mat A,MatFactorInfo *info,Mat *F)
112501b82886SBarry Smith {
112601b82886SBarry Smith   Mat_Plapack    *lu = (Mat_Plapack*)(*F)->spptr;
112701b82886SBarry Smith   PetscErrorCode ierr;
112801b82886SBarry Smith   PetscInt       M=A->rmap.N,m=A->rmap.n,rstart,rend;
112901b82886SBarry Smith   PetscInt       info_pla=0;
113001b82886SBarry Smith   PetscScalar    *array,one = 1.0;
113101b82886SBarry Smith 
113201b82886SBarry Smith   PetscFunctionBegin;
113362b4c0b3SBarry Smith   ierr = PLA_Obj_set_to_zero(lu->A);CHKERRQ(ierr);
113401b82886SBarry Smith 
113501b82886SBarry Smith   /* Copy A into lu->A */
113662b4c0b3SBarry Smith   ierr = PLA_API_begin();CHKERRQ(ierr);
113762b4c0b3SBarry Smith   ierr = PLA_Obj_API_open(lu->A);CHKERRQ(ierr);
113801b82886SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
113901b82886SBarry Smith   ierr = MatGetArray(A,&array);CHKERRQ(ierr);
114062b4c0b3SBarry Smith   ierr = PLA_API_axpy_matrix_to_global(m,M, &one,(void *)array,m,lu->A,rstart,0);CHKERRQ(ierr);
114101b82886SBarry Smith   ierr = MatRestoreArray(A,&array);CHKERRQ(ierr);
114262b4c0b3SBarry Smith   ierr = PLA_Obj_API_close(lu->A);CHKERRQ(ierr);
114362b4c0b3SBarry Smith   ierr = PLA_API_end();CHKERRQ(ierr);
114401b82886SBarry Smith 
114501b82886SBarry Smith   /* Factor P A -> L U overwriting lower triangular portion of A with L, upper, U */
114601b82886SBarry Smith   info_pla = PLA_LU(lu->A,lu->pivots);
114701b82886SBarry Smith   if (info_pla != 0)
114801b82886SBarry Smith     SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot encountered at row %d from PLA_LU()",info_pla);
114901b82886SBarry Smith 
115001b82886SBarry Smith   lu->CleanUpPlapack = PETSC_TRUE;
115101b82886SBarry Smith   lu->rstart         = rstart;
115201b82886SBarry Smith 
115301b82886SBarry Smith   (*F)->assembled    = PETSC_TRUE;  /* required by -ksp_view */
115401b82886SBarry Smith   PetscFunctionReturn(0);
115501b82886SBarry Smith }
115601b82886SBarry Smith 
115701b82886SBarry Smith #undef __FUNCT__
115801b82886SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_MPIDense"
115901b82886SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_MPIDense(Mat A,MatFactorInfo *info,Mat *F)
116001b82886SBarry Smith {
116101b82886SBarry Smith   Mat_Plapack    *lu = (Mat_Plapack*)(*F)->spptr;
116201b82886SBarry Smith   PetscErrorCode ierr;
116301b82886SBarry Smith   PetscInt       M=A->rmap.N,m=A->rmap.n,rstart,rend;
116401b82886SBarry Smith   PetscInt       info_pla=0;
116501b82886SBarry Smith   PetscScalar    *array,one = 1.0;
116601b82886SBarry Smith 
116701b82886SBarry Smith   PetscFunctionBegin;
1168e0fbc2eeSBarry Smith 
116901b82886SBarry Smith   /* Copy A into lu->A */
117062b4c0b3SBarry Smith   ierr = PLA_Obj_set_to_zero(lu->A);CHKERRQ(ierr);
117162b4c0b3SBarry Smith   ierr = PLA_API_begin();CHKERRQ(ierr);
117262b4c0b3SBarry Smith   ierr = PLA_Obj_API_open(lu->A);CHKERRQ(ierr);
117301b82886SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
117401b82886SBarry Smith   ierr = MatGetArray(A,&array);CHKERRQ(ierr);
117562b4c0b3SBarry Smith   ierr = PLA_API_axpy_matrix_to_global(m,M, &one,(void *)array,m,lu->A,rstart,0);CHKERRQ(ierr);
117601b82886SBarry Smith   ierr = MatRestoreArray(A,&array);CHKERRQ(ierr);
117762b4c0b3SBarry Smith   ierr = PLA_Obj_API_close(lu->A);CHKERRQ(ierr);
117862b4c0b3SBarry Smith   ierr = PLA_API_end();CHKERRQ(ierr);
117901b82886SBarry Smith 
118001b82886SBarry Smith   /* Factor P A -> Chol */
118101b82886SBarry Smith   info_pla = PLA_Chol(PLA_LOWER_TRIANGULAR,lu->A);
118201b82886SBarry Smith   if (info_pla != 0)
118301b82886SBarry Smith     SETERRQ1( PETSC_ERR_MAT_CH_ZRPVT,"Nonpositive definite matrix detected at row %d from PLA_Chol()",info_pla);
118401b82886SBarry Smith 
118501b82886SBarry Smith   lu->CleanUpPlapack = PETSC_TRUE;
118601b82886SBarry Smith   lu->rstart         = rstart;
118701b82886SBarry Smith 
118801b82886SBarry Smith   (*F)->assembled    = PETSC_TRUE;  /* required by -ksp_view */
118901b82886SBarry Smith   PetscFunctionReturn(0);
119001b82886SBarry Smith }
119101b82886SBarry Smith 
119201b82886SBarry Smith #undef __FUNCT__
119301b82886SBarry Smith #define __FUNCT__ "MatFactorSymbolic_MPIDense_Private"
119401b82886SBarry Smith PetscErrorCode MatFactorSymbolic_MPIDense_Private(Mat A,MatFactorInfo *info,Mat *F)
119501b82886SBarry Smith {
119601b82886SBarry Smith   Mat            B;
119701b82886SBarry Smith   Mat_Plapack    *lu;
119801b82886SBarry Smith   PetscErrorCode ierr;
119901b82886SBarry Smith   PetscInt       M=A->rmap.N,N=A->cmap.N;
120004fea9ffSBarry Smith   MPI_Comm       comm=((PetscObject)A)->comm;
120101b82886SBarry Smith   PetscMPIInt    size;
120201b82886SBarry Smith   PetscInt       ierror;
120301b82886SBarry Smith 
120401b82886SBarry Smith   PetscFunctionBegin;
120501b82886SBarry Smith   /* Create the factorization matrix */
120601b82886SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
120701b82886SBarry Smith   ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,M,N);CHKERRQ(ierr);
120801b82886SBarry Smith   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
120901b82886SBarry Smith 
121001b82886SBarry Smith   lu = (Mat_Plapack*)(B->spptr);
121101b82886SBarry Smith 
121201b82886SBarry Smith   /* Set default Plapack parameters */
121301b82886SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
121401b82886SBarry Smith   ierror = 0;
121501b82886SBarry Smith   lu->nb     = M/size;
121601b82886SBarry Smith   if (M - lu->nb*size) lu->nb++; /* without cyclic distribution */
121701b82886SBarry Smith 
121801b82886SBarry Smith   /* Set runtime options */
121901b82886SBarry Smith   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"PLAPACK Options","Mat");CHKERRQ(ierr);
122004fea9ffSBarry Smith   ierr = PetscOptionsInt("-mat_plapack_nprows","row dimension of 2D processor mesh","None",Plapack_nprows,&Plapack_nprows,PETSC_NULL);CHKERRQ(ierr);
122104fea9ffSBarry Smith   ierr = PetscOptionsInt("-mat_plapack_npcols","column dimension of 2D processor mesh","None",Plapack_npcols,&Plapack_npcols,PETSC_NULL);CHKERRQ(ierr);
122201b82886SBarry Smith 
122301b82886SBarry Smith   ierr = PetscOptionsInt("-mat_plapack_nb","block size of template vector","None",lu->nb,&lu->nb,PETSC_NULL);CHKERRQ(ierr);
122401b82886SBarry Smith   ierr = PetscOptionsInt("-mat_plapack_ckerror","error checking flag","None",ierror,&ierror,PETSC_NULL);CHKERRQ(ierr);
122501b82886SBarry Smith   if (ierror){
122662b4c0b3SBarry Smith     ierr = PLA_Set_error_checking(ierror,PETSC_TRUE,PETSC_TRUE,PETSC_FALSE );CHKERRQ(ierr);
122701b82886SBarry Smith   } else {
122862b4c0b3SBarry Smith     ierr = PLA_Set_error_checking(ierror,PETSC_FALSE,PETSC_FALSE,PETSC_FALSE );CHKERRQ(ierr);
122901b82886SBarry Smith   }
123001b82886SBarry Smith   lu->ierror = ierror;
123101b82886SBarry Smith 
123201b82886SBarry Smith   lu->nb_alg = 0;
123301b82886SBarry Smith   ierr = PetscOptionsInt("-mat_plapack_nb_alg","algorithmic block size","None",lu->nb_alg,&lu->nb_alg,PETSC_NULL);CHKERRQ(ierr);
123401b82886SBarry Smith   if (lu->nb_alg){
123562b4c0b3SBarry Smith     ierr = pla_Environ_set_nb_alg (PLA_OP_ALL_ALG,lu->nb_alg);CHKERRQ(ierr);
123601b82886SBarry Smith   }
123701b82886SBarry Smith   PetscOptionsEnd();
123801b82886SBarry Smith 
123901b82886SBarry Smith 
124001b82886SBarry Smith   /* Create object distribution template */
124101b82886SBarry Smith   lu->templ = NULL;
124262b4c0b3SBarry Smith   ierr = PLA_Temp_create(lu->nb, 0, &lu->templ);CHKERRQ(ierr);
124301b82886SBarry Smith 
124401b82886SBarry Smith   /* Use suggested nb_alg if it is not provided by user */
124501b82886SBarry Smith   if (lu->nb_alg == 0){
124662b4c0b3SBarry Smith     ierr = PLA_Environ_nb_alg(PLA_OP_PAN_PAN,lu->templ,&lu->nb_alg);CHKERRQ(ierr);
124762b4c0b3SBarry Smith     ierr = pla_Environ_set_nb_alg(PLA_OP_ALL_ALG,lu->nb_alg);CHKERRQ(ierr);
124801b82886SBarry Smith   }
124901b82886SBarry Smith 
125001b82886SBarry Smith   /* Set the datatype */
125101b82886SBarry Smith #if defined(PETSC_USE_COMPLEX)
125201b82886SBarry Smith   lu->datatype = MPI_DOUBLE_COMPLEX;
125301b82886SBarry Smith #else
125401b82886SBarry Smith   lu->datatype = MPI_DOUBLE;
125501b82886SBarry Smith #endif
125601b82886SBarry Smith 
1257e1156936SBarry Smith   ierr = PLA_Matrix_create(lu->datatype,M,M,lu->templ,PLA_ALIGN_FIRST,PLA_ALIGN_FIRST,&lu->A);CHKERRQ(ierr);
1258e1156936SBarry Smith 
125901b82886SBarry Smith   lu->pla_solved     = PETSC_FALSE; /* MatSolve_Plapack() is called yet */
126001b82886SBarry Smith   lu->CleanUpPlapack = PETSC_TRUE;
126101b82886SBarry Smith   *F                 = B;
126201b82886SBarry Smith   PetscFunctionReturn(0);
126301b82886SBarry Smith }
126401b82886SBarry Smith 
126501b82886SBarry Smith /* Note the Petsc r and c permutations are ignored */
126601b82886SBarry Smith #undef __FUNCT__
126701b82886SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_MPIDense"
126801b82886SBarry Smith PetscErrorCode MatLUFactorSymbolic_MPIDense(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
126901b82886SBarry Smith {
127001b82886SBarry Smith   PetscErrorCode ierr;
1271e1156936SBarry Smith   PetscInt       M = A->rmap.N;
1272e1156936SBarry Smith   Mat_Plapack    *lu;
127301b82886SBarry Smith 
127401b82886SBarry Smith   PetscFunctionBegin;
127501b82886SBarry Smith   ierr = MatFactorSymbolic_MPIDense_Private(A,info,F);CHKERRQ(ierr);
1276e1156936SBarry Smith   lu = (Mat_Plapack*)(*F)->spptr;
1277e1156936SBarry Smith   ierr = PLA_Mvector_create(MPI_INT,M,1,lu->templ,PLA_ALIGN_FIRST,&lu->pivots);CHKERRQ(ierr);
127801b82886SBarry Smith   (*F)->factor = FACTOR_LU;
127901b82886SBarry Smith   PetscFunctionReturn(0);
128001b82886SBarry Smith }
128101b82886SBarry Smith 
128201b82886SBarry Smith /* Note the Petsc perm permutation is ignored */
128301b82886SBarry Smith #undef __FUNCT__
128401b82886SBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_MPIDense"
128501b82886SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_MPIDense(Mat A,IS perm,MatFactorInfo *info,Mat *F)
128601b82886SBarry Smith {
128701b82886SBarry Smith   PetscErrorCode ierr;
128801b82886SBarry Smith   PetscTruth     issymmetric,set;
128901b82886SBarry Smith 
129001b82886SBarry Smith   PetscFunctionBegin;
129101b82886SBarry Smith   ierr = MatIsSymmetricKnown(A,&set,&issymmetric); CHKERRQ(ierr);
129201b82886SBarry Smith   if (!set || !issymmetric) SETERRQ(PETSC_ERR_USER,"Matrix must be set as MAT_SYMMETRIC for CholeskyFactor()");
129301b82886SBarry Smith   ierr = MatFactorSymbolic_MPIDense_Private(A,info,F);CHKERRQ(ierr);
129401b82886SBarry Smith   (*F)->factor = FACTOR_CHOLESKY;
129501b82886SBarry Smith   PetscFunctionReturn(0);
129601b82886SBarry Smith }
129701b82886SBarry Smith #endif
129801b82886SBarry Smith 
12998965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
130009dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
130109dc0095SBarry Smith        MatGetRow_MPIDense,
130209dc0095SBarry Smith        MatRestoreRow_MPIDense,
130309dc0095SBarry Smith        MatMult_MPIDense,
130497304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense,
13057c922b88SBarry Smith        MatMultTranspose_MPIDense,
13067c922b88SBarry Smith        MatMultTransposeAdd_MPIDense,
130701b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
130801b82886SBarry Smith        MatSolve_MPIDense,
130901b82886SBarry Smith #else
13108965ea79SLois Curfman McInnes        0,
131101b82886SBarry Smith #endif
131209dc0095SBarry Smith        0,
131309dc0095SBarry Smith        0,
131497304618SKris Buschelman /*10*/ 0,
131509dc0095SBarry Smith        0,
131609dc0095SBarry Smith        0,
131709dc0095SBarry Smith        0,
131809dc0095SBarry Smith        MatTranspose_MPIDense,
131997304618SKris Buschelman /*15*/ MatGetInfo_MPIDense,
13206e4ee0c6SHong Zhang        MatEqual_MPIDense,
132109dc0095SBarry Smith        MatGetDiagonal_MPIDense,
13225b2fa520SLois Curfman McInnes        MatDiagonalScale_MPIDense,
132309dc0095SBarry Smith        MatNorm_MPIDense,
132497304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense,
132509dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
132609dc0095SBarry Smith        0,
132709dc0095SBarry Smith        MatSetOption_MPIDense,
132809dc0095SBarry Smith        MatZeroEntries_MPIDense,
132997304618SKris Buschelman /*25*/ MatZeroRows_MPIDense,
133001b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
13310ad19415SHong Zhang        MatLUFactorSymbolic_MPIDense,
133201b82886SBarry Smith        MatLUFactorNumeric_MPIDense,
13330ad19415SHong Zhang        MatCholeskyFactorSymbolic_MPIDense,
133401b82886SBarry Smith        MatCholeskyFactorNumeric_MPIDense,
133501b82886SBarry Smith #else
1336919b68f7SBarry Smith        0,
133701b82886SBarry Smith        0,
133801b82886SBarry Smith        0,
133901b82886SBarry Smith        0,
134001b82886SBarry Smith #endif
134197304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense,
1342273d9f13SBarry Smith        0,
134309dc0095SBarry Smith        0,
134409dc0095SBarry Smith        MatGetArray_MPIDense,
134509dc0095SBarry Smith        MatRestoreArray_MPIDense,
134697304618SKris Buschelman /*35*/ MatDuplicate_MPIDense,
134709dc0095SBarry Smith        0,
134809dc0095SBarry Smith        0,
134909dc0095SBarry Smith        0,
135009dc0095SBarry Smith        0,
135197304618SKris Buschelman /*40*/ 0,
13522ce60cd0SSatish Balay        MatGetSubMatrices_MPIDense,
135309dc0095SBarry Smith        0,
135409dc0095SBarry Smith        MatGetValues_MPIDense,
135509dc0095SBarry Smith        0,
135697304618SKris Buschelman /*45*/ 0,
135709dc0095SBarry Smith        MatScale_MPIDense,
135809dc0095SBarry Smith        0,
135909dc0095SBarry Smith        0,
136009dc0095SBarry Smith        0,
1361521d7252SBarry Smith /*50*/ 0,
136209dc0095SBarry Smith        0,
136309dc0095SBarry Smith        0,
136409dc0095SBarry Smith        0,
136509dc0095SBarry Smith        0,
136697304618SKris Buschelman /*55*/ 0,
136709dc0095SBarry Smith        0,
136809dc0095SBarry Smith        0,
136909dc0095SBarry Smith        0,
137009dc0095SBarry Smith        0,
137197304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense,
1372b9b97703SBarry Smith        MatDestroy_MPIDense,
1373b9b97703SBarry Smith        MatView_MPIDense,
1374357abbc8SBarry Smith        0,
137597304618SKris Buschelman        0,
137697304618SKris Buschelman /*65*/ 0,
137797304618SKris Buschelman        0,
137897304618SKris Buschelman        0,
137997304618SKris Buschelman        0,
138097304618SKris Buschelman        0,
138197304618SKris Buschelman /*70*/ 0,
138297304618SKris Buschelman        0,
138397304618SKris Buschelman        0,
138497304618SKris Buschelman        0,
138597304618SKris Buschelman        0,
138697304618SKris Buschelman /*75*/ 0,
138797304618SKris Buschelman        0,
138897304618SKris Buschelman        0,
138997304618SKris Buschelman        0,
139097304618SKris Buschelman        0,
139197304618SKris Buschelman /*80*/ 0,
139297304618SKris Buschelman        0,
139397304618SKris Buschelman        0,
139497304618SKris Buschelman        0,
1395865e5f61SKris Buschelman /*84*/ MatLoad_MPIDense,
1396865e5f61SKris Buschelman        0,
1397865e5f61SKris Buschelman        0,
1398865e5f61SKris Buschelman        0,
1399865e5f61SKris Buschelman        0,
1400865e5f61SKris Buschelman        0,
1401865e5f61SKris Buschelman /*90*/ 0,
14024ae313f4SHong Zhang        MatMatMultSymbolic_MPIDense_MPIDense,
1403865e5f61SKris Buschelman        0,
1404865e5f61SKris Buschelman        0,
1405865e5f61SKris Buschelman        0,
1406865e5f61SKris Buschelman /*95*/ 0,
1407865e5f61SKris Buschelman        0,
1408865e5f61SKris Buschelman        0,
1409865e5f61SKris Buschelman        0};
14108965ea79SLois Curfman McInnes 
1411273d9f13SBarry Smith EXTERN_C_BEGIN
14124a2ae208SSatish Balay #undef __FUNCT__
1413a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense"
1414be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1415a23d5eceSKris Buschelman {
1416a23d5eceSKris Buschelman   Mat_MPIDense   *a;
1417dfbe8321SBarry Smith   PetscErrorCode ierr;
1418a23d5eceSKris Buschelman 
1419a23d5eceSKris Buschelman   PetscFunctionBegin;
1420a23d5eceSKris Buschelman   mat->preallocated = PETSC_TRUE;
1421a23d5eceSKris Buschelman   /* Note:  For now, when data is specified above, this assumes the user correctly
1422a23d5eceSKris Buschelman    allocates the local dense storage space.  We should add error checking. */
1423a23d5eceSKris Buschelman 
1424a23d5eceSKris Buschelman   a    = (Mat_MPIDense*)mat->data;
1425f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr);
1426899cda47SBarry Smith   ierr = MatSetSizes(a->A,mat->rmap.n,mat->cmap.N,mat->rmap.n,mat->cmap.N);CHKERRQ(ierr);
14275c5985e7SKris Buschelman   ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr);
14285c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr);
142952e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
1430a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1431a23d5eceSKris Buschelman }
1432a23d5eceSKris Buschelman EXTERN_C_END
1433a23d5eceSKris Buschelman 
14347878bbefSBarry Smith EXTERN_C_BEGIN
14357878bbefSBarry Smith #if defined(PETSC_HAVE_PLAPACK)
14367878bbefSBarry Smith #undef __FUNCT__
14377878bbefSBarry Smith #define __FUNCT__ "MatSetSolverType_MPIDense_PLAPACK"
14387878bbefSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatSetSolverType_MPIDense_PLAPACK(Mat mat,const char *type)
14397878bbefSBarry Smith {
14407878bbefSBarry Smith   PetscFunctionBegin;
14417878bbefSBarry Smith   /* dummy function so that -mat_solver_type plapack or MatSetSolverType(mat,"plapack"); silently work (since PLAPACK is on by default) */
14427878bbefSBarry Smith   PetscFunctionReturn(0);
14437878bbefSBarry Smith }
14447878bbefSBarry Smith #endif
14457a667e6fSSatish Balay EXTERN_C_END
14467878bbefSBarry Smith 
14470bad9183SKris Buschelman /*MC
1448fafad747SKris Buschelman    MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices.
14490bad9183SKris Buschelman 
14500bad9183SKris Buschelman    Options Database Keys:
14510bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions()
14520bad9183SKris Buschelman 
14530bad9183SKris Buschelman   Level: beginner
14540bad9183SKris Buschelman 
14557878bbefSBarry Smith   MATMPIDENSE matrices may use direct solvers (LU, Cholesky, and QR)
14567878bbefSBarry Smith   for parallel dense matrices via the external package PLAPACK, if PLAPACK is installed
14577878bbefSBarry Smith   (run config/configure.py with the option --download-plapack)
14587878bbefSBarry Smith 
14597878bbefSBarry Smith 
14607878bbefSBarry Smith   Options Database Keys:
14617878bbefSBarry Smith . -mat_plapack_nprows <n> - number of rows in processor partition
14627878bbefSBarry Smith . -mat_plapack_npcols <n> - number of columns in processor partition
14637878bbefSBarry Smith . -mat_plapack_nb <n> - block size of template vector
14647878bbefSBarry Smith . -mat_plapack_nb_alg <n> - algorithmic block size
14657878bbefSBarry Smith - -mat_plapack_ckerror <n> - error checking flag
14667878bbefSBarry Smith 
14677878bbefSBarry Smith .seealso: MatCreateMPIDense(), MATDENSE, MATSEQDENSE
14680bad9183SKris Buschelman M*/
14690bad9183SKris Buschelman 
1470a23d5eceSKris Buschelman EXTERN_C_BEGIN
1471a23d5eceSKris Buschelman #undef __FUNCT__
14724a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense"
1473be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPIDense(Mat mat)
1474273d9f13SBarry Smith {
1475273d9f13SBarry Smith   Mat_MPIDense   *a;
1476dfbe8321SBarry Smith   PetscErrorCode ierr;
147701b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
147801b82886SBarry Smith   Mat_Plapack    *lu;
147901b82886SBarry Smith #endif
1480273d9f13SBarry Smith 
1481273d9f13SBarry Smith   PetscFunctionBegin;
148238f2d2fdSLisandro Dalcin   ierr              = PetscNewLog(mat,Mat_MPIDense,&a);CHKERRQ(ierr);
1483b0a32e0cSBarry Smith   mat->data         = (void*)a;
1484273d9f13SBarry Smith   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1485273d9f13SBarry Smith   mat->factor       = 0;
1486273d9f13SBarry Smith   mat->mapping      = 0;
1487273d9f13SBarry Smith 
1488273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
14897adad957SLisandro Dalcin   ierr = MPI_Comm_rank(((PetscObject)mat)->comm,&a->rank);CHKERRQ(ierr);
14907adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)mat)->comm,&a->size);CHKERRQ(ierr);
1491273d9f13SBarry Smith 
1492899cda47SBarry Smith   mat->rmap.bs = mat->cmap.bs = 1;
14936148ca0dSBarry Smith   ierr = PetscMapSetUp(&mat->rmap);CHKERRQ(ierr);
14946148ca0dSBarry Smith   ierr = PetscMapSetUp(&mat->cmap);CHKERRQ(ierr);
1495899cda47SBarry Smith   a->nvec = mat->cmap.n;
1496273d9f13SBarry Smith 
1497273d9f13SBarry Smith   /* build cache for off array entries formed */
1498273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
14997adad957SLisandro Dalcin   ierr = MatStashCreate_Private(((PetscObject)mat)->comm,1,&mat->stash);CHKERRQ(ierr);
1500273d9f13SBarry Smith 
1501273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1502273d9f13SBarry Smith   a->lvec        = 0;
1503273d9f13SBarry Smith   a->Mvctx       = 0;
1504273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1505273d9f13SBarry Smith 
1506273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1507273d9f13SBarry Smith                                      "MatGetDiagonalBlock_MPIDense",
1508273d9f13SBarry Smith                                      MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1509a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C",
1510a23d5eceSKris Buschelman                                      "MatMPIDenseSetPreallocation_MPIDense",
1511a23d5eceSKris Buschelman                                      MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr);
15124ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",
15134ae313f4SHong Zhang                                      "MatMatMult_MPIAIJ_MPIDense",
15144ae313f4SHong Zhang                                       MatMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr);
15154ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",
15164ae313f4SHong Zhang                                      "MatMatMultSymbolic_MPIAIJ_MPIDense",
15174ae313f4SHong Zhang                                       MatMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr);
15184ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",
15194ae313f4SHong Zhang                                      "MatMatMultNumeric_MPIAIJ_MPIDense",
15204ae313f4SHong Zhang                                       MatMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr);
15217878bbefSBarry Smith #if defined(PETSC_HAVE_PLAPACK)
15227878bbefSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSetSolverType_mpidense_plapack_C",
15237878bbefSBarry Smith                                      "MatSetSolverType_MPIDense_PLAPACK",
15247878bbefSBarry Smith                                       MatSetSolverType_MPIDense_PLAPACK);CHKERRQ(ierr);
15257878bbefSBarry Smith #endif
152638aed534SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);CHKERRQ(ierr);
152701b82886SBarry Smith 
152801b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
152901b82886SBarry Smith   ierr = PetscNewLog(mat,Mat_Plapack,&lu);CHKERRQ(ierr);
153001b82886SBarry Smith   lu->CleanUpPlapack       = PETSC_FALSE;
153101b82886SBarry Smith   mat->spptr                 = (void*)lu;
153201b82886SBarry Smith #endif
1533273d9f13SBarry Smith   PetscFunctionReturn(0);
1534273d9f13SBarry Smith }
1535273d9f13SBarry Smith EXTERN_C_END
1536273d9f13SBarry Smith 
1537209238afSKris Buschelman /*MC
1538002d173eSKris Buschelman    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1539209238afSKris Buschelman 
1540209238afSKris Buschelman    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1541209238afSKris Buschelman    and MATMPIDENSE otherwise.
1542209238afSKris Buschelman 
1543209238afSKris Buschelman    Options Database Keys:
1544209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1545209238afSKris Buschelman 
1546209238afSKris Buschelman   Level: beginner
1547209238afSKris Buschelman 
154801b82886SBarry Smith 
1549209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1550209238afSKris Buschelman M*/
1551209238afSKris Buschelman 
1552209238afSKris Buschelman EXTERN_C_BEGIN
1553209238afSKris Buschelman #undef __FUNCT__
1554209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense"
1555be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Dense(Mat A)
1556dfbe8321SBarry Smith {
15576849ba73SBarry Smith   PetscErrorCode ierr;
155813f74950SBarry Smith   PetscMPIInt    size;
1559209238afSKris Buschelman 
1560209238afSKris Buschelman   PetscFunctionBegin;
15617adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1562209238afSKris Buschelman   if (size == 1) {
1563209238afSKris Buschelman     ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr);
1564209238afSKris Buschelman   } else {
1565209238afSKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
1566209238afSKris Buschelman   }
1567209238afSKris Buschelman   PetscFunctionReturn(0);
1568209238afSKris Buschelman }
1569209238afSKris Buschelman EXTERN_C_END
1570209238afSKris Buschelman 
15714a2ae208SSatish Balay #undef __FUNCT__
15724a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation"
1573273d9f13SBarry Smith /*@C
1574273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1575273d9f13SBarry Smith 
1576273d9f13SBarry Smith    Not collective
1577273d9f13SBarry Smith 
1578273d9f13SBarry Smith    Input Parameters:
1579273d9f13SBarry Smith .  A - the matrix
1580273d9f13SBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1581273d9f13SBarry Smith    to control all matrix memory allocation.
1582273d9f13SBarry Smith 
1583273d9f13SBarry Smith    Notes:
1584273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1585273d9f13SBarry Smith    storage by columns.
1586273d9f13SBarry Smith 
1587273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1588273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1589273d9f13SBarry Smith    set data=PETSC_NULL.
1590273d9f13SBarry Smith 
1591273d9f13SBarry Smith    Level: intermediate
1592273d9f13SBarry Smith 
1593273d9f13SBarry Smith .keywords: matrix,dense, parallel
1594273d9f13SBarry Smith 
1595273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1596273d9f13SBarry Smith @*/
1597be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data)
1598273d9f13SBarry Smith {
1599dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar *);
1600273d9f13SBarry Smith 
1601273d9f13SBarry Smith   PetscFunctionBegin;
1602565567f0SKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1603a23d5eceSKris Buschelman   if (f) {
1604a23d5eceSKris Buschelman     ierr = (*f)(mat,data);CHKERRQ(ierr);
1605a23d5eceSKris Buschelman   }
1606273d9f13SBarry Smith   PetscFunctionReturn(0);
1607273d9f13SBarry Smith }
1608273d9f13SBarry Smith 
16094a2ae208SSatish Balay #undef __FUNCT__
16104a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense"
16118965ea79SLois Curfman McInnes /*@C
161239ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
16138965ea79SLois Curfman McInnes 
1614db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1615db81eaa0SLois Curfman McInnes 
16168965ea79SLois Curfman McInnes    Input Parameters:
1617db81eaa0SLois Curfman McInnes +  comm - MPI communicator
16188965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1619db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
16208965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1621db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
16227f5ff6fdSBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1623dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
16248965ea79SLois Curfman McInnes 
16258965ea79SLois Curfman McInnes    Output Parameter:
1626477f1c0bSLois Curfman McInnes .  A - the matrix
16278965ea79SLois Curfman McInnes 
1628b259b22eSLois Curfman McInnes    Notes:
162939ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
163039ddd567SLois Curfman McInnes    storage by columns.
16318965ea79SLois Curfman McInnes 
163218f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
163318f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
16347f5ff6fdSBarry Smith    set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users).
163518f449edSLois Curfman McInnes 
16368965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
16378965ea79SLois Curfman McInnes    (possibly both).
16388965ea79SLois Curfman McInnes 
1639027ccd11SLois Curfman McInnes    Level: intermediate
1640027ccd11SLois Curfman McInnes 
164139ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
16428965ea79SLois Curfman McInnes 
164339ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
16448965ea79SLois Curfman McInnes @*/
1645be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
16468965ea79SLois Curfman McInnes {
16476849ba73SBarry Smith   PetscErrorCode ierr;
164813f74950SBarry Smith   PetscMPIInt    size;
16498965ea79SLois Curfman McInnes 
16503a40ed3dSBarry Smith   PetscFunctionBegin;
1651f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1652f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1653273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1654273d9f13SBarry Smith   if (size > 1) {
1655273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1656273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1657273d9f13SBarry Smith   } else {
1658273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1659273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
16608c469469SLois Curfman McInnes   }
16613a40ed3dSBarry Smith   PetscFunctionReturn(0);
16628965ea79SLois Curfman McInnes }
16638965ea79SLois Curfman McInnes 
16644a2ae208SSatish Balay #undef __FUNCT__
16654a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense"
16666849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
16678965ea79SLois Curfman McInnes {
16688965ea79SLois Curfman McInnes   Mat            mat;
16693501a2bdSLois Curfman McInnes   Mat_MPIDense   *a,*oldmat = (Mat_MPIDense*)A->data;
1670dfbe8321SBarry Smith   PetscErrorCode ierr;
16718965ea79SLois Curfman McInnes 
16723a40ed3dSBarry Smith   PetscFunctionBegin;
16738965ea79SLois Curfman McInnes   *newmat       = 0;
16747adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&mat);CHKERRQ(ierr);
1675899cda47SBarry Smith   ierr = MatSetSizes(mat,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr);
16767adad957SLisandro Dalcin   ierr = MatSetType(mat,((PetscObject)A)->type_name);CHKERRQ(ierr);
1677834f8fabSBarry Smith   a                 = (Mat_MPIDense*)mat->data;
1678e04c1aa4SHong Zhang   ierr              = PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
16793501a2bdSLois Curfman McInnes   mat->factor       = A->factor;
1680c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1681273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
16828965ea79SLois Curfman McInnes 
1683899cda47SBarry Smith   mat->rmap.rstart     = A->rmap.rstart;
1684899cda47SBarry Smith   mat->rmap.rend       = A->rmap.rend;
16858965ea79SLois Curfman McInnes   a->size              = oldmat->size;
16868965ea79SLois Curfman McInnes   a->rank              = oldmat->rank;
1687e0fa3b82SLois Curfman McInnes   mat->insertmode      = NOT_SET_VALUES;
1688b9b97703SBarry Smith   a->nvec              = oldmat->nvec;
16893782ba37SSatish Balay   a->donotstash        = oldmat->donotstash;
1690e04c1aa4SHong Zhang 
1691899cda47SBarry Smith   ierr = PetscMemcpy(mat->rmap.range,A->rmap.range,(a->size+1)*sizeof(PetscInt));CHKERRQ(ierr);
1692899cda47SBarry Smith   ierr = PetscMemcpy(mat->cmap.range,A->cmap.range,(a->size+1)*sizeof(PetscInt));CHKERRQ(ierr);
16937adad957SLisandro Dalcin   ierr = MatStashCreate_Private(((PetscObject)A)->comm,1,&mat->stash);CHKERRQ(ierr);
16948965ea79SLois Curfman McInnes 
1695329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
16965609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
169752e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
169801b82886SBarry Smith 
169901b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
170001b82886SBarry Smith   ierr = PetscMemcpy(mat->spptr,A->spptr,sizeof(Mat_Plapack));CHKERRQ(ierr);
170101b82886SBarry Smith #endif
17028965ea79SLois Curfman McInnes   *newmat = mat;
17033a40ed3dSBarry Smith   PetscFunctionReturn(0);
17048965ea79SLois Curfman McInnes }
17058965ea79SLois Curfman McInnes 
1706e090d566SSatish Balay #include "petscsys.h"
17078965ea79SLois Curfman McInnes 
17084a2ae208SSatish Balay #undef __FUNCT__
17094a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile"
1710f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N, MatType type,Mat *newmat)
171190ace30eSBarry Smith {
17126849ba73SBarry Smith   PetscErrorCode ierr;
171332dcc486SBarry Smith   PetscMPIInt    rank,size;
171413f74950SBarry Smith   PetscInt       *rowners,i,m,nz,j;
171587828ca2SBarry Smith   PetscScalar    *array,*vals,*vals_ptr;
171690ace30eSBarry Smith   MPI_Status     status;
171790ace30eSBarry Smith 
17183a40ed3dSBarry Smith   PetscFunctionBegin;
1719d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1720d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
172190ace30eSBarry Smith 
172290ace30eSBarry Smith   /* determine ownership of all rows */
172390ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
172413f74950SBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr);
172513f74950SBarry Smith   ierr       = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);CHKERRQ(ierr);
172690ace30eSBarry Smith   rowners[0] = 0;
172790ace30eSBarry Smith   for (i=2; i<=size; i++) {
172890ace30eSBarry Smith     rowners[i] += rowners[i-1];
172990ace30eSBarry Smith   }
173090ace30eSBarry Smith 
1731f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,newmat);CHKERRQ(ierr);
1732f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1733be5d1d56SKris Buschelman   ierr = MatSetType(*newmat,type);CHKERRQ(ierr);
1734878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
173590ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr);
173690ace30eSBarry Smith 
173790ace30eSBarry Smith   if (!rank) {
173887828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
173990ace30eSBarry Smith 
174090ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
17410752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
174290ace30eSBarry Smith 
174390ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
174490ace30eSBarry Smith     vals_ptr = vals;
174590ace30eSBarry Smith     for (i=0; i<m; i++) {
174690ace30eSBarry Smith       for (j=0; j<N; j++) {
174790ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
174890ace30eSBarry Smith       }
174990ace30eSBarry Smith     }
175090ace30eSBarry Smith 
175190ace30eSBarry Smith     /* read in other processors and ship out */
175290ace30eSBarry Smith     for (i=1; i<size; i++) {
175390ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
17540752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
17557adad957SLisandro Dalcin       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)(*newmat))->tag,comm);CHKERRQ(ierr);
175690ace30eSBarry Smith     }
17573a40ed3dSBarry Smith   } else {
175890ace30eSBarry Smith     /* receive numeric values */
175987828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
176090ace30eSBarry Smith 
176190ace30eSBarry Smith     /* receive message of values*/
17627adad957SLisandro Dalcin     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,((PetscObject)(*newmat))->tag,comm,&status);CHKERRQ(ierr);
176390ace30eSBarry Smith 
176490ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
176590ace30eSBarry Smith     vals_ptr = vals;
176690ace30eSBarry Smith     for (i=0; i<m; i++) {
176790ace30eSBarry Smith       for (j=0; j<N; j++) {
176890ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
176990ace30eSBarry Smith       }
177090ace30eSBarry Smith     }
177190ace30eSBarry Smith   }
1772606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
1773606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
17746d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17756d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17763a40ed3dSBarry Smith   PetscFunctionReturn(0);
177790ace30eSBarry Smith }
177890ace30eSBarry Smith 
17794a2ae208SSatish Balay #undef __FUNCT__
17804a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense"
1781f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense(PetscViewer viewer, MatType type,Mat *newmat)
17828965ea79SLois Curfman McInnes {
17838965ea79SLois Curfman McInnes   Mat            A;
178487828ca2SBarry Smith   PetscScalar    *vals,*svals;
178519bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
17868965ea79SLois Curfman McInnes   MPI_Status     status;
178713f74950SBarry Smith   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz;
178813f74950SBarry Smith   PetscInt       header[4],*rowlengths = 0,M,N,*cols;
178913f74950SBarry Smith   PetscInt       *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols;
179013f74950SBarry Smith   PetscInt       i,nz,j,rstart,rend;
179113f74950SBarry Smith   int            fd;
17926849ba73SBarry Smith   PetscErrorCode ierr;
17938965ea79SLois Curfman McInnes 
17943a40ed3dSBarry Smith   PetscFunctionBegin;
1795d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1796d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
17978965ea79SLois Curfman McInnes   if (!rank) {
1798b0a32e0cSBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
17990752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
1800552e946dSBarry Smith     if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
18018965ea79SLois Curfman McInnes   }
18028965ea79SLois Curfman McInnes 
180313f74950SBarry Smith   ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr);
180490ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
180590ace30eSBarry Smith 
180690ace30eSBarry Smith   /*
180790ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
180890ace30eSBarry Smith   */
180990ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
1810be5d1d56SKris Buschelman     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr);
18113a40ed3dSBarry Smith     PetscFunctionReturn(0);
181290ace30eSBarry Smith   }
181390ace30eSBarry Smith 
18148965ea79SLois Curfman McInnes   /* determine ownership of all rows */
1815e44c0bd4SBarry Smith   m          = PetscMPIIntCast(M/size + ((M % size) > rank));
1816e44c0bd4SBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(PetscMPIInt),&rowners);CHKERRQ(ierr);
1817ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
18188965ea79SLois Curfman McInnes   rowners[0] = 0;
18198965ea79SLois Curfman McInnes   for (i=2; i<=size; i++) {
18208965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
18218965ea79SLois Curfman McInnes   }
18228965ea79SLois Curfman McInnes   rstart = rowners[rank];
18238965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
18248965ea79SLois Curfman McInnes 
18258965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
182613f74950SBarry Smith   ierr    = PetscMalloc(2*(rend-rstart+1)*sizeof(PetscInt),&ourlens);CHKERRQ(ierr);
18278965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
18288965ea79SLois Curfman McInnes   if (!rank) {
182913f74950SBarry Smith     ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
18300752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
183113f74950SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr);
18328965ea79SLois Curfman McInnes     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
18331466f1e1SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
1834606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1835ca161407SBarry Smith   } else {
18361466f1e1SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
18378965ea79SLois Curfman McInnes   }
18388965ea79SLois Curfman McInnes 
18398965ea79SLois Curfman McInnes   if (!rank) {
18408965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
184113f74950SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr);
184213f74950SBarry Smith     ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr);
18438965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
18448965ea79SLois Curfman McInnes       for (j=rowners[i]; j< rowners[i+1]; j++) {
18458965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
18468965ea79SLois Curfman McInnes       }
18478965ea79SLois Curfman McInnes     }
1848606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
18498965ea79SLois Curfman McInnes 
18508965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
18518965ea79SLois Curfman McInnes     maxnz = 0;
18528965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
18530452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
18548965ea79SLois Curfman McInnes     }
185513f74950SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr);
18568965ea79SLois Curfman McInnes 
18578965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
18588965ea79SLois Curfman McInnes     nz = procsnz[0];
185913f74950SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
18600752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
18618965ea79SLois Curfman McInnes 
18628965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
18638965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
18648965ea79SLois Curfman McInnes       nz   = procsnz[i];
18650752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
186613f74950SBarry Smith       ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr);
18678965ea79SLois Curfman McInnes     }
1868606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
18693a40ed3dSBarry Smith   } else {
18708965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
18718965ea79SLois Curfman McInnes     nz = 0;
18728965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
18738965ea79SLois Curfman McInnes       nz += ourlens[i];
18748965ea79SLois Curfman McInnes     }
187513f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
18768965ea79SLois Curfman McInnes 
18778965ea79SLois Curfman McInnes     /* receive message of column indices*/
187813f74950SBarry Smith     ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr);
187913f74950SBarry Smith     ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr);
188029bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
18818965ea79SLois Curfman McInnes   }
18828965ea79SLois Curfman McInnes 
18838965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
188413f74950SBarry Smith   ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr);
18858965ea79SLois Curfman McInnes   jj = 0;
18868965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
18878965ea79SLois Curfman McInnes     for (j=0; j<ourlens[i]; j++) {
18888965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
18898965ea79SLois Curfman McInnes       jj++;
18908965ea79SLois Curfman McInnes     }
18918965ea79SLois Curfman McInnes   }
18928965ea79SLois Curfman McInnes 
18938965ea79SLois Curfman McInnes   /* create our matrix */
18948965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
18958965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
18968965ea79SLois Curfman McInnes   }
1897f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,newmat);CHKERRQ(ierr);
1898f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1899be5d1d56SKris Buschelman   ierr = MatSetType(*newmat,type);CHKERRQ(ierr);
1900878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
19018965ea79SLois Curfman McInnes   A = *newmat;
19028965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
19038965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
19048965ea79SLois Curfman McInnes   }
19058965ea79SLois Curfman McInnes 
19068965ea79SLois Curfman McInnes   if (!rank) {
190787828ca2SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
19088965ea79SLois Curfman McInnes 
19098965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
19108965ea79SLois Curfman McInnes     nz = procsnz[0];
19110752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
19128965ea79SLois Curfman McInnes 
19138965ea79SLois Curfman McInnes     /* insert into matrix */
19148965ea79SLois Curfman McInnes     jj      = rstart;
19158965ea79SLois Curfman McInnes     smycols = mycols;
19168965ea79SLois Curfman McInnes     svals   = vals;
19178965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
19188965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
19198965ea79SLois Curfman McInnes       smycols += ourlens[i];
19208965ea79SLois Curfman McInnes       svals   += ourlens[i];
19218965ea79SLois Curfman McInnes       jj++;
19228965ea79SLois Curfman McInnes     }
19238965ea79SLois Curfman McInnes 
19248965ea79SLois Curfman McInnes     /* read in other processors and ship out */
19258965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
19268965ea79SLois Curfman McInnes       nz   = procsnz[i];
19270752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
19287adad957SLisandro Dalcin       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)A)->tag,comm);CHKERRQ(ierr);
19298965ea79SLois Curfman McInnes     }
1930606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
19313a40ed3dSBarry Smith   } else {
19328965ea79SLois Curfman McInnes     /* receive numeric values */
193384d528b1SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
19348965ea79SLois Curfman McInnes 
19358965ea79SLois Curfman McInnes     /* receive message of values*/
19367adad957SLisandro Dalcin     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)A)->tag,comm,&status);CHKERRQ(ierr);
1937ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
193829bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
19398965ea79SLois Curfman McInnes 
19408965ea79SLois Curfman McInnes     /* insert into matrix */
19418965ea79SLois Curfman McInnes     jj      = rstart;
19428965ea79SLois Curfman McInnes     smycols = mycols;
19438965ea79SLois Curfman McInnes     svals   = vals;
19448965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
19458965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
19468965ea79SLois Curfman McInnes       smycols += ourlens[i];
19478965ea79SLois Curfman McInnes       svals   += ourlens[i];
19488965ea79SLois Curfman McInnes       jj++;
19498965ea79SLois Curfman McInnes     }
19508965ea79SLois Curfman McInnes   }
1951606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1952606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1953606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1954606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
19558965ea79SLois Curfman McInnes 
19566d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19576d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19583a40ed3dSBarry Smith   PetscFunctionReturn(0);
19598965ea79SLois Curfman McInnes }
196090ace30eSBarry Smith 
19616e4ee0c6SHong Zhang #undef __FUNCT__
19626e4ee0c6SHong Zhang #define __FUNCT__ "MatEqual_MPIDense"
19636e4ee0c6SHong Zhang PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscTruth *flag)
19646e4ee0c6SHong Zhang {
19656e4ee0c6SHong Zhang   Mat_MPIDense   *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data;
19666e4ee0c6SHong Zhang   Mat            a,b;
19676e4ee0c6SHong Zhang   PetscTruth     flg;
19686e4ee0c6SHong Zhang   PetscErrorCode ierr;
196990ace30eSBarry Smith 
19706e4ee0c6SHong Zhang   PetscFunctionBegin;
19716e4ee0c6SHong Zhang   a = matA->A;
19726e4ee0c6SHong Zhang   b = matB->A;
19736e4ee0c6SHong Zhang   ierr = MatEqual(a,b,&flg);CHKERRQ(ierr);
19747adad957SLisandro Dalcin   ierr = MPI_Allreduce(&flg,flag,1,MPI_INT,MPI_LAND,((PetscObject)A)->comm);CHKERRQ(ierr);
19756e4ee0c6SHong Zhang   PetscFunctionReturn(0);
19766e4ee0c6SHong Zhang }
197790ace30eSBarry Smith 
197809d27a7eSBarry Smith #if defined(PETSC_HAVE_PLAPACK)
197909d27a7eSBarry Smith 
198009d27a7eSBarry Smith #undef __FUNCT__
198109d27a7eSBarry Smith #define __FUNCT__ "PetscPLAPACKFinalizePackage"
198209d27a7eSBarry Smith /*@C
198309d27a7eSBarry Smith   PetscPLAPACKFinalizePackage - This function destroys everything in the Petsc interface to PLAPACK.
198409d27a7eSBarry Smith   Level: developer
198509d27a7eSBarry Smith 
198609d27a7eSBarry Smith .keywords: Petsc, destroy, package, PLAPACK
198709d27a7eSBarry Smith .seealso: PetscFinalize()
198809d27a7eSBarry Smith @*/
198909d27a7eSBarry Smith PetscErrorCode PETSC_DLLEXPORT PetscPLAPACKFinalizePackage(void)
199009d27a7eSBarry Smith {
199109d27a7eSBarry Smith   PetscErrorCode ierr;
199209d27a7eSBarry Smith 
199309d27a7eSBarry Smith   PetscFunctionBegin;
199409d27a7eSBarry Smith   ierr = PLA_Finalize();CHKERRQ(ierr);
199509d27a7eSBarry Smith   PetscFunctionReturn(0);
199609d27a7eSBarry Smith }
199709d27a7eSBarry Smith 
199809d27a7eSBarry Smith #undef __FUNCT__
199909d27a7eSBarry Smith #define __FUNCT__ "PetscPLAPACKInitializePackage"
200009d27a7eSBarry Smith /*@C
200109d27a7eSBarry Smith   PetscPLAPACKInitializePackage - This function initializes everything in the Petsc interface to PLAPACK. It is
200209d27a7eSBarry Smith   called from PetscDLLibraryRegister() when using dynamic libraries, and on the call to PetscInitialize()
200309d27a7eSBarry Smith   when using static libraries.
200409d27a7eSBarry Smith 
200509d27a7eSBarry Smith   Input Parameter:
200609d27a7eSBarry Smith   path - The dynamic library path, or PETSC_NULL
200709d27a7eSBarry Smith 
200809d27a7eSBarry Smith   Level: developer
200909d27a7eSBarry Smith 
201009d27a7eSBarry Smith .keywords: Petsc, initialize, package, PLAPACK
201109d27a7eSBarry Smith .seealso: PetscInitializePackage(), PetscInitialize()
201209d27a7eSBarry Smith @*/
201309d27a7eSBarry Smith PetscErrorCode PETSC_DLLEXPORT PetscPLAPACKInitializePackage(const char path[])
201409d27a7eSBarry Smith {
201504fea9ffSBarry Smith   MPI_Comm       comm = PETSC_COMM_WORLD;
2016*aeccfd6fSBarry Smith   PetscMPIInt    size,nb_alg;
201709d27a7eSBarry Smith   PetscErrorCode ierr;
2018*aeccfd6fSBarry Smith   PetscTruth     ierror;
201909d27a7eSBarry Smith 
202009d27a7eSBarry Smith   PetscFunctionBegin;
202109d27a7eSBarry Smith   if (!PLA_Initialized(PETSC_NULL)) {
202209d27a7eSBarry Smith 
202309d27a7eSBarry Smith     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
202404fea9ffSBarry Smith     Plapack_nprows = 1;
202504fea9ffSBarry Smith     Plapack_npcols = size;
202609d27a7eSBarry Smith 
2027*aeccfd6fSBarry Smith     ierr = PetscOptionsBegin(comm,PETSC_NULL,"PLAPACK Options","Mat");CHKERRQ(ierr);
2028*aeccfd6fSBarry Smith       ierr = PetscOptionsInt("-plapack_nprows","row dimension of 2D processor mesh","None",Plapack_nprows,&Plapack_nprows,PETSC_NULL);CHKERRQ(ierr);
2029*aeccfd6fSBarry Smith       ierr = PetscOptionsInt("-plapack_npcols","column dimension of 2D processor mesh","None",Plapack_npcols,&Plapack_npcols,PETSC_NULL);CHKERRQ(ierr);
2030*aeccfd6fSBarry Smith       ierr = PetscOptionsTruth("-plapack_ckerror","error checking flag","None",ierror,&ierror,PETSC_NULL);CHKERRQ(ierr);
2031*aeccfd6fSBarry Smith       if (ierror){
2032*aeccfd6fSBarry Smith 	ierr = PLA_Set_error_checking(ierror,PETSC_TRUE,PETSC_TRUE,PETSC_FALSE );CHKERRQ(ierr);
2033*aeccfd6fSBarry Smith       } else {
2034*aeccfd6fSBarry Smith 	ierr = PLA_Set_error_checking(ierror,PETSC_FALSE,PETSC_FALSE,PETSC_FALSE );CHKERRQ(ierr);
2035*aeccfd6fSBarry Smith       }
2036*aeccfd6fSBarry Smith 
2037*aeccfd6fSBarry Smith       nb_alg = 0;
2038*aeccfd6fSBarry Smith       ierr = PetscOptionsInt("-mat_plapack_nb_alg","algorithmic block size","None",nb_alg,&nb_alg,PETSC_NULL);CHKERRQ(ierr);
2039*aeccfd6fSBarry Smith       if (nb_alg) {
2040*aeccfd6fSBarry Smith         ierr = pla_Environ_set_nb_alg (PLA_OP_ALL_ALG,nb_alg);CHKERRQ(ierr);
2041*aeccfd6fSBarry Smith       }
2042*aeccfd6fSBarry Smith     PetscOptionsEnd();
2043*aeccfd6fSBarry Smith 
204404fea9ffSBarry Smith     ierr = PLA_Comm_1D_to_2D(comm,Plapack_nprows,Plapack_npcols,&Plapack_comm_2d);CHKERRQ(ierr);
204504fea9ffSBarry Smith     ierr = PLA_Init(Plapack_comm_2d);CHKERRQ(ierr);
204609d27a7eSBarry Smith     ierr = PetscRegisterFinalize(PetscPLAPACKFinalizePackage);CHKERRQ(ierr);
204709d27a7eSBarry Smith   }
204809d27a7eSBarry Smith   PetscFunctionReturn(0);
204909d27a7eSBarry Smith }
205090ace30eSBarry Smith 
205190ace30eSBarry Smith 
205209d27a7eSBarry Smith #endif
2053