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