xref: /petsc/src/mat/impls/dense/mpi/mpidense.c (revision eae6fb2ecf7051c084a928db0940ac509bd7cf3e)
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*/
9*eae6fb2eSBarry Smith #if defined(PETSC_HAVE_PLAPACK)
10*eae6fb2eSBarry Smith static PetscMPIInt Plapack_nprows,Plapack_npcols,Plapack_ierror,Plapack_nb_alg;
11*eae6fb2eSBarry Smith static MPI_Comm Plapack_comm_2d;
12*eae6fb2eSBarry Smith #endif
138965ea79SLois Curfman McInnes 
14ba8c8a56SBarry Smith #undef __FUNCT__
15ab92ecdeSBarry Smith #define __FUNCT__ "MatDenseGetLocalMatrix"
16ab92ecdeSBarry Smith /*@
17ab92ecdeSBarry Smith 
18ab92ecdeSBarry Smith       MatDenseGetLocalMatrix - For a MATMPIDENSE or MATSEQDENSE matrix returns the sequential
19ab92ecdeSBarry Smith               matrix that represents the operator. For sequential matrices it returns itself.
20ab92ecdeSBarry Smith 
21ab92ecdeSBarry Smith     Input Parameter:
22ab92ecdeSBarry Smith .      A - the Seq or MPI dense matrix
23ab92ecdeSBarry Smith 
24ab92ecdeSBarry Smith     Output Parameter:
25ab92ecdeSBarry Smith .      B - the inner matrix
26ab92ecdeSBarry Smith 
278e6c10adSSatish Balay     Level: intermediate
288e6c10adSSatish Balay 
29ab92ecdeSBarry Smith @*/
30ab92ecdeSBarry Smith PetscErrorCode MatDenseGetLocalMatrix(Mat A,Mat *B)
31ab92ecdeSBarry Smith {
32ab92ecdeSBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
33ab92ecdeSBarry Smith   PetscErrorCode ierr;
34ab92ecdeSBarry Smith   PetscTruth     flg;
35ab92ecdeSBarry Smith 
36ab92ecdeSBarry Smith   PetscFunctionBegin;
37ab92ecdeSBarry Smith   ierr = PetscTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr);
38ab92ecdeSBarry Smith   if (flg) {
39ab92ecdeSBarry Smith     *B = mat->A;
40ab92ecdeSBarry Smith   } else {
41ab92ecdeSBarry Smith     *B = A;
42ab92ecdeSBarry Smith   }
43ab92ecdeSBarry Smith   PetscFunctionReturn(0);
44ab92ecdeSBarry Smith }
45ab92ecdeSBarry Smith 
46ab92ecdeSBarry Smith #undef __FUNCT__
47ba8c8a56SBarry Smith #define __FUNCT__ "MatGetRow_MPIDense"
48ba8c8a56SBarry Smith PetscErrorCode MatGetRow_MPIDense(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
49ba8c8a56SBarry Smith {
50ba8c8a56SBarry Smith   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
51ba8c8a56SBarry Smith   PetscErrorCode ierr;
52899cda47SBarry Smith   PetscInt       lrow,rstart = A->rmap.rstart,rend = A->rmap.rend;
53ba8c8a56SBarry Smith 
54ba8c8a56SBarry Smith   PetscFunctionBegin;
55ba8c8a56SBarry Smith   if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows")
56ba8c8a56SBarry Smith   lrow = row - rstart;
57ba8c8a56SBarry Smith   ierr = MatGetRow(mat->A,lrow,nz,(const PetscInt **)idx,(const PetscScalar **)v);CHKERRQ(ierr);
58ba8c8a56SBarry Smith   PetscFunctionReturn(0);
59ba8c8a56SBarry Smith }
60ba8c8a56SBarry Smith 
61ba8c8a56SBarry Smith #undef __FUNCT__
62ba8c8a56SBarry Smith #define __FUNCT__ "MatRestoreRow_MPIDense"
63ba8c8a56SBarry Smith PetscErrorCode MatRestoreRow_MPIDense(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
64ba8c8a56SBarry Smith {
65ba8c8a56SBarry Smith   PetscErrorCode ierr;
66ba8c8a56SBarry Smith 
67ba8c8a56SBarry Smith   PetscFunctionBegin;
68ba8c8a56SBarry Smith   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
69ba8c8a56SBarry Smith   if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);}
70ba8c8a56SBarry Smith   PetscFunctionReturn(0);
71ba8c8a56SBarry Smith }
72ba8c8a56SBarry Smith 
730de54da6SSatish Balay EXTERN_C_BEGIN
744a2ae208SSatish Balay #undef __FUNCT__
754a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonalBlock_MPIDense"
76be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B)
770de54da6SSatish Balay {
780de54da6SSatish Balay   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
796849ba73SBarry Smith   PetscErrorCode ierr;
80899cda47SBarry Smith   PetscInt       m = A->rmap.n,rstart = A->rmap.rstart;
8187828ca2SBarry Smith   PetscScalar    *array;
820de54da6SSatish Balay   MPI_Comm       comm;
830de54da6SSatish Balay 
840de54da6SSatish Balay   PetscFunctionBegin;
85899cda47SBarry Smith   if (A->rmap.N != A->cmap.N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported.");
860de54da6SSatish Balay 
870de54da6SSatish Balay   /* The reuse aspect is not implemented efficiently */
880de54da6SSatish Balay   if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);}
890de54da6SSatish Balay 
900de54da6SSatish Balay   ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr);
910de54da6SSatish Balay   ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr);
92f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,B);CHKERRQ(ierr);
93f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*B,m,m,m,m);CHKERRQ(ierr);
947adad957SLisandro Dalcin   ierr = MatSetType(*B,((PetscObject)mdn->A)->type_name);CHKERRQ(ierr);
955c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr);
960de54da6SSatish Balay   ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr);
970de54da6SSatish Balay   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
980de54da6SSatish Balay   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
990de54da6SSatish Balay 
1000de54da6SSatish Balay   *iscopy = PETSC_TRUE;
1010de54da6SSatish Balay   PetscFunctionReturn(0);
1020de54da6SSatish Balay }
1030de54da6SSatish Balay EXTERN_C_END
1040de54da6SSatish Balay 
1054a2ae208SSatish Balay #undef __FUNCT__
1064a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_MPIDense"
10713f74950SBarry Smith PetscErrorCode MatSetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
1088965ea79SLois Curfman McInnes {
10939b7565bSBarry Smith   Mat_MPIDense   *A = (Mat_MPIDense*)mat->data;
110dfbe8321SBarry Smith   PetscErrorCode ierr;
111899cda47SBarry Smith   PetscInt       i,j,rstart = mat->rmap.rstart,rend = mat->rmap.rend,row;
112273d9f13SBarry Smith   PetscTruth     roworiented = A->roworiented;
1138965ea79SLois Curfman McInnes 
1143a40ed3dSBarry Smith   PetscFunctionBegin;
1158965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
1165ef9f2a5SBarry Smith     if (idxm[i] < 0) continue;
117899cda47SBarry Smith     if (idxm[i] >= mat->rmap.N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
1188965ea79SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
1198965ea79SLois Curfman McInnes       row = idxm[i] - rstart;
12039b7565bSBarry Smith       if (roworiented) {
12139b7565bSBarry Smith         ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr);
1223a40ed3dSBarry Smith       } else {
1238965ea79SLois Curfman McInnes         for (j=0; j<n; j++) {
1245ef9f2a5SBarry Smith           if (idxn[j] < 0) continue;
125899cda47SBarry Smith           if (idxn[j] >= mat->cmap.N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
12639b7565bSBarry Smith           ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr);
12739b7565bSBarry Smith         }
1288965ea79SLois Curfman McInnes       }
1293a40ed3dSBarry Smith     } else {
1303782ba37SSatish Balay       if (!A->donotstash) {
13139b7565bSBarry Smith         if (roworiented) {
1328798bf22SSatish Balay           ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr);
133d36fbae8SSatish Balay         } else {
1348798bf22SSatish Balay           ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr);
13539b7565bSBarry Smith         }
136b49de8d1SLois Curfman McInnes       }
137b49de8d1SLois Curfman McInnes     }
1383782ba37SSatish Balay   }
1393a40ed3dSBarry Smith   PetscFunctionReturn(0);
140b49de8d1SLois Curfman McInnes }
141b49de8d1SLois Curfman McInnes 
1424a2ae208SSatish Balay #undef __FUNCT__
1434a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_MPIDense"
14413f74950SBarry Smith PetscErrorCode MatGetValues_MPIDense(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
145b49de8d1SLois Curfman McInnes {
146b49de8d1SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
147dfbe8321SBarry Smith   PetscErrorCode ierr;
148899cda47SBarry Smith   PetscInt       i,j,rstart = mat->rmap.rstart,rend = mat->rmap.rend,row;
149b49de8d1SLois Curfman McInnes 
1503a40ed3dSBarry Smith   PetscFunctionBegin;
151b49de8d1SLois Curfman McInnes   for (i=0; i<m; i++) {
15297e567efSBarry Smith     if (idxm[i] < 0) continue; /* SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); */
153899cda47SBarry Smith     if (idxm[i] >= mat->rmap.N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
154b49de8d1SLois Curfman McInnes     if (idxm[i] >= rstart && idxm[i] < rend) {
155b49de8d1SLois Curfman McInnes       row = idxm[i] - rstart;
156b49de8d1SLois Curfman McInnes       for (j=0; j<n; j++) {
15797e567efSBarry Smith         if (idxn[j] < 0) continue; /* SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); */
158899cda47SBarry Smith         if (idxn[j] >= mat->cmap.N) {
15929bbc08cSBarry Smith           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
160a8c6a408SBarry Smith         }
161b49de8d1SLois Curfman McInnes         ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr);
162b49de8d1SLois Curfman McInnes       }
163a8c6a408SBarry Smith     } else {
16429bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"Only local values currently supported");
1658965ea79SLois Curfman McInnes     }
1668965ea79SLois Curfman McInnes   }
1673a40ed3dSBarry Smith   PetscFunctionReturn(0);
1688965ea79SLois Curfman McInnes }
1698965ea79SLois Curfman McInnes 
1704a2ae208SSatish Balay #undef __FUNCT__
1714a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_MPIDense"
172dfbe8321SBarry Smith PetscErrorCode MatGetArray_MPIDense(Mat A,PetscScalar *array[])
173ff14e315SSatish Balay {
174ff14e315SSatish Balay   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
175dfbe8321SBarry Smith   PetscErrorCode ierr;
176ff14e315SSatish Balay 
1773a40ed3dSBarry Smith   PetscFunctionBegin;
178ff14e315SSatish Balay   ierr = MatGetArray(a->A,array);CHKERRQ(ierr);
1793a40ed3dSBarry Smith   PetscFunctionReturn(0);
180ff14e315SSatish Balay }
181ff14e315SSatish Balay 
1824a2ae208SSatish Balay #undef __FUNCT__
1834a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_MPIDense"
18413f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
185ca3fa75bSLois Curfman McInnes {
186ca3fa75bSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data,*newmatd;
187ca3fa75bSLois Curfman McInnes   Mat_SeqDense   *lmat = (Mat_SeqDense*)mat->A->data;
1886849ba73SBarry Smith   PetscErrorCode ierr;
18913f74950SBarry Smith   PetscInt       i,j,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols;
19087828ca2SBarry Smith   PetscScalar    *av,*bv,*v = lmat->v;
191ca3fa75bSLois Curfman McInnes   Mat            newmat;
192ca3fa75bSLois Curfman McInnes 
193ca3fa75bSLois Curfman McInnes   PetscFunctionBegin;
194ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
195ca3fa75bSLois Curfman McInnes   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
196b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
197b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
198ca3fa75bSLois Curfman McInnes 
199ca3fa75bSLois Curfman McInnes   /* No parallel redistribution currently supported! Should really check each index set
2007eba5e9cSLois Curfman McInnes      to comfirm that it is OK.  ... Currently supports only submatrix same partitioning as
2017eba5e9cSLois Curfman McInnes      original matrix! */
202ca3fa75bSLois Curfman McInnes 
203ca3fa75bSLois Curfman McInnes   ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr);
2047eba5e9cSLois Curfman McInnes   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
205ca3fa75bSLois Curfman McInnes 
206ca3fa75bSLois Curfman McInnes   /* Check submatrix call */
207ca3fa75bSLois Curfman McInnes   if (scall == MAT_REUSE_MATRIX) {
20829bbc08cSBarry Smith     /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */
2097eba5e9cSLois Curfman McInnes     /* Really need to test rows and column sizes! */
210ca3fa75bSLois Curfman McInnes     newmat = *B;
211ca3fa75bSLois Curfman McInnes   } else {
212ca3fa75bSLois Curfman McInnes     /* Create and fill new matrix */
2137adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr);
214f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,cs,PETSC_DECIDE,ncols);CHKERRQ(ierr);
2157adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
216878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
217ca3fa75bSLois Curfman McInnes   }
218ca3fa75bSLois Curfman McInnes 
219ca3fa75bSLois Curfman McInnes   /* Now extract the data pointers and do the copy, column at a time */
220ca3fa75bSLois Curfman McInnes   newmatd = (Mat_MPIDense*)newmat->data;
221ca3fa75bSLois Curfman McInnes   bv      = ((Mat_SeqDense *)newmatd->A->data)->v;
222ca3fa75bSLois Curfman McInnes 
223ca3fa75bSLois Curfman McInnes   for (i=0; i<ncols; i++) {
224ca3fa75bSLois Curfman McInnes     av = v + nlrows*icol[i];
225ca3fa75bSLois Curfman McInnes     for (j=0; j<nrows; j++) {
2267eba5e9cSLois Curfman McInnes       *bv++ = av[irow[j] - rstart];
227ca3fa75bSLois Curfman McInnes     }
228ca3fa75bSLois Curfman McInnes   }
229ca3fa75bSLois Curfman McInnes 
230ca3fa75bSLois Curfman McInnes   /* Assemble the matrices so that the correct flags are set */
231ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
232ca3fa75bSLois Curfman McInnes   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
233ca3fa75bSLois Curfman McInnes 
234ca3fa75bSLois Curfman McInnes   /* Free work space */
235ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
236ca3fa75bSLois Curfman McInnes   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
237ca3fa75bSLois Curfman McInnes   *B = newmat;
238ca3fa75bSLois Curfman McInnes   PetscFunctionReturn(0);
239ca3fa75bSLois Curfman McInnes }
240ca3fa75bSLois Curfman McInnes 
2414a2ae208SSatish Balay #undef __FUNCT__
2424a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_MPIDense"
243dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_MPIDense(Mat A,PetscScalar *array[])
244ff14e315SSatish Balay {
2453a40ed3dSBarry Smith   PetscFunctionBegin;
2463a40ed3dSBarry Smith   PetscFunctionReturn(0);
247ff14e315SSatish Balay }
248ff14e315SSatish Balay 
2494a2ae208SSatish Balay #undef __FUNCT__
2504a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyBegin_MPIDense"
251dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode)
2528965ea79SLois Curfman McInnes {
25339ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
2547adad957SLisandro Dalcin   MPI_Comm       comm = ((PetscObject)mat)->comm;
255dfbe8321SBarry Smith   PetscErrorCode ierr;
25613f74950SBarry Smith   PetscInt       nstash,reallocs;
2578965ea79SLois Curfman McInnes   InsertMode     addv;
2588965ea79SLois Curfman McInnes 
2593a40ed3dSBarry Smith   PetscFunctionBegin;
2608965ea79SLois Curfman McInnes   /* make sure all processors are either in INSERTMODE or ADDMODE */
261ca161407SBarry Smith   ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr);
2627056b6fcSBarry Smith   if (addv == (ADD_VALUES|INSERT_VALUES)) {
26329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs");
2648965ea79SLois Curfman McInnes   }
265e0fa3b82SLois Curfman McInnes   mat->insertmode = addv; /* in case this processor had no cache */
2668965ea79SLois Curfman McInnes 
2671e2582c4SBarry Smith   ierr = MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap.range);CHKERRQ(ierr);
2688798bf22SSatish Balay   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
269ae15b995SBarry Smith   ierr = PetscInfo2(mdn->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);CHKERRQ(ierr);
2703a40ed3dSBarry Smith   PetscFunctionReturn(0);
2718965ea79SLois Curfman McInnes }
2728965ea79SLois Curfman McInnes 
2734a2ae208SSatish Balay #undef __FUNCT__
2744a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_MPIDense"
275dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode)
2768965ea79SLois Curfman McInnes {
27739ddd567SLois Curfman McInnes   Mat_MPIDense    *mdn=(Mat_MPIDense*)mat->data;
2786849ba73SBarry Smith   PetscErrorCode  ierr;
27913f74950SBarry Smith   PetscInt        i,*row,*col,flg,j,rstart,ncols;
28013f74950SBarry Smith   PetscMPIInt     n;
28187828ca2SBarry Smith   PetscScalar     *val;
282e0fa3b82SLois Curfman McInnes   InsertMode      addv=mat->insertmode;
2838965ea79SLois Curfman McInnes 
2843a40ed3dSBarry Smith   PetscFunctionBegin;
2858965ea79SLois Curfman McInnes   /*  wait on receives */
2867ef1d9bdSSatish Balay   while (1) {
2878798bf22SSatish Balay     ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
2887ef1d9bdSSatish Balay     if (!flg) break;
2898965ea79SLois Curfman McInnes 
2907ef1d9bdSSatish Balay     for (i=0; i<n;) {
2917ef1d9bdSSatish Balay       /* Now identify the consecutive vals belonging to the same row */
2927ef1d9bdSSatish Balay       for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
2937ef1d9bdSSatish Balay       if (j < n) ncols = j-i;
2947ef1d9bdSSatish Balay       else       ncols = n-i;
2957ef1d9bdSSatish Balay       /* Now assemble all these values with a single function call */
2967ef1d9bdSSatish Balay       ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
2977ef1d9bdSSatish Balay       i = j;
2988965ea79SLois Curfman McInnes     }
2997ef1d9bdSSatish Balay   }
3008798bf22SSatish Balay   ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
3018965ea79SLois Curfman McInnes 
30239ddd567SLois Curfman McInnes   ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr);
30339ddd567SLois Curfman McInnes   ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr);
3048965ea79SLois Curfman McInnes 
3056d4a8577SBarry Smith   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
30639ddd567SLois Curfman McInnes     ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
3078965ea79SLois Curfman McInnes   }
3083a40ed3dSBarry Smith   PetscFunctionReturn(0);
3098965ea79SLois Curfman McInnes }
3108965ea79SLois Curfman McInnes 
3114a2ae208SSatish Balay #undef __FUNCT__
3124a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_MPIDense"
313dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPIDense(Mat A)
3148965ea79SLois Curfman McInnes {
315dfbe8321SBarry Smith   PetscErrorCode ierr;
31639ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
3173a40ed3dSBarry Smith 
3183a40ed3dSBarry Smith   PetscFunctionBegin;
3193a40ed3dSBarry Smith   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
3203a40ed3dSBarry Smith   PetscFunctionReturn(0);
3218965ea79SLois Curfman McInnes }
3228965ea79SLois Curfman McInnes 
3238965ea79SLois Curfman McInnes /* the code does not do the diagonal entries correctly unless the
3248965ea79SLois Curfman McInnes    matrix is square and the column and row owerships are identical.
3258965ea79SLois Curfman McInnes    This is a BUG. The only way to fix it seems to be to access
3263501a2bdSLois Curfman McInnes    mdn->A and mdn->B directly and not through the MatZeroRows()
3278965ea79SLois Curfman McInnes    routine.
3288965ea79SLois Curfman McInnes */
3294a2ae208SSatish Balay #undef __FUNCT__
3304a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_MPIDense"
331f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_MPIDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
3328965ea79SLois Curfman McInnes {
33339ddd567SLois Curfman McInnes   Mat_MPIDense   *l = (Mat_MPIDense*)A->data;
3346849ba73SBarry Smith   PetscErrorCode ierr;
335899cda47SBarry Smith   PetscInt       i,*owners = A->rmap.range;
33613f74950SBarry Smith   PetscInt       *nprocs,j,idx,nsends;
33713f74950SBarry Smith   PetscInt       nmax,*svalues,*starts,*owner,nrecvs;
3387adad957SLisandro Dalcin   PetscInt       *rvalues,tag = ((PetscObject)A)->tag,count,base,slen,*source;
33913f74950SBarry Smith   PetscInt       *lens,*lrows,*values;
34013f74950SBarry Smith   PetscMPIInt    n,imdex,rank = l->rank,size = l->size;
3417adad957SLisandro Dalcin   MPI_Comm       comm = ((PetscObject)A)->comm;
3428965ea79SLois Curfman McInnes   MPI_Request    *send_waits,*recv_waits;
3438965ea79SLois Curfman McInnes   MPI_Status     recv_status,*send_status;
34435d8aa7fSBarry Smith   PetscTruth     found;
3458965ea79SLois Curfman McInnes 
3463a40ed3dSBarry Smith   PetscFunctionBegin;
3478965ea79SLois Curfman McInnes   /*  first count number of contributors to each processor */
34813f74950SBarry Smith   ierr  = PetscMalloc(2*size*sizeof(PetscInt),&nprocs);CHKERRQ(ierr);
34913f74950SBarry Smith   ierr  = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr);
35013f74950SBarry Smith   ierr  = PetscMalloc((N+1)*sizeof(PetscInt),&owner);CHKERRQ(ierr); /* see note*/
3518965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3528965ea79SLois Curfman McInnes     idx = rows[i];
35335d8aa7fSBarry Smith     found = PETSC_FALSE;
3548965ea79SLois Curfman McInnes     for (j=0; j<size; j++) {
3558965ea79SLois Curfman McInnes       if (idx >= owners[j] && idx < owners[j+1]) {
356c1dc657dSBarry Smith         nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
3578965ea79SLois Curfman McInnes       }
3588965ea79SLois Curfman McInnes     }
35929bbc08cSBarry Smith     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range");
3608965ea79SLois Curfman McInnes   }
361c1dc657dSBarry Smith   nsends = 0;  for (i=0; i<size; i++) { nsends += nprocs[2*i+1];}
3628965ea79SLois Curfman McInnes 
3638965ea79SLois Curfman McInnes   /* inform other processors of number of messages and max length*/
364c1dc657dSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
3658965ea79SLois Curfman McInnes 
3668965ea79SLois Curfman McInnes   /* post receives:   */
36713f74950SBarry Smith   ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(PetscInt),&rvalues);CHKERRQ(ierr);
368b0a32e0cSBarry Smith   ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr);
3698965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
37013f74950SBarry Smith     ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPIU_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr);
3718965ea79SLois Curfman McInnes   }
3728965ea79SLois Curfman McInnes 
3738965ea79SLois Curfman McInnes   /* do sends:
3748965ea79SLois Curfman McInnes       1) starts[i] gives the starting index in svalues for stuff going to
3758965ea79SLois Curfman McInnes          the ith processor
3768965ea79SLois Curfman McInnes   */
37713f74950SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&svalues);CHKERRQ(ierr);
378b0a32e0cSBarry Smith   ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr);
37913f74950SBarry Smith   ierr = PetscMalloc((size+1)*sizeof(PetscInt),&starts);CHKERRQ(ierr);
3808965ea79SLois Curfman McInnes   starts[0]  = 0;
381c1dc657dSBarry Smith   for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3828965ea79SLois Curfman McInnes   for (i=0; i<N; i++) {
3838965ea79SLois Curfman McInnes     svalues[starts[owner[i]]++] = rows[i];
3848965ea79SLois Curfman McInnes   }
3858965ea79SLois Curfman McInnes 
3868965ea79SLois Curfman McInnes   starts[0] = 0;
387c1dc657dSBarry Smith   for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
3888965ea79SLois Curfman McInnes   count = 0;
3898965ea79SLois Curfman McInnes   for (i=0; i<size; i++) {
390c1dc657dSBarry Smith     if (nprocs[2*i+1]) {
39113f74950SBarry Smith       ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPIU_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr);
3928965ea79SLois Curfman McInnes     }
3938965ea79SLois Curfman McInnes   }
394606d414cSSatish Balay   ierr = PetscFree(starts);CHKERRQ(ierr);
3958965ea79SLois Curfman McInnes 
3968965ea79SLois Curfman McInnes   base = owners[rank];
3978965ea79SLois Curfman McInnes 
3988965ea79SLois Curfman McInnes   /*  wait on receives */
39913f74950SBarry Smith   ierr   = PetscMalloc(2*(nrecvs+1)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
4008965ea79SLois Curfman McInnes   source = lens + nrecvs;
4018965ea79SLois Curfman McInnes   count  = nrecvs; slen = 0;
4028965ea79SLois Curfman McInnes   while (count) {
403ca161407SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
4048965ea79SLois Curfman McInnes     /* unpack receives into our local space */
40513f74950SBarry Smith     ierr = MPI_Get_count(&recv_status,MPIU_INT,&n);CHKERRQ(ierr);
4068965ea79SLois Curfman McInnes     source[imdex]  = recv_status.MPI_SOURCE;
4078965ea79SLois Curfman McInnes     lens[imdex]    = n;
4088965ea79SLois Curfman McInnes     slen += n;
4098965ea79SLois Curfman McInnes     count--;
4108965ea79SLois Curfman McInnes   }
411606d414cSSatish Balay   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
4128965ea79SLois Curfman McInnes 
4138965ea79SLois Curfman McInnes   /* move the data into the send scatter */
41413f74950SBarry Smith   ierr = PetscMalloc((slen+1)*sizeof(PetscInt),&lrows);CHKERRQ(ierr);
4158965ea79SLois Curfman McInnes   count = 0;
4168965ea79SLois Curfman McInnes   for (i=0; i<nrecvs; i++) {
4178965ea79SLois Curfman McInnes     values = rvalues + i*nmax;
4188965ea79SLois Curfman McInnes     for (j=0; j<lens[i]; j++) {
4198965ea79SLois Curfman McInnes       lrows[count++] = values[j] - base;
4208965ea79SLois Curfman McInnes     }
4218965ea79SLois Curfman McInnes   }
422606d414cSSatish Balay   ierr = PetscFree(rvalues);CHKERRQ(ierr);
423606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
424606d414cSSatish Balay   ierr = PetscFree(owner);CHKERRQ(ierr);
425606d414cSSatish Balay   ierr = PetscFree(nprocs);CHKERRQ(ierr);
4268965ea79SLois Curfman McInnes 
4278965ea79SLois Curfman McInnes   /* actually zap the local rows */
428f4df32b1SMatthew Knepley   ierr = MatZeroRows(l->A,slen,lrows,diag);CHKERRQ(ierr);
429606d414cSSatish Balay   ierr = PetscFree(lrows);CHKERRQ(ierr);
4308965ea79SLois Curfman McInnes 
4318965ea79SLois Curfman McInnes   /* wait on sends */
4328965ea79SLois Curfman McInnes   if (nsends) {
433b0a32e0cSBarry Smith     ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr);
434ca161407SBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
435606d414cSSatish Balay     ierr = PetscFree(send_status);CHKERRQ(ierr);
4368965ea79SLois Curfman McInnes   }
437606d414cSSatish Balay   ierr = PetscFree(send_waits);CHKERRQ(ierr);
438606d414cSSatish Balay   ierr = PetscFree(svalues);CHKERRQ(ierr);
4398965ea79SLois Curfman McInnes 
4403a40ed3dSBarry Smith   PetscFunctionReturn(0);
4418965ea79SLois Curfman McInnes }
4428965ea79SLois Curfman McInnes 
4434a2ae208SSatish Balay #undef __FUNCT__
4444a2ae208SSatish Balay #define __FUNCT__ "MatMult_MPIDense"
445dfbe8321SBarry Smith PetscErrorCode MatMult_MPIDense(Mat mat,Vec xx,Vec yy)
4468965ea79SLois Curfman McInnes {
44739ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
448dfbe8321SBarry Smith   PetscErrorCode ierr;
449c456f294SBarry Smith 
4503a40ed3dSBarry Smith   PetscFunctionBegin;
451ca9f406cSSatish Balay   ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
452ca9f406cSSatish Balay   ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
45344cd7ae7SLois Curfman McInnes   ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr);
4543a40ed3dSBarry Smith   PetscFunctionReturn(0);
4558965ea79SLois Curfman McInnes }
4568965ea79SLois Curfman McInnes 
4574a2ae208SSatish Balay #undef __FUNCT__
4584a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_MPIDense"
459dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz)
4608965ea79SLois Curfman McInnes {
46139ddd567SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
462dfbe8321SBarry Smith   PetscErrorCode ierr;
463c456f294SBarry Smith 
4643a40ed3dSBarry Smith   PetscFunctionBegin;
465ca9f406cSSatish Balay   ierr = VecScatterBegin(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
466ca9f406cSSatish Balay   ierr = VecScatterEnd(mdn->Mvctx,xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
46744cd7ae7SLois Curfman McInnes   ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr);
4683a40ed3dSBarry Smith   PetscFunctionReturn(0);
4698965ea79SLois Curfman McInnes }
4708965ea79SLois Curfman McInnes 
4714a2ae208SSatish Balay #undef __FUNCT__
4724a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_MPIDense"
473dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy)
474096963f5SLois Curfman McInnes {
475096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
476dfbe8321SBarry Smith   PetscErrorCode ierr;
47787828ca2SBarry Smith   PetscScalar    zero = 0.0;
478096963f5SLois Curfman McInnes 
4793a40ed3dSBarry Smith   PetscFunctionBegin;
4802dcb1b2aSMatthew Knepley   ierr = VecSet(yy,zero);CHKERRQ(ierr);
4817c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
482ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
483ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4843a40ed3dSBarry Smith   PetscFunctionReturn(0);
485096963f5SLois Curfman McInnes }
486096963f5SLois Curfman McInnes 
4874a2ae208SSatish Balay #undef __FUNCT__
4884a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_MPIDense"
489dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz)
490096963f5SLois Curfman McInnes {
491096963f5SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
492dfbe8321SBarry Smith   PetscErrorCode ierr;
493096963f5SLois Curfman McInnes 
4943a40ed3dSBarry Smith   PetscFunctionBegin;
4953501a2bdSLois Curfman McInnes   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
4967c922b88SBarry Smith   ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr);
497ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
498ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4993a40ed3dSBarry Smith   PetscFunctionReturn(0);
500096963f5SLois Curfman McInnes }
501096963f5SLois Curfman McInnes 
5024a2ae208SSatish Balay #undef __FUNCT__
5034a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_MPIDense"
504dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPIDense(Mat A,Vec v)
5058965ea79SLois Curfman McInnes {
50639ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
507096963f5SLois Curfman McInnes   Mat_SeqDense   *aloc = (Mat_SeqDense*)a->A->data;
508dfbe8321SBarry Smith   PetscErrorCode ierr;
509899cda47SBarry Smith   PetscInt       len,i,n,m = A->rmap.n,radd;
51087828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
511ed3cc1f0SBarry Smith 
5123a40ed3dSBarry Smith   PetscFunctionBegin;
5132dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
5141ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
515096963f5SLois Curfman McInnes   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
516899cda47SBarry Smith   if (n != A->rmap.N) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
517899cda47SBarry Smith   len  = PetscMin(a->A->rmap.n,a->A->cmap.n);
518899cda47SBarry Smith   radd = A->rmap.rstart*m;
51944cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
520096963f5SLois Curfman McInnes     x[i] = aloc->v[radd + i*m + i];
521096963f5SLois Curfman McInnes   }
5221ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
5233a40ed3dSBarry Smith   PetscFunctionReturn(0);
5248965ea79SLois Curfman McInnes }
5258965ea79SLois Curfman McInnes 
5264a2ae208SSatish Balay #undef __FUNCT__
5274a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIDense"
528dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIDense(Mat mat)
5298965ea79SLois Curfman McInnes {
5303501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)mat->data;
531dfbe8321SBarry Smith   PetscErrorCode ierr;
53201b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
53301b82886SBarry Smith   Mat_Plapack   *lu=(Mat_Plapack*)(mat->spptr);
53401b82886SBarry Smith #endif
535ed3cc1f0SBarry Smith 
5363a40ed3dSBarry Smith   PetscFunctionBegin;
53794d884c6SBarry Smith 
538aa482453SBarry Smith #if defined(PETSC_USE_LOG)
539899cda47SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->rmap.N,mat->cmap.N);
5408965ea79SLois Curfman McInnes #endif
5418798bf22SSatish Balay   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
5423501a2bdSLois Curfman McInnes   ierr = MatDestroy(mdn->A);CHKERRQ(ierr);
54305b42c5fSBarry Smith   if (mdn->lvec)   {ierr = VecDestroy(mdn->lvec);CHKERRQ(ierr);}
54405b42c5fSBarry Smith   if (mdn->Mvctx)  {ierr = VecScatterDestroy(mdn->Mvctx);CHKERRQ(ierr);}
54501b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
54601b82886SBarry Smith   if (lu->CleanUpPlapack) {
54762b4c0b3SBarry Smith     ierr = PLA_Obj_free(&lu->A);CHKERRQ(ierr);
54862b4c0b3SBarry Smith     ierr = PLA_Obj_free (&lu->pivots);CHKERRQ(ierr);
54962b4c0b3SBarry Smith     ierr = PLA_Temp_free(&lu->templ);CHKERRQ(ierr);
55009d27a7eSBarry Smith 
55101b82886SBarry Smith 
55201b82886SBarry Smith     ierr = ISDestroy(lu->is_pla);CHKERRQ(ierr);
55301b82886SBarry Smith     ierr = ISDestroy(lu->is_petsc);CHKERRQ(ierr);
55401b82886SBarry Smith     ierr = VecScatterDestroy(lu->ctx);CHKERRQ(ierr);
555622d7880SLois Curfman McInnes   }
55601b82886SBarry Smith #endif
55701b82886SBarry Smith 
558606d414cSSatish Balay   ierr = PetscFree(mdn);CHKERRQ(ierr);
559dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
560901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);CHKERRQ(ierr);
561901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
5624ae313f4SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C","",PETSC_NULL);CHKERRQ(ierr);
5634ae313f4SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C","",PETSC_NULL);CHKERRQ(ierr);
5644ae313f4SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C","",PETSC_NULL);CHKERRQ(ierr);
5653a40ed3dSBarry Smith   PetscFunctionReturn(0);
5668965ea79SLois Curfman McInnes }
56739ddd567SLois Curfman McInnes 
5684a2ae208SSatish Balay #undef __FUNCT__
5694a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_Binary"
5706849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_Binary(Mat mat,PetscViewer viewer)
5718965ea79SLois Curfman McInnes {
57239ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
573dfbe8321SBarry Smith   PetscErrorCode    ierr;
574aa05aa95SBarry Smith   PetscViewerFormat format;
575aa05aa95SBarry Smith   int               fd;
576aa05aa95SBarry Smith   PetscInt          header[4],mmax,N = mat->cmap.N,i,j,m,k;
577aa05aa95SBarry Smith   PetscMPIInt       rank,tag  = ((PetscObject)viewer)->tag,size;
578578230a0SSatish Balay   PetscScalar       *work,*v,*vv;
579aa05aa95SBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)mdn->A->data;
580aa05aa95SBarry Smith   MPI_Status        status;
5817056b6fcSBarry Smith 
5823a40ed3dSBarry Smith   PetscFunctionBegin;
58339ddd567SLois Curfman McInnes   if (mdn->size == 1) {
58439ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
585aa05aa95SBarry Smith   } else {
586aa05aa95SBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
5877adad957SLisandro Dalcin     ierr = MPI_Comm_rank(((PetscObject)mat)->comm,&rank);CHKERRQ(ierr);
5887adad957SLisandro Dalcin     ierr = MPI_Comm_size(((PetscObject)mat)->comm,&size);CHKERRQ(ierr);
589aa05aa95SBarry Smith 
590aa05aa95SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
591aa05aa95SBarry Smith     if (format == PETSC_VIEWER_BINARY_NATIVE) {
592aa05aa95SBarry Smith 
593aa05aa95SBarry Smith       if (!rank) {
594aa05aa95SBarry Smith 	/* store the matrix as a dense matrix */
595aa05aa95SBarry Smith 	header[0] = MAT_FILE_COOKIE;
596aa05aa95SBarry Smith 	header[1] = mat->rmap.N;
597aa05aa95SBarry Smith 	header[2] = N;
598aa05aa95SBarry Smith 	header[3] = MATRIX_BINARY_FORMAT_DENSE;
599aa05aa95SBarry Smith 	ierr = PetscBinaryWrite(fd,header,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
600aa05aa95SBarry Smith 
601aa05aa95SBarry Smith 	/* get largest work array needed for transposing array */
602aa05aa95SBarry Smith         mmax = mat->rmap.n;
603aa05aa95SBarry Smith         for (i=1; i<size; i++) {
604aa05aa95SBarry Smith           mmax = PetscMax(mmax,mat->rmap.range[i+1] - mat->rmap.range[i]);
6058965ea79SLois Curfman McInnes         }
606aa05aa95SBarry Smith 	ierr = PetscMalloc(mmax*N*sizeof(PetscScalar),&work);CHKERRQ(ierr);
607aa05aa95SBarry Smith 
608aa05aa95SBarry Smith 	/* write out local array, by rows */
609aa05aa95SBarry Smith         m    = mat->rmap.n;
610aa05aa95SBarry Smith 	v    = a->v;
611aa05aa95SBarry Smith         for (j=0; j<N; j++) {
612aa05aa95SBarry Smith           for (i=0; i<m; i++) {
613578230a0SSatish Balay 	    work[j + i*N] = *v++;
614aa05aa95SBarry Smith 	  }
615aa05aa95SBarry Smith 	}
616aa05aa95SBarry Smith 	ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
617aa05aa95SBarry Smith         /* get largest work array to receive messages from other processes, excludes process zero */
618aa05aa95SBarry Smith         mmax = 0;
619aa05aa95SBarry Smith         for (i=1; i<size; i++) {
620aa05aa95SBarry Smith           mmax = PetscMax(mmax,mat->rmap.range[i+1] - mat->rmap.range[i]);
621aa05aa95SBarry Smith         }
622578230a0SSatish Balay 	ierr = PetscMalloc(mmax*N*sizeof(PetscScalar),&vv);CHKERRQ(ierr);
623578230a0SSatish Balay         v = vv;
624aa05aa95SBarry Smith         for (k=1; k<size; k++) {
625aa05aa95SBarry Smith           m    = mat->rmap.range[k+1] - mat->rmap.range[k];
6267adad957SLisandro Dalcin           ierr = MPI_Recv(v,m*N,MPIU_SCALAR,k,tag,((PetscObject)mat)->comm,&status);CHKERRQ(ierr);
627aa05aa95SBarry Smith 
628aa05aa95SBarry Smith           for (j=0; j<N; j++) {
629aa05aa95SBarry Smith             for (i=0; i<m; i++) {
630578230a0SSatish Balay               work[j + i*N] = *v++;
631aa05aa95SBarry Smith 	    }
632aa05aa95SBarry Smith 	  }
633aa05aa95SBarry Smith 	  ierr = PetscBinaryWrite(fd,work,m*N,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
634aa05aa95SBarry Smith         }
635aa05aa95SBarry Smith         ierr = PetscFree(work);CHKERRQ(ierr);
636578230a0SSatish Balay         ierr = PetscFree(vv);CHKERRQ(ierr);
637aa05aa95SBarry Smith       } else {
6387adad957SLisandro Dalcin         ierr = MPI_Send(a->v,mat->rmap.n*mat->cmap.N,MPIU_SCALAR,0,tag,((PetscObject)mat)->comm);CHKERRQ(ierr);
639aa05aa95SBarry Smith       }
640aa05aa95SBarry Smith     }
641aa05aa95SBarry Smith   }
6423a40ed3dSBarry Smith   PetscFunctionReturn(0);
6438965ea79SLois Curfman McInnes }
6448965ea79SLois Curfman McInnes 
6454a2ae208SSatish Balay #undef __FUNCT__
6464a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket"
6476849ba73SBarry Smith static PetscErrorCode MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
6488965ea79SLois Curfman McInnes {
64939ddd567SLois Curfman McInnes   Mat_MPIDense      *mdn = (Mat_MPIDense*)mat->data;
650dfbe8321SBarry Smith   PetscErrorCode    ierr;
65132dcc486SBarry Smith   PetscMPIInt       size = mdn->size,rank = mdn->rank;
652b0a32e0cSBarry Smith   PetscViewerType   vtype;
65332077d6dSBarry Smith   PetscTruth        iascii,isdraw;
654b0a32e0cSBarry Smith   PetscViewer       sviewer;
655f3ef73ceSBarry Smith   PetscViewerFormat format;
65601b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
65701b82886SBarry Smith   Mat_Plapack       *lu=(Mat_Plapack*)(mat->spptr);
65801b82886SBarry Smith #endif
6598965ea79SLois Curfman McInnes 
6603a40ed3dSBarry Smith   PetscFunctionBegin;
66132077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
662fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
66332077d6dSBarry Smith   if (iascii) {
664b0a32e0cSBarry Smith     ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr);
665b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
666456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
6674e220ebcSLois Curfman McInnes       MatInfo info;
668888f2ed8SSatish Balay       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
669899cda47SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] local rows %D nz %D nz alloced %D mem %D \n",rank,mat->rmap.n,
67077431f27SBarry Smith                    (PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr);
671b0a32e0cSBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6723501a2bdSLois Curfman McInnes       ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr);
67301b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
67401b82886SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"PLAPACK run parameters:\n");CHKERRQ(ierr);
67504fea9ffSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Processor mesh: nprows %d, npcols %d\n",Plapack_nprows, Plapack_npcols);CHKERRQ(ierr);
67601b82886SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Distr. block size nb: %d \n",lu->nb);CHKERRQ(ierr);
677*eae6fb2eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Error checking: %d\n",Plapack_ierror);CHKERRQ(ierr);
678*eae6fb2eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Algorithmic block size: %d\n",Plapack_nb_alg);CHKERRQ(ierr);
67901b82886SBarry Smith #endif
6803a40ed3dSBarry Smith       PetscFunctionReturn(0);
681fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
6823a40ed3dSBarry Smith       PetscFunctionReturn(0);
6838965ea79SLois Curfman McInnes     }
684f1af5d2fSBarry Smith   } else if (isdraw) {
685b0a32e0cSBarry Smith     PetscDraw  draw;
686f1af5d2fSBarry Smith     PetscTruth isnull;
687f1af5d2fSBarry Smith 
688b0a32e0cSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
689b0a32e0cSBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
690f1af5d2fSBarry Smith     if (isnull) PetscFunctionReturn(0);
691f1af5d2fSBarry Smith   }
69277ed5343SBarry Smith 
6938965ea79SLois Curfman McInnes   if (size == 1) {
69439ddd567SLois Curfman McInnes     ierr = MatView(mdn->A,viewer);CHKERRQ(ierr);
6953a40ed3dSBarry Smith   } else {
6968965ea79SLois Curfman McInnes     /* assemble the entire matrix onto first processor. */
6978965ea79SLois Curfman McInnes     Mat         A;
698899cda47SBarry Smith     PetscInt    M = mat->rmap.N,N = mat->cmap.N,m,row,i,nz;
699ba8c8a56SBarry Smith     PetscInt    *cols;
700ba8c8a56SBarry Smith     PetscScalar *vals;
7018965ea79SLois Curfman McInnes 
7027adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)mat)->comm,&A);CHKERRQ(ierr);
7038965ea79SLois Curfman McInnes     if (!rank) {
704f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr);
7053a40ed3dSBarry Smith     } else {
706f69a0ea3SMatthew Knepley       ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr);
7078965ea79SLois Curfman McInnes     }
7087adad957SLisandro Dalcin     /* Since this is a temporary matrix, MATMPIDENSE instead of ((PetscObject)A)->type_name here is probably acceptable. */
709878740d9SKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
710878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL);
71152e6d16bSBarry Smith     ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr);
7128965ea79SLois Curfman McInnes 
71339ddd567SLois Curfman McInnes     /* Copy the matrix ... This isn't the most efficient means,
71439ddd567SLois Curfman McInnes        but it's quick for now */
71551022da4SBarry Smith     A->insertmode = INSERT_VALUES;
716899cda47SBarry Smith     row = mat->rmap.rstart; m = mdn->A->rmap.n;
7178965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
718ba8c8a56SBarry Smith       ierr = MatGetRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
719ba8c8a56SBarry Smith       ierr = MatSetValues_MPIDense(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
720ba8c8a56SBarry Smith       ierr = MatRestoreRow_MPIDense(mat,row,&nz,&cols,&vals);CHKERRQ(ierr);
72139ddd567SLois Curfman McInnes       row++;
7228965ea79SLois Curfman McInnes     }
7238965ea79SLois Curfman McInnes 
7246d4a8577SBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7256d4a8577SBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
726b0a32e0cSBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
727b9b97703SBarry Smith     if (!rank) {
7286831982aSBarry Smith       ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr);
7298965ea79SLois Curfman McInnes     }
730b0a32e0cSBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
731b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
7328965ea79SLois Curfman McInnes     ierr = MatDestroy(A);CHKERRQ(ierr);
7338965ea79SLois Curfman McInnes   }
7343a40ed3dSBarry Smith   PetscFunctionReturn(0);
7358965ea79SLois Curfman McInnes }
7368965ea79SLois Curfman McInnes 
7374a2ae208SSatish Balay #undef __FUNCT__
7384a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIDense"
739dfbe8321SBarry Smith PetscErrorCode MatView_MPIDense(Mat mat,PetscViewer viewer)
7408965ea79SLois Curfman McInnes {
741dfbe8321SBarry Smith   PetscErrorCode ierr;
74232077d6dSBarry Smith   PetscTruth     iascii,isbinary,isdraw,issocket;
7438965ea79SLois Curfman McInnes 
744433994e6SBarry Smith   PetscFunctionBegin;
7450f5bd95cSBarry Smith 
74632077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
747fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
748b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
749fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
7500f5bd95cSBarry Smith 
75132077d6dSBarry Smith   if (iascii || issocket || isdraw) {
752f1af5d2fSBarry Smith     ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
7530f5bd95cSBarry Smith   } else if (isbinary) {
7543a40ed3dSBarry Smith     ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr);
7555cd90555SBarry Smith   } else {
756958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name);
7578965ea79SLois Curfman McInnes   }
7583a40ed3dSBarry Smith   PetscFunctionReturn(0);
7598965ea79SLois Curfman McInnes }
7608965ea79SLois Curfman McInnes 
7614a2ae208SSatish Balay #undef __FUNCT__
7624a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_MPIDense"
763dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info)
7648965ea79SLois Curfman McInnes {
7653501a2bdSLois Curfman McInnes   Mat_MPIDense   *mat = (Mat_MPIDense*)A->data;
7663501a2bdSLois Curfman McInnes   Mat            mdn = mat->A;
767dfbe8321SBarry Smith   PetscErrorCode ierr;
768329f5518SBarry Smith   PetscReal      isend[5],irecv[5];
7698965ea79SLois Curfman McInnes 
7703a40ed3dSBarry Smith   PetscFunctionBegin;
771899cda47SBarry Smith   info->rows_global    = (double)A->rmap.N;
772899cda47SBarry Smith   info->columns_global = (double)A->cmap.N;
773899cda47SBarry Smith   info->rows_local     = (double)A->rmap.n;
774899cda47SBarry Smith   info->columns_local  = (double)A->cmap.N;
7754e220ebcSLois Curfman McInnes   info->block_size     = 1.0;
7764e220ebcSLois Curfman McInnes   ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr);
7774e220ebcSLois Curfman McInnes   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
7784e220ebcSLois Curfman McInnes   isend[3] = info->memory;  isend[4] = info->mallocs;
7798965ea79SLois Curfman McInnes   if (flag == MAT_LOCAL) {
7804e220ebcSLois Curfman McInnes     info->nz_used      = isend[0];
7814e220ebcSLois Curfman McInnes     info->nz_allocated = isend[1];
7824e220ebcSLois Curfman McInnes     info->nz_unneeded  = isend[2];
7834e220ebcSLois Curfman McInnes     info->memory       = isend[3];
7844e220ebcSLois Curfman McInnes     info->mallocs      = isend[4];
7858965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_MAX) {
7867adad957SLisandro Dalcin     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,((PetscObject)A)->comm);CHKERRQ(ierr);
7874e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7884e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7894e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7904e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7914e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7928965ea79SLois Curfman McInnes   } else if (flag == MAT_GLOBAL_SUM) {
7937adad957SLisandro Dalcin     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
7944e220ebcSLois Curfman McInnes     info->nz_used      = irecv[0];
7954e220ebcSLois Curfman McInnes     info->nz_allocated = irecv[1];
7964e220ebcSLois Curfman McInnes     info->nz_unneeded  = irecv[2];
7974e220ebcSLois Curfman McInnes     info->memory       = irecv[3];
7984e220ebcSLois Curfman McInnes     info->mallocs      = irecv[4];
7998965ea79SLois Curfman McInnes   }
8004e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
8014e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
8024e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
8033a40ed3dSBarry Smith   PetscFunctionReturn(0);
8048965ea79SLois Curfman McInnes }
8058965ea79SLois Curfman McInnes 
8064a2ae208SSatish Balay #undef __FUNCT__
8074a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIDense"
8084e0d8c25SBarry Smith PetscErrorCode MatSetOption_MPIDense(Mat A,MatOption op,PetscTruth flg)
8098965ea79SLois Curfman McInnes {
81039ddd567SLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
811dfbe8321SBarry Smith   PetscErrorCode ierr;
8128965ea79SLois Curfman McInnes 
8133a40ed3dSBarry Smith   PetscFunctionBegin;
81412c028f9SKris Buschelman   switch (op) {
815512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
81612c028f9SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
81712c028f9SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
8184e0d8c25SBarry Smith     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
81912c028f9SKris Buschelman     break;
82012c028f9SKris Buschelman   case MAT_ROW_ORIENTED:
8214e0d8c25SBarry Smith     a->roworiented = flg;
8224e0d8c25SBarry Smith     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
82312c028f9SKris Buschelman     break;
8244e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
82512c028f9SKris Buschelman   case MAT_USE_HASH_TABLE:
826290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
82712c028f9SKris Buschelman     break;
82812c028f9SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
8294e0d8c25SBarry Smith     a->donotstash = flg;
83012c028f9SKris Buschelman     break;
83177e54ba9SKris Buschelman   case MAT_SYMMETRIC:
83277e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
8339a4540c5SBarry Smith   case MAT_HERMITIAN:
8349a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
835600fe468SBarry Smith   case MAT_IGNORE_LOWER_TRIANGULAR:
836290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
83777e54ba9SKris Buschelman     break;
83812c028f9SKris Buschelman   default:
839600fe468SBarry Smith     SETERRQ1(PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
8403a40ed3dSBarry Smith   }
8413a40ed3dSBarry Smith   PetscFunctionReturn(0);
8428965ea79SLois Curfman McInnes }
8438965ea79SLois Curfman McInnes 
8448965ea79SLois Curfman McInnes 
8454a2ae208SSatish Balay #undef __FUNCT__
8464a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_MPIDense"
847dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr)
8485b2fa520SLois Curfman McInnes {
8495b2fa520SLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8505b2fa520SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
85187828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
852dfbe8321SBarry Smith   PetscErrorCode ierr;
853899cda47SBarry Smith   PetscInt       i,j,s2a,s3a,s2,s3,m=mdn->A->rmap.n,n=mdn->A->cmap.n;
8545b2fa520SLois Curfman McInnes 
8555b2fa520SLois Curfman McInnes   PetscFunctionBegin;
85672d926a5SLois Curfman McInnes   ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr);
8575b2fa520SLois Curfman McInnes   if (ll) {
85872d926a5SLois Curfman McInnes     ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr);
859175be7b4SMatthew Knepley     if (s2a != s2) SETERRQ2(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size, %d != %d.", s2a, s2);
8601ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
8615b2fa520SLois Curfman McInnes     for (i=0; i<m; i++) {
8625b2fa520SLois Curfman McInnes       x = l[i];
8635b2fa520SLois Curfman McInnes       v = mat->v + i;
8645b2fa520SLois Curfman McInnes       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
8655b2fa520SLois Curfman McInnes     }
8661ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
867efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8685b2fa520SLois Curfman McInnes   }
8695b2fa520SLois Curfman McInnes   if (rr) {
870175be7b4SMatthew Knepley     ierr = VecGetLocalSize(rr,&s3a);CHKERRQ(ierr);
871175be7b4SMatthew Knepley     if (s3a != s3) SETERRQ2(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size, %d != %d.", s3a, s3);
872ca9f406cSSatish Balay     ierr = VecScatterBegin(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
873ca9f406cSSatish Balay     ierr = VecScatterEnd(mdn->Mvctx,rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8741ebc52fbSHong Zhang     ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr);
8755b2fa520SLois Curfman McInnes     for (i=0; i<n; i++) {
8765b2fa520SLois Curfman McInnes       x = r[i];
8775b2fa520SLois Curfman McInnes       v = mat->v + i*m;
8785b2fa520SLois Curfman McInnes       for (j=0; j<m; j++) { (*v++) *= x;}
8795b2fa520SLois Curfman McInnes     }
8801ebc52fbSHong Zhang     ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr);
881efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
8825b2fa520SLois Curfman McInnes   }
8835b2fa520SLois Curfman McInnes   PetscFunctionReturn(0);
8845b2fa520SLois Curfman McInnes }
8855b2fa520SLois Curfman McInnes 
8864a2ae208SSatish Balay #undef __FUNCT__
8874a2ae208SSatish Balay #define __FUNCT__ "MatNorm_MPIDense"
888dfbe8321SBarry Smith PetscErrorCode MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm)
889096963f5SLois Curfman McInnes {
8903501a2bdSLois Curfman McInnes   Mat_MPIDense   *mdn = (Mat_MPIDense*)A->data;
8913501a2bdSLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)mdn->A->data;
892dfbe8321SBarry Smith   PetscErrorCode ierr;
89313f74950SBarry Smith   PetscInt       i,j;
894329f5518SBarry Smith   PetscReal      sum = 0.0;
89587828ca2SBarry Smith   PetscScalar    *v = mat->v;
8963501a2bdSLois Curfman McInnes 
8973a40ed3dSBarry Smith   PetscFunctionBegin;
8983501a2bdSLois Curfman McInnes   if (mdn->size == 1) {
899064f8208SBarry Smith     ierr =  MatNorm(mdn->A,type,nrm);CHKERRQ(ierr);
9003501a2bdSLois Curfman McInnes   } else {
9013501a2bdSLois Curfman McInnes     if (type == NORM_FROBENIUS) {
902899cda47SBarry Smith       for (i=0; i<mdn->A->cmap.n*mdn->A->rmap.n; i++) {
903aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
904329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
9053501a2bdSLois Curfman McInnes #else
9063501a2bdSLois Curfman McInnes         sum += (*v)*(*v); v++;
9073501a2bdSLois Curfman McInnes #endif
9083501a2bdSLois Curfman McInnes       }
9097adad957SLisandro Dalcin       ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
910064f8208SBarry Smith       *nrm = sqrt(*nrm);
911899cda47SBarry Smith       ierr = PetscLogFlops(2*mdn->A->cmap.n*mdn->A->rmap.n);CHKERRQ(ierr);
9123a40ed3dSBarry Smith     } else if (type == NORM_1) {
913329f5518SBarry Smith       PetscReal *tmp,*tmp2;
914899cda47SBarry Smith       ierr = PetscMalloc(2*A->cmap.N*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
915899cda47SBarry Smith       tmp2 = tmp + A->cmap.N;
916899cda47SBarry Smith       ierr = PetscMemzero(tmp,2*A->cmap.N*sizeof(PetscReal));CHKERRQ(ierr);
917064f8208SBarry Smith       *nrm = 0.0;
9183501a2bdSLois Curfman McInnes       v = mat->v;
919899cda47SBarry Smith       for (j=0; j<mdn->A->cmap.n; j++) {
920899cda47SBarry Smith         for (i=0; i<mdn->A->rmap.n; i++) {
92167e560aaSBarry Smith           tmp[j] += PetscAbsScalar(*v);  v++;
9223501a2bdSLois Curfman McInnes         }
9233501a2bdSLois Curfman McInnes       }
9247adad957SLisandro Dalcin       ierr = MPI_Allreduce(tmp,tmp2,A->cmap.N,MPIU_REAL,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
925899cda47SBarry Smith       for (j=0; j<A->cmap.N; j++) {
926064f8208SBarry Smith         if (tmp2[j] > *nrm) *nrm = tmp2[j];
9273501a2bdSLois Curfman McInnes       }
928606d414cSSatish Balay       ierr = PetscFree(tmp);CHKERRQ(ierr);
929899cda47SBarry Smith       ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr);
9303a40ed3dSBarry Smith     } else if (type == NORM_INFINITY) { /* max row norm */
931329f5518SBarry Smith       PetscReal ntemp;
9323501a2bdSLois Curfman McInnes       ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr);
9337adad957SLisandro Dalcin       ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,((PetscObject)A)->comm);CHKERRQ(ierr);
9343a40ed3dSBarry Smith     } else {
93529bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"No support for two norm");
9363501a2bdSLois Curfman McInnes     }
9373501a2bdSLois Curfman McInnes   }
9383a40ed3dSBarry Smith   PetscFunctionReturn(0);
9393501a2bdSLois Curfman McInnes }
9403501a2bdSLois Curfman McInnes 
9414a2ae208SSatish Balay #undef __FUNCT__
9424a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_MPIDense"
943fc4dec0aSBarry Smith PetscErrorCode MatTranspose_MPIDense(Mat A,MatReuse reuse,Mat *matout)
9443501a2bdSLois Curfman McInnes {
9453501a2bdSLois Curfman McInnes   Mat_MPIDense   *a = (Mat_MPIDense*)A->data;
9463501a2bdSLois Curfman McInnes   Mat_SeqDense   *Aloc = (Mat_SeqDense*)a->A->data;
9473501a2bdSLois Curfman McInnes   Mat            B;
948899cda47SBarry Smith   PetscInt       M = A->rmap.N,N = A->cmap.N,m,n,*rwork,rstart = A->rmap.rstart;
9496849ba73SBarry Smith   PetscErrorCode ierr;
95013f74950SBarry Smith   PetscInt       j,i;
95187828ca2SBarry Smith   PetscScalar    *v;
9523501a2bdSLois Curfman McInnes 
9533a40ed3dSBarry Smith   PetscFunctionBegin;
954fc4dec0aSBarry Smith   if (A == *matout && M != N) {
95529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place");
9567056b6fcSBarry Smith   }
957fc4dec0aSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || A == *matout) {
9587adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
959f69a0ea3SMatthew Knepley     ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,M);CHKERRQ(ierr);
9607adad957SLisandro Dalcin     ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
961878740d9SKris Buschelman     ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr);
962fc4dec0aSBarry Smith   } else {
963fc4dec0aSBarry Smith     B = *matout;
964fc4dec0aSBarry Smith   }
9653501a2bdSLois Curfman McInnes 
966899cda47SBarry Smith   m = a->A->rmap.n; n = a->A->cmap.n; v = Aloc->v;
9671acff37aSSatish Balay   ierr = PetscMalloc(m*sizeof(PetscInt),&rwork);CHKERRQ(ierr);
9683501a2bdSLois Curfman McInnes   for (i=0; i<m; i++) rwork[i] = rstart + i;
9691acff37aSSatish Balay   for (j=0; j<n; j++) {
9703501a2bdSLois Curfman McInnes     ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr);
9713501a2bdSLois Curfman McInnes     v   += m;
9723501a2bdSLois Curfman McInnes   }
973606d414cSSatish Balay   ierr = PetscFree(rwork);CHKERRQ(ierr);
9746d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9756d4a8577SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
976fc4dec0aSBarry Smith   if (*matout != A) {
9773501a2bdSLois Curfman McInnes     *matout = B;
9783501a2bdSLois Curfman McInnes   } else {
979273d9f13SBarry Smith     ierr = MatHeaderCopy(A,B);CHKERRQ(ierr);
9803501a2bdSLois Curfman McInnes   }
9813a40ed3dSBarry Smith   PetscFunctionReturn(0);
982096963f5SLois Curfman McInnes }
983096963f5SLois Curfman McInnes 
984d9eff348SSatish Balay #include "petscblaslapack.h"
9854a2ae208SSatish Balay #undef __FUNCT__
9864a2ae208SSatish Balay #define __FUNCT__ "MatScale_MPIDense"
987f4df32b1SMatthew Knepley PetscErrorCode MatScale_MPIDense(Mat inA,PetscScalar alpha)
98844cd7ae7SLois Curfman McInnes {
98944cd7ae7SLois Curfman McInnes   Mat_MPIDense   *A = (Mat_MPIDense*)inA->data;
99044cd7ae7SLois Curfman McInnes   Mat_SeqDense   *a = (Mat_SeqDense*)A->A->data;
991f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
992efee365bSSatish Balay   PetscErrorCode ierr;
9930805154bSBarry Smith   PetscBLASInt   one = 1,nz = PetscBLASIntCast(inA->rmap.n*inA->cmap.N);
99444cd7ae7SLois Curfman McInnes 
9953a40ed3dSBarry Smith   PetscFunctionBegin;
996f4df32b1SMatthew Knepley   BLASscal_(&nz,&oalpha,a->v,&one);
997efee365bSSatish Balay   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
9983a40ed3dSBarry Smith   PetscFunctionReturn(0);
99944cd7ae7SLois Curfman McInnes }
100044cd7ae7SLois Curfman McInnes 
10016849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *);
10028965ea79SLois Curfman McInnes 
10034a2ae208SSatish Balay #undef __FUNCT__
10044a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_MPIDense"
1005dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_MPIDense(Mat A)
1006273d9f13SBarry Smith {
1007dfbe8321SBarry Smith   PetscErrorCode ierr;
1008273d9f13SBarry Smith 
1009273d9f13SBarry Smith   PetscFunctionBegin;
1010273d9f13SBarry Smith   ierr =  MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr);
1011273d9f13SBarry Smith   PetscFunctionReturn(0);
1012273d9f13SBarry Smith }
1013273d9f13SBarry Smith 
10144ae313f4SHong Zhang #undef __FUNCT__
10154ae313f4SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIDense_MPIDense"
10164ae313f4SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
10174ae313f4SHong Zhang {
10184ae313f4SHong Zhang   PetscErrorCode ierr;
10194ae313f4SHong Zhang   PetscInt       m=A->rmap.n,n=B->cmap.n;
10204ae313f4SHong Zhang   Mat            Cmat;
10214ae313f4SHong Zhang 
10224ae313f4SHong Zhang   PetscFunctionBegin;
10234ae313f4SHong 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);
10247adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)B)->comm,&Cmat);CHKERRQ(ierr);
10254ae313f4SHong Zhang   ierr = MatSetSizes(Cmat,m,n,A->rmap.N,B->cmap.N);CHKERRQ(ierr);
10264ae313f4SHong Zhang   ierr = MatSetType(Cmat,MATMPIDENSE);CHKERRQ(ierr);
10274ae313f4SHong Zhang   ierr = MatMPIDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr);
10282a6744ebSBarry Smith   ierr = MatAssemblyBegin(Cmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10290a71e23eSMatthew Knepley   ierr = MatAssemblyEnd(Cmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10304ae313f4SHong Zhang   *C = Cmat;
10314ae313f4SHong Zhang   PetscFunctionReturn(0);
10324ae313f4SHong Zhang }
10334ae313f4SHong Zhang 
103401b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
1035*eae6fb2eSBarry Smith 
1036*eae6fb2eSBarry Smith #undef __FUNCT__
1037*eae6fb2eSBarry Smith #define __FUNCT__ "MatMPIDenseCreatePlapack"
1038*eae6fb2eSBarry Smith PetscErrorCode MatMPIDenseCreatePlapack(Mat A)
1039*eae6fb2eSBarry Smith {
1040*eae6fb2eSBarry Smith   Mat_Plapack    *lu;
1041*eae6fb2eSBarry Smith   PetscErrorCode ierr;
1042*eae6fb2eSBarry Smith   PetscInt       M=A->rmap.N;
1043*eae6fb2eSBarry Smith   MPI_Comm       comm=((PetscObject)A)->comm;
1044*eae6fb2eSBarry Smith   PetscMPIInt    size;
1045*eae6fb2eSBarry Smith 
1046*eae6fb2eSBarry Smith   PetscFunctionBegin;
1047*eae6fb2eSBarry Smith   lu = (Mat_Plapack*)(A->spptr);
1048*eae6fb2eSBarry Smith 
1049*eae6fb2eSBarry Smith   /* Set default Plapack parameters */
1050*eae6fb2eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1051*eae6fb2eSBarry Smith   lu->nb     = M/size;
1052*eae6fb2eSBarry Smith   if (M - lu->nb*size) lu->nb++; /* without cyclic distribution */
1053*eae6fb2eSBarry Smith 
1054*eae6fb2eSBarry Smith   /* Set runtime options */
1055*eae6fb2eSBarry Smith   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"PLAPACK Options","Mat");CHKERRQ(ierr);
1056*eae6fb2eSBarry Smith     ierr = PetscOptionsInt("-mat_plapack_nb","block size of template vector","None",lu->nb,&lu->nb,PETSC_NULL);CHKERRQ(ierr);
1057*eae6fb2eSBarry Smith   PetscOptionsEnd();
1058*eae6fb2eSBarry Smith 
1059*eae6fb2eSBarry Smith   /* Create object distribution template */
1060*eae6fb2eSBarry Smith   lu->templ = NULL;
1061*eae6fb2eSBarry Smith   ierr = PLA_Temp_create(lu->nb, 0, &lu->templ);CHKERRQ(ierr);
1062*eae6fb2eSBarry Smith 
1063*eae6fb2eSBarry Smith   /* Set the datatype */
1064*eae6fb2eSBarry Smith #if defined(PETSC_USE_COMPLEX)
1065*eae6fb2eSBarry Smith   lu->datatype = MPI_DOUBLE_COMPLEX;
1066*eae6fb2eSBarry Smith #else
1067*eae6fb2eSBarry Smith   lu->datatype = MPI_DOUBLE;
1068*eae6fb2eSBarry Smith #endif
1069*eae6fb2eSBarry Smith 
1070*eae6fb2eSBarry Smith   ierr = PLA_Matrix_create(lu->datatype,M,M,lu->templ,PLA_ALIGN_FIRST,PLA_ALIGN_FIRST,&lu->A);CHKERRQ(ierr);
1071*eae6fb2eSBarry Smith 
1072*eae6fb2eSBarry Smith   lu->pla_solved     = PETSC_FALSE; /* MatSolve_Plapack() is called yet */
1073*eae6fb2eSBarry Smith   lu->CleanUpPlapack = PETSC_TRUE;
1074*eae6fb2eSBarry Smith   PetscFunctionReturn(0);
1075*eae6fb2eSBarry Smith }
1076*eae6fb2eSBarry Smith 
1077*eae6fb2eSBarry Smith #undef __FUNCT__
1078*eae6fb2eSBarry Smith #define __FUNCT__ "MatMPIDenseCopyToPlapack"
1079*eae6fb2eSBarry Smith PetscErrorCode MatMPIDenseCopyToPlapack(Mat A,Mat F)
1080*eae6fb2eSBarry Smith {
1081*eae6fb2eSBarry Smith   Mat_Plapack    *lu = (Mat_Plapack*)(F)->spptr;
1082*eae6fb2eSBarry Smith   PetscErrorCode ierr;
1083*eae6fb2eSBarry Smith   PetscInt       M=A->rmap.N,m=A->rmap.n,rstart,rend;
1084*eae6fb2eSBarry Smith   PetscScalar    *array;
1085*eae6fb2eSBarry Smith   PetscReal      one = 1.0;
1086*eae6fb2eSBarry Smith 
1087*eae6fb2eSBarry Smith   PetscFunctionBegin;
1088*eae6fb2eSBarry Smith   /* Copy A into lu->A */
1089*eae6fb2eSBarry Smith   ierr = PLA_Obj_set_to_zero(lu->A);CHKERRQ(ierr);
1090*eae6fb2eSBarry Smith   ierr = PLA_API_begin();CHKERRQ(ierr);
1091*eae6fb2eSBarry Smith   ierr = PLA_Obj_API_open(lu->A);CHKERRQ(ierr);
1092*eae6fb2eSBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
1093*eae6fb2eSBarry Smith   ierr = MatGetArray(A,&array);CHKERRQ(ierr);
1094*eae6fb2eSBarry Smith   ierr = PLA_API_axpy_matrix_to_global(m,M, &one,(void *)array,m,lu->A,rstart,0);CHKERRQ(ierr);
1095*eae6fb2eSBarry Smith   ierr = MatRestoreArray(A,&array);CHKERRQ(ierr);
1096*eae6fb2eSBarry Smith   ierr = PLA_Obj_API_close(lu->A);CHKERRQ(ierr);
1097*eae6fb2eSBarry Smith   ierr = PLA_API_end();CHKERRQ(ierr);
1098*eae6fb2eSBarry Smith   lu->rstart         = rstart;
1099*eae6fb2eSBarry Smith   PetscFunctionReturn(0);
1100*eae6fb2eSBarry Smith }
1101*eae6fb2eSBarry Smith 
110201b82886SBarry Smith #undef __FUNCT__
110301b82886SBarry Smith #define __FUNCT__ "MatSolve_MPIDense"
110401b82886SBarry Smith PetscErrorCode MatSolve_MPIDense(Mat A,Vec b,Vec x)
110501b82886SBarry Smith {
110604fea9ffSBarry Smith   MPI_Comm       comm = ((PetscObject)A)->comm,dummy_comm;
110701b82886SBarry Smith   Mat_Plapack    *lu = (Mat_Plapack*)A->spptr;
110801b82886SBarry Smith   PetscErrorCode ierr;
110901b82886SBarry Smith   PetscInt       M=A->rmap.N,m=A->rmap.n,rstart,i,j,*idx_pla,*idx_petsc,loc_m,loc_stride;
111001b82886SBarry Smith   PetscScalar    *array;
111101b82886SBarry Smith   PetscReal      one = 1.0;
1112aeccfd6fSBarry Smith   PetscMPIInt    size,r_rank,r_nproc,c_rank,c_nproc;;
111301b82886SBarry Smith   PLA_Obj        v_pla = NULL;
111401b82886SBarry Smith   PetscScalar    *loc_buf;
111501b82886SBarry Smith   Vec            loc_x;
111601b82886SBarry Smith 
111701b82886SBarry Smith   PetscFunctionBegin;
111801b82886SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
111901b82886SBarry Smith 
112001b82886SBarry Smith   /* Create PLAPACK vector objects, then copy b into PLAPACK b */
112162b4c0b3SBarry Smith   ierr = PLA_Mvector_create(lu->datatype,M,1,lu->templ,PLA_ALIGN_FIRST,&v_pla);CHKERRQ(ierr);
112262b4c0b3SBarry Smith   ierr = PLA_Obj_set_to_zero(v_pla);CHKERRQ(ierr);
112301b82886SBarry Smith 
112401b82886SBarry Smith   /* Copy b into rhs_pla */
112562b4c0b3SBarry Smith   ierr = PLA_API_begin();CHKERRQ(ierr);
112662b4c0b3SBarry Smith   ierr = PLA_Obj_API_open(v_pla);CHKERRQ(ierr);
112701b82886SBarry Smith   ierr = VecGetArray(b,&array);CHKERRQ(ierr);
112862b4c0b3SBarry Smith   ierr = PLA_API_axpy_vector_to_global(m,&one,(void *)array,1,v_pla,lu->rstart);CHKERRQ(ierr);
112901b82886SBarry Smith   ierr = VecRestoreArray(b,&array);CHKERRQ(ierr);
113062b4c0b3SBarry Smith   ierr = PLA_Obj_API_close(v_pla);CHKERRQ(ierr);
113162b4c0b3SBarry Smith   ierr = PLA_API_end();CHKERRQ(ierr);
113201b82886SBarry Smith 
113301b82886SBarry Smith   if (A->factor == FACTOR_LU){
113401b82886SBarry Smith     /* Apply the permutations to the right hand sides */
113562b4c0b3SBarry Smith     ierr = PLA_Apply_pivots_to_rows (v_pla,lu->pivots);CHKERRQ(ierr);
113601b82886SBarry Smith 
113701b82886SBarry Smith     /* Solve L y = b, overwriting b with y */
113862b4c0b3SBarry Smith     ierr = PLA_Trsv( PLA_LOWER_TRIANGULAR,PLA_NO_TRANSPOSE,PLA_UNIT_DIAG,lu->A,v_pla );CHKERRQ(ierr);
113901b82886SBarry Smith 
114001b82886SBarry Smith     /* Solve U x = y (=b), overwriting b with x */
114162b4c0b3SBarry Smith     ierr = PLA_Trsv( PLA_UPPER_TRIANGULAR,PLA_NO_TRANSPOSE,PLA_NONUNIT_DIAG,lu->A,v_pla );CHKERRQ(ierr);
114201b82886SBarry Smith   } else { /* FACTOR_CHOLESKY */
114362b4c0b3SBarry Smith     ierr = PLA_Trsv( PLA_LOWER_TRIANGULAR,PLA_NO_TRANSPOSE,PLA_NONUNIT_DIAG,lu->A,v_pla);CHKERRQ(ierr);
114462b4c0b3SBarry 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);
114501b82886SBarry Smith   }
114601b82886SBarry Smith 
114701b82886SBarry Smith   /* Copy PLAPACK x into Petsc vector x  */
114862b4c0b3SBarry Smith   ierr = PLA_Obj_local_length(v_pla, &loc_m);CHKERRQ(ierr);
114962b4c0b3SBarry Smith   ierr = PLA_Obj_local_buffer(v_pla, (void**)&loc_buf);CHKERRQ(ierr);
115062b4c0b3SBarry Smith   ierr = PLA_Obj_local_stride(v_pla, &loc_stride);CHKERRQ(ierr);
115101b82886SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,loc_m*loc_stride,loc_buf,&loc_x);CHKERRQ(ierr);
115201b82886SBarry Smith   if (!lu->pla_solved){
115301b82886SBarry Smith 
115404fea9ffSBarry Smith     ierr = PLA_Temp_comm_row_info(lu->templ,&dummy_comm,&r_rank,&r_nproc);CHKERRQ(ierr);
115504fea9ffSBarry Smith     ierr = PLA_Temp_comm_col_info(lu->templ,&dummy_comm,&c_rank,&c_nproc);CHKERRQ(ierr);
115601b82886SBarry Smith     /* printf(" [%d] rank: %d %d, nproc: %d %d\n",rank,r_rank,c_rank,r_nproc,c_nproc); */
115701b82886SBarry Smith 
115801b82886SBarry Smith     /* Create IS and cts for VecScatterring */
115962b4c0b3SBarry Smith     ierr = PLA_Obj_local_length(v_pla, &loc_m);CHKERRQ(ierr);
116062b4c0b3SBarry Smith     ierr = PLA_Obj_local_stride(v_pla, &loc_stride);CHKERRQ(ierr);
116101b82886SBarry Smith     ierr = PetscMalloc((2*loc_m+1)*sizeof(PetscInt),&idx_pla);CHKERRQ(ierr);
116201b82886SBarry Smith     idx_petsc = idx_pla + loc_m;
116301b82886SBarry Smith 
116401b82886SBarry Smith     rstart = (r_rank*c_nproc+c_rank)*lu->nb;
116501b82886SBarry Smith     for (i=0; i<loc_m; i+=lu->nb){
116601b82886SBarry Smith       j = 0;
116701b82886SBarry Smith       while (j < lu->nb && i+j < loc_m){
116801b82886SBarry Smith         idx_petsc[i+j] = rstart + j; j++;
116901b82886SBarry Smith       }
117001b82886SBarry Smith       rstart += size*lu->nb;
117101b82886SBarry Smith     }
117201b82886SBarry Smith 
117301b82886SBarry Smith     for (i=0; i<loc_m; i++) idx_pla[i] = i*loc_stride;
117401b82886SBarry Smith 
117501b82886SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,loc_m,idx_pla,&lu->is_pla);CHKERRQ(ierr);
117601b82886SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,loc_m,idx_petsc,&lu->is_petsc);CHKERRQ(ierr);
117701b82886SBarry Smith     ierr = PetscFree(idx_pla);CHKERRQ(ierr);
117801b82886SBarry Smith     ierr = VecScatterCreate(loc_x,lu->is_pla,x,lu->is_petsc,&lu->ctx);CHKERRQ(ierr);
117901b82886SBarry Smith   }
118001b82886SBarry Smith   ierr = VecScatterBegin(lu->ctx,loc_x,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
118101b82886SBarry Smith   ierr = VecScatterEnd(lu->ctx,loc_x,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
118201b82886SBarry Smith 
118301b82886SBarry Smith   /* Free data */
118401b82886SBarry Smith   ierr = VecDestroy(loc_x);CHKERRQ(ierr);
118562b4c0b3SBarry Smith   ierr = PLA_Obj_free(&v_pla);CHKERRQ(ierr);
118601b82886SBarry Smith 
118701b82886SBarry Smith   lu->pla_solved = PETSC_TRUE;
118801b82886SBarry Smith   PetscFunctionReturn(0);
118901b82886SBarry Smith }
119001b82886SBarry Smith 
119101b82886SBarry Smith #undef __FUNCT__
119201b82886SBarry Smith #define __FUNCT__ "MatLUFactorNumeric_MPIDense"
119301b82886SBarry Smith PetscErrorCode MatLUFactorNumeric_MPIDense(Mat A,MatFactorInfo *info,Mat *F)
119401b82886SBarry Smith {
119501b82886SBarry Smith   Mat_Plapack    *lu = (Mat_Plapack*)(*F)->spptr;
119601b82886SBarry Smith   PetscErrorCode ierr;
1197*eae6fb2eSBarry Smith   PetscInt       info_pla;
119801b82886SBarry Smith 
119901b82886SBarry Smith   PetscFunctionBegin;
120062b4c0b3SBarry Smith   ierr = PLA_Obj_set_to_zero(lu->A);CHKERRQ(ierr);
120101b82886SBarry Smith 
1202*eae6fb2eSBarry Smith   /* Copy A into F->lu->A */
1203*eae6fb2eSBarry Smith   ierr = MatMPIDenseCopyToPlapack(A,*F);CHKERRQ(ierr);
120401b82886SBarry Smith 
120501b82886SBarry Smith   /* Factor P A -> L U overwriting lower triangular portion of A with L, upper, U */
120601b82886SBarry Smith   info_pla = PLA_LU(lu->A,lu->pivots);
1207*eae6fb2eSBarry Smith   if (info_pla) SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot encountered at row %d from PLA_LU()",info_pla);
120801b82886SBarry Smith 
1209*eae6fb2eSBarry Smith   (*F)->assembled    = PETSC_TRUE;
121001b82886SBarry Smith   PetscFunctionReturn(0);
121101b82886SBarry Smith }
121201b82886SBarry Smith 
121301b82886SBarry Smith #undef __FUNCT__
121401b82886SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_MPIDense"
121501b82886SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_MPIDense(Mat A,MatFactorInfo *info,Mat *F)
121601b82886SBarry Smith {
121701b82886SBarry Smith   Mat_Plapack    *lu = (Mat_Plapack*)(*F)->spptr;
121801b82886SBarry Smith   PetscErrorCode ierr;
1219*eae6fb2eSBarry Smith   PetscInt       info_pla;
122001b82886SBarry Smith 
122101b82886SBarry Smith   PetscFunctionBegin;
1222e0fbc2eeSBarry Smith 
1223*eae6fb2eSBarry Smith   /* Copy A into F->lu->A */
1224*eae6fb2eSBarry Smith   ierr = MatMPIDenseCopyToPlapack(A,*F);CHKERRQ(ierr);
122501b82886SBarry Smith 
122601b82886SBarry Smith   /* Factor P A -> Chol */
122701b82886SBarry Smith   info_pla = PLA_Chol(PLA_LOWER_TRIANGULAR,lu->A);
1228*eae6fb2eSBarry Smith   if (info_pla) SETERRQ1( PETSC_ERR_MAT_CH_ZRPVT,"Nonpositive definite matrix detected at row %d from PLA_Chol()",info_pla);
122901b82886SBarry Smith 
1230*eae6fb2eSBarry Smith   (*F)->assembled = PETSC_TRUE;
123101b82886SBarry Smith   PetscFunctionReturn(0);
123201b82886SBarry Smith }
123301b82886SBarry Smith 
1234*eae6fb2eSBarry Smith 
123501b82886SBarry Smith #undef __FUNCT__
123601b82886SBarry Smith #define __FUNCT__ "MatFactorSymbolic_MPIDense_Private"
123701b82886SBarry Smith PetscErrorCode MatFactorSymbolic_MPIDense_Private(Mat A,MatFactorInfo *info,Mat *F)
123801b82886SBarry Smith {
123901b82886SBarry Smith   PetscErrorCode ierr;
124001b82886SBarry Smith 
124101b82886SBarry Smith   PetscFunctionBegin;
124201b82886SBarry Smith   /* Create the factorization matrix */
1243*eae6fb2eSBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,F);CHKERRQ(ierr);
1244*eae6fb2eSBarry Smith   ierr = MatSetSizes(*F,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr);
1245*eae6fb2eSBarry Smith   ierr = MatSetType(*F,((PetscObject)A)->type_name);CHKERRQ(ierr);
1246*eae6fb2eSBarry Smith   ierr = MatMPIDenseCreatePlapack(*F);CHKERRQ(ierr);
124701b82886SBarry Smith   PetscFunctionReturn(0);
124801b82886SBarry Smith }
124901b82886SBarry Smith 
125001b82886SBarry Smith /* Note the Petsc r and c permutations are ignored */
125101b82886SBarry Smith #undef __FUNCT__
125201b82886SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_MPIDense"
125301b82886SBarry Smith PetscErrorCode MatLUFactorSymbolic_MPIDense(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
125401b82886SBarry Smith {
125501b82886SBarry Smith   PetscErrorCode ierr;
1256e1156936SBarry Smith   PetscInt       M = A->rmap.N;
1257e1156936SBarry Smith   Mat_Plapack    *lu;
125801b82886SBarry Smith 
125901b82886SBarry Smith   PetscFunctionBegin;
126001b82886SBarry Smith   ierr = MatFactorSymbolic_MPIDense_Private(A,info,F);CHKERRQ(ierr);
1261e1156936SBarry Smith   lu = (Mat_Plapack*)(*F)->spptr;
1262e1156936SBarry Smith   ierr = PLA_Mvector_create(MPI_INT,M,1,lu->templ,PLA_ALIGN_FIRST,&lu->pivots);CHKERRQ(ierr);
126301b82886SBarry Smith   (*F)->factor = FACTOR_LU;
126401b82886SBarry Smith   PetscFunctionReturn(0);
126501b82886SBarry Smith }
126601b82886SBarry Smith 
126701b82886SBarry Smith /* Note the Petsc perm permutation is ignored */
126801b82886SBarry Smith #undef __FUNCT__
126901b82886SBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_MPIDense"
127001b82886SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_MPIDense(Mat A,IS perm,MatFactorInfo *info,Mat *F)
127101b82886SBarry Smith {
127201b82886SBarry Smith   PetscErrorCode ierr;
127301b82886SBarry Smith   PetscTruth     issymmetric,set;
127401b82886SBarry Smith 
127501b82886SBarry Smith   PetscFunctionBegin;
127601b82886SBarry Smith   ierr = MatIsSymmetricKnown(A,&set,&issymmetric); CHKERRQ(ierr);
127701b82886SBarry Smith   if (!set || !issymmetric) SETERRQ(PETSC_ERR_USER,"Matrix must be set as MAT_SYMMETRIC for CholeskyFactor()");
127801b82886SBarry Smith   ierr = MatFactorSymbolic_MPIDense_Private(A,info,F);CHKERRQ(ierr);
127901b82886SBarry Smith   (*F)->factor = FACTOR_CHOLESKY;
128001b82886SBarry Smith   PetscFunctionReturn(0);
128101b82886SBarry Smith }
128201b82886SBarry Smith #endif
128301b82886SBarry Smith 
12848965ea79SLois Curfman McInnes /* -------------------------------------------------------------------*/
128509dc0095SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIDense,
128609dc0095SBarry Smith        MatGetRow_MPIDense,
128709dc0095SBarry Smith        MatRestoreRow_MPIDense,
128809dc0095SBarry Smith        MatMult_MPIDense,
128997304618SKris Buschelman /* 4*/ MatMultAdd_MPIDense,
12907c922b88SBarry Smith        MatMultTranspose_MPIDense,
12917c922b88SBarry Smith        MatMultTransposeAdd_MPIDense,
129201b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
129301b82886SBarry Smith        MatSolve_MPIDense,
129401b82886SBarry Smith #else
12958965ea79SLois Curfman McInnes        0,
129601b82886SBarry Smith #endif
129709dc0095SBarry Smith        0,
129809dc0095SBarry Smith        0,
129997304618SKris Buschelman /*10*/ 0,
130009dc0095SBarry Smith        0,
130109dc0095SBarry Smith        0,
130209dc0095SBarry Smith        0,
130309dc0095SBarry Smith        MatTranspose_MPIDense,
130497304618SKris Buschelman /*15*/ MatGetInfo_MPIDense,
13056e4ee0c6SHong Zhang        MatEqual_MPIDense,
130609dc0095SBarry Smith        MatGetDiagonal_MPIDense,
13075b2fa520SLois Curfman McInnes        MatDiagonalScale_MPIDense,
130809dc0095SBarry Smith        MatNorm_MPIDense,
130997304618SKris Buschelman /*20*/ MatAssemblyBegin_MPIDense,
131009dc0095SBarry Smith        MatAssemblyEnd_MPIDense,
131109dc0095SBarry Smith        0,
131209dc0095SBarry Smith        MatSetOption_MPIDense,
131309dc0095SBarry Smith        MatZeroEntries_MPIDense,
131497304618SKris Buschelman /*25*/ MatZeroRows_MPIDense,
131501b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
13160ad19415SHong Zhang        MatLUFactorSymbolic_MPIDense,
131701b82886SBarry Smith        MatLUFactorNumeric_MPIDense,
13180ad19415SHong Zhang        MatCholeskyFactorSymbolic_MPIDense,
131901b82886SBarry Smith        MatCholeskyFactorNumeric_MPIDense,
132001b82886SBarry Smith #else
1321919b68f7SBarry Smith        0,
132201b82886SBarry Smith        0,
132301b82886SBarry Smith        0,
132401b82886SBarry Smith        0,
132501b82886SBarry Smith #endif
132697304618SKris Buschelman /*30*/ MatSetUpPreallocation_MPIDense,
1327273d9f13SBarry Smith        0,
132809dc0095SBarry Smith        0,
132909dc0095SBarry Smith        MatGetArray_MPIDense,
133009dc0095SBarry Smith        MatRestoreArray_MPIDense,
133197304618SKris Buschelman /*35*/ MatDuplicate_MPIDense,
133209dc0095SBarry Smith        0,
133309dc0095SBarry Smith        0,
133409dc0095SBarry Smith        0,
133509dc0095SBarry Smith        0,
133697304618SKris Buschelman /*40*/ 0,
13372ce60cd0SSatish Balay        MatGetSubMatrices_MPIDense,
133809dc0095SBarry Smith        0,
133909dc0095SBarry Smith        MatGetValues_MPIDense,
134009dc0095SBarry Smith        0,
134197304618SKris Buschelman /*45*/ 0,
134209dc0095SBarry Smith        MatScale_MPIDense,
134309dc0095SBarry Smith        0,
134409dc0095SBarry Smith        0,
134509dc0095SBarry Smith        0,
1346521d7252SBarry Smith /*50*/ 0,
134709dc0095SBarry Smith        0,
134809dc0095SBarry Smith        0,
134909dc0095SBarry Smith        0,
135009dc0095SBarry Smith        0,
135197304618SKris Buschelman /*55*/ 0,
135209dc0095SBarry Smith        0,
135309dc0095SBarry Smith        0,
135409dc0095SBarry Smith        0,
135509dc0095SBarry Smith        0,
135697304618SKris Buschelman /*60*/ MatGetSubMatrix_MPIDense,
1357b9b97703SBarry Smith        MatDestroy_MPIDense,
1358b9b97703SBarry Smith        MatView_MPIDense,
1359357abbc8SBarry Smith        0,
136097304618SKris Buschelman        0,
136197304618SKris Buschelman /*65*/ 0,
136297304618SKris Buschelman        0,
136397304618SKris Buschelman        0,
136497304618SKris Buschelman        0,
136597304618SKris Buschelman        0,
136697304618SKris Buschelman /*70*/ 0,
136797304618SKris Buschelman        0,
136897304618SKris Buschelman        0,
136997304618SKris Buschelman        0,
137097304618SKris Buschelman        0,
137197304618SKris Buschelman /*75*/ 0,
137297304618SKris Buschelman        0,
137397304618SKris Buschelman        0,
137497304618SKris Buschelman        0,
137597304618SKris Buschelman        0,
137697304618SKris Buschelman /*80*/ 0,
137797304618SKris Buschelman        0,
137897304618SKris Buschelman        0,
137997304618SKris Buschelman        0,
1380865e5f61SKris Buschelman /*84*/ MatLoad_MPIDense,
1381865e5f61SKris Buschelman        0,
1382865e5f61SKris Buschelman        0,
1383865e5f61SKris Buschelman        0,
1384865e5f61SKris Buschelman        0,
1385865e5f61SKris Buschelman        0,
1386865e5f61SKris Buschelman /*90*/ 0,
13874ae313f4SHong Zhang        MatMatMultSymbolic_MPIDense_MPIDense,
1388865e5f61SKris Buschelman        0,
1389865e5f61SKris Buschelman        0,
1390865e5f61SKris Buschelman        0,
1391865e5f61SKris Buschelman /*95*/ 0,
1392865e5f61SKris Buschelman        0,
1393865e5f61SKris Buschelman        0,
1394865e5f61SKris Buschelman        0};
13958965ea79SLois Curfman McInnes 
1396273d9f13SBarry Smith EXTERN_C_BEGIN
13974a2ae208SSatish Balay #undef __FUNCT__
1398a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense"
1399be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data)
1400a23d5eceSKris Buschelman {
1401a23d5eceSKris Buschelman   Mat_MPIDense   *a;
1402dfbe8321SBarry Smith   PetscErrorCode ierr;
1403a23d5eceSKris Buschelman 
1404a23d5eceSKris Buschelman   PetscFunctionBegin;
1405a23d5eceSKris Buschelman   mat->preallocated = PETSC_TRUE;
1406a23d5eceSKris Buschelman   /* Note:  For now, when data is specified above, this assumes the user correctly
1407a23d5eceSKris Buschelman    allocates the local dense storage space.  We should add error checking. */
1408a23d5eceSKris Buschelman 
1409a23d5eceSKris Buschelman   a    = (Mat_MPIDense*)mat->data;
1410f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&a->A);CHKERRQ(ierr);
1411899cda47SBarry Smith   ierr = MatSetSizes(a->A,mat->rmap.n,mat->cmap.N,mat->rmap.n,mat->cmap.N);CHKERRQ(ierr);
14125c5985e7SKris Buschelman   ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr);
14135c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr);
141452e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
1415a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1416a23d5eceSKris Buschelman }
1417a23d5eceSKris Buschelman EXTERN_C_END
1418a23d5eceSKris Buschelman 
14197878bbefSBarry Smith EXTERN_C_BEGIN
14207878bbefSBarry Smith #if defined(PETSC_HAVE_PLAPACK)
14217878bbefSBarry Smith #undef __FUNCT__
14227878bbefSBarry Smith #define __FUNCT__ "MatSetSolverType_MPIDense_PLAPACK"
14237878bbefSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatSetSolverType_MPIDense_PLAPACK(Mat mat,const char *type)
14247878bbefSBarry Smith {
14257878bbefSBarry Smith   PetscFunctionBegin;
14267878bbefSBarry Smith   /* dummy function so that -mat_solver_type plapack or MatSetSolverType(mat,"plapack"); silently work (since PLAPACK is on by default) */
14277878bbefSBarry Smith   PetscFunctionReturn(0);
14287878bbefSBarry Smith }
14297878bbefSBarry Smith #endif
14307a667e6fSSatish Balay EXTERN_C_END
14317878bbefSBarry Smith 
14320bad9183SKris Buschelman /*MC
1433fafad747SKris Buschelman    MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices.
14340bad9183SKris Buschelman 
14350bad9183SKris Buschelman    Options Database Keys:
14360bad9183SKris Buschelman . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions()
14370bad9183SKris Buschelman 
14380bad9183SKris Buschelman   Level: beginner
14390bad9183SKris Buschelman 
14407878bbefSBarry Smith   MATMPIDENSE matrices may use direct solvers (LU, Cholesky, and QR)
14417878bbefSBarry Smith   for parallel dense matrices via the external package PLAPACK, if PLAPACK is installed
14427878bbefSBarry Smith   (run config/configure.py with the option --download-plapack)
14437878bbefSBarry Smith 
14447878bbefSBarry Smith 
14457878bbefSBarry Smith   Options Database Keys:
14467878bbefSBarry Smith . -mat_plapack_nprows <n> - number of rows in processor partition
14477878bbefSBarry Smith . -mat_plapack_npcols <n> - number of columns in processor partition
14487878bbefSBarry Smith . -mat_plapack_nb <n> - block size of template vector
14497878bbefSBarry Smith . -mat_plapack_nb_alg <n> - algorithmic block size
14507878bbefSBarry Smith - -mat_plapack_ckerror <n> - error checking flag
14517878bbefSBarry Smith 
14527878bbefSBarry Smith .seealso: MatCreateMPIDense(), MATDENSE, MATSEQDENSE
14530bad9183SKris Buschelman M*/
14540bad9183SKris Buschelman 
1455a23d5eceSKris Buschelman EXTERN_C_BEGIN
1456a23d5eceSKris Buschelman #undef __FUNCT__
14574a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIDense"
1458be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPIDense(Mat mat)
1459273d9f13SBarry Smith {
1460273d9f13SBarry Smith   Mat_MPIDense   *a;
1461dfbe8321SBarry Smith   PetscErrorCode ierr;
146201b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
146301b82886SBarry Smith   Mat_Plapack    *lu;
146401b82886SBarry Smith #endif
1465273d9f13SBarry Smith 
1466273d9f13SBarry Smith   PetscFunctionBegin;
146738f2d2fdSLisandro Dalcin   ierr              = PetscNewLog(mat,Mat_MPIDense,&a);CHKERRQ(ierr);
1468b0a32e0cSBarry Smith   mat->data         = (void*)a;
1469273d9f13SBarry Smith   ierr              = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1470273d9f13SBarry Smith   mat->factor       = 0;
1471273d9f13SBarry Smith   mat->mapping      = 0;
1472273d9f13SBarry Smith 
1473273d9f13SBarry Smith   mat->insertmode = NOT_SET_VALUES;
14747adad957SLisandro Dalcin   ierr = MPI_Comm_rank(((PetscObject)mat)->comm,&a->rank);CHKERRQ(ierr);
14757adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)mat)->comm,&a->size);CHKERRQ(ierr);
1476273d9f13SBarry Smith 
1477899cda47SBarry Smith   mat->rmap.bs = mat->cmap.bs = 1;
14786148ca0dSBarry Smith   ierr = PetscMapSetUp(&mat->rmap);CHKERRQ(ierr);
14796148ca0dSBarry Smith   ierr = PetscMapSetUp(&mat->cmap);CHKERRQ(ierr);
1480899cda47SBarry Smith   a->nvec = mat->cmap.n;
1481273d9f13SBarry Smith 
1482273d9f13SBarry Smith   /* build cache for off array entries formed */
1483273d9f13SBarry Smith   a->donotstash = PETSC_FALSE;
14847adad957SLisandro Dalcin   ierr = MatStashCreate_Private(((PetscObject)mat)->comm,1,&mat->stash);CHKERRQ(ierr);
1485273d9f13SBarry Smith 
1486273d9f13SBarry Smith   /* stuff used for matrix vector multiply */
1487273d9f13SBarry Smith   a->lvec        = 0;
1488273d9f13SBarry Smith   a->Mvctx       = 0;
1489273d9f13SBarry Smith   a->roworiented = PETSC_TRUE;
1490273d9f13SBarry Smith 
1491273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C",
1492273d9f13SBarry Smith                                      "MatGetDiagonalBlock_MPIDense",
1493273d9f13SBarry Smith                                      MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr);
1494a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C",
1495a23d5eceSKris Buschelman                                      "MatMPIDenseSetPreallocation_MPIDense",
1496a23d5eceSKris Buschelman                                      MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr);
14974ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_mpiaij_mpidense_C",
14984ae313f4SHong Zhang                                      "MatMatMult_MPIAIJ_MPIDense",
14994ae313f4SHong Zhang                                       MatMatMult_MPIAIJ_MPIDense);CHKERRQ(ierr);
15004ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_mpiaij_mpidense_C",
15014ae313f4SHong Zhang                                      "MatMatMultSymbolic_MPIAIJ_MPIDense",
15024ae313f4SHong Zhang                                       MatMatMultSymbolic_MPIAIJ_MPIDense);CHKERRQ(ierr);
15034ae313f4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_mpiaij_mpidense_C",
15044ae313f4SHong Zhang                                      "MatMatMultNumeric_MPIAIJ_MPIDense",
15054ae313f4SHong Zhang                                       MatMatMultNumeric_MPIAIJ_MPIDense);CHKERRQ(ierr);
15067878bbefSBarry Smith #if defined(PETSC_HAVE_PLAPACK)
15077878bbefSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSetSolverType_mpidense_plapack_C",
15087878bbefSBarry Smith                                      "MatSetSolverType_MPIDense_PLAPACK",
15097878bbefSBarry Smith                                       MatSetSolverType_MPIDense_PLAPACK);CHKERRQ(ierr);
15107878bbefSBarry Smith #endif
151138aed534SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)mat,MATMPIDENSE);CHKERRQ(ierr);
151201b82886SBarry Smith 
151301b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
151401b82886SBarry Smith   ierr = PetscNewLog(mat,Mat_Plapack,&lu);CHKERRQ(ierr);
151501b82886SBarry Smith   lu->CleanUpPlapack       = PETSC_FALSE;
151601b82886SBarry Smith   mat->spptr                 = (void*)lu;
151701b82886SBarry Smith #endif
1518273d9f13SBarry Smith   PetscFunctionReturn(0);
1519273d9f13SBarry Smith }
1520273d9f13SBarry Smith EXTERN_C_END
1521273d9f13SBarry Smith 
1522209238afSKris Buschelman /*MC
1523002d173eSKris Buschelman    MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices.
1524209238afSKris Buschelman 
1525209238afSKris Buschelman    This matrix type is identical to MATSEQDENSE when constructed with a single process communicator,
1526209238afSKris Buschelman    and MATMPIDENSE otherwise.
1527209238afSKris Buschelman 
1528209238afSKris Buschelman    Options Database Keys:
1529209238afSKris Buschelman . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions()
1530209238afSKris Buschelman 
1531209238afSKris Buschelman   Level: beginner
1532209238afSKris Buschelman 
153301b82886SBarry Smith 
1534209238afSKris Buschelman .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE
1535209238afSKris Buschelman M*/
1536209238afSKris Buschelman 
1537209238afSKris Buschelman EXTERN_C_BEGIN
1538209238afSKris Buschelman #undef __FUNCT__
1539209238afSKris Buschelman #define __FUNCT__ "MatCreate_Dense"
1540be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_Dense(Mat A)
1541dfbe8321SBarry Smith {
15426849ba73SBarry Smith   PetscErrorCode ierr;
154313f74950SBarry Smith   PetscMPIInt    size;
1544209238afSKris Buschelman 
1545209238afSKris Buschelman   PetscFunctionBegin;
15467adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1547209238afSKris Buschelman   if (size == 1) {
1548209238afSKris Buschelman     ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr);
1549209238afSKris Buschelman   } else {
1550209238afSKris Buschelman     ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
1551209238afSKris Buschelman   }
1552209238afSKris Buschelman   PetscFunctionReturn(0);
1553209238afSKris Buschelman }
1554209238afSKris Buschelman EXTERN_C_END
1555209238afSKris Buschelman 
15564a2ae208SSatish Balay #undef __FUNCT__
15574a2ae208SSatish Balay #define __FUNCT__ "MatMPIDenseSetPreallocation"
1558273d9f13SBarry Smith /*@C
1559273d9f13SBarry Smith    MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries
1560273d9f13SBarry Smith 
1561273d9f13SBarry Smith    Not collective
1562273d9f13SBarry Smith 
1563273d9f13SBarry Smith    Input Parameters:
1564273d9f13SBarry Smith .  A - the matrix
1565273d9f13SBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1566273d9f13SBarry Smith    to control all matrix memory allocation.
1567273d9f13SBarry Smith 
1568273d9f13SBarry Smith    Notes:
1569273d9f13SBarry Smith    The dense format is fully compatible with standard Fortran 77
1570273d9f13SBarry Smith    storage by columns.
1571273d9f13SBarry Smith 
1572273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1573273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1574273d9f13SBarry Smith    set data=PETSC_NULL.
1575273d9f13SBarry Smith 
1576273d9f13SBarry Smith    Level: intermediate
1577273d9f13SBarry Smith 
1578273d9f13SBarry Smith .keywords: matrix,dense, parallel
1579273d9f13SBarry Smith 
1580273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
1581273d9f13SBarry Smith @*/
1582be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data)
1583273d9f13SBarry Smith {
1584dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar *);
1585273d9f13SBarry Smith 
1586273d9f13SBarry Smith   PetscFunctionBegin;
1587565567f0SKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1588a23d5eceSKris Buschelman   if (f) {
1589a23d5eceSKris Buschelman     ierr = (*f)(mat,data);CHKERRQ(ierr);
1590a23d5eceSKris Buschelman   }
1591273d9f13SBarry Smith   PetscFunctionReturn(0);
1592273d9f13SBarry Smith }
1593273d9f13SBarry Smith 
15944a2ae208SSatish Balay #undef __FUNCT__
15954a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIDense"
15968965ea79SLois Curfman McInnes /*@C
159739ddd567SLois Curfman McInnes    MatCreateMPIDense - Creates a sparse parallel matrix in dense format.
15988965ea79SLois Curfman McInnes 
1599db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1600db81eaa0SLois Curfman McInnes 
16018965ea79SLois Curfman McInnes    Input Parameters:
1602db81eaa0SLois Curfman McInnes +  comm - MPI communicator
16038965ea79SLois Curfman McInnes .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1604db81eaa0SLois Curfman McInnes .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
16058965ea79SLois Curfman McInnes .  M - number of global rows (or PETSC_DECIDE to have calculated if m is given)
1606db81eaa0SLois Curfman McInnes .  N - number of global columns (or PETSC_DECIDE to have calculated if n is given)
16077f5ff6fdSBarry Smith -  data - optional location of matrix data.  Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc
1608dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
16098965ea79SLois Curfman McInnes 
16108965ea79SLois Curfman McInnes    Output Parameter:
1611477f1c0bSLois Curfman McInnes .  A - the matrix
16128965ea79SLois Curfman McInnes 
1613b259b22eSLois Curfman McInnes    Notes:
161439ddd567SLois Curfman McInnes    The dense format is fully compatible with standard Fortran 77
161539ddd567SLois Curfman McInnes    storage by columns.
16168965ea79SLois Curfman McInnes 
161718f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
161818f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
16197f5ff6fdSBarry Smith    set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users).
162018f449edSLois Curfman McInnes 
16218965ea79SLois Curfman McInnes    The user MUST specify either the local or global matrix dimensions
16228965ea79SLois Curfman McInnes    (possibly both).
16238965ea79SLois Curfman McInnes 
1624027ccd11SLois Curfman McInnes    Level: intermediate
1625027ccd11SLois Curfman McInnes 
162639ddd567SLois Curfman McInnes .keywords: matrix,dense, parallel
16278965ea79SLois Curfman McInnes 
162839ddd567SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues()
16298965ea79SLois Curfman McInnes @*/
1630be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar *data,Mat *A)
16318965ea79SLois Curfman McInnes {
16326849ba73SBarry Smith   PetscErrorCode ierr;
163313f74950SBarry Smith   PetscMPIInt    size;
16348965ea79SLois Curfman McInnes 
16353a40ed3dSBarry Smith   PetscFunctionBegin;
1636f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1637f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1638273d9f13SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1639273d9f13SBarry Smith   if (size > 1) {
1640273d9f13SBarry Smith     ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr);
1641273d9f13SBarry Smith     ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1642273d9f13SBarry Smith   } else {
1643273d9f13SBarry Smith     ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1644273d9f13SBarry Smith     ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
16458c469469SLois Curfman McInnes   }
16463a40ed3dSBarry Smith   PetscFunctionReturn(0);
16478965ea79SLois Curfman McInnes }
16488965ea79SLois Curfman McInnes 
16494a2ae208SSatish Balay #undef __FUNCT__
16504a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_MPIDense"
16516849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
16528965ea79SLois Curfman McInnes {
16538965ea79SLois Curfman McInnes   Mat            mat;
16543501a2bdSLois Curfman McInnes   Mat_MPIDense   *a,*oldmat = (Mat_MPIDense*)A->data;
1655dfbe8321SBarry Smith   PetscErrorCode ierr;
16568965ea79SLois Curfman McInnes 
16573a40ed3dSBarry Smith   PetscFunctionBegin;
16588965ea79SLois Curfman McInnes   *newmat       = 0;
16597adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&mat);CHKERRQ(ierr);
1660899cda47SBarry Smith   ierr = MatSetSizes(mat,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr);
16617adad957SLisandro Dalcin   ierr = MatSetType(mat,((PetscObject)A)->type_name);CHKERRQ(ierr);
1662834f8fabSBarry Smith   a                 = (Mat_MPIDense*)mat->data;
1663e04c1aa4SHong Zhang   ierr              = PetscMemcpy(mat->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
16643501a2bdSLois Curfman McInnes   mat->factor       = A->factor;
1665c456f294SBarry Smith   mat->assembled    = PETSC_TRUE;
1666273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
16678965ea79SLois Curfman McInnes 
1668899cda47SBarry Smith   mat->rmap.rstart     = A->rmap.rstart;
1669899cda47SBarry Smith   mat->rmap.rend       = A->rmap.rend;
16708965ea79SLois Curfman McInnes   a->size              = oldmat->size;
16718965ea79SLois Curfman McInnes   a->rank              = oldmat->rank;
1672e0fa3b82SLois Curfman McInnes   mat->insertmode      = NOT_SET_VALUES;
1673b9b97703SBarry Smith   a->nvec              = oldmat->nvec;
16743782ba37SSatish Balay   a->donotstash        = oldmat->donotstash;
1675e04c1aa4SHong Zhang 
1676899cda47SBarry Smith   ierr = PetscMemcpy(mat->rmap.range,A->rmap.range,(a->size+1)*sizeof(PetscInt));CHKERRQ(ierr);
1677899cda47SBarry Smith   ierr = PetscMemcpy(mat->cmap.range,A->cmap.range,(a->size+1)*sizeof(PetscInt));CHKERRQ(ierr);
16787adad957SLisandro Dalcin   ierr = MatStashCreate_Private(((PetscObject)A)->comm,1,&mat->stash);CHKERRQ(ierr);
16798965ea79SLois Curfman McInnes 
1680329f5518SBarry Smith   ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr);
16815609ef8eSBarry Smith   ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
168252e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
168301b82886SBarry Smith 
168401b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
168501b82886SBarry Smith   ierr = PetscMemcpy(mat->spptr,A->spptr,sizeof(Mat_Plapack));CHKERRQ(ierr);
168601b82886SBarry Smith #endif
16878965ea79SLois Curfman McInnes   *newmat = mat;
16883a40ed3dSBarry Smith   PetscFunctionReturn(0);
16898965ea79SLois Curfman McInnes }
16908965ea79SLois Curfman McInnes 
1691e090d566SSatish Balay #include "petscsys.h"
16928965ea79SLois Curfman McInnes 
16934a2ae208SSatish Balay #undef __FUNCT__
16944a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense_DenseInFile"
1695f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense_DenseInFile(MPI_Comm comm,PetscInt fd,PetscInt M,PetscInt N, MatType type,Mat *newmat)
169690ace30eSBarry Smith {
16976849ba73SBarry Smith   PetscErrorCode ierr;
169832dcc486SBarry Smith   PetscMPIInt    rank,size;
169913f74950SBarry Smith   PetscInt       *rowners,i,m,nz,j;
170087828ca2SBarry Smith   PetscScalar    *array,*vals,*vals_ptr;
170190ace30eSBarry Smith   MPI_Status     status;
170290ace30eSBarry Smith 
17033a40ed3dSBarry Smith   PetscFunctionBegin;
1704d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1705d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
170690ace30eSBarry Smith 
170790ace30eSBarry Smith   /* determine ownership of all rows */
170890ace30eSBarry Smith   m          = M/size + ((M % size) > rank);
170913f74950SBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(PetscInt),&rowners);CHKERRQ(ierr);
171013f74950SBarry Smith   ierr       = MPI_Allgather(&m,1,MPIU_INT,rowners+1,1,MPIU_INT,comm);CHKERRQ(ierr);
171190ace30eSBarry Smith   rowners[0] = 0;
171290ace30eSBarry Smith   for (i=2; i<=size; i++) {
171390ace30eSBarry Smith     rowners[i] += rowners[i-1];
171490ace30eSBarry Smith   }
171590ace30eSBarry Smith 
1716f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,newmat);CHKERRQ(ierr);
1717f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1718be5d1d56SKris Buschelman   ierr = MatSetType(*newmat,type);CHKERRQ(ierr);
1719878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
172090ace30eSBarry Smith   ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr);
172190ace30eSBarry Smith 
172290ace30eSBarry Smith   if (!rank) {
172387828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
172490ace30eSBarry Smith 
172590ace30eSBarry Smith     /* read in my part of the matrix numerical values  */
17260752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr);
172790ace30eSBarry Smith 
172890ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
172990ace30eSBarry Smith     vals_ptr = vals;
173090ace30eSBarry Smith     for (i=0; i<m; i++) {
173190ace30eSBarry Smith       for (j=0; j<N; j++) {
173290ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
173390ace30eSBarry Smith       }
173490ace30eSBarry Smith     }
173590ace30eSBarry Smith 
173690ace30eSBarry Smith     /* read in other processors and ship out */
173790ace30eSBarry Smith     for (i=1; i<size; i++) {
173890ace30eSBarry Smith       nz   = (rowners[i+1] - rowners[i])*N;
17390752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
17407adad957SLisandro Dalcin       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)(*newmat))->tag,comm);CHKERRQ(ierr);
174190ace30eSBarry Smith     }
17423a40ed3dSBarry Smith   } else {
174390ace30eSBarry Smith     /* receive numeric values */
174487828ca2SBarry Smith     ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
174590ace30eSBarry Smith 
174690ace30eSBarry Smith     /* receive message of values*/
17477adad957SLisandro Dalcin     ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,((PetscObject)(*newmat))->tag,comm,&status);CHKERRQ(ierr);
174890ace30eSBarry Smith 
174990ace30eSBarry Smith     /* insert into matrix-by row (this is why cannot directly read into array */
175090ace30eSBarry Smith     vals_ptr = vals;
175190ace30eSBarry Smith     for (i=0; i<m; i++) {
175290ace30eSBarry Smith       for (j=0; j<N; j++) {
175390ace30eSBarry Smith         array[i + j*m] = *vals_ptr++;
175490ace30eSBarry Smith       }
175590ace30eSBarry Smith     }
175690ace30eSBarry Smith   }
1757606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
1758606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
17596d4a8577SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17606d4a8577SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17613a40ed3dSBarry Smith   PetscFunctionReturn(0);
176290ace30eSBarry Smith }
176390ace30eSBarry Smith 
17644a2ae208SSatish Balay #undef __FUNCT__
17654a2ae208SSatish Balay #define __FUNCT__ "MatLoad_MPIDense"
1766f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_MPIDense(PetscViewer viewer, MatType type,Mat *newmat)
17678965ea79SLois Curfman McInnes {
17688965ea79SLois Curfman McInnes   Mat            A;
176987828ca2SBarry Smith   PetscScalar    *vals,*svals;
177019bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
17718965ea79SLois Curfman McInnes   MPI_Status     status;
177213f74950SBarry Smith   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*rowners,*sndcounts,m,maxnz;
177313f74950SBarry Smith   PetscInt       header[4],*rowlengths = 0,M,N,*cols;
177413f74950SBarry Smith   PetscInt       *ourlens,*procsnz = 0,*offlens,jj,*mycols,*smycols;
177513f74950SBarry Smith   PetscInt       i,nz,j,rstart,rend;
177613f74950SBarry Smith   int            fd;
17776849ba73SBarry Smith   PetscErrorCode ierr;
17788965ea79SLois Curfman McInnes 
17793a40ed3dSBarry Smith   PetscFunctionBegin;
1780d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1781d132466eSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
17828965ea79SLois Curfman McInnes   if (!rank) {
1783b0a32e0cSBarry Smith     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
17840752156aSBarry Smith     ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr);
1785552e946dSBarry Smith     if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
17868965ea79SLois Curfman McInnes   }
17878965ea79SLois Curfman McInnes 
178813f74950SBarry Smith   ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr);
178990ace30eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
179090ace30eSBarry Smith 
179190ace30eSBarry Smith   /*
179290ace30eSBarry Smith        Handle case where matrix is stored on disk as a dense matrix
179390ace30eSBarry Smith   */
179490ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) {
1795be5d1d56SKris Buschelman     ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr);
17963a40ed3dSBarry Smith     PetscFunctionReturn(0);
179790ace30eSBarry Smith   }
179890ace30eSBarry Smith 
17998965ea79SLois Curfman McInnes   /* determine ownership of all rows */
1800e44c0bd4SBarry Smith   m          = PetscMPIIntCast(M/size + ((M % size) > rank));
1801e44c0bd4SBarry Smith   ierr       = PetscMalloc((size+2)*sizeof(PetscMPIInt),&rowners);CHKERRQ(ierr);
1802ca161407SBarry Smith   ierr       = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
18038965ea79SLois Curfman McInnes   rowners[0] = 0;
18048965ea79SLois Curfman McInnes   for (i=2; i<=size; i++) {
18058965ea79SLois Curfman McInnes     rowners[i] += rowners[i-1];
18068965ea79SLois Curfman McInnes   }
18078965ea79SLois Curfman McInnes   rstart = rowners[rank];
18088965ea79SLois Curfman McInnes   rend   = rowners[rank+1];
18098965ea79SLois Curfman McInnes 
18108965ea79SLois Curfman McInnes   /* distribute row lengths to all processors */
181113f74950SBarry Smith   ierr    = PetscMalloc(2*(rend-rstart+1)*sizeof(PetscInt),&ourlens);CHKERRQ(ierr);
18128965ea79SLois Curfman McInnes   offlens = ourlens + (rend-rstart);
18138965ea79SLois Curfman McInnes   if (!rank) {
181413f74950SBarry Smith     ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
18150752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
181613f74950SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr);
18178965ea79SLois Curfman McInnes     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
18181466f1e1SBarry Smith     ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
1819606d414cSSatish Balay     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
1820ca161407SBarry Smith   } else {
18211466f1e1SBarry Smith     ierr = MPI_Scatterv(0,0,0,MPIU_INT,ourlens,rend-rstart,MPIU_INT,0,comm);CHKERRQ(ierr);
18228965ea79SLois Curfman McInnes   }
18238965ea79SLois Curfman McInnes 
18248965ea79SLois Curfman McInnes   if (!rank) {
18258965ea79SLois Curfman McInnes     /* calculate the number of nonzeros on each processor */
182613f74950SBarry Smith     ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr);
182713f74950SBarry Smith     ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr);
18288965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
18298965ea79SLois Curfman McInnes       for (j=rowners[i]; j< rowners[i+1]; j++) {
18308965ea79SLois Curfman McInnes         procsnz[i] += rowlengths[j];
18318965ea79SLois Curfman McInnes       }
18328965ea79SLois Curfman McInnes     }
1833606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
18348965ea79SLois Curfman McInnes 
18358965ea79SLois Curfman McInnes     /* determine max buffer needed and allocate it */
18368965ea79SLois Curfman McInnes     maxnz = 0;
18378965ea79SLois Curfman McInnes     for (i=0; i<size; i++) {
18380452661fSBarry Smith       maxnz = PetscMax(maxnz,procsnz[i]);
18398965ea79SLois Curfman McInnes     }
184013f74950SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr);
18418965ea79SLois Curfman McInnes 
18428965ea79SLois Curfman McInnes     /* read in my part of the matrix column indices  */
18438965ea79SLois Curfman McInnes     nz = procsnz[0];
184413f74950SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
18450752156aSBarry Smith     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
18468965ea79SLois Curfman McInnes 
18478965ea79SLois Curfman McInnes     /* read in every one elses and ship off */
18488965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
18498965ea79SLois Curfman McInnes       nz   = procsnz[i];
18500752156aSBarry Smith       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
185113f74950SBarry Smith       ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr);
18528965ea79SLois Curfman McInnes     }
1853606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
18543a40ed3dSBarry Smith   } else {
18558965ea79SLois Curfman McInnes     /* determine buffer space needed for message */
18568965ea79SLois Curfman McInnes     nz = 0;
18578965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
18588965ea79SLois Curfman McInnes       nz += ourlens[i];
18598965ea79SLois Curfman McInnes     }
186013f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&mycols);CHKERRQ(ierr);
18618965ea79SLois Curfman McInnes 
18628965ea79SLois Curfman McInnes     /* receive message of column indices*/
186313f74950SBarry Smith     ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr);
186413f74950SBarry Smith     ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr);
186529bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
18668965ea79SLois Curfman McInnes   }
18678965ea79SLois Curfman McInnes 
18688965ea79SLois Curfman McInnes   /* loop over local rows, determining number of off diagonal entries */
186913f74950SBarry Smith   ierr = PetscMemzero(offlens,m*sizeof(PetscInt));CHKERRQ(ierr);
18708965ea79SLois Curfman McInnes   jj = 0;
18718965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
18728965ea79SLois Curfman McInnes     for (j=0; j<ourlens[i]; j++) {
18738965ea79SLois Curfman McInnes       if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++;
18748965ea79SLois Curfman McInnes       jj++;
18758965ea79SLois Curfman McInnes     }
18768965ea79SLois Curfman McInnes   }
18778965ea79SLois Curfman McInnes 
18788965ea79SLois Curfman McInnes   /* create our matrix */
18798965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
18808965ea79SLois Curfman McInnes     ourlens[i] -= offlens[i];
18818965ea79SLois Curfman McInnes   }
1882f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,newmat);CHKERRQ(ierr);
1883f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1884be5d1d56SKris Buschelman   ierr = MatSetType(*newmat,type);CHKERRQ(ierr);
1885878740d9SKris Buschelman   ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr);
18868965ea79SLois Curfman McInnes   A = *newmat;
18878965ea79SLois Curfman McInnes   for (i=0; i<m; i++) {
18888965ea79SLois Curfman McInnes     ourlens[i] += offlens[i];
18898965ea79SLois Curfman McInnes   }
18908965ea79SLois Curfman McInnes 
18918965ea79SLois Curfman McInnes   if (!rank) {
189287828ca2SBarry Smith     ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
18938965ea79SLois Curfman McInnes 
18948965ea79SLois Curfman McInnes     /* read in my part of the matrix numerical values  */
18958965ea79SLois Curfman McInnes     nz = procsnz[0];
18960752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
18978965ea79SLois Curfman McInnes 
18988965ea79SLois Curfman McInnes     /* insert into matrix */
18998965ea79SLois Curfman McInnes     jj      = rstart;
19008965ea79SLois Curfman McInnes     smycols = mycols;
19018965ea79SLois Curfman McInnes     svals   = vals;
19028965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
19038965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
19048965ea79SLois Curfman McInnes       smycols += ourlens[i];
19058965ea79SLois Curfman McInnes       svals   += ourlens[i];
19068965ea79SLois Curfman McInnes       jj++;
19078965ea79SLois Curfman McInnes     }
19088965ea79SLois Curfman McInnes 
19098965ea79SLois Curfman McInnes     /* read in other processors and ship out */
19108965ea79SLois Curfman McInnes     for (i=1; i<size; i++) {
19118965ea79SLois Curfman McInnes       nz   = procsnz[i];
19120752156aSBarry Smith       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
19137adad957SLisandro Dalcin       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)A)->tag,comm);CHKERRQ(ierr);
19148965ea79SLois Curfman McInnes     }
1915606d414cSSatish Balay     ierr = PetscFree(procsnz);CHKERRQ(ierr);
19163a40ed3dSBarry Smith   } else {
19178965ea79SLois Curfman McInnes     /* receive numeric values */
191884d528b1SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
19198965ea79SLois Curfman McInnes 
19208965ea79SLois Curfman McInnes     /* receive message of values*/
19217adad957SLisandro Dalcin     ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)A)->tag,comm,&status);CHKERRQ(ierr);
1922ca161407SBarry Smith     ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
192329bbc08cSBarry Smith     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
19248965ea79SLois Curfman McInnes 
19258965ea79SLois Curfman McInnes     /* insert into matrix */
19268965ea79SLois Curfman McInnes     jj      = rstart;
19278965ea79SLois Curfman McInnes     smycols = mycols;
19288965ea79SLois Curfman McInnes     svals   = vals;
19298965ea79SLois Curfman McInnes     for (i=0; i<m; i++) {
19308965ea79SLois Curfman McInnes       ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr);
19318965ea79SLois Curfman McInnes       smycols += ourlens[i];
19328965ea79SLois Curfman McInnes       svals   += ourlens[i];
19338965ea79SLois Curfman McInnes       jj++;
19348965ea79SLois Curfman McInnes     }
19358965ea79SLois Curfman McInnes   }
1936606d414cSSatish Balay   ierr = PetscFree(ourlens);CHKERRQ(ierr);
1937606d414cSSatish Balay   ierr = PetscFree(vals);CHKERRQ(ierr);
1938606d414cSSatish Balay   ierr = PetscFree(mycols);CHKERRQ(ierr);
1939606d414cSSatish Balay   ierr = PetscFree(rowners);CHKERRQ(ierr);
19408965ea79SLois Curfman McInnes 
19416d4a8577SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19426d4a8577SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19433a40ed3dSBarry Smith   PetscFunctionReturn(0);
19448965ea79SLois Curfman McInnes }
194590ace30eSBarry Smith 
19466e4ee0c6SHong Zhang #undef __FUNCT__
19476e4ee0c6SHong Zhang #define __FUNCT__ "MatEqual_MPIDense"
19486e4ee0c6SHong Zhang PetscErrorCode MatEqual_MPIDense(Mat A,Mat B,PetscTruth *flag)
19496e4ee0c6SHong Zhang {
19506e4ee0c6SHong Zhang   Mat_MPIDense   *matB = (Mat_MPIDense*)B->data,*matA = (Mat_MPIDense*)A->data;
19516e4ee0c6SHong Zhang   Mat            a,b;
19526e4ee0c6SHong Zhang   PetscTruth     flg;
19536e4ee0c6SHong Zhang   PetscErrorCode ierr;
195490ace30eSBarry Smith 
19556e4ee0c6SHong Zhang   PetscFunctionBegin;
19566e4ee0c6SHong Zhang   a = matA->A;
19576e4ee0c6SHong Zhang   b = matB->A;
19586e4ee0c6SHong Zhang   ierr = MatEqual(a,b,&flg);CHKERRQ(ierr);
19597adad957SLisandro Dalcin   ierr = MPI_Allreduce(&flg,flag,1,MPI_INT,MPI_LAND,((PetscObject)A)->comm);CHKERRQ(ierr);
19606e4ee0c6SHong Zhang   PetscFunctionReturn(0);
19616e4ee0c6SHong Zhang }
196290ace30eSBarry Smith 
196309d27a7eSBarry Smith #if defined(PETSC_HAVE_PLAPACK)
196409d27a7eSBarry Smith 
196509d27a7eSBarry Smith #undef __FUNCT__
196609d27a7eSBarry Smith #define __FUNCT__ "PetscPLAPACKFinalizePackage"
196709d27a7eSBarry Smith /*@C
196809d27a7eSBarry Smith   PetscPLAPACKFinalizePackage - This function destroys everything in the Petsc interface to PLAPACK.
196909d27a7eSBarry Smith   Level: developer
197009d27a7eSBarry Smith 
197109d27a7eSBarry Smith .keywords: Petsc, destroy, package, PLAPACK
197209d27a7eSBarry Smith .seealso: PetscFinalize()
197309d27a7eSBarry Smith @*/
197409d27a7eSBarry Smith PetscErrorCode PETSC_DLLEXPORT PetscPLAPACKFinalizePackage(void)
197509d27a7eSBarry Smith {
197609d27a7eSBarry Smith   PetscErrorCode ierr;
197709d27a7eSBarry Smith 
197809d27a7eSBarry Smith   PetscFunctionBegin;
197909d27a7eSBarry Smith   ierr = PLA_Finalize();CHKERRQ(ierr);
198009d27a7eSBarry Smith   PetscFunctionReturn(0);
198109d27a7eSBarry Smith }
198209d27a7eSBarry Smith 
198309d27a7eSBarry Smith #undef __FUNCT__
198409d27a7eSBarry Smith #define __FUNCT__ "PetscPLAPACKInitializePackage"
198509d27a7eSBarry Smith /*@C
198609d27a7eSBarry Smith   PetscPLAPACKInitializePackage - This function initializes everything in the Petsc interface to PLAPACK. It is
198709d27a7eSBarry Smith   called from PetscDLLibraryRegister() when using dynamic libraries, and on the call to PetscInitialize()
198809d27a7eSBarry Smith   when using static libraries.
198909d27a7eSBarry Smith 
199009d27a7eSBarry Smith   Input Parameter:
199109d27a7eSBarry Smith   path - The dynamic library path, or PETSC_NULL
199209d27a7eSBarry Smith 
199309d27a7eSBarry Smith   Level: developer
199409d27a7eSBarry Smith 
199509d27a7eSBarry Smith .keywords: Petsc, initialize, package, PLAPACK
199609d27a7eSBarry Smith .seealso: PetscInitializePackage(), PetscInitialize()
199709d27a7eSBarry Smith @*/
199809d27a7eSBarry Smith PetscErrorCode PETSC_DLLEXPORT PetscPLAPACKInitializePackage(const char path[])
199909d27a7eSBarry Smith {
200004fea9ffSBarry Smith   MPI_Comm       comm = PETSC_COMM_WORLD;
2001*eae6fb2eSBarry Smith   PetscMPIInt    size;
200209d27a7eSBarry Smith   PetscErrorCode ierr;
200309d27a7eSBarry Smith 
200409d27a7eSBarry Smith   PetscFunctionBegin;
200509d27a7eSBarry Smith   if (!PLA_Initialized(PETSC_NULL)) {
200609d27a7eSBarry Smith 
200709d27a7eSBarry Smith     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
200804fea9ffSBarry Smith     Plapack_nprows = 1;
200904fea9ffSBarry Smith     Plapack_npcols = size;
201009d27a7eSBarry Smith 
2011aeccfd6fSBarry Smith     ierr = PetscOptionsBegin(comm,PETSC_NULL,"PLAPACK Options","Mat");CHKERRQ(ierr);
2012aeccfd6fSBarry Smith       ierr = PetscOptionsInt("-plapack_nprows","row dimension of 2D processor mesh","None",Plapack_nprows,&Plapack_nprows,PETSC_NULL);CHKERRQ(ierr);
2013aeccfd6fSBarry Smith       ierr = PetscOptionsInt("-plapack_npcols","column dimension of 2D processor mesh","None",Plapack_npcols,&Plapack_npcols,PETSC_NULL);CHKERRQ(ierr);
2014*eae6fb2eSBarry Smith       ierr = PetscOptionsInt("-plapack_ckerror","error checking flag","None",Plapack_ierror,&Plapack_ierror,PETSC_NULL);CHKERRQ(ierr);
2015*eae6fb2eSBarry Smith       if (Plapack_ierror){
2016*eae6fb2eSBarry Smith 	ierr = PLA_Set_error_checking(Plapack_ierror,PETSC_TRUE,PETSC_TRUE,PETSC_FALSE );CHKERRQ(ierr);
2017aeccfd6fSBarry Smith       } else {
2018*eae6fb2eSBarry Smith 	ierr = PLA_Set_error_checking(Plapack_ierror,PETSC_FALSE,PETSC_FALSE,PETSC_FALSE );CHKERRQ(ierr);
2019aeccfd6fSBarry Smith       }
2020aeccfd6fSBarry Smith 
2021*eae6fb2eSBarry Smith       Plapack_nb_alg = 0;
2022*eae6fb2eSBarry Smith       ierr = PetscOptionsInt("-mat_plapack_nb_alg","algorithmic block size","None",Plapack_nb_alg,&Plapack_nb_alg,PETSC_NULL);CHKERRQ(ierr);
2023*eae6fb2eSBarry Smith       if (Plapack_nb_alg) {
2024*eae6fb2eSBarry Smith         ierr = pla_Environ_set_nb_alg (PLA_OP_ALL_ALG,Plapack_nb_alg);CHKERRQ(ierr);
2025aeccfd6fSBarry Smith       }
2026aeccfd6fSBarry Smith     PetscOptionsEnd();
2027aeccfd6fSBarry Smith 
202804fea9ffSBarry Smith     ierr = PLA_Comm_1D_to_2D(comm,Plapack_nprows,Plapack_npcols,&Plapack_comm_2d);CHKERRQ(ierr);
202904fea9ffSBarry Smith     ierr = PLA_Init(Plapack_comm_2d);CHKERRQ(ierr);
203009d27a7eSBarry Smith     ierr = PetscRegisterFinalize(PetscPLAPACKFinalizePackage);CHKERRQ(ierr);
203109d27a7eSBarry Smith   }
203209d27a7eSBarry Smith   PetscFunctionReturn(0);
203309d27a7eSBarry Smith }
203490ace30eSBarry Smith 
203590ace30eSBarry Smith 
203609d27a7eSBarry Smith #endif
2037