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