xref: /petsc/src/mat/impls/elemental/matelem.cxx (revision eb1ec7c1db909b5b8e5e7107245f4f2eaa6849e7)
1db31f6deSJed Brown #include <../src/mat/impls/elemental/matelemimpl.h> /*I "petscmat.h" I*/
2db31f6deSJed Brown 
35e9f5b67SHong Zhang /*
45e9f5b67SHong Zhang     The variable Petsc_Elemental_keyval is used to indicate an MPI attribute that
55e9f5b67SHong Zhang   is attached to a communicator, in this case the attribute is a Mat_Elemental_Grid
65e9f5b67SHong Zhang */
75e9f5b67SHong Zhang static PetscMPIInt Petsc_Elemental_keyval = MPI_KEYVAL_INVALID;
85e9f5b67SHong Zhang 
9db31f6deSJed Brown /*@C
10db31f6deSJed Brown    PetscElementalInitializePackage - Initialize Elemental package
11db31f6deSJed Brown 
12db31f6deSJed Brown    Logically Collective
13db31f6deSJed Brown 
14db31f6deSJed Brown    Level: developer
15db31f6deSJed Brown 
16db31f6deSJed Brown .seealso: MATELEMENTAL, PetscElementalFinalizePackage()
17db31f6deSJed Brown @*/
18607a6623SBarry Smith PetscErrorCode PetscElementalInitializePackage(void)
19db31f6deSJed Brown {
20db31f6deSJed Brown   PetscErrorCode ierr;
21db31f6deSJed Brown 
22db31f6deSJed Brown   PetscFunctionBegin;
23e8146892SSatish Balay   if (El::Initialized()) PetscFunctionReturn(0);
24e05a9d59SSatish Balay   El::Initialize();   /* called by the 1st call of MatCreate_Elemental */
25db31f6deSJed Brown   ierr = PetscRegisterFinalize(PetscElementalFinalizePackage);CHKERRQ(ierr);
26db31f6deSJed Brown   PetscFunctionReturn(0);
27db31f6deSJed Brown }
28db31f6deSJed Brown 
29db31f6deSJed Brown /*@C
30db31f6deSJed Brown    PetscElementalFinalizePackage - Finalize Elemental package
31db31f6deSJed Brown 
32db31f6deSJed Brown    Logically Collective
33db31f6deSJed Brown 
34db31f6deSJed Brown    Level: developer
35db31f6deSJed Brown 
36db31f6deSJed Brown .seealso: MATELEMENTAL, PetscElementalInitializePackage()
37db31f6deSJed Brown @*/
38db31f6deSJed Brown PetscErrorCode PetscElementalFinalizePackage(void)
39db31f6deSJed Brown {
40db31f6deSJed Brown   PetscFunctionBegin;
41e8146892SSatish Balay   El::Finalize();  /* called by PetscFinalize() */
42db31f6deSJed Brown   PetscFunctionReturn(0);
43db31f6deSJed Brown }
44db31f6deSJed Brown 
45db31f6deSJed Brown static PetscErrorCode MatView_Elemental(Mat A,PetscViewer viewer)
46db31f6deSJed Brown {
47db31f6deSJed Brown   PetscErrorCode ierr;
48db31f6deSJed Brown   Mat_Elemental  *a = (Mat_Elemental*)A->data;
49db31f6deSJed Brown   PetscBool      iascii;
50db31f6deSJed Brown 
51db31f6deSJed Brown   PetscFunctionBegin;
52db31f6deSJed Brown   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
53db31f6deSJed Brown   if (iascii) {
54db31f6deSJed Brown     PetscViewerFormat format;
55db31f6deSJed Brown     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
56db31f6deSJed Brown     if (format == PETSC_VIEWER_ASCII_INFO) {
5779673f7bSHong Zhang       /* call elemental viewing function */
582d8adcc7SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Elemental run parameters:\n");CHKERRQ(ierr);
59ed667823SXuan Zhou       ierr = PetscViewerASCIIPrintf(viewer,"  allocated entries=%d\n",(*a->emat).AllocatedMemory());CHKERRQ(ierr);
60ed667823SXuan Zhou       ierr = PetscViewerASCIIPrintf(viewer,"  grid height=%d, grid width=%d\n",(*a->emat).Grid().Height(),(*a->emat).Grid().Width());CHKERRQ(ierr);
614fe7bbcaSHong Zhang       if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
6279673f7bSHong Zhang         /* call elemental viewing function */
63ce94432eSBarry Smith         ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer),"test matview_elemental 2\n");CHKERRQ(ierr);
644fe7bbcaSHong Zhang       }
6579673f7bSHong Zhang 
66db31f6deSJed Brown     } else if (format == PETSC_VIEWER_DEFAULT) {
67db31f6deSJed Brown       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
68e8146892SSatish Balay       El::Print( *a->emat, "Elemental matrix (cyclic ordering)" );
69db31f6deSJed Brown       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
70834d3fecSHong Zhang       if (A->factortype == MAT_FACTOR_NONE){
7161119200SXuan Zhou         Mat Adense;
72ce94432eSBarry Smith         ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer),"Elemental matrix (explicit ordering)\n");CHKERRQ(ierr);
7361119200SXuan Zhou         ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense);CHKERRQ(ierr);
7461119200SXuan Zhou         ierr = MatView(Adense,viewer);CHKERRQ(ierr);
7561119200SXuan Zhou         ierr = MatDestroy(&Adense);CHKERRQ(ierr);
76834d3fecSHong Zhang       }
77ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Format");
78d2daa67eSHong Zhang   } else {
795cb544a0SHong Zhang     /* convert to dense format and call MatView() */
8061119200SXuan Zhou     Mat Adense;
81ce94432eSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer),"Elemental matrix (explicit ordering)\n");CHKERRQ(ierr);
8261119200SXuan Zhou     ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense);CHKERRQ(ierr);
8361119200SXuan Zhou     ierr = MatView(Adense,viewer);CHKERRQ(ierr);
8461119200SXuan Zhou     ierr = MatDestroy(&Adense);CHKERRQ(ierr);
85d2daa67eSHong Zhang   }
86db31f6deSJed Brown   PetscFunctionReturn(0);
87db31f6deSJed Brown }
88db31f6deSJed Brown 
8915767789SHong Zhang static PetscErrorCode MatGetInfo_Elemental(Mat A,MatInfoType flag,MatInfo *info)
90180a43e4SHong Zhang {
9115767789SHong Zhang   Mat_Elemental  *a = (Mat_Elemental*)A->data;
9215767789SHong Zhang 
93180a43e4SHong Zhang   PetscFunctionBegin;
945cb544a0SHong Zhang   info->block_size = 1.0;
9515767789SHong Zhang 
9615767789SHong Zhang   if (flag == MAT_LOCAL) {
973966268fSBarry Smith     info->nz_allocated   = (*a->emat).AllocatedMemory(); /* locally allocated */
9815767789SHong Zhang     info->nz_used        = info->nz_allocated;
9915767789SHong Zhang   } else if (flag == MAT_GLOBAL_MAX) {
100b2566f29SBarry Smith     //ierr = MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)matin));CHKERRQ(ierr);
10115767789SHong Zhang     /* see MatGetInfo_MPIAIJ() for getting global info->nz_allocated! */
10215767789SHong Zhang     //SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_GLOBAL_MAX not written yet");
10315767789SHong Zhang   } else if (flag == MAT_GLOBAL_SUM) {
10415767789SHong Zhang     //SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_GLOBAL_SUM not written yet");
1053966268fSBarry Smith     info->nz_allocated   = (*a->emat).AllocatedMemory(); /* locally allocated */
10615767789SHong Zhang     info->nz_used        = info->nz_allocated; /* assume Elemental does accurate allocation */
107b2566f29SBarry Smith     //ierr = MPIU_Allreduce(isend,irecv,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
10815767789SHong Zhang     //PetscPrintf(PETSC_COMM_SELF,"    ... [%d] locally allocated %g\n",rank,info->nz_allocated);
10915767789SHong Zhang   }
11015767789SHong Zhang 
11115767789SHong Zhang   info->nz_unneeded       = 0.0;
1123966268fSBarry Smith   info->assemblies        = A->num_ass;
11315767789SHong Zhang   info->mallocs           = 0;
11415767789SHong Zhang   info->memory            = ((PetscObject)A)->mem;
11515767789SHong Zhang   info->fill_ratio_given  = 0; /* determined by Elemental */
11615767789SHong Zhang   info->fill_ratio_needed = 0;
11715767789SHong Zhang   info->factor_mallocs    = 0;
118180a43e4SHong Zhang   PetscFunctionReturn(0);
119180a43e4SHong Zhang }
120180a43e4SHong Zhang 
12132d7a744SHong Zhang PetscErrorCode MatSetOption_Elemental(Mat A,MatOption op,PetscBool flg)
12232d7a744SHong Zhang {
12332d7a744SHong Zhang   Mat_Elemental  *a = (Mat_Elemental*)A->data;
12432d7a744SHong Zhang 
12532d7a744SHong Zhang   PetscFunctionBegin;
12632d7a744SHong Zhang   switch (op) {
12732d7a744SHong Zhang   case MAT_NEW_NONZERO_LOCATIONS:
12832d7a744SHong Zhang   case MAT_NEW_NONZERO_LOCATION_ERR:
12932d7a744SHong Zhang   case MAT_NEW_NONZERO_ALLOCATION_ERR:
130891a0571SHong Zhang   case MAT_SYMMETRIC:
131071fcb05SBarry Smith   case MAT_SORTED_FULL:
132*eb1ec7c1SStefano Zampini   case MAT_HERMITIAN:
133891a0571SHong Zhang     break;
13432d7a744SHong Zhang   case MAT_ROW_ORIENTED:
13532d7a744SHong Zhang     a->roworiented = flg;
13632d7a744SHong Zhang     break;
13732d7a744SHong Zhang   default:
13832d7a744SHong Zhang     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
13932d7a744SHong Zhang   }
14032d7a744SHong Zhang   PetscFunctionReturn(0);
14132d7a744SHong Zhang }
14232d7a744SHong Zhang 
143e6dea9dbSXuan Zhou static PetscErrorCode MatSetValues_Elemental(Mat A,PetscInt nr,const PetscInt *rows,PetscInt nc,const PetscInt *cols,const PetscScalar *vals,InsertMode imode)
144db31f6deSJed Brown {
145db31f6deSJed Brown   Mat_Elemental  *a = (Mat_Elemental*)A->data;
146fbbcb730SJack Poulson   PetscInt       i,j,rrank,ridx,crank,cidx,erow,ecol,numQueues=0;
147db31f6deSJed Brown 
148db31f6deSJed Brown   PetscFunctionBegin;
149fbbcb730SJack Poulson   // TODO: Initialize matrix to all zeros?
150fbbcb730SJack Poulson 
151fbbcb730SJack Poulson   // Count the number of queues from this process
15232d7a744SHong Zhang   if (a->roworiented) {
153db31f6deSJed Brown     for (i=0; i<nr; i++) {
154db31f6deSJed Brown       if (rows[i] < 0) continue;
155db31f6deSJed Brown       P2RO(A,0,rows[i],&rrank,&ridx);
156db31f6deSJed Brown       RO2E(A,0,rrank,ridx,&erow);
157ce94432eSBarry Smith       if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Incorrect row translation");
158db31f6deSJed Brown       for (j=0; j<nc; j++) {
159db31f6deSJed Brown         if (cols[j] < 0) continue;
160db31f6deSJed Brown         P2RO(A,1,cols[j],&crank,&cidx);
161db31f6deSJed Brown         RO2E(A,1,crank,cidx,&ecol);
162ce94432eSBarry Smith         if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Incorrect col translation");
163fbbcb730SJack Poulson         if (!a->emat->IsLocal(erow,ecol) ){ /* off-proc entry */
164fbbcb730SJack Poulson           /* printf("Will later remotely update (%d,%d)\n",erow,ecol); */
165aae2c449SHong Zhang           if (imode != ADD_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only ADD_VALUES to off-processor entry is supported");
166fbbcb730SJack Poulson           ++numQueues;
167aae2c449SHong Zhang           continue;
168ed36708cSHong Zhang         }
169fbbcb730SJack Poulson         /* printf("Locally updating (%d,%d)\n",erow,ecol); */
170db31f6deSJed Brown         switch (imode) {
171fbbcb730SJack Poulson         case INSERT_VALUES: a->emat->Set(erow,ecol,(PetscElemScalar)vals[i*nc+j]); break;
172fbbcb730SJack Poulson         case ADD_VALUES: a->emat->Update(erow,ecol,(PetscElemScalar)vals[i*nc+j]); break;
173ce94432eSBarry Smith         default: SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for InsertMode %d",(int)imode);
174db31f6deSJed Brown         }
175db31f6deSJed Brown       }
176db31f6deSJed Brown     }
177fbbcb730SJack Poulson 
178fbbcb730SJack Poulson     /* printf("numQueues=%d\n",numQueues); */
179fbbcb730SJack Poulson     a->emat->Reserve( numQueues );
180fbbcb730SJack Poulson     for (i=0; i<nr; i++) {
181fbbcb730SJack Poulson       if (rows[i] < 0) continue;
182fbbcb730SJack Poulson       P2RO(A,0,rows[i],&rrank,&ridx);
183fbbcb730SJack Poulson       RO2E(A,0,rrank,ridx,&erow);
184fbbcb730SJack Poulson       for (j=0; j<nc; j++) {
185fbbcb730SJack Poulson         if (cols[j] < 0) continue;
186fbbcb730SJack Poulson         P2RO(A,1,cols[j],&crank,&cidx);
187fbbcb730SJack Poulson         RO2E(A,1,crank,cidx,&ecol);
188fbbcb730SJack Poulson         if ( !a->emat->IsLocal(erow,ecol) ) { /*off-proc entry*/
189fbbcb730SJack Poulson           /* printf("Queueing remotely update of (%d,%d)\n",erow,ecol); */
190fbbcb730SJack Poulson           a->emat->QueueUpdate( erow, ecol, vals[i*nc+j] );
191fbbcb730SJack Poulson         }
192fbbcb730SJack Poulson       }
193fbbcb730SJack Poulson     }
19432d7a744SHong Zhang   } else { /* columnoriented */
19532d7a744SHong Zhang     for (j=0; j<nc; j++) {
19632d7a744SHong Zhang       if (cols[j] < 0) continue;
19732d7a744SHong Zhang       P2RO(A,1,cols[j],&crank,&cidx);
19832d7a744SHong Zhang       RO2E(A,1,crank,cidx,&ecol);
19932d7a744SHong Zhang       if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Incorrect col translation");
20032d7a744SHong Zhang       for (i=0; i<nr; i++) {
20132d7a744SHong Zhang         if (rows[i] < 0) continue;
20232d7a744SHong Zhang         P2RO(A,0,rows[i],&rrank,&ridx);
20332d7a744SHong Zhang         RO2E(A,0,rrank,ridx,&erow);
20432d7a744SHong Zhang         if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Incorrect row translation");
20532d7a744SHong Zhang         if (!a->emat->IsLocal(erow,ecol) ){ /* off-proc entry */
20632d7a744SHong Zhang           /* printf("Will later remotely update (%d,%d)\n",erow,ecol); */
20732d7a744SHong Zhang           if (imode != ADD_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only ADD_VALUES to off-processor entry is supported");
20832d7a744SHong Zhang           ++numQueues;
20932d7a744SHong Zhang           continue;
21032d7a744SHong Zhang         }
21132d7a744SHong Zhang         /* printf("Locally updating (%d,%d)\n",erow,ecol); */
21232d7a744SHong Zhang         switch (imode) {
21332d7a744SHong Zhang         case INSERT_VALUES: a->emat->Set(erow,ecol,(PetscElemScalar)vals[i+j*nr]); break;
21432d7a744SHong Zhang         case ADD_VALUES: a->emat->Update(erow,ecol,(PetscElemScalar)vals[i+j*nr]); break;
21532d7a744SHong Zhang         default: SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for InsertMode %d",(int)imode);
21632d7a744SHong Zhang         }
21732d7a744SHong Zhang       }
21832d7a744SHong Zhang     }
21932d7a744SHong Zhang 
22032d7a744SHong Zhang     /* printf("numQueues=%d\n",numQueues); */
22132d7a744SHong Zhang     a->emat->Reserve( numQueues );
22232d7a744SHong Zhang     for (j=0; j<nc; j++) {
22332d7a744SHong Zhang       if (cols[j] < 0) continue;
22432d7a744SHong Zhang       P2RO(A,1,cols[j],&crank,&cidx);
22532d7a744SHong Zhang       RO2E(A,1,crank,cidx,&ecol);
22632d7a744SHong Zhang 
22732d7a744SHong Zhang       for (i=0; i<nr; i++) {
22832d7a744SHong Zhang         if (rows[i] < 0) continue;
22932d7a744SHong Zhang         P2RO(A,0,rows[i],&rrank,&ridx);
23032d7a744SHong Zhang         RO2E(A,0,rrank,ridx,&erow);
23132d7a744SHong Zhang         if ( !a->emat->IsLocal(erow,ecol) ) { /*off-proc entry*/
23232d7a744SHong Zhang           /* printf("Queueing remotely update of (%d,%d)\n",erow,ecol); */
23332d7a744SHong Zhang           a->emat->QueueUpdate( erow, ecol, vals[i+j*nr] );
23432d7a744SHong Zhang         }
23532d7a744SHong Zhang       }
23632d7a744SHong Zhang     }
23732d7a744SHong Zhang   }
238db31f6deSJed Brown   PetscFunctionReturn(0);
239db31f6deSJed Brown }
240db31f6deSJed Brown 
241db31f6deSJed Brown static PetscErrorCode MatMult_Elemental(Mat A,Vec X,Vec Y)
242db31f6deSJed Brown {
243db31f6deSJed Brown   Mat_Elemental         *a = (Mat_Elemental*)A->data;
244db31f6deSJed Brown   PetscErrorCode        ierr;
245e6dea9dbSXuan Zhou   const PetscElemScalar *x;
246e6dea9dbSXuan Zhou   PetscElemScalar       *y;
247df311e6cSXuan Zhou   PetscElemScalar       one = 1,zero = 0;
248db31f6deSJed Brown 
249db31f6deSJed Brown   PetscFunctionBegin;
250e6dea9dbSXuan Zhou   ierr = VecGetArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr);
251e6dea9dbSXuan Zhou   ierr = VecGetArray(Y,(PetscScalar **)&y);CHKERRQ(ierr);
252db31f6deSJed Brown   { /* Scoping so that constructor is called before pointer is returned */
253e8146892SSatish Balay     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ye;
2540c18141cSBarry Smith     xe.LockedAttach(A->cmap->N,1,*a->grid,0,0,x,A->cmap->n);
2550c18141cSBarry Smith     ye.Attach(A->rmap->N,1,*a->grid,0,0,y,A->rmap->n);
256e8146892SSatish Balay     El::Gemv(El::NORMAL,one,*a->emat,xe,zero,ye);
257db31f6deSJed Brown   }
258e6dea9dbSXuan Zhou   ierr = VecRestoreArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr);
259e6dea9dbSXuan Zhou   ierr = VecRestoreArray(Y,(PetscScalar **)&y);CHKERRQ(ierr);
260db31f6deSJed Brown   PetscFunctionReturn(0);
261db31f6deSJed Brown }
262db31f6deSJed Brown 
2639426833fSXuan Zhou static PetscErrorCode MatMultTranspose_Elemental(Mat A,Vec X,Vec Y)
2649426833fSXuan Zhou {
2659426833fSXuan Zhou   Mat_Elemental         *a = (Mat_Elemental*)A->data;
2669426833fSXuan Zhou   PetscErrorCode        ierr;
267df311e6cSXuan Zhou   const PetscElemScalar *x;
268df311e6cSXuan Zhou   PetscElemScalar       *y;
269e6dea9dbSXuan Zhou   PetscElemScalar       one = 1,zero = 0;
2709426833fSXuan Zhou 
2719426833fSXuan Zhou   PetscFunctionBegin;
272e6dea9dbSXuan Zhou   ierr = VecGetArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr);
273e6dea9dbSXuan Zhou   ierr = VecGetArray(Y,(PetscScalar **)&y);CHKERRQ(ierr);
2749426833fSXuan Zhou   { /* Scoping so that constructor is called before pointer is returned */
275e8146892SSatish Balay     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ye;
2760c18141cSBarry Smith     xe.LockedAttach(A->rmap->N,1,*a->grid,0,0,x,A->rmap->n);
2770c18141cSBarry Smith     ye.Attach(A->cmap->N,1,*a->grid,0,0,y,A->cmap->n);
278e8146892SSatish Balay     El::Gemv(El::TRANSPOSE,one,*a->emat,xe,zero,ye);
2799426833fSXuan Zhou   }
280e6dea9dbSXuan Zhou   ierr = VecRestoreArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr);
281e6dea9dbSXuan Zhou   ierr = VecRestoreArray(Y,(PetscScalar **)&y);CHKERRQ(ierr);
2829426833fSXuan Zhou   PetscFunctionReturn(0);
2839426833fSXuan Zhou }
2849426833fSXuan Zhou 
285db31f6deSJed Brown static PetscErrorCode MatMultAdd_Elemental(Mat A,Vec X,Vec Y,Vec Z)
286db31f6deSJed Brown {
287db31f6deSJed Brown   Mat_Elemental         *a = (Mat_Elemental*)A->data;
288db31f6deSJed Brown   PetscErrorCode        ierr;
289df311e6cSXuan Zhou   const PetscElemScalar *x;
290df311e6cSXuan Zhou   PetscElemScalar       *z;
291e6dea9dbSXuan Zhou   PetscElemScalar       one = 1;
292db31f6deSJed Brown 
293db31f6deSJed Brown   PetscFunctionBegin;
294db31f6deSJed Brown   if (Y != Z) {ierr = VecCopy(Y,Z);CHKERRQ(ierr);}
295e6dea9dbSXuan Zhou   ierr = VecGetArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr);
296e6dea9dbSXuan Zhou   ierr = VecGetArray(Z,(PetscScalar **)&z);CHKERRQ(ierr);
297db31f6deSJed Brown   { /* Scoping so that constructor is called before pointer is returned */
298e8146892SSatish Balay     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ze;
2990c18141cSBarry Smith     xe.LockedAttach(A->cmap->N,1,*a->grid,0,0,x,A->cmap->n);
3000c18141cSBarry Smith     ze.Attach(A->rmap->N,1,*a->grid,0,0,z,A->rmap->n);
301e8146892SSatish Balay     El::Gemv(El::NORMAL,one,*a->emat,xe,one,ze);
302db31f6deSJed Brown   }
303e6dea9dbSXuan Zhou   ierr = VecRestoreArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr);
304e6dea9dbSXuan Zhou   ierr = VecRestoreArray(Z,(PetscScalar **)&z);CHKERRQ(ierr);
305db31f6deSJed Brown   PetscFunctionReturn(0);
306db31f6deSJed Brown }
307db31f6deSJed Brown 
308e883f9d5SXuan Zhou static PetscErrorCode MatMultTransposeAdd_Elemental(Mat A,Vec X,Vec Y,Vec Z)
309e883f9d5SXuan Zhou {
310e883f9d5SXuan Zhou   Mat_Elemental         *a = (Mat_Elemental*)A->data;
311e883f9d5SXuan Zhou   PetscErrorCode        ierr;
312df311e6cSXuan Zhou   const PetscElemScalar *x;
313df311e6cSXuan Zhou   PetscElemScalar       *z;
314e6dea9dbSXuan Zhou   PetscElemScalar       one = 1;
315e883f9d5SXuan Zhou 
316e883f9d5SXuan Zhou   PetscFunctionBegin;
317e883f9d5SXuan Zhou   if (Y != Z) {ierr = VecCopy(Y,Z);CHKERRQ(ierr);}
318e6dea9dbSXuan Zhou   ierr = VecGetArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr);
319e6dea9dbSXuan Zhou   ierr = VecGetArray(Z,(PetscScalar **)&z);CHKERRQ(ierr);
320e883f9d5SXuan Zhou   { /* Scoping so that constructor is called before pointer is returned */
321e8146892SSatish Balay     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ze;
3220c18141cSBarry Smith     xe.LockedAttach(A->rmap->N,1,*a->grid,0,0,x,A->rmap->n);
3230c18141cSBarry Smith     ze.Attach(A->cmap->N,1,*a->grid,0,0,z,A->cmap->n);
324e8146892SSatish Balay     El::Gemv(El::TRANSPOSE,one,*a->emat,xe,one,ze);
325e883f9d5SXuan Zhou   }
326e6dea9dbSXuan Zhou   ierr = VecRestoreArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr);
327e6dea9dbSXuan Zhou   ierr = VecRestoreArray(Z,(PetscScalar **)&z);CHKERRQ(ierr);
328e883f9d5SXuan Zhou   PetscFunctionReturn(0);
329e883f9d5SXuan Zhou }
330e883f9d5SXuan Zhou 
3319a9e8502SHong Zhang static PetscErrorCode MatMatMultNumeric_Elemental(Mat A,Mat B,Mat C)
332c1d1b975SXuan Zhou {
333c1d1b975SXuan Zhou   Mat_Elemental    *a = (Mat_Elemental*)A->data;
334c1d1b975SXuan Zhou   Mat_Elemental    *b = (Mat_Elemental*)B->data;
3359a9e8502SHong Zhang   Mat_Elemental    *c = (Mat_Elemental*)C->data;
336e6dea9dbSXuan Zhou   PetscElemScalar  one = 1,zero = 0;
337c1d1b975SXuan Zhou 
338c1d1b975SXuan Zhou   PetscFunctionBegin;
339aae2c449SHong Zhang   { /* Scoping so that constructor is called before pointer is returned */
340e8146892SSatish Balay     El::Gemm(El::NORMAL,El::NORMAL,one,*a->emat,*b->emat,zero,*c->emat);
341aae2c449SHong Zhang   }
3429a9e8502SHong Zhang   C->assembled = PETSC_TRUE;
3439a9e8502SHong Zhang   PetscFunctionReturn(0);
3449a9e8502SHong Zhang }
3459a9e8502SHong Zhang 
3469a9e8502SHong Zhang static PetscErrorCode MatMatMultSymbolic_Elemental(Mat A,Mat B,PetscReal fill,Mat *C)
3479a9e8502SHong Zhang {
3489a9e8502SHong Zhang   PetscErrorCode ierr;
3499a9e8502SHong Zhang   Mat            Ce;
350ce94432eSBarry Smith   MPI_Comm       comm;
3519a9e8502SHong Zhang 
3529a9e8502SHong Zhang   PetscFunctionBegin;
353ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
3549a9e8502SHong Zhang   ierr = MatCreate(comm,&Ce);CHKERRQ(ierr);
3559a9e8502SHong Zhang   ierr = MatSetSizes(Ce,A->rmap->n,B->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
3569a9e8502SHong Zhang   ierr = MatSetType(Ce,MATELEMENTAL);CHKERRQ(ierr);
3579a9e8502SHong Zhang   ierr = MatSetUp(Ce);CHKERRQ(ierr);
3589a9e8502SHong Zhang   *C = Ce;
3599a9e8502SHong Zhang   PetscFunctionReturn(0);
3609a9e8502SHong Zhang }
3619a9e8502SHong Zhang 
3629a9e8502SHong Zhang static PetscErrorCode MatMatMult_Elemental(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3639a9e8502SHong Zhang {
3649a9e8502SHong Zhang   PetscErrorCode ierr;
3659a9e8502SHong Zhang 
3669a9e8502SHong Zhang   PetscFunctionBegin;
3679a9e8502SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
3683ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3699a9e8502SHong Zhang     ierr = MatMatMultSymbolic_Elemental(A,B,1.0,C);CHKERRQ(ierr);
3703ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3719a9e8502SHong Zhang   }
3723ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3739a9e8502SHong Zhang   ierr = MatMatMultNumeric_Elemental(A,B,*C);CHKERRQ(ierr);
3743ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
375c1d1b975SXuan Zhou   PetscFunctionReturn(0);
376c1d1b975SXuan Zhou }
377c1d1b975SXuan Zhou 
378df311e6cSXuan Zhou static PetscErrorCode MatMatTransposeMultNumeric_Elemental(Mat A,Mat B,Mat C)
379df311e6cSXuan Zhou {
380df311e6cSXuan Zhou   Mat_Elemental      *a = (Mat_Elemental*)A->data;
381df311e6cSXuan Zhou   Mat_Elemental      *b = (Mat_Elemental*)B->data;
382df311e6cSXuan Zhou   Mat_Elemental      *c = (Mat_Elemental*)C->data;
383e6dea9dbSXuan Zhou   PetscElemScalar    one = 1,zero = 0;
384df311e6cSXuan Zhou 
385df311e6cSXuan Zhou   PetscFunctionBegin;
386df311e6cSXuan Zhou   { /* Scoping so that constructor is called before pointer is returned */
387e8146892SSatish Balay     El::Gemm(El::NORMAL,El::TRANSPOSE,one,*a->emat,*b->emat,zero,*c->emat);
388df311e6cSXuan Zhou   }
389df311e6cSXuan Zhou   C->assembled = PETSC_TRUE;
390df311e6cSXuan Zhou   PetscFunctionReturn(0);
391df311e6cSXuan Zhou }
392df311e6cSXuan Zhou 
393df311e6cSXuan Zhou static PetscErrorCode MatMatTransposeMultSymbolic_Elemental(Mat A,Mat B,PetscReal fill,Mat *C)
394df311e6cSXuan Zhou {
395df311e6cSXuan Zhou   PetscErrorCode ierr;
396df311e6cSXuan Zhou   Mat            Ce;
397ce94432eSBarry Smith   MPI_Comm       comm;
398df311e6cSXuan Zhou 
399df311e6cSXuan Zhou   PetscFunctionBegin;
400ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
401df311e6cSXuan Zhou   ierr = MatCreate(comm,&Ce);CHKERRQ(ierr);
402df311e6cSXuan Zhou   ierr = MatSetSizes(Ce,A->rmap->n,B->rmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
403df311e6cSXuan Zhou   ierr = MatSetType(Ce,MATELEMENTAL);CHKERRQ(ierr);
404df311e6cSXuan Zhou   ierr = MatSetUp(Ce);CHKERRQ(ierr);
405df311e6cSXuan Zhou   *C = Ce;
406df311e6cSXuan Zhou   PetscFunctionReturn(0);
407df311e6cSXuan Zhou }
408df311e6cSXuan Zhou 
409df311e6cSXuan Zhou static PetscErrorCode MatMatTransposeMult_Elemental(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
410df311e6cSXuan Zhou {
411df311e6cSXuan Zhou   PetscErrorCode ierr;
412df311e6cSXuan Zhou 
413df311e6cSXuan Zhou   PetscFunctionBegin;
414df311e6cSXuan Zhou   if (scall == MAT_INITIAL_MATRIX){
415df311e6cSXuan Zhou     ierr = PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
416df311e6cSXuan Zhou     ierr = MatMatMultSymbolic_Elemental(A,B,1.0,C);CHKERRQ(ierr);
417df311e6cSXuan Zhou     ierr = PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
418df311e6cSXuan Zhou   }
419df311e6cSXuan Zhou   ierr = PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
420df311e6cSXuan Zhou   ierr = MatMatTransposeMultNumeric_Elemental(A,B,*C);CHKERRQ(ierr);
421df311e6cSXuan Zhou   ierr = PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
422df311e6cSXuan Zhou   PetscFunctionReturn(0);
423df311e6cSXuan Zhou }
424df311e6cSXuan Zhou 
425a9d89745SXuan Zhou static PetscErrorCode MatGetDiagonal_Elemental(Mat A,Vec D)
42661119200SXuan Zhou {
427a9d89745SXuan Zhou   PetscInt        i,nrows,ncols,nD,rrank,ridx,crank,cidx;
428a9d89745SXuan Zhou   Mat_Elemental   *a = (Mat_Elemental*)A->data;
42961119200SXuan Zhou   PetscErrorCode  ierr;
430a9d89745SXuan Zhou   PetscElemScalar v;
431ce94432eSBarry Smith   MPI_Comm        comm;
43261119200SXuan Zhou 
43361119200SXuan Zhou   PetscFunctionBegin;
434ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
435a9d89745SXuan Zhou   ierr = MatGetSize(A,&nrows,&ncols);CHKERRQ(ierr);
436a9d89745SXuan Zhou   nD = nrows>ncols ? ncols : nrows;
437a9d89745SXuan Zhou   for (i=0; i<nD; i++) {
438a9d89745SXuan Zhou     PetscInt erow,ecol;
439a9d89745SXuan Zhou     P2RO(A,0,i,&rrank,&ridx);
440a9d89745SXuan Zhou     RO2E(A,0,rrank,ridx,&erow);
441a9d89745SXuan Zhou     if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect row translation");
442a9d89745SXuan Zhou     P2RO(A,1,i,&crank,&cidx);
443a9d89745SXuan Zhou     RO2E(A,1,crank,cidx,&ecol);
444a9d89745SXuan Zhou     if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect col translation");
445a9d89745SXuan Zhou     v = a->emat->Get(erow,ecol);
4467c8b904fSSatish Balay     ierr = VecSetValues(D,1,&i,(PetscScalar*)&v,INSERT_VALUES);CHKERRQ(ierr);
447ade3cc5eSXuan Zhou   }
448a9d89745SXuan Zhou   ierr = VecAssemblyBegin(D);CHKERRQ(ierr);
449a9d89745SXuan Zhou   ierr = VecAssemblyEnd(D);CHKERRQ(ierr);
45061119200SXuan Zhou   PetscFunctionReturn(0);
45161119200SXuan Zhou }
45261119200SXuan Zhou 
453ade3cc5eSXuan Zhou static PetscErrorCode MatDiagonalScale_Elemental(Mat X,Vec L,Vec R)
454ade3cc5eSXuan Zhou {
455ade3cc5eSXuan Zhou   Mat_Elemental         *x = (Mat_Elemental*)X->data;
456ade3cc5eSXuan Zhou   const PetscElemScalar *d;
457ade3cc5eSXuan Zhou   PetscErrorCode        ierr;
458ade3cc5eSXuan Zhou 
459ade3cc5eSXuan Zhou   PetscFunctionBegin;
4609065cd98SJed Brown   if (R) {
461ade3cc5eSXuan Zhou     ierr = VecGetArrayRead(R,(const PetscScalar **)&d);CHKERRQ(ierr);
462e8146892SSatish Balay     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> de;
4630c18141cSBarry Smith     de.LockedAttach(X->cmap->N,1,*x->grid,0,0,d,X->cmap->n);
464e8146892SSatish Balay     El::DiagonalScale(El::RIGHT,El::NORMAL,de,*x->emat);
465ade3cc5eSXuan Zhou     ierr = VecRestoreArrayRead(R,(const PetscScalar **)&d);CHKERRQ(ierr);
4669065cd98SJed Brown   }
4679065cd98SJed Brown   if (L) {
468ade3cc5eSXuan Zhou     ierr = VecGetArrayRead(L,(const PetscScalar **)&d);CHKERRQ(ierr);
469e8146892SSatish Balay     El::DistMatrix<PetscElemScalar,El::VC,El::STAR> de;
4700c18141cSBarry Smith     de.LockedAttach(X->rmap->N,1,*x->grid,0,0,d,X->rmap->n);
471e8146892SSatish Balay     El::DiagonalScale(El::LEFT,El::NORMAL,de,*x->emat);
472ade3cc5eSXuan Zhou     ierr = VecRestoreArrayRead(L,(const PetscScalar **)&d);CHKERRQ(ierr);
473ade3cc5eSXuan Zhou   }
474ade3cc5eSXuan Zhou   PetscFunctionReturn(0);
475ade3cc5eSXuan Zhou }
476ade3cc5eSXuan Zhou 
47734fa46e1SFrancesco Ballarin static PetscErrorCode MatMissingDiagonal_Elemental(Mat A,PetscBool *missing,PetscInt *d)
47834fa46e1SFrancesco Ballarin {
47934fa46e1SFrancesco Ballarin   PetscFunctionBegin;
48034fa46e1SFrancesco Ballarin   *missing = PETSC_FALSE;
48134fa46e1SFrancesco Ballarin   PetscFunctionReturn(0);
48234fa46e1SFrancesco Ballarin }
48334fa46e1SFrancesco Ballarin 
484e6dea9dbSXuan Zhou static PetscErrorCode MatScale_Elemental(Mat X,PetscScalar a)
48565b78793SXuan Zhou {
48665b78793SXuan Zhou   Mat_Elemental  *x = (Mat_Elemental*)X->data;
48765b78793SXuan Zhou 
48865b78793SXuan Zhou   PetscFunctionBegin;
489e8146892SSatish Balay   El::Scale((PetscElemScalar)a,*x->emat);
49065b78793SXuan Zhou   PetscFunctionReturn(0);
49165b78793SXuan Zhou }
49265b78793SXuan Zhou 
493f82baa17SHong Zhang /*
494f82baa17SHong Zhang   MatAXPY - Computes Y = a*X + Y.
495f82baa17SHong Zhang */
496e6dea9dbSXuan Zhou static PetscErrorCode MatAXPY_Elemental(Mat Y,PetscScalar a,Mat X,MatStructure str)
497e09a3074SHong Zhang {
498e09a3074SHong Zhang   Mat_Elemental  *x = (Mat_Elemental*)X->data;
499e09a3074SHong Zhang   Mat_Elemental  *y = (Mat_Elemental*)Y->data;
50068446cf8SJed Brown   PetscErrorCode ierr;
501e09a3074SHong Zhang 
502e09a3074SHong Zhang   PetscFunctionBegin;
503e8146892SSatish Balay   El::Axpy((PetscElemScalar)a,*x->emat,*y->emat);
50421f6c9c4SJed Brown   ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);
505e09a3074SHong Zhang   PetscFunctionReturn(0);
506e09a3074SHong Zhang }
507e09a3074SHong Zhang 
508d6223691SXuan Zhou static PetscErrorCode MatCopy_Elemental(Mat A,Mat B,MatStructure str)
509d6223691SXuan Zhou {
510d6223691SXuan Zhou   Mat_Elemental *a=(Mat_Elemental*)A->data;
511d6223691SXuan Zhou   Mat_Elemental *b=(Mat_Elemental*)B->data;
512a4ee3bb6SSatish Balay   PetscErrorCode ierr;
513d6223691SXuan Zhou 
514d6223691SXuan Zhou   PetscFunctionBegin;
515e8146892SSatish Balay   El::Copy(*a->emat,*b->emat);
516cdc753b6SBarry Smith   ierr = PetscObjectStateIncrease((PetscObject)B);CHKERRQ(ierr);
517d6223691SXuan Zhou   PetscFunctionReturn(0);
518d6223691SXuan Zhou }
519d6223691SXuan Zhou 
520df311e6cSXuan Zhou static PetscErrorCode MatDuplicate_Elemental(Mat A,MatDuplicateOption op,Mat *B)
521df311e6cSXuan Zhou {
522df311e6cSXuan Zhou   Mat            Be;
523ce94432eSBarry Smith   MPI_Comm       comm;
524df311e6cSXuan Zhou   Mat_Elemental  *a=(Mat_Elemental*)A->data;
525df311e6cSXuan Zhou   PetscErrorCode ierr;
526df311e6cSXuan Zhou 
527df311e6cSXuan Zhou   PetscFunctionBegin;
528ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
529df311e6cSXuan Zhou   ierr = MatCreate(comm,&Be);CHKERRQ(ierr);
530df311e6cSXuan Zhou   ierr = MatSetSizes(Be,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
531df311e6cSXuan Zhou   ierr = MatSetType(Be,MATELEMENTAL);CHKERRQ(ierr);
532df311e6cSXuan Zhou   ierr = MatSetUp(Be);CHKERRQ(ierr);
533df311e6cSXuan Zhou   *B = Be;
534df311e6cSXuan Zhou   if (op == MAT_COPY_VALUES) {
535df311e6cSXuan Zhou     Mat_Elemental *b=(Mat_Elemental*)Be->data;
536e8146892SSatish Balay     El::Copy(*a->emat,*b->emat);
537df311e6cSXuan Zhou   }
538df311e6cSXuan Zhou   Be->assembled = PETSC_TRUE;
539df311e6cSXuan Zhou   PetscFunctionReturn(0);
540df311e6cSXuan Zhou }
541df311e6cSXuan Zhou 
542d6223691SXuan Zhou static PetscErrorCode MatTranspose_Elemental(Mat A,MatReuse reuse,Mat *B)
543d6223691SXuan Zhou {
544ec8cb81fSBarry Smith   Mat            Be = *B;
5455262d616SXuan Zhou   PetscErrorCode ierr;
546ce94432eSBarry Smith   MPI_Comm       comm;
5474fe7bbcaSHong Zhang   Mat_Elemental  *a = (Mat_Elemental*)A->data, *b;
548d6223691SXuan Zhou 
549d6223691SXuan Zhou   PetscFunctionBegin;
550ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
5515cb544a0SHong Zhang   /* Only out-of-place supported */
5522847e3fdSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Only out-of-place supported");
5535262d616SXuan Zhou   if (reuse == MAT_INITIAL_MATRIX) {
5545262d616SXuan Zhou     ierr = MatCreate(comm,&Be);CHKERRQ(ierr);
5555262d616SXuan Zhou     ierr = MatSetSizes(Be,A->cmap->n,A->rmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
5565262d616SXuan Zhou     ierr = MatSetType(Be,MATELEMENTAL);CHKERRQ(ierr);
5575262d616SXuan Zhou     ierr = MatSetUp(Be);CHKERRQ(ierr);
5585262d616SXuan Zhou     *B = Be;
5595262d616SXuan Zhou   }
5604fe7bbcaSHong Zhang   b = (Mat_Elemental*)Be->data;
561e8146892SSatish Balay   El::Transpose(*a->emat,*b->emat);
5625262d616SXuan Zhou   Be->assembled = PETSC_TRUE;
563d6223691SXuan Zhou   PetscFunctionReturn(0);
564d6223691SXuan Zhou }
565d6223691SXuan Zhou 
566dfcb0403SXuan Zhou static PetscErrorCode MatConjugate_Elemental(Mat A)
567dfcb0403SXuan Zhou {
568dfcb0403SXuan Zhou   Mat_Elemental  *a = (Mat_Elemental*)A->data;
569dfcb0403SXuan Zhou 
570dfcb0403SXuan Zhou   PetscFunctionBegin;
571e8146892SSatish Balay   El::Conjugate(*a->emat);
572dfcb0403SXuan Zhou   PetscFunctionReturn(0);
573dfcb0403SXuan Zhou }
574dfcb0403SXuan Zhou 
5754a29722dSXuan Zhou static PetscErrorCode MatHermitianTranspose_Elemental(Mat A,MatReuse reuse,Mat *B)
5764a29722dSXuan Zhou {
577ec8cb81fSBarry Smith   Mat            Be = *B;
5784a29722dSXuan Zhou   PetscErrorCode ierr;
579ce94432eSBarry Smith   MPI_Comm       comm;
5804a29722dSXuan Zhou   Mat_Elemental  *a = (Mat_Elemental*)A->data, *b;
5814a29722dSXuan Zhou 
5824a29722dSXuan Zhou   PetscFunctionBegin;
583ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
5845cb544a0SHong Zhang   /* Only out-of-place supported */
5854a29722dSXuan Zhou   if (reuse == MAT_INITIAL_MATRIX){
5864a29722dSXuan Zhou     ierr = MatCreate(comm,&Be);CHKERRQ(ierr);
5874a29722dSXuan Zhou     ierr = MatSetSizes(Be,A->cmap->n,A->rmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
5884a29722dSXuan Zhou     ierr = MatSetType(Be,MATELEMENTAL);CHKERRQ(ierr);
5894a29722dSXuan Zhou     ierr = MatSetUp(Be);CHKERRQ(ierr);
5904a29722dSXuan Zhou     *B = Be;
5914a29722dSXuan Zhou   }
5924a29722dSXuan Zhou   b = (Mat_Elemental*)Be->data;
593e8146892SSatish Balay   El::Adjoint(*a->emat,*b->emat);
5944a29722dSXuan Zhou   Be->assembled = PETSC_TRUE;
5954a29722dSXuan Zhou   PetscFunctionReturn(0);
5964a29722dSXuan Zhou }
5974a29722dSXuan Zhou 
5981f881ff8SXuan Zhou static PetscErrorCode MatSolve_Elemental(Mat A,Vec B,Vec X)
5991f881ff8SXuan Zhou {
6001f881ff8SXuan Zhou   Mat_Elemental     *a = (Mat_Elemental*)A->data;
6011f881ff8SXuan Zhou   PetscErrorCode    ierr;
602df311e6cSXuan Zhou   PetscElemScalar   *x;
603ea394475SSatish Balay   PetscInt          pivoting = a->pivoting;
6041f881ff8SXuan Zhou 
6051f881ff8SXuan Zhou   PetscFunctionBegin;
60645cf121fSXuan Zhou   ierr = VecCopy(B,X);CHKERRQ(ierr);
607e6dea9dbSXuan Zhou   ierr = VecGetArray(X,(PetscScalar **)&x);CHKERRQ(ierr);
608ea394475SSatish Balay 
609e8146892SSatish Balay   El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe;
6100c18141cSBarry Smith   xe.Attach(A->rmap->N,1,*a->grid,0,0,x,A->rmap->n);
611e8146892SSatish Balay   El::DistMatrix<PetscElemScalar,El::MC,El::MR> xer(xe);
612fc54b460SXuan Zhou   switch (A->factortype) {
613fc54b460SXuan Zhou   case MAT_FACTOR_LU:
614ea394475SSatish Balay     if (pivoting == 0) {
615e8146892SSatish Balay       El::lu::SolveAfter(El::NORMAL,*a->emat,xer);
616ea394475SSatish Balay     } else if (pivoting == 1) {
617ea394475SSatish Balay       El::lu::SolveAfter(El::NORMAL,*a->emat,*a->P,xer);
618ea394475SSatish Balay     } else { /* pivoting == 2 */
619ea394475SSatish Balay       El::lu::SolveAfter(El::NORMAL,*a->emat,*a->P,*a->Q,xer);
6201f881ff8SXuan Zhou     }
621fc54b460SXuan Zhou     break;
622fc54b460SXuan Zhou   case MAT_FACTOR_CHOLESKY:
623e8146892SSatish Balay     El::cholesky::SolveAfter(El::UPPER,El::NORMAL,*a->emat,xer);
624fc54b460SXuan Zhou     break;
625fc54b460SXuan Zhou   default:
6264fe7bbcaSHong Zhang     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unfactored Matrix or Unsupported MatFactorType");
627fc54b460SXuan Zhou     break;
6281f881ff8SXuan Zhou   }
629ea394475SSatish Balay   El::Copy(xer,xe);
630ea394475SSatish Balay 
631e6dea9dbSXuan Zhou   ierr = VecRestoreArray(X,(PetscScalar **)&x);CHKERRQ(ierr);
6321f881ff8SXuan Zhou   PetscFunctionReturn(0);
6331f881ff8SXuan Zhou }
6341f881ff8SXuan Zhou 
635df311e6cSXuan Zhou static PetscErrorCode MatSolveAdd_Elemental(Mat A,Vec B,Vec Y,Vec X)
636df311e6cSXuan Zhou {
637df311e6cSXuan Zhou   PetscErrorCode    ierr;
638df311e6cSXuan Zhou 
639df311e6cSXuan Zhou   PetscFunctionBegin;
640df311e6cSXuan Zhou   ierr = MatSolve_Elemental(A,B,X);CHKERRQ(ierr);
6413d7f40dbSXuan Zhou   ierr = VecAXPY(X,1,Y);CHKERRQ(ierr);
642df311e6cSXuan Zhou   PetscFunctionReturn(0);
643df311e6cSXuan Zhou }
644df311e6cSXuan Zhou 
645ae844d54SHong Zhang static PetscErrorCode MatMatSolve_Elemental(Mat A,Mat B,Mat X)
646ae844d54SHong Zhang {
6471f0e42cfSHong Zhang   Mat_Elemental *a=(Mat_Elemental*)A->data;
648d6223691SXuan Zhou   Mat_Elemental *b=(Mat_Elemental*)B->data;
6491f0e42cfSHong Zhang   Mat_Elemental *x=(Mat_Elemental*)X->data;
650ea394475SSatish Balay   PetscInt      pivoting = a->pivoting;
6511f0e42cfSHong Zhang 
652ae844d54SHong Zhang   PetscFunctionBegin;
653e8146892SSatish Balay   El::Copy(*b->emat,*x->emat);
654fc54b460SXuan Zhou   switch (A->factortype) {
655fc54b460SXuan Zhou   case MAT_FACTOR_LU:
656ea394475SSatish Balay     if (pivoting == 0) {
657e8146892SSatish Balay       El::lu::SolveAfter(El::NORMAL,*a->emat,*x->emat);
658ea394475SSatish Balay     } else if (pivoting == 1) {
659ea394475SSatish Balay       El::lu::SolveAfter(El::NORMAL,*a->emat,*a->P,*x->emat);
660ea394475SSatish Balay     } else {
661ea394475SSatish Balay       El::lu::SolveAfter(El::NORMAL,*a->emat,*a->P,*a->Q,*x->emat);
662d6223691SXuan Zhou     }
663fc54b460SXuan Zhou     break;
664fc54b460SXuan Zhou   case MAT_FACTOR_CHOLESKY:
665e8146892SSatish Balay     El::cholesky::SolveAfter(El::UPPER,El::NORMAL,*a->emat,*x->emat);
666fc54b460SXuan Zhou     break;
667fc54b460SXuan Zhou   default:
6684fe7bbcaSHong Zhang     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unfactored Matrix or Unsupported MatFactorType");
669fc54b460SXuan Zhou     break;
670fc54b460SXuan Zhou   }
671ae844d54SHong Zhang   PetscFunctionReturn(0);
672ae844d54SHong Zhang }
673ae844d54SHong Zhang 
674ae844d54SHong Zhang static PetscErrorCode MatLUFactor_Elemental(Mat A,IS row,IS col,const MatFactorInfo *info)
675ae844d54SHong Zhang {
6767c920d81SXuan Zhou   Mat_Elemental  *a = (Mat_Elemental*)A->data;
677124af5d2SSatish Balay   PetscErrorCode ierr;
678ea394475SSatish Balay   PetscInt       pivoting = a->pivoting;
6797c920d81SXuan Zhou 
680ae844d54SHong Zhang   PetscFunctionBegin;
681ea394475SSatish Balay   if (pivoting == 0) {
682e8146892SSatish Balay     El::LU(*a->emat);
683ea394475SSatish Balay   } else if (pivoting == 1) {
684ea394475SSatish Balay     El::LU(*a->emat,*a->P);
685ea394475SSatish Balay   } else {
686ea394475SSatish Balay     El::LU(*a->emat,*a->P,*a->Q);
687d6223691SXuan Zhou   }
6881293973dSHong Zhang   A->factortype = MAT_FACTOR_LU;
689834d3fecSHong Zhang   A->assembled  = PETSC_TRUE;
690f6224b95SHong Zhang 
691f6224b95SHong Zhang   ierr = PetscFree(A->solvertype);CHKERRQ(ierr);
692f6224b95SHong Zhang   ierr = PetscStrallocpy(MATSOLVERELEMENTAL,&A->solvertype);CHKERRQ(ierr);
693ae844d54SHong Zhang   PetscFunctionReturn(0);
694ae844d54SHong Zhang }
695ae844d54SHong Zhang 
696d7c3f9d8SHong Zhang static PetscErrorCode  MatLUFactorNumeric_Elemental(Mat F,Mat A,const MatFactorInfo *info)
697d7c3f9d8SHong Zhang {
698d7c3f9d8SHong Zhang   PetscErrorCode ierr;
699d7c3f9d8SHong Zhang 
700d7c3f9d8SHong Zhang   PetscFunctionBegin;
701d7c3f9d8SHong Zhang   ierr = MatCopy(A,F,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
702d7c3f9d8SHong Zhang   ierr = MatLUFactor_Elemental(F,0,0,info);CHKERRQ(ierr);
703d7c3f9d8SHong Zhang   PetscFunctionReturn(0);
704d7c3f9d8SHong Zhang }
705d7c3f9d8SHong Zhang 
706d7c3f9d8SHong Zhang static PetscErrorCode  MatLUFactorSymbolic_Elemental(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
707d7c3f9d8SHong Zhang {
708d7c3f9d8SHong Zhang   PetscFunctionBegin;
709d7c3f9d8SHong Zhang   /* F is create and allocated by MatGetFactor_elemental_petsc(), skip this routine. */
710d7c3f9d8SHong Zhang   PetscFunctionReturn(0);
711d7c3f9d8SHong Zhang }
712d7c3f9d8SHong Zhang 
71345cf121fSXuan Zhou static PetscErrorCode MatCholeskyFactor_Elemental(Mat A,IS perm,const MatFactorInfo *info)
71445cf121fSXuan Zhou {
71545cf121fSXuan Zhou   Mat_Elemental  *a = (Mat_Elemental*)A->data;
716e8146892SSatish Balay   El::DistMatrix<PetscElemScalar,El::MC,El::STAR> d;
717124af5d2SSatish Balay   PetscErrorCode ierr;
71845cf121fSXuan Zhou 
71945cf121fSXuan Zhou   PetscFunctionBegin;
720e8146892SSatish Balay   El::Cholesky(El::UPPER,*a->emat);
721fc54b460SXuan Zhou   A->factortype = MAT_FACTOR_CHOLESKY;
722834d3fecSHong Zhang   A->assembled  = PETSC_TRUE;
723f6224b95SHong Zhang 
724f6224b95SHong Zhang   ierr = PetscFree(A->solvertype);CHKERRQ(ierr);
725f6224b95SHong Zhang   ierr = PetscStrallocpy(MATSOLVERELEMENTAL,&A->solvertype);CHKERRQ(ierr);
72645cf121fSXuan Zhou   PetscFunctionReturn(0);
72745cf121fSXuan Zhou }
72845cf121fSXuan Zhou 
72979673f7bSHong Zhang static PetscErrorCode MatCholeskyFactorNumeric_Elemental(Mat F,Mat A,const MatFactorInfo *info)
73079673f7bSHong Zhang {
731cb76c1d8SXuan Zhou   PetscErrorCode ierr;
732cb76c1d8SXuan Zhou 
733cb76c1d8SXuan Zhou   PetscFunctionBegin;
734cb76c1d8SXuan Zhou   ierr = MatCopy(A,F,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
735cb76c1d8SXuan Zhou   ierr = MatCholeskyFactor_Elemental(F,0,info);CHKERRQ(ierr);
736cb76c1d8SXuan Zhou   PetscFunctionReturn(0);
73779673f7bSHong Zhang }
738cb76c1d8SXuan Zhou 
73979673f7bSHong Zhang static PetscErrorCode MatCholeskyFactorSymbolic_Elemental(Mat F,Mat A,IS perm,const MatFactorInfo *info)
74079673f7bSHong Zhang {
74179673f7bSHong Zhang   PetscFunctionBegin;
74279673f7bSHong Zhang   /* F is create and allocated by MatGetFactor_elemental_petsc(), skip this routine. */
74379673f7bSHong Zhang   PetscFunctionReturn(0);
74479673f7bSHong Zhang }
74579673f7bSHong Zhang 
746ea799195SBarry Smith PetscErrorCode MatFactorGetSolverType_elemental_elemental(Mat A,MatSolverType *type)
74715767789SHong Zhang {
74815767789SHong Zhang   PetscFunctionBegin;
74915767789SHong Zhang   *type = MATSOLVERELEMENTAL;
75015767789SHong Zhang   PetscFunctionReturn(0);
75115767789SHong Zhang }
75215767789SHong Zhang 
75315767789SHong Zhang static PetscErrorCode MatGetFactor_elemental_elemental(Mat A,MatFactorType ftype,Mat *F)
7541293973dSHong Zhang {
7551293973dSHong Zhang   Mat            B;
7561293973dSHong Zhang   PetscErrorCode ierr;
7571293973dSHong Zhang 
7581293973dSHong Zhang   PetscFunctionBegin;
7591293973dSHong Zhang   /* Create the factorization matrix */
760ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
7611293973dSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
7621293973dSHong Zhang   ierr = MatSetType(B,MATELEMENTAL);CHKERRQ(ierr);
7631293973dSHong Zhang   ierr = MatSetUp(B);CHKERRQ(ierr);
7641293973dSHong Zhang   B->factortype = ftype;
76500c67f3bSHong Zhang   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
76600c67f3bSHong Zhang   ierr = PetscStrallocpy(MATSOLVERELEMENTAL,&B->solvertype);CHKERRQ(ierr);
76700c67f3bSHong Zhang 
7683ca39a21SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_elemental_elemental);CHKERRQ(ierr);
7691293973dSHong Zhang   *F            = B;
7701293973dSHong Zhang   PetscFunctionReturn(0);
7711293973dSHong Zhang }
7721293973dSHong Zhang 
7733ca39a21SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_Elemental(void)
77442c9c57cSBarry Smith {
77518713533SBarry Smith   PetscErrorCode ierr;
77618713533SBarry Smith 
77742c9c57cSBarry Smith   PetscFunctionBegin;
7783ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERELEMENTAL,MATELEMENTAL,        MAT_FACTOR_LU,MatGetFactor_elemental_elemental);CHKERRQ(ierr);
7793ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERELEMENTAL,MATELEMENTAL,        MAT_FACTOR_CHOLESKY,MatGetFactor_elemental_elemental);CHKERRQ(ierr);
78042c9c57cSBarry Smith   PetscFunctionReturn(0);
78142c9c57cSBarry Smith }
78242c9c57cSBarry Smith 
7831f881ff8SXuan Zhou static PetscErrorCode MatNorm_Elemental(Mat A,NormType type,PetscReal *nrm)
7841f881ff8SXuan Zhou {
7851f881ff8SXuan Zhou   Mat_Elemental *a=(Mat_Elemental*)A->data;
7861f881ff8SXuan Zhou 
7871f881ff8SXuan Zhou   PetscFunctionBegin;
7881f881ff8SXuan Zhou   switch (type){
7891f881ff8SXuan Zhou   case NORM_1:
790e8146892SSatish Balay     *nrm = El::OneNorm(*a->emat);
7911f881ff8SXuan Zhou     break;
7921f881ff8SXuan Zhou   case NORM_FROBENIUS:
793e8146892SSatish Balay     *nrm = El::FrobeniusNorm(*a->emat);
7941f881ff8SXuan Zhou     break;
7951f881ff8SXuan Zhou   case NORM_INFINITY:
796e8146892SSatish Balay     *nrm = El::InfinityNorm(*a->emat);
7971f881ff8SXuan Zhou     break;
7981f881ff8SXuan Zhou   default:
7991f881ff8SXuan Zhou     printf("Error: unsupported norm type!\n");
8001f881ff8SXuan Zhou   }
8011f881ff8SXuan Zhou   PetscFunctionReturn(0);
8021f881ff8SXuan Zhou }
8031f881ff8SXuan Zhou 
8045262d616SXuan Zhou static PetscErrorCode MatZeroEntries_Elemental(Mat A)
8055262d616SXuan Zhou {
8065262d616SXuan Zhou   Mat_Elemental *a=(Mat_Elemental*)A->data;
8075262d616SXuan Zhou 
8085262d616SXuan Zhou   PetscFunctionBegin;
809e8146892SSatish Balay   El::Zero(*a->emat);
8105262d616SXuan Zhou   PetscFunctionReturn(0);
8115262d616SXuan Zhou }
8125262d616SXuan Zhou 
813db31f6deSJed Brown static PetscErrorCode MatGetOwnershipIS_Elemental(Mat A,IS *rows,IS *cols)
814db31f6deSJed Brown {
815db31f6deSJed Brown   Mat_Elemental  *a = (Mat_Elemental*)A->data;
816db31f6deSJed Brown   PetscErrorCode ierr;
817db31f6deSJed Brown   PetscInt       i,m,shift,stride,*idx;
818db31f6deSJed Brown 
819db31f6deSJed Brown   PetscFunctionBegin;
820db31f6deSJed Brown   if (rows) {
821db31f6deSJed Brown     m = a->emat->LocalHeight();
822db31f6deSJed Brown     shift = a->emat->ColShift();
823db31f6deSJed Brown     stride = a->emat->ColStride();
824785e854fSJed Brown     ierr = PetscMalloc1(m,&idx);CHKERRQ(ierr);
825db31f6deSJed Brown     for (i=0; i<m; i++) {
826db31f6deSJed Brown       PetscInt rank,offset;
827db31f6deSJed Brown       E2RO(A,0,shift+i*stride,&rank,&offset);
828db31f6deSJed Brown       RO2P(A,0,rank,offset,&idx[i]);
829db31f6deSJed Brown     }
830db31f6deSJed Brown     ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx,PETSC_OWN_POINTER,rows);CHKERRQ(ierr);
831db31f6deSJed Brown   }
832db31f6deSJed Brown   if (cols) {
833db31f6deSJed Brown     m = a->emat->LocalWidth();
834db31f6deSJed Brown     shift = a->emat->RowShift();
835db31f6deSJed Brown     stride = a->emat->RowStride();
836785e854fSJed Brown     ierr = PetscMalloc1(m,&idx);CHKERRQ(ierr);
837db31f6deSJed Brown     for (i=0; i<m; i++) {
838db31f6deSJed Brown       PetscInt rank,offset;
839db31f6deSJed Brown       E2RO(A,1,shift+i*stride,&rank,&offset);
840db31f6deSJed Brown       RO2P(A,1,rank,offset,&idx[i]);
841db31f6deSJed Brown     }
842db31f6deSJed Brown     ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx,PETSC_OWN_POINTER,cols);CHKERRQ(ierr);
843db31f6deSJed Brown   }
844db31f6deSJed Brown   PetscFunctionReturn(0);
845db31f6deSJed Brown }
846db31f6deSJed Brown 
84719fd82e9SBarry Smith static PetscErrorCode MatConvert_Elemental_Dense(Mat A,MatType newtype,MatReuse reuse,Mat *B)
848af295397SXuan Zhou {
8492ef0cf24SXuan Zhou   Mat                Bmpi;
850af295397SXuan Zhou   Mat_Elemental      *a = (Mat_Elemental*)A->data;
851ce94432eSBarry Smith   MPI_Comm           comm;
8522ef0cf24SXuan Zhou   PetscErrorCode     ierr;
853fa9e26c8SHong Zhang   IS                 isrows,iscols;
854fa9e26c8SHong Zhang   PetscInt           rrank,ridx,crank,cidx,nrows,ncols,i,j,erow,ecol,elrow,elcol;
855fa9e26c8SHong Zhang   const PetscInt     *rows,*cols;
856df311e6cSXuan Zhou   PetscElemScalar    v;
857fa9e26c8SHong Zhang   const El::Grid     &grid = a->emat->Grid();
858af295397SXuan Zhou 
859af295397SXuan Zhou   PetscFunctionBegin;
860ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
861a000e53bSHong Zhang 
862de0a22f0SHong Zhang   if (reuse == MAT_REUSE_MATRIX) {
863de0a22f0SHong Zhang     Bmpi = *B;
864de0a22f0SHong Zhang   } else {
865af295397SXuan Zhou     ierr = MatCreate(comm,&Bmpi);CHKERRQ(ierr);
866af295397SXuan Zhou     ierr = MatSetSizes(Bmpi,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
8672ef0cf24SXuan Zhou     ierr = MatSetType(Bmpi,MATDENSE);CHKERRQ(ierr);
868af295397SXuan Zhou     ierr = MatSetUp(Bmpi);CHKERRQ(ierr);
869de0a22f0SHong Zhang   }
870fa9e26c8SHong Zhang 
871fa9e26c8SHong Zhang   /* Get local entries of A */
872fa9e26c8SHong Zhang   ierr = MatGetOwnershipIS(A,&isrows,&iscols);CHKERRQ(ierr);
873fa9e26c8SHong Zhang   ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr);
874fa9e26c8SHong Zhang   ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr);
875fa9e26c8SHong Zhang   ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr);
876fa9e26c8SHong Zhang   ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr);
877fa9e26c8SHong Zhang 
87857299f2fSHong Zhang   if (a->roworiented) {
8792ef0cf24SXuan Zhou     for (i=0; i<nrows; i++) {
880fa9e26c8SHong Zhang       P2RO(A,0,rows[i],&rrank,&ridx); /* convert indices between PETSc <-> (Rank,Offset) <-> Elemental */
8812ef0cf24SXuan Zhou       RO2E(A,0,rrank,ridx,&erow);
8822ef0cf24SXuan Zhou       if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect row translation");
8832ef0cf24SXuan Zhou       for (j=0; j<ncols; j++) {
884fa9e26c8SHong Zhang         P2RO(A,1,cols[j],&crank,&cidx);
8852ef0cf24SXuan Zhou         RO2E(A,1,crank,cidx,&ecol);
8862ef0cf24SXuan Zhou         if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect col translation");
887fa9e26c8SHong Zhang 
888fa9e26c8SHong Zhang         elrow = erow / grid.MCSize(); /* Elemental local row index */
889fa9e26c8SHong Zhang         elcol = ecol / grid.MRSize(); /* Elemental local column index */
890fa9e26c8SHong Zhang         v = a->emat->GetLocal(elrow,elcol);
891fa9e26c8SHong Zhang         ierr = MatSetValues(Bmpi,1,&rows[i],1,&cols[j],(PetscScalar *)&v,INSERT_VALUES);CHKERRQ(ierr);
8922ef0cf24SXuan Zhou       }
8932ef0cf24SXuan Zhou     }
89457299f2fSHong Zhang   } else { /* column-oriented */
89557299f2fSHong Zhang     for (j=0; j<ncols; j++) {
89657299f2fSHong Zhang       P2RO(A,1,cols[j],&crank,&cidx);
89757299f2fSHong Zhang       RO2E(A,1,crank,cidx,&ecol);
89857299f2fSHong Zhang       if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect col translation");
89957299f2fSHong Zhang       for (i=0; i<nrows; i++) {
90057299f2fSHong Zhang         P2RO(A,0,rows[i],&rrank,&ridx); /* convert indices between PETSc <-> (Rank,Offset) <-> Elemental */
90157299f2fSHong Zhang         RO2E(A,0,rrank,ridx,&erow);
90257299f2fSHong Zhang         if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect row translation");
90357299f2fSHong Zhang 
90457299f2fSHong Zhang         elrow = erow / grid.MCSize(); /* Elemental local row index */
90557299f2fSHong Zhang         elcol = ecol / grid.MRSize(); /* Elemental local column index */
90657299f2fSHong Zhang         v = a->emat->GetLocal(elrow,elcol);
90757299f2fSHong Zhang         ierr = MatSetValues(Bmpi,1,&rows[i],1,&cols[j],(PetscScalar *)&v,INSERT_VALUES);CHKERRQ(ierr);
90857299f2fSHong Zhang       }
90957299f2fSHong Zhang     }
91057299f2fSHong Zhang   }
911af295397SXuan Zhou   ierr = MatAssemblyBegin(Bmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
912af295397SXuan Zhou   ierr = MatAssemblyEnd(Bmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
913511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
91428be2f97SBarry Smith     ierr = MatHeaderReplace(A,&Bmpi);CHKERRQ(ierr);
915c4ad791aSXuan Zhou   } else {
916c4ad791aSXuan Zhou     *B = Bmpi;
917c4ad791aSXuan Zhou   }
918fa9e26c8SHong Zhang   ierr = ISDestroy(&isrows);CHKERRQ(ierr);
919fa9e26c8SHong Zhang   ierr = ISDestroy(&iscols);CHKERRQ(ierr);
920af295397SXuan Zhou   PetscFunctionReturn(0);
921af295397SXuan Zhou }
922af295397SXuan Zhou 
923cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
924af8000cdSHong Zhang {
925af8000cdSHong Zhang   Mat               mat_elemental;
926af8000cdSHong Zhang   PetscErrorCode    ierr;
927af8000cdSHong Zhang   PetscInt          M=A->rmap->N,N=A->cmap->N,row,ncols;
928af8000cdSHong Zhang   const PetscInt    *cols;
929af8000cdSHong Zhang   const PetscScalar *vals;
930af8000cdSHong Zhang 
931af8000cdSHong Zhang   PetscFunctionBegin;
932bdeae278SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
933bdeae278SStefano Zampini     mat_elemental = *newmat;
934bdeae278SStefano Zampini     ierr = MatZeroEntries(mat_elemental);CHKERRQ(ierr);
935bdeae278SStefano Zampini   } else {
936af8000cdSHong Zhang     ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr);
937af8000cdSHong Zhang     ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
938af8000cdSHong Zhang     ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr);
939af8000cdSHong Zhang     ierr = MatSetUp(mat_elemental);CHKERRQ(ierr);
940bdeae278SStefano Zampini   }
941af8000cdSHong Zhang   for (row=0; row<M; row++) {
942af8000cdSHong Zhang     ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr);
943bdeae278SStefano Zampini     /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
944af8000cdSHong Zhang     ierr = MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
945af8000cdSHong Zhang     ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr);
946af8000cdSHong Zhang   }
947af8000cdSHong Zhang   ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
948af8000cdSHong Zhang   ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
949af8000cdSHong Zhang 
950511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
95128be2f97SBarry Smith     ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr);
952af8000cdSHong Zhang   } else {
953af8000cdSHong Zhang     *newmat = mat_elemental;
954af8000cdSHong Zhang   }
955af8000cdSHong Zhang   PetscFunctionReturn(0);
956af8000cdSHong Zhang }
957af8000cdSHong Zhang 
958cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
9595d7652ecSHong Zhang {
9605d7652ecSHong Zhang   Mat               mat_elemental;
9615d7652ecSHong Zhang   PetscErrorCode    ierr;
9625d7652ecSHong Zhang   PetscInt          row,ncols,rstart=A->rmap->rstart,rend=A->rmap->rend,j;
9635d7652ecSHong Zhang   const PetscInt    *cols;
9645d7652ecSHong Zhang   const PetscScalar *vals;
9655d7652ecSHong Zhang 
9665d7652ecSHong Zhang   PetscFunctionBegin;
967bdeae278SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
968bdeae278SStefano Zampini     mat_elemental = *newmat;
969bdeae278SStefano Zampini     ierr = MatZeroEntries(mat_elemental);CHKERRQ(ierr);
970bdeae278SStefano Zampini   } else {
9715d7652ecSHong Zhang     ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr);
9725d7652ecSHong Zhang     ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
9735d7652ecSHong Zhang     ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr);
9745d7652ecSHong Zhang     ierr = MatSetUp(mat_elemental);CHKERRQ(ierr);
975bdeae278SStefano Zampini   }
9765d7652ecSHong Zhang   for (row=rstart; row<rend; row++) {
9775d7652ecSHong Zhang     ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr);
9785d7652ecSHong Zhang     for (j=0; j<ncols; j++) {
979bdeae278SStefano Zampini       /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
9805d7652ecSHong Zhang       ierr = MatSetValues(mat_elemental,1,&row,1,&cols[j],&vals[j],ADD_VALUES);CHKERRQ(ierr);
9815d7652ecSHong Zhang     }
9825d7652ecSHong Zhang     ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr);
9835d7652ecSHong Zhang   }
9845d7652ecSHong Zhang   ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9855d7652ecSHong Zhang   ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9865d7652ecSHong Zhang 
987511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
98828be2f97SBarry Smith     ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr);
9895d7652ecSHong Zhang   } else {
9905d7652ecSHong Zhang     *newmat = mat_elemental;
9915d7652ecSHong Zhang   }
9925d7652ecSHong Zhang   PetscFunctionReturn(0);
9935d7652ecSHong Zhang }
9945d7652ecSHong Zhang 
995cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
9966214f412SHong Zhang {
9976214f412SHong Zhang   Mat               mat_elemental;
9986214f412SHong Zhang   PetscErrorCode    ierr;
9996214f412SHong Zhang   PetscInt          M=A->rmap->N,N=A->cmap->N,row,ncols,j;
10006214f412SHong Zhang   const PetscInt    *cols;
10016214f412SHong Zhang   const PetscScalar *vals;
10026214f412SHong Zhang 
10036214f412SHong Zhang   PetscFunctionBegin;
1004bdeae278SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
1005bdeae278SStefano Zampini     mat_elemental = *newmat;
1006bdeae278SStefano Zampini     ierr = MatZeroEntries(mat_elemental);CHKERRQ(ierr);
1007bdeae278SStefano Zampini   } else {
10086214f412SHong Zhang     ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr);
10096214f412SHong Zhang     ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
10106214f412SHong Zhang     ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr);
10116214f412SHong Zhang     ierr = MatSetUp(mat_elemental);CHKERRQ(ierr);
1012bdeae278SStefano Zampini   }
10136214f412SHong Zhang   ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr);
10146214f412SHong Zhang   for (row=0; row<M; row++) {
10156214f412SHong Zhang     ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr);
1016bdeae278SStefano Zampini     /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
10176214f412SHong Zhang     ierr = MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
10186214f412SHong Zhang     for (j=0; j<ncols; j++) { /* lower triangular part */
1019*eb1ec7c1SStefano Zampini       PetscScalar v;
10206214f412SHong Zhang       if (cols[j] == row) continue;
1021*eb1ec7c1SStefano Zampini       v    = A->hermitian ? PetscConj(vals[j]) : vals[j];
1022*eb1ec7c1SStefano Zampini       ierr = MatSetValues(mat_elemental,1,&cols[j],1,&row,&v,ADD_VALUES);CHKERRQ(ierr);
10236214f412SHong Zhang     }
10246214f412SHong Zhang     ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr);
10256214f412SHong Zhang   }
10266214f412SHong Zhang   ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr);
10276214f412SHong Zhang   ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10286214f412SHong Zhang   ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10296214f412SHong Zhang 
1030511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
103128be2f97SBarry Smith     ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr);
10326214f412SHong Zhang   } else {
10336214f412SHong Zhang     *newmat = mat_elemental;
10346214f412SHong Zhang   }
10356214f412SHong Zhang   PetscFunctionReturn(0);
10366214f412SHong Zhang }
10376214f412SHong Zhang 
1038cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
10396214f412SHong Zhang {
10406214f412SHong Zhang   Mat               mat_elemental;
10416214f412SHong Zhang   PetscErrorCode    ierr;
10426214f412SHong Zhang   PetscInt          M=A->rmap->N,N=A->cmap->N,row,ncols,j,rstart=A->rmap->rstart,rend=A->rmap->rend;
10436214f412SHong Zhang   const PetscInt    *cols;
10446214f412SHong Zhang   const PetscScalar *vals;
10456214f412SHong Zhang 
10466214f412SHong Zhang   PetscFunctionBegin;
1047bdeae278SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
1048bdeae278SStefano Zampini     mat_elemental = *newmat;
1049bdeae278SStefano Zampini     ierr = MatZeroEntries(mat_elemental);CHKERRQ(ierr);
1050bdeae278SStefano Zampini   } else {
10516214f412SHong Zhang     ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr);
10526214f412SHong Zhang     ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
10536214f412SHong Zhang     ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr);
10546214f412SHong Zhang     ierr = MatSetUp(mat_elemental);CHKERRQ(ierr);
1055bdeae278SStefano Zampini   }
10566214f412SHong Zhang   ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr);
10576214f412SHong Zhang   for (row=rstart; row<rend; row++) {
10586214f412SHong Zhang     ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr);
1059bdeae278SStefano Zampini     /* PETSc-Elemental interface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
10606214f412SHong Zhang     ierr = MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
10616214f412SHong Zhang     for (j=0; j<ncols; j++) { /* lower triangular part */
1062*eb1ec7c1SStefano Zampini       PetscScalar v;
10636214f412SHong Zhang       if (cols[j] == row) continue;
1064*eb1ec7c1SStefano Zampini       v    = A->hermitian ? PetscConj(vals[j]) : vals[j];
1065*eb1ec7c1SStefano Zampini       ierr = MatSetValues(mat_elemental,1,&cols[j],1,&row,&v,ADD_VALUES);CHKERRQ(ierr);
10666214f412SHong Zhang     }
10676214f412SHong Zhang     ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr);
10686214f412SHong Zhang   }
10696214f412SHong Zhang   ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr);
10706214f412SHong Zhang   ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10716214f412SHong Zhang   ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10726214f412SHong Zhang 
1073511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
107428be2f97SBarry Smith     ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr);
10756214f412SHong Zhang   } else {
10766214f412SHong Zhang     *newmat = mat_elemental;
10776214f412SHong Zhang   }
10786214f412SHong Zhang   PetscFunctionReturn(0);
10796214f412SHong Zhang }
10806214f412SHong Zhang 
1081db31f6deSJed Brown static PetscErrorCode MatDestroy_Elemental(Mat A)
1082db31f6deSJed Brown {
1083db31f6deSJed Brown   Mat_Elemental      *a = (Mat_Elemental*)A->data;
1084db31f6deSJed Brown   PetscErrorCode     ierr;
10855e9f5b67SHong Zhang   Mat_Elemental_Grid *commgrid;
10865e9f5b67SHong Zhang   PetscBool          flg;
10875e9f5b67SHong Zhang   MPI_Comm           icomm;
1088db31f6deSJed Brown 
1089db31f6deSJed Brown   PetscFunctionBegin;
1090db31f6deSJed Brown   delete a->emat;
1091ea394475SSatish Balay   delete a->P;
1092ea394475SSatish Balay   delete a->Q;
10935e9f5b67SHong Zhang 
1094e8146892SSatish Balay   El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A));
10950c18141cSBarry Smith   ierr = PetscCommDuplicate(cxxcomm.comm,&icomm,NULL);CHKERRQ(ierr);
109647435625SJed Brown   ierr = MPI_Comm_get_attr(icomm,Petsc_Elemental_keyval,(void**)&commgrid,(int*)&flg);CHKERRQ(ierr);
10975e9f5b67SHong Zhang   if (--commgrid->grid_refct == 0) {
10985e9f5b67SHong Zhang     delete commgrid->grid;
10995e9f5b67SHong Zhang     ierr = PetscFree(commgrid);CHKERRQ(ierr);
110047435625SJed Brown     ierr = MPI_Comm_free_keyval(&Petsc_Elemental_keyval);CHKERRQ(ierr);
11015e9f5b67SHong Zhang   }
11025e9f5b67SHong Zhang   ierr = PetscCommDestroy(&icomm);CHKERRQ(ierr);
1103bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatGetOwnershipIS_C",NULL);CHKERRQ(ierr);
11043ca39a21SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);CHKERRQ(ierr);
1105db31f6deSJed Brown   ierr = PetscFree(A->data);CHKERRQ(ierr);
1106db31f6deSJed Brown   PetscFunctionReturn(0);
1107db31f6deSJed Brown }
1108db31f6deSJed Brown 
1109db31f6deSJed Brown PetscErrorCode MatSetUp_Elemental(Mat A)
1110db31f6deSJed Brown {
1111db31f6deSJed Brown   Mat_Elemental  *a = (Mat_Elemental*)A->data;
1112db31f6deSJed Brown   PetscErrorCode ierr;
1113f19600a0SHong Zhang   MPI_Comm       comm;
1114db31f6deSJed Brown   PetscMPIInt    rsize,csize;
1115f19600a0SHong Zhang   PetscInt       n;
1116db31f6deSJed Brown 
1117db31f6deSJed Brown   PetscFunctionBegin;
1118db31f6deSJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1119db31f6deSJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1120db31f6deSJed Brown 
1121f19600a0SHong Zhang   /* Check if local row and clomun sizes are equally distributed.
1122f19600a0SHong Zhang      Jed: Elemental uses "element" cyclic ordering so the sizes need to match that
1123f19600a0SHong Zhang      exactly.  The strategy in MatElemental is for PETSc to implicitly permute to block ordering (like would be returned by
1124f19600a0SHong Zhang      PetscSplitOwnership(comm,&n,&N)), at which point Elemental matrices can act on PETSc vectors without redistributing the vectors. */
1125f19600a0SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1126f19600a0SHong Zhang   n = PETSC_DECIDE;
1127f19600a0SHong Zhang   ierr = PetscSplitOwnership(comm,&n,&A->rmap->N);CHKERRQ(ierr);
1128f19600a0SHong Zhang   if (n != A->rmap->n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local row size %D of ELEMENTAL matrix must be equally distributed",A->rmap->n);
1129f19600a0SHong Zhang 
1130f19600a0SHong Zhang   n = PETSC_DECIDE;
1131f19600a0SHong Zhang   ierr = PetscSplitOwnership(comm,&n,&A->cmap->N);CHKERRQ(ierr);
1132f19600a0SHong Zhang   if (n != A->cmap->n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local column size %D of ELEMENTAL matrix must be equally distributed",A->cmap->n);
1133f19600a0SHong Zhang 
1134efb79153SJack Poulson   a->emat->Resize(A->rmap->N,A->cmap->N);CHKERRQ(ierr);
1135e8146892SSatish Balay   El::Zero(*a->emat);
1136db31f6deSJed Brown 
1137db31f6deSJed Brown   ierr = MPI_Comm_size(A->rmap->comm,&rsize);CHKERRQ(ierr);
1138db31f6deSJed Brown   ierr = MPI_Comm_size(A->cmap->comm,&csize);CHKERRQ(ierr);
1139ce94432eSBarry Smith   if (csize != rsize) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Cannot use row and column communicators of different sizes");
1140db31f6deSJed Brown   a->commsize = rsize;
1141db31f6deSJed Brown   a->mr[0] = A->rmap->N % rsize; if (!a->mr[0]) a->mr[0] = rsize;
1142db31f6deSJed Brown   a->mr[1] = A->cmap->N % csize; if (!a->mr[1]) a->mr[1] = csize;
1143db31f6deSJed Brown   a->m[0]  = A->rmap->N / rsize + (a->mr[0] != rsize);
1144db31f6deSJed Brown   a->m[1]  = A->cmap->N / csize + (a->mr[1] != csize);
1145db31f6deSJed Brown   PetscFunctionReturn(0);
1146db31f6deSJed Brown }
1147db31f6deSJed Brown 
1148aae2c449SHong Zhang PetscErrorCode MatAssemblyBegin_Elemental(Mat A, MatAssemblyType type)
1149aae2c449SHong Zhang {
1150aae2c449SHong Zhang   Mat_Elemental  *a = (Mat_Elemental*)A->data;
1151aae2c449SHong Zhang 
1152aae2c449SHong Zhang   PetscFunctionBegin;
1153fbbcb730SJack Poulson   /* printf("Calling ProcessQueues\n"); */
1154fbbcb730SJack Poulson   a->emat->ProcessQueues();
1155fbbcb730SJack Poulson   /* printf("Finished ProcessQueues\n"); */
1156aae2c449SHong Zhang   PetscFunctionReturn(0);
1157aae2c449SHong Zhang }
1158aae2c449SHong Zhang 
1159aae2c449SHong Zhang PetscErrorCode MatAssemblyEnd_Elemental(Mat A, MatAssemblyType type)
1160aae2c449SHong Zhang {
1161aae2c449SHong Zhang   PetscFunctionBegin;
1162aae2c449SHong Zhang   /* Currently does nothing */
1163aae2c449SHong Zhang   PetscFunctionReturn(0);
1164aae2c449SHong Zhang }
1165aae2c449SHong Zhang 
11667b30e0f1SHong Zhang PetscErrorCode MatLoad_Elemental(Mat newMat, PetscViewer viewer)
11677b30e0f1SHong Zhang {
11687b30e0f1SHong Zhang   PetscErrorCode ierr;
11697b30e0f1SHong Zhang   Mat            Adense,Ae;
11707b30e0f1SHong Zhang   MPI_Comm       comm;
11717b30e0f1SHong Zhang 
11727b30e0f1SHong Zhang   PetscFunctionBegin;
11737b30e0f1SHong Zhang   ierr = PetscObjectGetComm((PetscObject)newMat,&comm);CHKERRQ(ierr);
11747b30e0f1SHong Zhang   ierr = MatCreate(comm,&Adense);CHKERRQ(ierr);
11757b30e0f1SHong Zhang   ierr = MatSetType(Adense,MATDENSE);CHKERRQ(ierr);
11767b30e0f1SHong Zhang   ierr = MatLoad(Adense,viewer);CHKERRQ(ierr);
11777b30e0f1SHong Zhang   ierr = MatConvert(Adense, MATELEMENTAL, MAT_INITIAL_MATRIX,&Ae);CHKERRQ(ierr);
11787b30e0f1SHong Zhang   ierr = MatDestroy(&Adense);CHKERRQ(ierr);
117928be2f97SBarry Smith   ierr = MatHeaderReplace(newMat,&Ae);CHKERRQ(ierr);
11807b30e0f1SHong Zhang   PetscFunctionReturn(0);
11817b30e0f1SHong Zhang }
11827b30e0f1SHong Zhang 
118340d92e34SHong Zhang /* -------------------------------------------------------------------*/
118440d92e34SHong Zhang static struct _MatOps MatOps_Values = {
118540d92e34SHong Zhang        MatSetValues_Elemental,
118640d92e34SHong Zhang        0,
118740d92e34SHong Zhang        0,
118840d92e34SHong Zhang        MatMult_Elemental,
118940d92e34SHong Zhang /* 4*/ MatMultAdd_Elemental,
11909426833fSXuan Zhou        MatMultTranspose_Elemental,
1191e883f9d5SXuan Zhou        MatMultTransposeAdd_Elemental,
119240d92e34SHong Zhang        MatSolve_Elemental,
1193df311e6cSXuan Zhou        MatSolveAdd_Elemental,
11942ef15aa2SHong Zhang        0,
11952ef15aa2SHong Zhang /*10*/ 0,
119640d92e34SHong Zhang        MatLUFactor_Elemental,
119740d92e34SHong Zhang        MatCholeskyFactor_Elemental,
119840d92e34SHong Zhang        0,
119940d92e34SHong Zhang        MatTranspose_Elemental,
120040d92e34SHong Zhang /*15*/ MatGetInfo_Elemental,
120140d92e34SHong Zhang        0,
120261119200SXuan Zhou        MatGetDiagonal_Elemental,
1203ade3cc5eSXuan Zhou        MatDiagonalScale_Elemental,
120440d92e34SHong Zhang        MatNorm_Elemental,
120540d92e34SHong Zhang /*20*/ MatAssemblyBegin_Elemental,
120640d92e34SHong Zhang        MatAssemblyEnd_Elemental,
120732d7a744SHong Zhang        MatSetOption_Elemental,
120840d92e34SHong Zhang        MatZeroEntries_Elemental,
120940d92e34SHong Zhang /*24*/ 0,
121040d92e34SHong Zhang        MatLUFactorSymbolic_Elemental,
121140d92e34SHong Zhang        MatLUFactorNumeric_Elemental,
121240d92e34SHong Zhang        MatCholeskyFactorSymbolic_Elemental,
121340d92e34SHong Zhang        MatCholeskyFactorNumeric_Elemental,
121440d92e34SHong Zhang /*29*/ MatSetUp_Elemental,
121540d92e34SHong Zhang        0,
121640d92e34SHong Zhang        0,
121740d92e34SHong Zhang        0,
121840d92e34SHong Zhang        0,
1219df311e6cSXuan Zhou /*34*/ MatDuplicate_Elemental,
122040d92e34SHong Zhang        0,
122140d92e34SHong Zhang        0,
122240d92e34SHong Zhang        0,
122340d92e34SHong Zhang        0,
122440d92e34SHong Zhang /*39*/ MatAXPY_Elemental,
122540d92e34SHong Zhang        0,
122640d92e34SHong Zhang        0,
122740d92e34SHong Zhang        0,
122840d92e34SHong Zhang        MatCopy_Elemental,
122940d92e34SHong Zhang /*44*/ 0,
123040d92e34SHong Zhang        MatScale_Elemental,
12317d68702bSBarry Smith        MatShift_Basic,
123240d92e34SHong Zhang        0,
123340d92e34SHong Zhang        0,
123440d92e34SHong Zhang /*49*/ 0,
123540d92e34SHong Zhang        0,
123640d92e34SHong Zhang        0,
123740d92e34SHong Zhang        0,
123840d92e34SHong Zhang        0,
123940d92e34SHong Zhang /*54*/ 0,
124040d92e34SHong Zhang        0,
124140d92e34SHong Zhang        0,
124240d92e34SHong Zhang        0,
124340d92e34SHong Zhang        0,
124440d92e34SHong Zhang /*59*/ 0,
124540d92e34SHong Zhang        MatDestroy_Elemental,
124640d92e34SHong Zhang        MatView_Elemental,
124740d92e34SHong Zhang        0,
124840d92e34SHong Zhang        0,
124940d92e34SHong Zhang /*64*/ 0,
125040d92e34SHong Zhang        0,
125140d92e34SHong Zhang        0,
125240d92e34SHong Zhang        0,
125340d92e34SHong Zhang        0,
125440d92e34SHong Zhang /*69*/ 0,
125540d92e34SHong Zhang        0,
12562ef0cf24SXuan Zhou        MatConvert_Elemental_Dense,
125740d92e34SHong Zhang        0,
125840d92e34SHong Zhang        0,
125940d92e34SHong Zhang /*74*/ 0,
126040d92e34SHong Zhang        0,
126140d92e34SHong Zhang        0,
126240d92e34SHong Zhang        0,
126340d92e34SHong Zhang        0,
126440d92e34SHong Zhang /*79*/ 0,
126540d92e34SHong Zhang        0,
126640d92e34SHong Zhang        0,
126740d92e34SHong Zhang        0,
12687b30e0f1SHong Zhang        MatLoad_Elemental,
126940d92e34SHong Zhang /*84*/ 0,
127040d92e34SHong Zhang        0,
127140d92e34SHong Zhang        0,
127240d92e34SHong Zhang        0,
127340d92e34SHong Zhang        0,
127440d92e34SHong Zhang /*89*/ MatMatMult_Elemental,
127540d92e34SHong Zhang        MatMatMultSymbolic_Elemental,
127640d92e34SHong Zhang        MatMatMultNumeric_Elemental,
127740d92e34SHong Zhang        0,
127840d92e34SHong Zhang        0,
127940d92e34SHong Zhang /*94*/ 0,
1280df311e6cSXuan Zhou        MatMatTransposeMult_Elemental,
1281df311e6cSXuan Zhou        MatMatTransposeMultSymbolic_Elemental,
1282df311e6cSXuan Zhou        MatMatTransposeMultNumeric_Elemental,
128340d92e34SHong Zhang        0,
128440d92e34SHong Zhang /*99*/ 0,
128540d92e34SHong Zhang        0,
128640d92e34SHong Zhang        0,
1287dfcb0403SXuan Zhou        MatConjugate_Elemental,
128840d92e34SHong Zhang        0,
128940d92e34SHong Zhang /*104*/0,
129040d92e34SHong Zhang        0,
129140d92e34SHong Zhang        0,
129240d92e34SHong Zhang        0,
129340d92e34SHong Zhang        0,
129440d92e34SHong Zhang /*109*/MatMatSolve_Elemental,
129540d92e34SHong Zhang        0,
129640d92e34SHong Zhang        0,
129740d92e34SHong Zhang        0,
129834fa46e1SFrancesco Ballarin        MatMissingDiagonal_Elemental,
129940d92e34SHong Zhang /*114*/0,
130040d92e34SHong Zhang        0,
130140d92e34SHong Zhang        0,
130240d92e34SHong Zhang        0,
130340d92e34SHong Zhang        0,
130440d92e34SHong Zhang /*119*/0,
13054a29722dSXuan Zhou        MatHermitianTranspose_Elemental,
130640d92e34SHong Zhang        0,
130740d92e34SHong Zhang        0,
130840d92e34SHong Zhang        0,
130940d92e34SHong Zhang /*124*/0,
131040d92e34SHong Zhang        0,
131140d92e34SHong Zhang        0,
131240d92e34SHong Zhang        0,
131340d92e34SHong Zhang        0,
131440d92e34SHong Zhang /*129*/0,
131540d92e34SHong Zhang        0,
131640d92e34SHong Zhang        0,
131740d92e34SHong Zhang        0,
131840d92e34SHong Zhang        0,
131940d92e34SHong Zhang /*134*/0,
132040d92e34SHong Zhang        0,
132140d92e34SHong Zhang        0,
132240d92e34SHong Zhang        0,
132340d92e34SHong Zhang        0
132440d92e34SHong Zhang };
132540d92e34SHong Zhang 
1326ed36708cSHong Zhang /*MC
1327ed36708cSHong Zhang    MATELEMENTAL = "elemental" - A matrix type for dense matrices using the Elemental package
1328ed36708cSHong Zhang 
1329c2b89b5dSBarry Smith   Use ./configure --download-elemental to install PETSc to use Elemental
1330c2b89b5dSBarry Smith 
13313ca39a21SBarry Smith   Use -pc_type lu -pc_factor_mat_solver_type elemental to use this direct solver
1332c2b89b5dSBarry Smith 
1333ed36708cSHong Zhang    Options Database Keys:
13345cc86fc1SJed Brown + -mat_type elemental - sets the matrix type to "elemental" during a call to MatSetFromOptions()
13355cc86fc1SJed Brown - -mat_elemental_grid_height - sets Grid Height for 2D cyclic ordering of internal matrix
1336ed36708cSHong Zhang 
1337ed36708cSHong Zhang   Level: beginner
1338ed36708cSHong Zhang 
13395cb544a0SHong Zhang .seealso: MATDENSE
1340ed36708cSHong Zhang M*/
13414a29722dSXuan Zhou 
13428cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Elemental(Mat A)
1343db31f6deSJed Brown {
1344db31f6deSJed Brown   Mat_Elemental      *a;
1345db31f6deSJed Brown   PetscErrorCode     ierr;
13465682a260SJack Poulson   PetscBool          flg,flg1;
13475e9f5b67SHong Zhang   Mat_Elemental_Grid *commgrid;
13485e9f5b67SHong Zhang   MPI_Comm           icomm;
13495682a260SJack Poulson   PetscInt           optv1;
1350db31f6deSJed Brown 
1351db31f6deSJed Brown   PetscFunctionBegin;
1352607a6623SBarry Smith   ierr = PetscElementalInitializePackage();CHKERRQ(ierr);
135340d92e34SHong Zhang   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
135440d92e34SHong Zhang   A->insertmode = NOT_SET_VALUES;
1355db31f6deSJed Brown 
1356b00a9115SJed Brown   ierr = PetscNewLog(A,&a);CHKERRQ(ierr);
1357db31f6deSJed Brown   A->data = (void*)a;
1358db31f6deSJed Brown 
1359db31f6deSJed Brown   /* Set up the elemental matrix */
1360e8146892SSatish Balay   El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A));
13615e9f5b67SHong Zhang 
13625e9f5b67SHong Zhang   /* Grid needs to be shared between multiple Mats on the same communicator, implement by attribute caching on the MPI_Comm */
13635e9f5b67SHong Zhang   if (Petsc_Elemental_keyval == MPI_KEYVAL_INVALID) {
136412801b39SBarry Smith     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Elemental_keyval,(void*)0);CHKERRQ(ierr);
13655e9f5b67SHong Zhang   }
13660c18141cSBarry Smith   ierr = PetscCommDuplicate(cxxcomm.comm,&icomm,NULL);CHKERRQ(ierr);
136747435625SJed Brown   ierr = MPI_Comm_get_attr(icomm,Petsc_Elemental_keyval,(void**)&commgrid,(int*)&flg);CHKERRQ(ierr);
13685e9f5b67SHong Zhang   if (!flg) {
1369b00a9115SJed Brown     ierr = PetscNewLog(A,&commgrid);CHKERRQ(ierr);
13705cb544a0SHong Zhang 
1371ce94432eSBarry Smith     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"Elemental Options","Mat");CHKERRQ(ierr);
1372e8146892SSatish Balay     /* displayed default grid sizes (CommSize,1) are set by us arbitrarily until El::Grid() is called */
1373e8146892SSatish Balay     ierr = PetscOptionsInt("-mat_elemental_grid_height","Grid Height","None",El::mpi::Size(cxxcomm),&optv1,&flg1);CHKERRQ(ierr);
13745682a260SJack Poulson     if (flg1) {
1375e8146892SSatish Balay       if (El::mpi::Size(cxxcomm) % optv1 != 0) {
1376e8146892SSatish Balay         SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Grid Height %D must evenly divide CommSize %D",optv1,(PetscInt)El::mpi::Size(cxxcomm));
1377ed667823SXuan Zhou       }
1378e8146892SSatish Balay       commgrid->grid = new El::Grid(cxxcomm,optv1); /* use user-provided grid height */
13792ef0cf24SXuan Zhou     } else {
1380e8146892SSatish Balay       commgrid->grid = new El::Grid(cxxcomm); /* use Elemental default grid sizes */
1381da0640a4SHong Zhang       /* printf("new commgrid->grid = %p\n",commgrid->grid);  -- memory leak revealed by valgrind? */
1382ed667823SXuan Zhou     }
13835e9f5b67SHong Zhang     commgrid->grid_refct = 1;
138447435625SJed Brown     ierr = MPI_Comm_set_attr(icomm,Petsc_Elemental_keyval,(void*)commgrid);CHKERRQ(ierr);
1385ea394475SSatish Balay 
1386ea394475SSatish Balay     a->pivoting    = 1;
1387ea394475SSatish Balay     ierr = PetscOptionsInt("-mat_elemental_pivoting","Pivoting","None",a->pivoting,&a->pivoting,NULL);CHKERRQ(ierr);
1388ea394475SSatish Balay 
13892a808120SBarry Smith     ierr = PetscOptionsEnd();CHKERRQ(ierr);
13905e9f5b67SHong Zhang   } else {
13915e9f5b67SHong Zhang     commgrid->grid_refct++;
13925e9f5b67SHong Zhang   }
13935e9f5b67SHong Zhang   ierr = PetscCommDestroy(&icomm);CHKERRQ(ierr);
13945e9f5b67SHong Zhang   a->grid        = commgrid->grid;
1395e8146892SSatish Balay   a->emat        = new El::DistMatrix<PetscElemScalar>(*a->grid);
139632d7a744SHong Zhang   a->roworiented = PETSC_TRUE;
1397db31f6deSJed Brown 
1398bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatGetOwnershipIS_C",MatGetOwnershipIS_Elemental);CHKERRQ(ierr);
1399db31f6deSJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATELEMENTAL);CHKERRQ(ierr);
1400db31f6deSJed Brown   PetscFunctionReturn(0);
1401db31f6deSJed Brown }
1402