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