1 #include <../src/mat/impls/elemental/matelemimpl.h> /*I "petscmat.h" I*/ 2 3 /* 4 The variable Petsc_Elemental_keyval is used to indicate an MPI attribute that 5 is attached to a communicator, in this case the attribute is a Mat_Elemental_Grid 6 */ 7 static PetscMPIInt Petsc_Elemental_keyval = MPI_KEYVAL_INVALID; 8 9 #undef __FUNCT__ 10 #define __FUNCT__ "PetscElementalInitializePackage" 11 /*@C 12 PetscElementalInitializePackage - Initialize Elemental package 13 14 Logically Collective 15 16 Level: developer 17 18 .seealso: MATELEMENTAL, PetscElementalFinalizePackage() 19 @*/ 20 PetscErrorCode PetscElementalInitializePackage(void) 21 { 22 PetscErrorCode ierr; 23 24 PetscFunctionBegin; 25 if (El::Initialized()) PetscFunctionReturn(0); 26 { /* We have already initialized MPI, so this song and dance is just to pass these variables (which won't be used by Elemental) through the interface that needs references */ 27 int zero = 0; 28 char **nothing = 0; 29 El::Initialize(zero,nothing); /* called by the 1st call of MatCreate_Elemental */ 30 } 31 ierr = PetscRegisterFinalize(PetscElementalFinalizePackage);CHKERRQ(ierr); 32 PetscFunctionReturn(0); 33 } 34 35 #undef __FUNCT__ 36 #define __FUNCT__ "PetscElementalFinalizePackage" 37 /*@C 38 PetscElementalFinalizePackage - Finalize Elemental package 39 40 Logically Collective 41 42 Level: developer 43 44 .seealso: MATELEMENTAL, PetscElementalInitializePackage() 45 @*/ 46 PetscErrorCode PetscElementalFinalizePackage(void) 47 { 48 PetscFunctionBegin; 49 El::Finalize(); /* called by PetscFinalize() */ 50 PetscFunctionReturn(0); 51 } 52 53 #undef __FUNCT__ 54 #define __FUNCT__ "MatView_Elemental" 55 static PetscErrorCode MatView_Elemental(Mat A,PetscViewer viewer) 56 { 57 PetscErrorCode ierr; 58 Mat_Elemental *a = (Mat_Elemental*)A->data; 59 PetscBool iascii; 60 61 PetscFunctionBegin; 62 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 63 if (iascii) { 64 PetscViewerFormat format; 65 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 66 if (format == PETSC_VIEWER_ASCII_INFO) { 67 /* call elemental viewing function */ 68 ierr = PetscViewerASCIIPrintf(viewer,"Elemental run parameters:\n");CHKERRQ(ierr); 69 ierr = PetscViewerASCIIPrintf(viewer," allocated entries=%d\n",(*a->emat).AllocatedMemory());CHKERRQ(ierr); 70 ierr = PetscViewerASCIIPrintf(viewer," grid height=%d, grid width=%d\n",(*a->emat).Grid().Height(),(*a->emat).Grid().Width());CHKERRQ(ierr); 71 if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 72 /* call elemental viewing function */ 73 ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer),"test matview_elemental 2\n");CHKERRQ(ierr); 74 } 75 76 } else if (format == PETSC_VIEWER_DEFAULT) { 77 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 78 El::Print( *a->emat, "Elemental matrix (cyclic ordering)" ); 79 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 80 if (A->factortype == MAT_FACTOR_NONE){ 81 Mat Adense; 82 ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer),"Elemental matrix (explicit ordering)\n");CHKERRQ(ierr); 83 ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense);CHKERRQ(ierr); 84 ierr = MatView(Adense,viewer);CHKERRQ(ierr); 85 ierr = MatDestroy(&Adense);CHKERRQ(ierr); 86 } 87 } else SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Format"); 88 } else { 89 /* convert to dense format and call MatView() */ 90 Mat Adense; 91 ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer),"Elemental matrix (explicit ordering)\n");CHKERRQ(ierr); 92 ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense);CHKERRQ(ierr); 93 ierr = MatView(Adense,viewer);CHKERRQ(ierr); 94 ierr = MatDestroy(&Adense);CHKERRQ(ierr); 95 } 96 PetscFunctionReturn(0); 97 } 98 99 #undef __FUNCT__ 100 #define __FUNCT__ "MatGetInfo_Elemental" 101 static PetscErrorCode MatGetInfo_Elemental(Mat A,MatInfoType flag,MatInfo *info) 102 { 103 Mat_Elemental *a = (Mat_Elemental*)A->data; 104 PetscMPIInt rank; 105 106 PetscFunctionBegin; 107 MPI_Comm_rank(PetscObjectComm((PetscObject)A),&rank); 108 109 /* if (!rank) printf(" .........MatGetInfo_Elemental ...\n"); */ 110 info->block_size = 1.0; 111 112 if (flag == MAT_LOCAL) { 113 info->nz_allocated = (double)(*a->emat).AllocatedMemory(); /* locally allocated */ 114 info->nz_used = info->nz_allocated; 115 } else if (flag == MAT_GLOBAL_MAX) { 116 //ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)matin));CHKERRQ(ierr); 117 /* see MatGetInfo_MPIAIJ() for getting global info->nz_allocated! */ 118 //SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_GLOBAL_MAX not written yet"); 119 } else if (flag == MAT_GLOBAL_SUM) { 120 //SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_GLOBAL_SUM not written yet"); 121 info->nz_allocated = (double)(*a->emat).AllocatedMemory(); /* locally allocated */ 122 info->nz_used = info->nz_allocated; /* assume Elemental does accurate allocation */ 123 //ierr = MPI_Allreduce(isend,irecv,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 124 //PetscPrintf(PETSC_COMM_SELF," ... [%d] locally allocated %g\n",rank,info->nz_allocated); 125 } 126 127 info->nz_unneeded = 0.0; 128 info->assemblies = (double)A->num_ass; 129 info->mallocs = 0; 130 info->memory = ((PetscObject)A)->mem; 131 info->fill_ratio_given = 0; /* determined by Elemental */ 132 info->fill_ratio_needed = 0; 133 info->factor_mallocs = 0; 134 PetscFunctionReturn(0); 135 } 136 137 #undef __FUNCT__ 138 #define __FUNCT__ "MatSetValues_Elemental" 139 static PetscErrorCode MatSetValues_Elemental(Mat A,PetscInt nr,const PetscInt *rows,PetscInt nc,const PetscInt *cols,const PetscScalar *vals,InsertMode imode) 140 { 141 PetscErrorCode ierr; 142 Mat_Elemental *a = (Mat_Elemental*)A->data; 143 PetscMPIInt rank; 144 PetscInt i,j,rrank,ridx,crank,cidx; 145 146 PetscFunctionBegin; 147 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A),&rank);CHKERRQ(ierr); 148 149 const El::Grid &grid = a->emat->Grid(); 150 for (i=0; i<nr; i++) { 151 PetscInt erow,ecol,elrow,elcol; 152 if (rows[i] < 0) continue; 153 P2RO(A,0,rows[i],&rrank,&ridx); 154 RO2E(A,0,rrank,ridx,&erow); 155 if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Incorrect row translation"); 156 for (j=0; j<nc; j++) { 157 if (cols[j] < 0) continue; 158 P2RO(A,1,cols[j],&crank,&cidx); 159 RO2E(A,1,crank,cidx,&ecol); 160 if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Incorrect col translation"); 161 if (erow % grid.MCSize() != grid.MCRank() || ecol % grid.MRSize() != grid.MRRank()){ /* off-proc entry */ 162 if (imode != ADD_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only ADD_VALUES to off-processor entry is supported"); 163 /* PetscPrintf(PETSC_COMM_SELF,"[%D] add off-proc entry (%D,%D, %g) (%D %D)\n",rank,rows[i],cols[j],*(vals+i*nc),erow,ecol); */ 164 a->esubmat->Set(0,0, (PetscElemScalar)vals[i*nc+j]); 165 a->interface->Axpy(1.0,*(a->esubmat),erow,ecol); 166 continue; 167 } 168 elrow = erow / grid.MCSize(); 169 elcol = ecol / grid.MRSize(); 170 switch (imode) { 171 case INSERT_VALUES: a->emat->SetLocal(elrow,elcol,(PetscElemScalar)vals[i*nc+j]); break; 172 case ADD_VALUES: a->emat->UpdateLocal(elrow,elcol,(PetscElemScalar)vals[i*nc+j]); break; 173 default: SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for InsertMode %d",(int)imode); 174 } 175 } 176 } 177 PetscFunctionReturn(0); 178 } 179 180 #undef __FUNCT__ 181 #define __FUNCT__ "MatMult_Elemental" 182 static PetscErrorCode MatMult_Elemental(Mat A,Vec X,Vec Y) 183 { 184 Mat_Elemental *a = (Mat_Elemental*)A->data; 185 PetscErrorCode ierr; 186 const PetscElemScalar *x; 187 PetscElemScalar *y; 188 PetscElemScalar one = 1,zero = 0; 189 190 PetscFunctionBegin; 191 ierr = VecGetArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr); 192 ierr = VecGetArray(Y,(PetscScalar **)&y);CHKERRQ(ierr); 193 { /* Scoping so that constructor is called before pointer is returned */ 194 El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ye; 195 xe.LockedAttach(A->cmap->N,1,*a->grid,0,0,x,A->cmap->n); 196 ye.Attach(A->rmap->N,1,*a->grid,0,0,y,A->rmap->n); 197 El::Gemv(El::NORMAL,one,*a->emat,xe,zero,ye); 198 } 199 ierr = VecRestoreArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr); 200 ierr = VecRestoreArray(Y,(PetscScalar **)&y);CHKERRQ(ierr); 201 PetscFunctionReturn(0); 202 } 203 204 #undef __FUNCT__ 205 #define __FUNCT__ "MatMultTranspose_Elemental" 206 static PetscErrorCode MatMultTranspose_Elemental(Mat A,Vec X,Vec Y) 207 { 208 Mat_Elemental *a = (Mat_Elemental*)A->data; 209 PetscErrorCode ierr; 210 const PetscElemScalar *x; 211 PetscElemScalar *y; 212 PetscElemScalar one = 1,zero = 0; 213 214 PetscFunctionBegin; 215 ierr = VecGetArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr); 216 ierr = VecGetArray(Y,(PetscScalar **)&y);CHKERRQ(ierr); 217 { /* Scoping so that constructor is called before pointer is returned */ 218 El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ye; 219 xe.LockedAttach(A->rmap->N,1,*a->grid,0,0,x,A->rmap->n); 220 ye.Attach(A->cmap->N,1,*a->grid,0,0,y,A->cmap->n); 221 El::Gemv(El::TRANSPOSE,one,*a->emat,xe,zero,ye); 222 } 223 ierr = VecRestoreArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr); 224 ierr = VecRestoreArray(Y,(PetscScalar **)&y);CHKERRQ(ierr); 225 PetscFunctionReturn(0); 226 } 227 228 #undef __FUNCT__ 229 #define __FUNCT__ "MatMultAdd_Elemental" 230 static PetscErrorCode MatMultAdd_Elemental(Mat A,Vec X,Vec Y,Vec Z) 231 { 232 Mat_Elemental *a = (Mat_Elemental*)A->data; 233 PetscErrorCode ierr; 234 const PetscElemScalar *x; 235 PetscElemScalar *z; 236 PetscElemScalar one = 1; 237 238 PetscFunctionBegin; 239 if (Y != Z) {ierr = VecCopy(Y,Z);CHKERRQ(ierr);} 240 ierr = VecGetArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr); 241 ierr = VecGetArray(Z,(PetscScalar **)&z);CHKERRQ(ierr); 242 { /* Scoping so that constructor is called before pointer is returned */ 243 El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ze; 244 xe.LockedAttach(A->cmap->N,1,*a->grid,0,0,x,A->cmap->n); 245 ze.Attach(A->rmap->N,1,*a->grid,0,0,z,A->rmap->n); 246 El::Gemv(El::NORMAL,one,*a->emat,xe,one,ze); 247 } 248 ierr = VecRestoreArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr); 249 ierr = VecRestoreArray(Z,(PetscScalar **)&z);CHKERRQ(ierr); 250 PetscFunctionReturn(0); 251 } 252 253 #undef __FUNCT__ 254 #define __FUNCT__ "MatMultTransposeAdd_Elemental" 255 static PetscErrorCode MatMultTransposeAdd_Elemental(Mat A,Vec X,Vec Y,Vec Z) 256 { 257 Mat_Elemental *a = (Mat_Elemental*)A->data; 258 PetscErrorCode ierr; 259 const PetscElemScalar *x; 260 PetscElemScalar *z; 261 PetscElemScalar one = 1; 262 263 PetscFunctionBegin; 264 if (Y != Z) {ierr = VecCopy(Y,Z);CHKERRQ(ierr);} 265 ierr = VecGetArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr); 266 ierr = VecGetArray(Z,(PetscScalar **)&z);CHKERRQ(ierr); 267 { /* Scoping so that constructor is called before pointer is returned */ 268 El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe, ze; 269 xe.LockedAttach(A->rmap->N,1,*a->grid,0,0,x,A->rmap->n); 270 ze.Attach(A->cmap->N,1,*a->grid,0,0,z,A->cmap->n); 271 El::Gemv(El::TRANSPOSE,one,*a->emat,xe,one,ze); 272 } 273 ierr = VecRestoreArrayRead(X,(const PetscScalar **)&x);CHKERRQ(ierr); 274 ierr = VecRestoreArray(Z,(PetscScalar **)&z);CHKERRQ(ierr); 275 PetscFunctionReturn(0); 276 } 277 278 #undef __FUNCT__ 279 #define __FUNCT__ "MatMatMultNumeric_Elemental" 280 static PetscErrorCode MatMatMultNumeric_Elemental(Mat A,Mat B,Mat C) 281 { 282 Mat_Elemental *a = (Mat_Elemental*)A->data; 283 Mat_Elemental *b = (Mat_Elemental*)B->data; 284 Mat_Elemental *c = (Mat_Elemental*)C->data; 285 PetscElemScalar one = 1,zero = 0; 286 287 PetscFunctionBegin; 288 { /* Scoping so that constructor is called before pointer is returned */ 289 El::Gemm(El::NORMAL,El::NORMAL,one,*a->emat,*b->emat,zero,*c->emat); 290 } 291 C->assembled = PETSC_TRUE; 292 PetscFunctionReturn(0); 293 } 294 295 #undef __FUNCT__ 296 #define __FUNCT__ "MatMatMultSymbolic_Elemental" 297 static PetscErrorCode MatMatMultSymbolic_Elemental(Mat A,Mat B,PetscReal fill,Mat *C) 298 { 299 PetscErrorCode ierr; 300 Mat Ce; 301 MPI_Comm comm; 302 303 PetscFunctionBegin; 304 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 305 ierr = MatCreate(comm,&Ce);CHKERRQ(ierr); 306 ierr = MatSetSizes(Ce,A->rmap->n,B->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 307 ierr = MatSetType(Ce,MATELEMENTAL);CHKERRQ(ierr); 308 ierr = MatSetUp(Ce);CHKERRQ(ierr); 309 *C = Ce; 310 PetscFunctionReturn(0); 311 } 312 313 #undef __FUNCT__ 314 #define __FUNCT__ "MatMatMult_Elemental" 315 static PetscErrorCode MatMatMult_Elemental(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 316 { 317 PetscErrorCode ierr; 318 319 PetscFunctionBegin; 320 if (scall == MAT_INITIAL_MATRIX){ 321 ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 322 ierr = MatMatMultSymbolic_Elemental(A,B,1.0,C);CHKERRQ(ierr); 323 ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 324 } 325 ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 326 ierr = MatMatMultNumeric_Elemental(A,B,*C);CHKERRQ(ierr); 327 ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 328 PetscFunctionReturn(0); 329 } 330 331 #undef __FUNCT__ 332 #define __FUNCT__ "MatMatTransposeMultNumeric_Elemental" 333 static PetscErrorCode MatMatTransposeMultNumeric_Elemental(Mat A,Mat B,Mat C) 334 { 335 Mat_Elemental *a = (Mat_Elemental*)A->data; 336 Mat_Elemental *b = (Mat_Elemental*)B->data; 337 Mat_Elemental *c = (Mat_Elemental*)C->data; 338 PetscElemScalar one = 1,zero = 0; 339 340 PetscFunctionBegin; 341 { /* Scoping so that constructor is called before pointer is returned */ 342 El::Gemm(El::NORMAL,El::TRANSPOSE,one,*a->emat,*b->emat,zero,*c->emat); 343 } 344 C->assembled = PETSC_TRUE; 345 PetscFunctionReturn(0); 346 } 347 348 #undef __FUNCT__ 349 #define __FUNCT__ "MatMatTransposeMultSymbolic_Elemental" 350 static PetscErrorCode MatMatTransposeMultSymbolic_Elemental(Mat A,Mat B,PetscReal fill,Mat *C) 351 { 352 PetscErrorCode ierr; 353 Mat Ce; 354 MPI_Comm comm; 355 356 PetscFunctionBegin; 357 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 358 ierr = MatCreate(comm,&Ce);CHKERRQ(ierr); 359 ierr = MatSetSizes(Ce,A->rmap->n,B->rmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 360 ierr = MatSetType(Ce,MATELEMENTAL);CHKERRQ(ierr); 361 ierr = MatSetUp(Ce);CHKERRQ(ierr); 362 *C = Ce; 363 PetscFunctionReturn(0); 364 } 365 366 #undef __FUNCT__ 367 #define __FUNCT__ "MatMatTransposeMult_Elemental" 368 static PetscErrorCode MatMatTransposeMult_Elemental(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 369 { 370 PetscErrorCode ierr; 371 372 PetscFunctionBegin; 373 if (scall == MAT_INITIAL_MATRIX){ 374 ierr = PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr); 375 ierr = MatMatMultSymbolic_Elemental(A,B,1.0,C);CHKERRQ(ierr); 376 ierr = PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr); 377 } 378 ierr = PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr); 379 ierr = MatMatTransposeMultNumeric_Elemental(A,B,*C);CHKERRQ(ierr); 380 ierr = PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr); 381 PetscFunctionReturn(0); 382 } 383 384 #undef __FUNCT__ 385 #define __FUNCT__ "MatGetDiagonal_Elemental" 386 static PetscErrorCode MatGetDiagonal_Elemental(Mat A,Vec D) 387 { 388 PetscInt i,nrows,ncols,nD,rrank,ridx,crank,cidx; 389 Mat_Elemental *a = (Mat_Elemental*)A->data; 390 PetscErrorCode ierr; 391 PetscElemScalar v; 392 MPI_Comm comm; 393 394 PetscFunctionBegin; 395 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 396 ierr = MatGetSize(A,&nrows,&ncols);CHKERRQ(ierr); 397 nD = nrows>ncols ? ncols : nrows; 398 for (i=0; i<nD; i++) { 399 PetscInt erow,ecol; 400 P2RO(A,0,i,&rrank,&ridx); 401 RO2E(A,0,rrank,ridx,&erow); 402 if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect row translation"); 403 P2RO(A,1,i,&crank,&cidx); 404 RO2E(A,1,crank,cidx,&ecol); 405 if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect col translation"); 406 v = a->emat->Get(erow,ecol); 407 ierr = VecSetValues(D,1,&i,(PetscScalar*)&v,INSERT_VALUES);CHKERRQ(ierr); 408 } 409 ierr = VecAssemblyBegin(D);CHKERRQ(ierr); 410 ierr = VecAssemblyEnd(D);CHKERRQ(ierr); 411 PetscFunctionReturn(0); 412 } 413 414 #undef __FUNCT__ 415 #define __FUNCT__ "MatDiagonalScale_Elemental" 416 static PetscErrorCode MatDiagonalScale_Elemental(Mat X,Vec L,Vec R) 417 { 418 Mat_Elemental *x = (Mat_Elemental*)X->data; 419 const PetscElemScalar *d; 420 PetscErrorCode ierr; 421 422 PetscFunctionBegin; 423 if (R) { 424 ierr = VecGetArrayRead(R,(const PetscScalar **)&d);CHKERRQ(ierr); 425 El::DistMatrix<PetscElemScalar,El::VC,El::STAR> de; 426 de.LockedAttach(X->cmap->N,1,*x->grid,0,0,d,X->cmap->n); 427 El::DiagonalScale(El::RIGHT,El::NORMAL,de,*x->emat); 428 ierr = VecRestoreArrayRead(R,(const PetscScalar **)&d);CHKERRQ(ierr); 429 } 430 if (L) { 431 ierr = VecGetArrayRead(L,(const PetscScalar **)&d);CHKERRQ(ierr); 432 El::DistMatrix<PetscElemScalar,El::VC,El::STAR> de; 433 de.LockedAttach(X->rmap->N,1,*x->grid,0,0,d,X->rmap->n); 434 El::DiagonalScale(El::LEFT,El::NORMAL,de,*x->emat); 435 ierr = VecRestoreArrayRead(L,(const PetscScalar **)&d);CHKERRQ(ierr); 436 } 437 PetscFunctionReturn(0); 438 } 439 440 #undef __FUNCT__ 441 #define __FUNCT__ "MatScale_Elemental" 442 static PetscErrorCode MatScale_Elemental(Mat X,PetscScalar a) 443 { 444 Mat_Elemental *x = (Mat_Elemental*)X->data; 445 446 PetscFunctionBegin; 447 El::Scale((PetscElemScalar)a,*x->emat); 448 PetscFunctionReturn(0); 449 } 450 451 /* 452 MatAXPY - Computes Y = a*X + Y. 453 */ 454 #undef __FUNCT__ 455 #define __FUNCT__ "MatAXPY_Elemental" 456 static PetscErrorCode MatAXPY_Elemental(Mat Y,PetscScalar a,Mat X,MatStructure str) 457 { 458 Mat_Elemental *x = (Mat_Elemental*)X->data; 459 Mat_Elemental *y = (Mat_Elemental*)Y->data; 460 PetscErrorCode ierr; 461 462 PetscFunctionBegin; 463 El::Axpy((PetscElemScalar)a,*x->emat,*y->emat); 464 ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr); 465 PetscFunctionReturn(0); 466 } 467 468 #undef __FUNCT__ 469 #define __FUNCT__ "MatCopy_Elemental" 470 static PetscErrorCode MatCopy_Elemental(Mat A,Mat B,MatStructure str) 471 { 472 Mat_Elemental *a=(Mat_Elemental*)A->data; 473 Mat_Elemental *b=(Mat_Elemental*)B->data; 474 475 PetscFunctionBegin; 476 El::Copy(*a->emat,*b->emat); 477 PetscFunctionReturn(0); 478 } 479 480 #undef __FUNCT__ 481 #define __FUNCT__ "MatDuplicate_Elemental" 482 static PetscErrorCode MatDuplicate_Elemental(Mat A,MatDuplicateOption op,Mat *B) 483 { 484 Mat Be; 485 MPI_Comm comm; 486 Mat_Elemental *a=(Mat_Elemental*)A->data; 487 PetscErrorCode ierr; 488 489 PetscFunctionBegin; 490 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 491 ierr = MatCreate(comm,&Be);CHKERRQ(ierr); 492 ierr = MatSetSizes(Be,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 493 ierr = MatSetType(Be,MATELEMENTAL);CHKERRQ(ierr); 494 ierr = MatSetUp(Be);CHKERRQ(ierr); 495 *B = Be; 496 if (op == MAT_COPY_VALUES) { 497 Mat_Elemental *b=(Mat_Elemental*)Be->data; 498 El::Copy(*a->emat,*b->emat); 499 } 500 Be->assembled = PETSC_TRUE; 501 PetscFunctionReturn(0); 502 } 503 504 #undef __FUNCT__ 505 #define __FUNCT__ "MatTranspose_Elemental" 506 static PetscErrorCode MatTranspose_Elemental(Mat A,MatReuse reuse,Mat *B) 507 { 508 Mat Be = *B; 509 PetscErrorCode ierr; 510 MPI_Comm comm; 511 Mat_Elemental *a = (Mat_Elemental*)A->data, *b; 512 513 PetscFunctionBegin; 514 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 515 /* Only out-of-place supported */ 516 if (reuse == MAT_INITIAL_MATRIX){ 517 ierr = MatCreate(comm,&Be);CHKERRQ(ierr); 518 ierr = MatSetSizes(Be,A->cmap->n,A->rmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 519 ierr = MatSetType(Be,MATELEMENTAL);CHKERRQ(ierr); 520 ierr = MatSetUp(Be);CHKERRQ(ierr); 521 *B = Be; 522 } 523 b = (Mat_Elemental*)Be->data; 524 El::Transpose(*a->emat,*b->emat); 525 Be->assembled = PETSC_TRUE; 526 PetscFunctionReturn(0); 527 } 528 529 #undef __FUNCT__ 530 #define __FUNCT__ "MatConjugate_Elemental" 531 static PetscErrorCode MatConjugate_Elemental(Mat A) 532 { 533 Mat_Elemental *a = (Mat_Elemental*)A->data; 534 535 PetscFunctionBegin; 536 El::Conjugate(*a->emat); 537 PetscFunctionReturn(0); 538 } 539 540 #undef __FUNCT__ 541 #define __FUNCT__ "MatHermitianTranspose_Elemental" 542 static PetscErrorCode MatHermitianTranspose_Elemental(Mat A,MatReuse reuse,Mat *B) 543 { 544 Mat Be = *B; 545 PetscErrorCode ierr; 546 MPI_Comm comm; 547 Mat_Elemental *a = (Mat_Elemental*)A->data, *b; 548 549 PetscFunctionBegin; 550 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 551 /* Only out-of-place supported */ 552 if (reuse == MAT_INITIAL_MATRIX){ 553 ierr = MatCreate(comm,&Be);CHKERRQ(ierr); 554 ierr = MatSetSizes(Be,A->cmap->n,A->rmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 555 ierr = MatSetType(Be,MATELEMENTAL);CHKERRQ(ierr); 556 ierr = MatSetUp(Be);CHKERRQ(ierr); 557 *B = Be; 558 } 559 b = (Mat_Elemental*)Be->data; 560 El::Adjoint(*a->emat,*b->emat); 561 Be->assembled = PETSC_TRUE; 562 PetscFunctionReturn(0); 563 } 564 565 #undef __FUNCT__ 566 #define __FUNCT__ "MatSolve_Elemental" 567 static PetscErrorCode MatSolve_Elemental(Mat A,Vec B,Vec X) 568 { 569 Mat_Elemental *a = (Mat_Elemental*)A->data; 570 PetscErrorCode ierr; 571 PetscElemScalar *x; 572 573 PetscFunctionBegin; 574 ierr = VecCopy(B,X);CHKERRQ(ierr); 575 ierr = VecGetArray(X,(PetscScalar **)&x);CHKERRQ(ierr); 576 El::DistMatrix<PetscElemScalar,El::VC,El::STAR> xe; 577 xe.Attach(A->rmap->N,1,*a->grid,0,0,x,A->rmap->n); 578 El::DistMatrix<PetscElemScalar,El::MC,El::MR> xer(xe); 579 switch (A->factortype) { 580 case MAT_FACTOR_LU: 581 if ((*a->pivot).AllocatedMemory()) { 582 El::lu::SolveAfter(El::NORMAL,*a->emat,*a->pivot,xer); 583 El::Copy(xer,xe); 584 } else { 585 El::lu::SolveAfter(El::NORMAL,*a->emat,xer); 586 El::Copy(xer,xe); 587 } 588 break; 589 case MAT_FACTOR_CHOLESKY: 590 El::cholesky::SolveAfter(El::UPPER,El::NORMAL,*a->emat,xer); 591 El::Copy(xer,xe); 592 break; 593 default: 594 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unfactored Matrix or Unsupported MatFactorType"); 595 break; 596 } 597 ierr = VecRestoreArray(X,(PetscScalar **)&x);CHKERRQ(ierr); 598 PetscFunctionReturn(0); 599 } 600 601 #undef __FUNCT__ 602 #define __FUNCT__ "MatSolveAdd_Elemental" 603 static PetscErrorCode MatSolveAdd_Elemental(Mat A,Vec B,Vec Y,Vec X) 604 { 605 PetscErrorCode ierr; 606 607 PetscFunctionBegin; 608 ierr = MatSolve_Elemental(A,B,X);CHKERRQ(ierr); 609 ierr = VecAXPY(X,1,Y);CHKERRQ(ierr); 610 PetscFunctionReturn(0); 611 } 612 613 #undef __FUNCT__ 614 #define __FUNCT__ "MatMatSolve_Elemental" 615 static PetscErrorCode MatMatSolve_Elemental(Mat A,Mat B,Mat X) 616 { 617 Mat_Elemental *a=(Mat_Elemental*)A->data; 618 Mat_Elemental *b=(Mat_Elemental*)B->data; 619 Mat_Elemental *x=(Mat_Elemental*)X->data; 620 621 PetscFunctionBegin; 622 El::Copy(*b->emat,*x->emat); 623 switch (A->factortype) { 624 case MAT_FACTOR_LU: 625 if ((*a->pivot).AllocatedMemory()) { 626 El::lu::SolveAfter(El::NORMAL,*a->emat,*a->pivot,*x->emat); 627 } else { 628 El::lu::SolveAfter(El::NORMAL,*a->emat,*x->emat); 629 } 630 break; 631 case MAT_FACTOR_CHOLESKY: 632 El::cholesky::SolveAfter(El::UPPER,El::NORMAL,*a->emat,*x->emat); 633 break; 634 default: 635 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unfactored Matrix or Unsupported MatFactorType"); 636 break; 637 } 638 PetscFunctionReturn(0); 639 } 640 641 #undef __FUNCT__ 642 #define __FUNCT__ "MatLUFactor_Elemental" 643 static PetscErrorCode MatLUFactor_Elemental(Mat A,IS row,IS col,const MatFactorInfo *info) 644 { 645 Mat_Elemental *a = (Mat_Elemental*)A->data; 646 647 PetscFunctionBegin; 648 if (info->dtcol){ 649 El::LU(*a->emat,*a->pivot); 650 } else { 651 El::LU(*a->emat); 652 } 653 A->factortype = MAT_FACTOR_LU; 654 A->assembled = PETSC_TRUE; 655 PetscFunctionReturn(0); 656 } 657 658 #undef __FUNCT__ 659 #define __FUNCT__ "MatLUFactorNumeric_Elemental" 660 static PetscErrorCode MatLUFactorNumeric_Elemental(Mat F,Mat A,const MatFactorInfo *info) 661 { 662 PetscErrorCode ierr; 663 664 PetscFunctionBegin; 665 ierr = MatCopy(A,F,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 666 ierr = MatLUFactor_Elemental(F,0,0,info);CHKERRQ(ierr); 667 PetscFunctionReturn(0); 668 } 669 670 #undef __FUNCT__ 671 #define __FUNCT__ "MatLUFactorSymbolic_Elemental" 672 static PetscErrorCode MatLUFactorSymbolic_Elemental(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 673 { 674 PetscFunctionBegin; 675 /* F is create and allocated by MatGetFactor_elemental_petsc(), skip this routine. */ 676 PetscFunctionReturn(0); 677 } 678 679 #undef __FUNCT__ 680 #define __FUNCT__ "MatCholeskyFactor_Elemental" 681 static PetscErrorCode MatCholeskyFactor_Elemental(Mat A,IS perm,const MatFactorInfo *info) 682 { 683 Mat_Elemental *a = (Mat_Elemental*)A->data; 684 El::DistMatrix<PetscElemScalar,El::MC,El::STAR> d; 685 686 PetscFunctionBegin; 687 El::Cholesky(El::UPPER,*a->emat); 688 A->factortype = MAT_FACTOR_CHOLESKY; 689 A->assembled = PETSC_TRUE; 690 PetscFunctionReturn(0); 691 } 692 693 #undef __FUNCT__ 694 #define __FUNCT__ "MatCholeskyFactorNumeric_Elemental" 695 static PetscErrorCode MatCholeskyFactorNumeric_Elemental(Mat F,Mat A,const MatFactorInfo *info) 696 { 697 PetscErrorCode ierr; 698 699 PetscFunctionBegin; 700 ierr = MatCopy(A,F,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 701 ierr = MatCholeskyFactor_Elemental(F,0,info);CHKERRQ(ierr); 702 PetscFunctionReturn(0); 703 } 704 705 #undef __FUNCT__ 706 #define __FUNCT__ "MatCholeskyFactorSymbolic_Elemental" 707 static PetscErrorCode MatCholeskyFactorSymbolic_Elemental(Mat F,Mat A,IS perm,const MatFactorInfo *info) 708 { 709 PetscFunctionBegin; 710 /* F is create and allocated by MatGetFactor_elemental_petsc(), skip this routine. */ 711 PetscFunctionReturn(0); 712 } 713 714 #undef __FUNCT__ 715 #define __FUNCT__ "MatFactorGetSolverPackage_elemental_elemental" 716 PetscErrorCode MatFactorGetSolverPackage_elemental_elemental(Mat A,const MatSolverPackage *type) 717 { 718 PetscFunctionBegin; 719 *type = MATSOLVERELEMENTAL; 720 PetscFunctionReturn(0); 721 } 722 723 #undef __FUNCT__ 724 #define __FUNCT__ "MatGetFactor_elemental_elemental" 725 static PetscErrorCode MatGetFactor_elemental_elemental(Mat A,MatFactorType ftype,Mat *F) 726 { 727 Mat B; 728 PetscErrorCode ierr; 729 730 PetscFunctionBegin; 731 /* Create the factorization matrix */ 732 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 733 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 734 ierr = MatSetType(B,MATELEMENTAL);CHKERRQ(ierr); 735 ierr = MatSetUp(B);CHKERRQ(ierr); 736 B->factortype = ftype; 737 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_elemental_elemental);CHKERRQ(ierr); 738 *F = B; 739 PetscFunctionReturn(0); 740 } 741 742 #undef __FUNCT__ 743 #define __FUNCT__ "MatSolverPackageRegister_Elemental" 744 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_Elemental(void) 745 { 746 PetscErrorCode ierr; 747 748 PetscFunctionBegin; 749 ierr = MatSolverPackageRegister(MATSOLVERELEMENTAL,MATELEMENTAL, MAT_FACTOR_LU,MatGetFactor_elemental_elemental);CHKERRQ(ierr); 750 ierr = MatSolverPackageRegister(MATSOLVERELEMENTAL,MATELEMENTAL, MAT_FACTOR_CHOLESKY,MatGetFactor_elemental_elemental);CHKERRQ(ierr); 751 PetscFunctionReturn(0); 752 } 753 754 #undef __FUNCT__ 755 #define __FUNCT__ "MatNorm_Elemental" 756 static PetscErrorCode MatNorm_Elemental(Mat A,NormType type,PetscReal *nrm) 757 { 758 Mat_Elemental *a=(Mat_Elemental*)A->data; 759 760 PetscFunctionBegin; 761 switch (type){ 762 case NORM_1: 763 *nrm = El::OneNorm(*a->emat); 764 break; 765 case NORM_FROBENIUS: 766 *nrm = El::FrobeniusNorm(*a->emat); 767 break; 768 case NORM_INFINITY: 769 *nrm = El::InfinityNorm(*a->emat); 770 break; 771 default: 772 printf("Error: unsupported norm type!\n"); 773 } 774 PetscFunctionReturn(0); 775 } 776 777 #undef __FUNCT__ 778 #define __FUNCT__ "MatZeroEntries_Elemental" 779 static PetscErrorCode MatZeroEntries_Elemental(Mat A) 780 { 781 Mat_Elemental *a=(Mat_Elemental*)A->data; 782 783 PetscFunctionBegin; 784 El::Zero(*a->emat); 785 PetscFunctionReturn(0); 786 } 787 788 #undef __FUNCT__ 789 #define __FUNCT__ "MatGetOwnershipIS_Elemental" 790 static PetscErrorCode MatGetOwnershipIS_Elemental(Mat A,IS *rows,IS *cols) 791 { 792 Mat_Elemental *a = (Mat_Elemental*)A->data; 793 PetscErrorCode ierr; 794 PetscInt i,m,shift,stride,*idx; 795 796 PetscFunctionBegin; 797 if (rows) { 798 m = a->emat->LocalHeight(); 799 shift = a->emat->ColShift(); 800 stride = a->emat->ColStride(); 801 ierr = PetscMalloc1(m,&idx);CHKERRQ(ierr); 802 for (i=0; i<m; i++) { 803 PetscInt rank,offset; 804 E2RO(A,0,shift+i*stride,&rank,&offset); 805 RO2P(A,0,rank,offset,&idx[i]); 806 } 807 ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx,PETSC_OWN_POINTER,rows);CHKERRQ(ierr); 808 } 809 if (cols) { 810 m = a->emat->LocalWidth(); 811 shift = a->emat->RowShift(); 812 stride = a->emat->RowStride(); 813 ierr = PetscMalloc1(m,&idx);CHKERRQ(ierr); 814 for (i=0; i<m; i++) { 815 PetscInt rank,offset; 816 E2RO(A,1,shift+i*stride,&rank,&offset); 817 RO2P(A,1,rank,offset,&idx[i]); 818 } 819 ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx,PETSC_OWN_POINTER,cols);CHKERRQ(ierr); 820 } 821 PetscFunctionReturn(0); 822 } 823 824 #undef __FUNCT__ 825 #define __FUNCT__ "MatConvert_Elemental_Dense" 826 static PetscErrorCode MatConvert_Elemental_Dense(Mat A,MatType newtype,MatReuse reuse,Mat *B) 827 { 828 Mat Bmpi; 829 Mat_Elemental *a = (Mat_Elemental*)A->data; 830 MPI_Comm comm; 831 PetscErrorCode ierr; 832 PetscInt rrank,ridx,crank,cidx,nrows,ncols,i,j; 833 PetscElemScalar v; 834 PetscBool s1,s2,s3; 835 836 PetscFunctionBegin; 837 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 838 ierr = PetscStrcmp(newtype,MATDENSE,&s1);CHKERRQ(ierr); 839 ierr = PetscStrcmp(newtype,MATSEQDENSE,&s2);CHKERRQ(ierr); 840 ierr = PetscStrcmp(newtype,MATMPIDENSE,&s3);CHKERRQ(ierr); 841 if (!s1 && !s2 && !s3) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported New MatType: must be MATDENSE, MATSEQDENSE or MATMPIDENSE"); 842 ierr = MatCreate(comm,&Bmpi);CHKERRQ(ierr); 843 ierr = MatSetSizes(Bmpi,A->rmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 844 ierr = MatSetType(Bmpi,MATDENSE);CHKERRQ(ierr); 845 ierr = MatSetUp(Bmpi);CHKERRQ(ierr); 846 ierr = MatGetSize(A,&nrows,&ncols);CHKERRQ(ierr); 847 for (i=0; i<nrows; i++) { 848 PetscInt erow,ecol; 849 P2RO(A,0,i,&rrank,&ridx); 850 RO2E(A,0,rrank,ridx,&erow); 851 if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect row translation"); 852 for (j=0; j<ncols; j++) { 853 P2RO(A,1,j,&crank,&cidx); 854 RO2E(A,1,crank,cidx,&ecol); 855 if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect col translation"); 856 v = a->emat->Get(erow,ecol); 857 ierr = MatSetValues(Bmpi,1,&i,1,&j,(PetscScalar *)&v,INSERT_VALUES);CHKERRQ(ierr); 858 } 859 } 860 ierr = MatAssemblyBegin(Bmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 861 ierr = MatAssemblyEnd(Bmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 862 if (reuse == MAT_REUSE_MATRIX) { 863 ierr = MatHeaderReplace(A,Bmpi);CHKERRQ(ierr); 864 } else { 865 *B = Bmpi; 866 } 867 PetscFunctionReturn(0); 868 } 869 870 #undef __FUNCT__ 871 #define __FUNCT__ "MatConvert_SeqAIJ_Elemental" 872 PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 873 { 874 Mat mat_elemental; 875 PetscErrorCode ierr; 876 PetscInt M=A->rmap->N,N=A->cmap->N,row,ncols; 877 const PetscInt *cols; 878 const PetscScalar *vals; 879 880 PetscFunctionBegin; 881 ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr); 882 ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 883 ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr); 884 ierr = MatSetUp(mat_elemental);CHKERRQ(ierr); 885 for (row=0; row<M; row++) { 886 ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 887 /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */ 888 ierr = MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 889 ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 890 } 891 ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 892 ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 893 894 if (reuse == MAT_REUSE_MATRIX) { 895 ierr = MatHeaderReplace(A,mat_elemental);CHKERRQ(ierr); 896 } else { 897 *newmat = mat_elemental; 898 } 899 PetscFunctionReturn(0); 900 } 901 902 #undef __FUNCT__ 903 #define __FUNCT__ "MatConvert_MPIAIJ_Elemental" 904 PETSC_EXTERN PetscErrorCode MatConvert_MPIAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 905 { 906 Mat mat_elemental; 907 PetscErrorCode ierr; 908 PetscInt row,ncols,rstart=A->rmap->rstart,rend=A->rmap->rend,j; 909 const PetscInt *cols; 910 const PetscScalar *vals; 911 912 PetscFunctionBegin; 913 ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr); 914 ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 915 ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr); 916 ierr = MatSetUp(mat_elemental);CHKERRQ(ierr); 917 for (row=rstart; row<rend; row++) { 918 ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 919 for (j=0; j<ncols; j++) { 920 /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */ 921 ierr = MatSetValues(mat_elemental,1,&row,1,&cols[j],&vals[j],ADD_VALUES);CHKERRQ(ierr); 922 } 923 ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 924 } 925 ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 926 ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 927 928 if (reuse == MAT_REUSE_MATRIX) { 929 ierr = MatHeaderReplace(A,mat_elemental);CHKERRQ(ierr); 930 } else { 931 *newmat = mat_elemental; 932 } 933 PetscFunctionReturn(0); 934 } 935 936 #undef __FUNCT__ 937 #define __FUNCT__ "MatConvert_SeqSBAIJ_Elemental" 938 PETSC_EXTERN PetscErrorCode MatConvert_SeqSBAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 939 { 940 Mat mat_elemental; 941 PetscErrorCode ierr; 942 PetscInt M=A->rmap->N,N=A->cmap->N,row,ncols,j; 943 const PetscInt *cols; 944 const PetscScalar *vals; 945 946 PetscFunctionBegin; 947 ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr); 948 ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 949 ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr); 950 ierr = MatSetUp(mat_elemental);CHKERRQ(ierr); 951 ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr); 952 for (row=0; row<M; row++) { 953 ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 954 /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */ 955 ierr = MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 956 for (j=0; j<ncols; j++) { /* lower triangular part */ 957 if (cols[j] == row) continue; 958 ierr = MatSetValues(mat_elemental,1,&cols[j],1,&row,&vals[j],ADD_VALUES);CHKERRQ(ierr); 959 } 960 ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 961 } 962 ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr); 963 ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 964 ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 965 966 if (reuse == MAT_REUSE_MATRIX) { 967 ierr = MatHeaderReplace(A,mat_elemental);CHKERRQ(ierr); 968 } else { 969 *newmat = mat_elemental; 970 } 971 PetscFunctionReturn(0); 972 } 973 974 #undef __FUNCT__ 975 #define __FUNCT__ "MatConvert_MPISBAIJ_Elemental" 976 PETSC_EXTERN PetscErrorCode MatConvert_MPISBAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 977 { 978 Mat mat_elemental; 979 PetscErrorCode ierr; 980 PetscInt M=A->rmap->N,N=A->cmap->N,row,ncols,j,rstart=A->rmap->rstart,rend=A->rmap->rend; 981 const PetscInt *cols; 982 const PetscScalar *vals; 983 984 PetscFunctionBegin; 985 ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr); 986 ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 987 ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr); 988 ierr = MatSetUp(mat_elemental);CHKERRQ(ierr); 989 ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr); 990 for (row=rstart; row<rend; row++) { 991 ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 992 /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */ 993 ierr = MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 994 for (j=0; j<ncols; j++) { /* lower triangular part */ 995 if (cols[j] == row) continue; 996 ierr = MatSetValues(mat_elemental,1,&cols[j],1,&row,&vals[j],ADD_VALUES);CHKERRQ(ierr); 997 } 998 ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 999 } 1000 ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr); 1001 ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1002 ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1003 1004 if (reuse == MAT_REUSE_MATRIX) { 1005 ierr = MatHeaderReplace(A,mat_elemental);CHKERRQ(ierr); 1006 } else { 1007 *newmat = mat_elemental; 1008 } 1009 PetscFunctionReturn(0); 1010 } 1011 1012 #undef __FUNCT__ 1013 #define __FUNCT__ "MatDestroy_Elemental" 1014 static PetscErrorCode MatDestroy_Elemental(Mat A) 1015 { 1016 Mat_Elemental *a = (Mat_Elemental*)A->data; 1017 PetscErrorCode ierr; 1018 Mat_Elemental_Grid *commgrid; 1019 PetscBool flg; 1020 MPI_Comm icomm; 1021 1022 PetscFunctionBegin; 1023 a->interface->Detach(); 1024 delete a->interface; 1025 delete a->esubmat; 1026 delete a->emat; 1027 delete a->pivot; 1028 1029 El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A)); 1030 ierr = PetscCommDuplicate(cxxcomm.comm,&icomm,NULL);CHKERRQ(ierr); 1031 ierr = MPI_Attr_get(icomm,Petsc_Elemental_keyval,(void**)&commgrid,(int*)&flg);CHKERRQ(ierr); 1032 /* printf("commgrid->grid_refct = %d, grid=%p\n",commgrid->grid_refct,commgrid->grid); -- memory leak revealed by valgrind? */ 1033 if (--commgrid->grid_refct == 0) { 1034 delete commgrid->grid; 1035 ierr = PetscFree(commgrid);CHKERRQ(ierr); 1036 } 1037 ierr = PetscCommDestroy(&icomm);CHKERRQ(ierr); 1038 ierr = PetscObjectComposeFunction((PetscObject)A,"MatGetOwnershipIS_C",NULL);CHKERRQ(ierr); 1039 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr); 1040 ierr = PetscObjectComposeFunction((PetscObject)A,"MatElementalHermitianGenDefEig_C",NULL);CHKERRQ(ierr); 1041 ierr = PetscFree(A->data);CHKERRQ(ierr); 1042 PetscFunctionReturn(0); 1043 } 1044 1045 #undef __FUNCT__ 1046 #define __FUNCT__ "MatSetUp_Elemental" 1047 PetscErrorCode MatSetUp_Elemental(Mat A) 1048 { 1049 Mat_Elemental *a = (Mat_Elemental*)A->data; 1050 PetscErrorCode ierr; 1051 PetscMPIInt rsize,csize; 1052 1053 PetscFunctionBegin; 1054 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1055 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1056 1057 a->emat->Resize(A->rmap->N,A->cmap->N);CHKERRQ(ierr); 1058 El::Zero(*a->emat); 1059 1060 ierr = MPI_Comm_size(A->rmap->comm,&rsize);CHKERRQ(ierr); 1061 ierr = MPI_Comm_size(A->cmap->comm,&csize);CHKERRQ(ierr); 1062 if (csize != rsize) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Cannot use row and column communicators of different sizes"); 1063 a->commsize = rsize; 1064 a->mr[0] = A->rmap->N % rsize; if (!a->mr[0]) a->mr[0] = rsize; 1065 a->mr[1] = A->cmap->N % csize; if (!a->mr[1]) a->mr[1] = csize; 1066 a->m[0] = A->rmap->N / rsize + (a->mr[0] != rsize); 1067 a->m[1] = A->cmap->N / csize + (a->mr[1] != csize); 1068 PetscFunctionReturn(0); 1069 } 1070 1071 #undef __FUNCT__ 1072 #define __FUNCT__ "MatAssemblyBegin_Elemental" 1073 PetscErrorCode MatAssemblyBegin_Elemental(Mat A, MatAssemblyType type) 1074 { 1075 Mat_Elemental *a = (Mat_Elemental*)A->data; 1076 1077 PetscFunctionBegin; 1078 a->interface->Detach(); 1079 a->interface->Attach(El::LOCAL_TO_GLOBAL,*(a->emat)); 1080 PetscFunctionReturn(0); 1081 } 1082 1083 #undef __FUNCT__ 1084 #define __FUNCT__ "MatAssemblyEnd_Elemental" 1085 PetscErrorCode MatAssemblyEnd_Elemental(Mat A, MatAssemblyType type) 1086 { 1087 PetscFunctionBegin; 1088 /* Currently does nothing */ 1089 PetscFunctionReturn(0); 1090 } 1091 1092 #undef __FUNCT__ 1093 #define __FUNCT__ "MatLoad_Elemental" 1094 PetscErrorCode MatLoad_Elemental(Mat newMat, PetscViewer viewer) 1095 { 1096 PetscErrorCode ierr; 1097 Mat Adense,Ae; 1098 MPI_Comm comm; 1099 1100 PetscFunctionBegin; 1101 ierr = PetscObjectGetComm((PetscObject)newMat,&comm);CHKERRQ(ierr); 1102 ierr = MatCreate(comm,&Adense);CHKERRQ(ierr); 1103 ierr = MatSetType(Adense,MATDENSE);CHKERRQ(ierr); 1104 ierr = MatLoad(Adense,viewer);CHKERRQ(ierr); 1105 ierr = MatConvert(Adense, MATELEMENTAL, MAT_INITIAL_MATRIX,&Ae);CHKERRQ(ierr); 1106 ierr = MatDestroy(&Adense);CHKERRQ(ierr); 1107 ierr = MatHeaderReplace(newMat,Ae);CHKERRQ(ierr); 1108 PetscFunctionReturn(0); 1109 } 1110 1111 #undef __FUNCT__ 1112 #define __FUNCT__ "MatElementalHermitianGenDefEig_Elemental" 1113 PetscErrorCode MatElementalHermitianGenDefEig_Elemental(El::Pencil eigtype,El::UpperOrLower uplo,Mat A,Mat B,Mat *evals,Mat *evec,El::SortType sort,El::HermitianEigSubset<PetscElemScalar> subset,const El::HermitianEigCtrl<PetscElemScalar> ctrl) 1114 { 1115 PetscErrorCode ierr; 1116 Mat_Elemental *a=(Mat_Elemental*)A->data,*b=(Mat_Elemental*)B->data,*x; 1117 MPI_Comm comm; 1118 Mat EVAL; 1119 Mat_Elemental *e; 1120 1121 PetscFunctionBegin; 1122 /* Compute eigenvalues and eigenvectors */ 1123 El::DistMatrix<PetscElemScalar,El::VR,El::STAR> w( *a->grid ); /* holding eigenvalues */ 1124 El::DistMatrix<PetscElemScalar> X( *a->grid ); /* holding eigenvectors */ 1125 El::HermitianGenDefEig(eigtype,uplo,*a->emat,*b->emat,w,X,sort,subset,ctrl); 1126 /* El::Print(w, "Eigenvalues"); */ 1127 1128 /* Wrap w and X into PETSc's MATMATELEMENTAL matrices */ 1129 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1130 ierr = MatCreate(comm,evec);CHKERRQ(ierr); 1131 ierr = MatSetSizes(*evec,PETSC_DECIDE,PETSC_DECIDE,X.Height(),X.Width());CHKERRQ(ierr); 1132 ierr = MatSetType(*evec,MATELEMENTAL);CHKERRQ(ierr); 1133 ierr = MatSetFromOptions(*evec);CHKERRQ(ierr); 1134 ierr = MatSetUp(*evec);CHKERRQ(ierr); 1135 ierr = MatAssemblyBegin(*evec,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1136 ierr = MatAssemblyEnd(*evec,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1137 1138 x = (Mat_Elemental*)(*evec)->data; 1139 //delete x->emat; //-- memory leak??? 1140 *x->emat = X; 1141 1142 ierr = MatCreate(comm,&EVAL);CHKERRQ(ierr); 1143 ierr = MatSetSizes(EVAL,PETSC_DECIDE,PETSC_DECIDE,w.Height(),w.Width());CHKERRQ(ierr); 1144 ierr = MatSetType(EVAL,MATELEMENTAL);CHKERRQ(ierr); 1145 ierr = MatSetFromOptions(EVAL);CHKERRQ(ierr); 1146 ierr = MatSetUp(EVAL);CHKERRQ(ierr); 1147 ierr = MatAssemblyBegin(EVAL,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1148 ierr = MatAssemblyEnd(EVAL,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1149 e = (Mat_Elemental*)EVAL->data; 1150 *e->emat = w; //-- memory leak??? 1151 *evals = EVAL; 1152 1153 #if defined(MV) 1154 /* Test correctness norm = || - A*X + B*X*w || */ 1155 { 1156 PetscElemScalar alpha,beta; 1157 El::DistMatrix<PetscElemScalar> Y(*a->grid); //tmp matrix 1158 alpha = 1.0; beta=0.0; 1159 El::Gemm(El::NORMAL,El::NORMAL,alpha,*b->emat,X,beta,Y); //Y = B*X 1160 El::DiagonalScale(El::RIGHT,El::NORMAL, w, Y); //Y = Y*w 1161 alpha = -1.0; beta=1.0; 1162 El::Gemm(El::NORMAL,El::NORMAL,alpha,*a->emat,X,beta,Y); //Y = - A*X + B*X*w 1163 1164 PetscElemScalar norm = El::FrobeniusNorm(Y); 1165 if ((*a->grid).Rank()==0) printf(" norm (- A*X + B*X*w) = %g\n",norm); 1166 } 1167 1168 { 1169 PetscMPIInt rank; 1170 ierr = MPI_Comm_rank(comm,&rank); 1171 printf("w: [%d] [%d, %d %d] %d; X: %d %d\n",rank,w.DistRank(),w.ColRank(),w.RowRank(),w.LocalHeight(),X.LocalHeight(),X.LocalWidth()); 1172 } 1173 #endif 1174 PetscFunctionReturn(0); 1175 } 1176 1177 #undef __FUNCT__ 1178 #define __FUNCT__ "MatElementalHermitianGenDefEig" 1179 /*@ 1180 MatElementalHermitianGenDefEig - Compute the set of eigenvalues of the Hermitian-definite matrix pencil determined by the subset structure 1181 1182 Logically Collective on Mat 1183 1184 Level: beginner 1185 1186 References: Elemental Users' Guide 1187 @*/ 1188 PetscErrorCode MatElementalHermitianGenDefEig(El::Pencil type,El::UpperOrLower uplo,Mat A,Mat B,Mat *evals,Mat *evec,El::SortType sort,El::HermitianEigSubset<PetscElemScalar> subset,const El::HermitianEigCtrl<PetscElemScalar> ctrl) 1189 { 1190 PetscErrorCode ierr; 1191 1192 PetscFunctionBegin; 1193 ierr = PetscTryMethod(A,"MatElementalHermitianGenDefEig_C",(El::Pencil,El::UpperOrLower,Mat,Mat,Mat*,Mat*,El::SortType,El::HermitianEigSubset<PetscElemScalar>,const El::HermitianEigCtrl<PetscElemScalar>),(type,uplo,A,B,evals,evec,sort,subset,ctrl));CHKERRQ(ierr); 1194 PetscFunctionReturn(0); 1195 } 1196 1197 /* -------------------------------------------------------------------*/ 1198 static struct _MatOps MatOps_Values = { 1199 MatSetValues_Elemental, 1200 0, 1201 0, 1202 MatMult_Elemental, 1203 /* 4*/ MatMultAdd_Elemental, 1204 MatMultTranspose_Elemental, 1205 MatMultTransposeAdd_Elemental, 1206 MatSolve_Elemental, 1207 MatSolveAdd_Elemental, 1208 0, 1209 /*10*/ 0, 1210 MatLUFactor_Elemental, 1211 MatCholeskyFactor_Elemental, 1212 0, 1213 MatTranspose_Elemental, 1214 /*15*/ MatGetInfo_Elemental, 1215 0, 1216 MatGetDiagonal_Elemental, 1217 MatDiagonalScale_Elemental, 1218 MatNorm_Elemental, 1219 /*20*/ MatAssemblyBegin_Elemental, 1220 MatAssemblyEnd_Elemental, 1221 0, 1222 MatZeroEntries_Elemental, 1223 /*24*/ 0, 1224 MatLUFactorSymbolic_Elemental, 1225 MatLUFactorNumeric_Elemental, 1226 MatCholeskyFactorSymbolic_Elemental, 1227 MatCholeskyFactorNumeric_Elemental, 1228 /*29*/ MatSetUp_Elemental, 1229 0, 1230 0, 1231 0, 1232 0, 1233 /*34*/ MatDuplicate_Elemental, 1234 0, 1235 0, 1236 0, 1237 0, 1238 /*39*/ MatAXPY_Elemental, 1239 0, 1240 0, 1241 0, 1242 MatCopy_Elemental, 1243 /*44*/ 0, 1244 MatScale_Elemental, 1245 MatShift_Basic, 1246 0, 1247 0, 1248 /*49*/ 0, 1249 0, 1250 0, 1251 0, 1252 0, 1253 /*54*/ 0, 1254 0, 1255 0, 1256 0, 1257 0, 1258 /*59*/ 0, 1259 MatDestroy_Elemental, 1260 MatView_Elemental, 1261 0, 1262 0, 1263 /*64*/ 0, 1264 0, 1265 0, 1266 0, 1267 0, 1268 /*69*/ 0, 1269 0, 1270 MatConvert_Elemental_Dense, 1271 0, 1272 0, 1273 /*74*/ 0, 1274 0, 1275 0, 1276 0, 1277 0, 1278 /*79*/ 0, 1279 0, 1280 0, 1281 0, 1282 MatLoad_Elemental, 1283 /*84*/ 0, 1284 0, 1285 0, 1286 0, 1287 0, 1288 /*89*/ MatMatMult_Elemental, 1289 MatMatMultSymbolic_Elemental, 1290 MatMatMultNumeric_Elemental, 1291 0, 1292 0, 1293 /*94*/ 0, 1294 MatMatTransposeMult_Elemental, 1295 MatMatTransposeMultSymbolic_Elemental, 1296 MatMatTransposeMultNumeric_Elemental, 1297 0, 1298 /*99*/ 0, 1299 0, 1300 0, 1301 MatConjugate_Elemental, 1302 0, 1303 /*104*/0, 1304 0, 1305 0, 1306 0, 1307 0, 1308 /*109*/MatMatSolve_Elemental, 1309 0, 1310 0, 1311 0, 1312 0, 1313 /*114*/0, 1314 0, 1315 0, 1316 0, 1317 0, 1318 /*119*/0, 1319 MatHermitianTranspose_Elemental, 1320 0, 1321 0, 1322 0, 1323 /*124*/0, 1324 0, 1325 0, 1326 0, 1327 0, 1328 /*129*/0, 1329 0, 1330 0, 1331 0, 1332 0, 1333 /*134*/0, 1334 0, 1335 0, 1336 0, 1337 0 1338 }; 1339 1340 /*MC 1341 MATELEMENTAL = "elemental" - A matrix type for dense matrices using the Elemental package 1342 1343 Options Database Keys: 1344 + -mat_type elemental - sets the matrix type to "elemental" during a call to MatSetFromOptions() 1345 - -mat_elemental_grid_height - sets Grid Height for 2D cyclic ordering of internal matrix 1346 1347 Level: beginner 1348 1349 .seealso: MATDENSE 1350 M*/ 1351 1352 #undef __FUNCT__ 1353 #define __FUNCT__ "MatCreate_Elemental" 1354 PETSC_EXTERN PetscErrorCode MatCreate_Elemental(Mat A) 1355 { 1356 Mat_Elemental *a; 1357 PetscErrorCode ierr; 1358 PetscBool flg,flg1; 1359 Mat_Elemental_Grid *commgrid; 1360 MPI_Comm icomm; 1361 PetscInt optv1; 1362 1363 PetscFunctionBegin; 1364 ierr = PetscElementalInitializePackage();CHKERRQ(ierr); 1365 ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1366 A->insertmode = NOT_SET_VALUES; 1367 1368 ierr = PetscNewLog(A,&a);CHKERRQ(ierr); 1369 A->data = (void*)a; 1370 1371 /* Set up the elemental matrix */ 1372 El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A)); 1373 1374 /* Grid needs to be shared between multiple Mats on the same communicator, implement by attribute caching on the MPI_Comm */ 1375 if (Petsc_Elemental_keyval == MPI_KEYVAL_INVALID) { 1376 ierr = MPI_Keyval_create(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Elemental_keyval,(void*)0); 1377 } 1378 ierr = PetscCommDuplicate(cxxcomm.comm,&icomm,NULL);CHKERRQ(ierr); 1379 ierr = MPI_Attr_get(icomm,Petsc_Elemental_keyval,(void**)&commgrid,(int*)&flg);CHKERRQ(ierr); 1380 if (!flg) { 1381 ierr = PetscNewLog(A,&commgrid);CHKERRQ(ierr); 1382 1383 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"Elemental Options","Mat");CHKERRQ(ierr); 1384 /* displayed default grid sizes (CommSize,1) are set by us arbitrarily until El::Grid() is called */ 1385 ierr = PetscOptionsInt("-mat_elemental_grid_height","Grid Height","None",El::mpi::Size(cxxcomm),&optv1,&flg1);CHKERRQ(ierr); 1386 if (flg1) { 1387 if (El::mpi::Size(cxxcomm) % optv1 != 0) { 1388 SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Grid Height %D must evenly divide CommSize %D",optv1,(PetscInt)El::mpi::Size(cxxcomm)); 1389 } 1390 commgrid->grid = new El::Grid(cxxcomm,optv1); /* use user-provided grid height */ 1391 } else { 1392 commgrid->grid = new El::Grid(cxxcomm); /* use Elemental default grid sizes */ 1393 /* printf("new commgrid->grid = %p\n",commgrid->grid); -- memory leak revealed by valgrind? */ 1394 } 1395 commgrid->grid_refct = 1; 1396 ierr = MPI_Attr_put(icomm,Petsc_Elemental_keyval,(void*)commgrid);CHKERRQ(ierr); 1397 PetscOptionsEnd(); 1398 } else { 1399 commgrid->grid_refct++; 1400 } 1401 ierr = PetscCommDestroy(&icomm);CHKERRQ(ierr); 1402 a->grid = commgrid->grid; 1403 a->emat = new El::DistMatrix<PetscElemScalar>(*a->grid); 1404 a->esubmat = new El::Matrix<PetscElemScalar>(1,1); 1405 a->interface = new El::AxpyInterface<PetscElemScalar>; 1406 a->pivot = new El::DistMatrix<PetscInt,El::VC,El::STAR>; 1407 1408 /* build cache for off array entries formed */ 1409 a->interface->Attach(El::LOCAL_TO_GLOBAL,*(a->emat)); 1410 1411 ierr = PetscObjectComposeFunction((PetscObject)A,"MatGetOwnershipIS_C",MatGetOwnershipIS_Elemental);CHKERRQ(ierr); 1412 ierr = PetscObjectComposeFunction((PetscObject)A,"MatElementalHermitianGenDefEig_C",MatElementalHermitianGenDefEig_Elemental);CHKERRQ(ierr); 1413 1414 ierr = PetscObjectChangeTypeName((PetscObject)A,MATELEMENTAL);CHKERRQ(ierr); 1415 PetscFunctionReturn(0); 1416 } 1417