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