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