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 for (i=0; i<nrows; i++) { 927 P2RO(A,0,rows[i],&rrank,&ridx); /* convert indices between PETSc <-> (Rank,Offset) <-> Elemental */ 928 RO2E(A,0,rrank,ridx,&erow); 929 if (rrank < 0 || ridx < 0 || erow < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect row translation"); 930 for (j=0; j<ncols; j++) { 931 P2RO(A,1,cols[j],&crank,&cidx); 932 RO2E(A,1,crank,cidx,&ecol); 933 if (crank < 0 || cidx < 0 || ecol < 0) SETERRQ(comm,PETSC_ERR_PLIB,"Incorrect col translation"); 934 935 elrow = erow / grid.MCSize(); /* Elemental local row index */ 936 elcol = ecol / grid.MRSize(); /* Elemental local column index */ 937 v = a->emat->GetLocal(elrow,elcol); 938 ierr = MatSetValues(Bmpi,1,&rows[i],1,&cols[j],(PetscScalar *)&v,INSERT_VALUES);CHKERRQ(ierr); 939 } 940 } 941 ierr = MatAssemblyBegin(Bmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 942 ierr = MatAssemblyEnd(Bmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 943 if (reuse == MAT_INPLACE_MATRIX) { 944 ierr = MatHeaderReplace(A,&Bmpi);CHKERRQ(ierr); 945 } else { 946 *B = Bmpi; 947 } 948 ierr = ISDestroy(&isrows);CHKERRQ(ierr); 949 ierr = ISDestroy(&iscols);CHKERRQ(ierr); 950 PetscFunctionReturn(0); 951 } 952 953 #undef __FUNCT__ 954 #define __FUNCT__ "MatConvert_SeqAIJ_Elemental" 955 PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 956 { 957 Mat mat_elemental; 958 PetscErrorCode ierr; 959 PetscInt M=A->rmap->N,N=A->cmap->N,row,ncols; 960 const PetscInt *cols; 961 const PetscScalar *vals; 962 963 PetscFunctionBegin; 964 ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr); 965 ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 966 ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr); 967 ierr = MatSetUp(mat_elemental);CHKERRQ(ierr); 968 for (row=0; row<M; row++) { 969 ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 970 /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */ 971 ierr = MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 972 ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 973 } 974 ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 975 ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 976 977 if (reuse == MAT_INPLACE_MATRIX) { 978 ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr); 979 } else { 980 *newmat = mat_elemental; 981 } 982 PetscFunctionReturn(0); 983 } 984 985 #undef __FUNCT__ 986 #define __FUNCT__ "MatConvert_MPIAIJ_Elemental" 987 PETSC_EXTERN PetscErrorCode MatConvert_MPIAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 988 { 989 Mat mat_elemental; 990 PetscErrorCode ierr; 991 PetscInt row,ncols,rstart=A->rmap->rstart,rend=A->rmap->rend,j; 992 const PetscInt *cols; 993 const PetscScalar *vals; 994 995 PetscFunctionBegin; 996 ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr); 997 ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 998 ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr); 999 ierr = MatSetUp(mat_elemental);CHKERRQ(ierr); 1000 for (row=rstart; row<rend; row++) { 1001 ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 1002 for (j=0; j<ncols; j++) { 1003 /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */ 1004 ierr = MatSetValues(mat_elemental,1,&row,1,&cols[j],&vals[j],ADD_VALUES);CHKERRQ(ierr); 1005 } 1006 ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 1007 } 1008 ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1009 ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1010 1011 if (reuse == MAT_INPLACE_MATRIX) { 1012 ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr); 1013 } else { 1014 *newmat = mat_elemental; 1015 } 1016 PetscFunctionReturn(0); 1017 } 1018 1019 #undef __FUNCT__ 1020 #define __FUNCT__ "MatConvert_SeqSBAIJ_Elemental" 1021 PETSC_EXTERN PetscErrorCode MatConvert_SeqSBAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 1022 { 1023 Mat mat_elemental; 1024 PetscErrorCode ierr; 1025 PetscInt M=A->rmap->N,N=A->cmap->N,row,ncols,j; 1026 const PetscInt *cols; 1027 const PetscScalar *vals; 1028 1029 PetscFunctionBegin; 1030 ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr); 1031 ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1032 ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr); 1033 ierr = MatSetUp(mat_elemental);CHKERRQ(ierr); 1034 ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr); 1035 for (row=0; row<M; row++) { 1036 ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 1037 /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */ 1038 ierr = MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 1039 for (j=0; j<ncols; j++) { /* lower triangular part */ 1040 if (cols[j] == row) continue; 1041 ierr = MatSetValues(mat_elemental,1,&cols[j],1,&row,&vals[j],ADD_VALUES);CHKERRQ(ierr); 1042 } 1043 ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 1044 } 1045 ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr); 1046 ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1047 ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1048 1049 if (reuse == MAT_INPLACE_MATRIX) { 1050 ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr); 1051 } else { 1052 *newmat = mat_elemental; 1053 } 1054 PetscFunctionReturn(0); 1055 } 1056 1057 #undef __FUNCT__ 1058 #define __FUNCT__ "MatConvert_MPISBAIJ_Elemental" 1059 PETSC_EXTERN PetscErrorCode MatConvert_MPISBAIJ_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 1060 { 1061 Mat mat_elemental; 1062 PetscErrorCode ierr; 1063 PetscInt M=A->rmap->N,N=A->cmap->N,row,ncols,j,rstart=A->rmap->rstart,rend=A->rmap->rend; 1064 const PetscInt *cols; 1065 const PetscScalar *vals; 1066 1067 PetscFunctionBegin; 1068 ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr); 1069 ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1070 ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr); 1071 ierr = MatSetUp(mat_elemental);CHKERRQ(ierr); 1072 ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr); 1073 for (row=rstart; row<rend; row++) { 1074 ierr = MatGetRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 1075 /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */ 1076 ierr = MatSetValues(mat_elemental,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 1077 for (j=0; j<ncols; j++) { /* lower triangular part */ 1078 if (cols[j] == row) continue; 1079 ierr = MatSetValues(mat_elemental,1,&cols[j],1,&row,&vals[j],ADD_VALUES);CHKERRQ(ierr); 1080 } 1081 ierr = MatRestoreRow(A,row,&ncols,&cols,&vals);CHKERRQ(ierr); 1082 } 1083 ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr); 1084 ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1085 ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1086 1087 if (reuse == MAT_INPLACE_MATRIX) { 1088 ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr); 1089 } else { 1090 *newmat = mat_elemental; 1091 } 1092 PetscFunctionReturn(0); 1093 } 1094 1095 #undef __FUNCT__ 1096 #define __FUNCT__ "MatDestroy_Elemental" 1097 static PetscErrorCode MatDestroy_Elemental(Mat A) 1098 { 1099 Mat_Elemental *a = (Mat_Elemental*)A->data; 1100 PetscErrorCode ierr; 1101 Mat_Elemental_Grid *commgrid; 1102 PetscBool flg; 1103 MPI_Comm icomm; 1104 1105 PetscFunctionBegin; 1106 delete a->emat; 1107 delete a->pivot; 1108 1109 El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A)); 1110 ierr = PetscCommDuplicate(cxxcomm.comm,&icomm,NULL);CHKERRQ(ierr); 1111 ierr = MPI_Attr_get(icomm,Petsc_Elemental_keyval,(void**)&commgrid,(int*)&flg);CHKERRQ(ierr); 1112 if (--commgrid->grid_refct == 0) { 1113 delete commgrid->grid; 1114 ierr = PetscFree(commgrid);CHKERRQ(ierr); 1115 ierr = MPI_Keyval_free(&Petsc_Elemental_keyval);CHKERRQ(ierr); 1116 } 1117 ierr = PetscCommDestroy(&icomm);CHKERRQ(ierr); 1118 ierr = PetscObjectComposeFunction((PetscObject)A,"MatGetOwnershipIS_C",NULL);CHKERRQ(ierr); 1119 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr); 1120 ierr = PetscObjectComposeFunction((PetscObject)A,"MatElementalHermitianGenDefEig_C",NULL);CHKERRQ(ierr); 1121 ierr = PetscFree(A->data);CHKERRQ(ierr); 1122 PetscFunctionReturn(0); 1123 } 1124 1125 #undef __FUNCT__ 1126 #define __FUNCT__ "MatSetUp_Elemental" 1127 PetscErrorCode MatSetUp_Elemental(Mat A) 1128 { 1129 Mat_Elemental *a = (Mat_Elemental*)A->data; 1130 PetscErrorCode ierr; 1131 PetscMPIInt rsize,csize; 1132 1133 PetscFunctionBegin; 1134 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1135 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1136 1137 a->emat->Resize(A->rmap->N,A->cmap->N);CHKERRQ(ierr); 1138 El::Zero(*a->emat); 1139 1140 ierr = MPI_Comm_size(A->rmap->comm,&rsize);CHKERRQ(ierr); 1141 ierr = MPI_Comm_size(A->cmap->comm,&csize);CHKERRQ(ierr); 1142 if (csize != rsize) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Cannot use row and column communicators of different sizes"); 1143 a->commsize = rsize; 1144 a->mr[0] = A->rmap->N % rsize; if (!a->mr[0]) a->mr[0] = rsize; 1145 a->mr[1] = A->cmap->N % csize; if (!a->mr[1]) a->mr[1] = csize; 1146 a->m[0] = A->rmap->N / rsize + (a->mr[0] != rsize); 1147 a->m[1] = A->cmap->N / csize + (a->mr[1] != csize); 1148 PetscFunctionReturn(0); 1149 } 1150 1151 #undef __FUNCT__ 1152 #define __FUNCT__ "MatAssemblyBegin_Elemental" 1153 PetscErrorCode MatAssemblyBegin_Elemental(Mat A, MatAssemblyType type) 1154 { 1155 Mat_Elemental *a = (Mat_Elemental*)A->data; 1156 1157 PetscFunctionBegin; 1158 /* printf("Calling ProcessQueues\n"); */ 1159 a->emat->ProcessQueues(); 1160 /* printf("Finished ProcessQueues\n"); */ 1161 PetscFunctionReturn(0); 1162 } 1163 1164 #undef __FUNCT__ 1165 #define __FUNCT__ "MatAssemblyEnd_Elemental" 1166 PetscErrorCode MatAssemblyEnd_Elemental(Mat A, MatAssemblyType type) 1167 { 1168 PetscFunctionBegin; 1169 /* Currently does nothing */ 1170 PetscFunctionReturn(0); 1171 } 1172 1173 #undef __FUNCT__ 1174 #define __FUNCT__ "MatLoad_Elemental" 1175 PetscErrorCode MatLoad_Elemental(Mat newMat, PetscViewer viewer) 1176 { 1177 PetscErrorCode ierr; 1178 Mat Adense,Ae; 1179 MPI_Comm comm; 1180 1181 PetscFunctionBegin; 1182 ierr = PetscObjectGetComm((PetscObject)newMat,&comm);CHKERRQ(ierr); 1183 ierr = MatCreate(comm,&Adense);CHKERRQ(ierr); 1184 ierr = MatSetType(Adense,MATDENSE);CHKERRQ(ierr); 1185 ierr = MatLoad(Adense,viewer);CHKERRQ(ierr); 1186 ierr = MatConvert(Adense, MATELEMENTAL, MAT_INITIAL_MATRIX,&Ae);CHKERRQ(ierr); 1187 ierr = MatDestroy(&Adense);CHKERRQ(ierr); 1188 ierr = MatHeaderReplace(newMat,&Ae);CHKERRQ(ierr); 1189 PetscFunctionReturn(0); 1190 } 1191 1192 #undef __FUNCT__ 1193 #define __FUNCT__ "MatElementalHermitianGenDefEig_Elemental" 1194 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) 1195 { 1196 PetscErrorCode ierr; 1197 Mat_Elemental *a=(Mat_Elemental*)A->data,*b=(Mat_Elemental*)B->data,*x; 1198 MPI_Comm comm; 1199 Mat EVAL; 1200 Mat_Elemental *e; 1201 1202 PetscFunctionBegin; 1203 /* Compute eigenvalues and eigenvectors */ 1204 El::DistMatrix<PetscElemScalar,El::VR,El::STAR> w( *a->grid ); /* holding eigenvalues */ 1205 El::DistMatrix<PetscElemScalar> X( *a->grid ); /* holding eigenvectors */ 1206 El::HermitianGenDefEig(eigtype,uplo,*a->emat,*b->emat,w,X,sort,subset,ctrl); 1207 /* El::Print(w, "Eigenvalues"); */ 1208 1209 /* Wrap w and X into PETSc's MATMATELEMENTAL matrices */ 1210 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1211 ierr = MatCreate(comm,evec);CHKERRQ(ierr); 1212 ierr = MatSetSizes(*evec,PETSC_DECIDE,PETSC_DECIDE,X.Height(),X.Width());CHKERRQ(ierr); 1213 ierr = MatSetType(*evec,MATELEMENTAL);CHKERRQ(ierr); 1214 ierr = MatSetFromOptions(*evec);CHKERRQ(ierr); 1215 ierr = MatSetUp(*evec);CHKERRQ(ierr); 1216 ierr = MatAssemblyBegin(*evec,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1217 ierr = MatAssemblyEnd(*evec,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1218 1219 x = (Mat_Elemental*)(*evec)->data; 1220 //delete x->emat; //-- memory leak??? 1221 *x->emat = X; 1222 1223 ierr = MatCreate(comm,&EVAL);CHKERRQ(ierr); 1224 ierr = MatSetSizes(EVAL,PETSC_DECIDE,PETSC_DECIDE,w.Height(),w.Width());CHKERRQ(ierr); 1225 ierr = MatSetType(EVAL,MATELEMENTAL);CHKERRQ(ierr); 1226 ierr = MatSetFromOptions(EVAL);CHKERRQ(ierr); 1227 ierr = MatSetUp(EVAL);CHKERRQ(ierr); 1228 ierr = MatAssemblyBegin(EVAL,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1229 ierr = MatAssemblyEnd(EVAL,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1230 e = (Mat_Elemental*)EVAL->data; 1231 *e->emat = w; //-- memory leak??? 1232 *evals = EVAL; 1233 1234 #if defined(MV) 1235 /* Test correctness norm = || - A*X + B*X*w || */ 1236 { 1237 PetscElemScalar alpha,beta; 1238 El::DistMatrix<PetscElemScalar> Y(*a->grid); //tmp matrix 1239 alpha = 1.0; beta=0.0; 1240 El::Gemm(El::NORMAL,El::NORMAL,alpha,*b->emat,X,beta,Y); //Y = B*X 1241 El::DiagonalScale(El::RIGHT,El::NORMAL, w, Y); //Y = Y*w 1242 alpha = -1.0; beta=1.0; 1243 El::Gemm(El::NORMAL,El::NORMAL,alpha,*a->emat,X,beta,Y); //Y = - A*X + B*X*w 1244 1245 PetscElemScalar norm = El::FrobeniusNorm(Y); 1246 if ((*a->grid).Rank()==0) printf(" norm (- A*X + B*X*w) = %g\n",norm); 1247 } 1248 1249 { 1250 PetscMPIInt rank; 1251 ierr = MPI_Comm_rank(comm,&rank); 1252 printf("w: [%d] [%d, %d %d] %d; X: %d %d\n",rank,w.DistRank(),w.ColRank(),w.RowRank(),w.LocalHeight(),X.LocalHeight(),X.LocalWidth()); 1253 } 1254 #endif 1255 PetscFunctionReturn(0); 1256 } 1257 1258 #undef __FUNCT__ 1259 #define __FUNCT__ "MatElementalHermitianGenDefEig" 1260 /*@ 1261 MatElementalHermitianGenDefEig - Compute the set of eigenvalues of the Hermitian-definite matrix pencil determined by the subset structure 1262 1263 Logically Collective on Mat 1264 1265 Level: beginner 1266 1267 References: Elemental Users' Guide 1268 @*/ 1269 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) 1270 { 1271 PetscErrorCode ierr; 1272 1273 PetscFunctionBegin; 1274 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); 1275 PetscFunctionReturn(0); 1276 } 1277 1278 /* -------------------------------------------------------------------*/ 1279 static struct _MatOps MatOps_Values = { 1280 MatSetValues_Elemental, 1281 0, 1282 0, 1283 MatMult_Elemental, 1284 /* 4*/ MatMultAdd_Elemental, 1285 MatMultTranspose_Elemental, 1286 MatMultTransposeAdd_Elemental, 1287 MatSolve_Elemental, 1288 MatSolveAdd_Elemental, 1289 0, 1290 /*10*/ 0, 1291 MatLUFactor_Elemental, 1292 MatCholeskyFactor_Elemental, 1293 0, 1294 MatTranspose_Elemental, 1295 /*15*/ MatGetInfo_Elemental, 1296 0, 1297 MatGetDiagonal_Elemental, 1298 MatDiagonalScale_Elemental, 1299 MatNorm_Elemental, 1300 /*20*/ MatAssemblyBegin_Elemental, 1301 MatAssemblyEnd_Elemental, 1302 MatSetOption_Elemental, 1303 MatZeroEntries_Elemental, 1304 /*24*/ 0, 1305 MatLUFactorSymbolic_Elemental, 1306 MatLUFactorNumeric_Elemental, 1307 MatCholeskyFactorSymbolic_Elemental, 1308 MatCholeskyFactorNumeric_Elemental, 1309 /*29*/ MatSetUp_Elemental, 1310 0, 1311 0, 1312 0, 1313 0, 1314 /*34*/ MatDuplicate_Elemental, 1315 0, 1316 0, 1317 0, 1318 0, 1319 /*39*/ MatAXPY_Elemental, 1320 0, 1321 0, 1322 0, 1323 MatCopy_Elemental, 1324 /*44*/ 0, 1325 MatScale_Elemental, 1326 MatShift_Basic, 1327 0, 1328 0, 1329 /*49*/ 0, 1330 0, 1331 0, 1332 0, 1333 0, 1334 /*54*/ 0, 1335 0, 1336 0, 1337 0, 1338 0, 1339 /*59*/ 0, 1340 MatDestroy_Elemental, 1341 MatView_Elemental, 1342 0, 1343 0, 1344 /*64*/ 0, 1345 0, 1346 0, 1347 0, 1348 0, 1349 /*69*/ 0, 1350 0, 1351 MatConvert_Elemental_Dense, 1352 0, 1353 0, 1354 /*74*/ 0, 1355 0, 1356 0, 1357 0, 1358 0, 1359 /*79*/ 0, 1360 0, 1361 0, 1362 0, 1363 MatLoad_Elemental, 1364 /*84*/ 0, 1365 0, 1366 0, 1367 0, 1368 0, 1369 /*89*/ MatMatMult_Elemental, 1370 MatMatMultSymbolic_Elemental, 1371 MatMatMultNumeric_Elemental, 1372 0, 1373 0, 1374 /*94*/ 0, 1375 MatMatTransposeMult_Elemental, 1376 MatMatTransposeMultSymbolic_Elemental, 1377 MatMatTransposeMultNumeric_Elemental, 1378 0, 1379 /*99*/ 0, 1380 0, 1381 0, 1382 MatConjugate_Elemental, 1383 0, 1384 /*104*/0, 1385 0, 1386 0, 1387 0, 1388 0, 1389 /*109*/MatMatSolve_Elemental, 1390 0, 1391 0, 1392 0, 1393 0, 1394 /*114*/0, 1395 0, 1396 0, 1397 0, 1398 0, 1399 /*119*/0, 1400 MatHermitianTranspose_Elemental, 1401 0, 1402 0, 1403 0, 1404 /*124*/0, 1405 0, 1406 0, 1407 0, 1408 0, 1409 /*129*/0, 1410 0, 1411 0, 1412 0, 1413 0, 1414 /*134*/0, 1415 0, 1416 0, 1417 0, 1418 0 1419 }; 1420 1421 /*MC 1422 MATELEMENTAL = "elemental" - A matrix type for dense matrices using the Elemental package 1423 1424 Use ./configure --download-elemental to install PETSc to use Elemental 1425 1426 Use -pc_type lu -pc_factor_mat_solver_package elemental to us this direct solver 1427 1428 Options Database Keys: 1429 + -mat_type elemental - sets the matrix type to "elemental" during a call to MatSetFromOptions() 1430 - -mat_elemental_grid_height - sets Grid Height for 2D cyclic ordering of internal matrix 1431 1432 Level: beginner 1433 1434 .seealso: MATDENSE 1435 M*/ 1436 1437 #undef __FUNCT__ 1438 #define __FUNCT__ "MatCreate_Elemental" 1439 PETSC_EXTERN PetscErrorCode MatCreate_Elemental(Mat A) 1440 { 1441 Mat_Elemental *a; 1442 PetscErrorCode ierr; 1443 PetscBool flg,flg1; 1444 Mat_Elemental_Grid *commgrid; 1445 MPI_Comm icomm; 1446 PetscInt optv1; 1447 1448 PetscFunctionBegin; 1449 ierr = PetscElementalInitializePackage();CHKERRQ(ierr); 1450 ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1451 A->insertmode = NOT_SET_VALUES; 1452 1453 ierr = PetscNewLog(A,&a);CHKERRQ(ierr); 1454 A->data = (void*)a; 1455 1456 /* Set up the elemental matrix */ 1457 El::mpi::Comm cxxcomm(PetscObjectComm((PetscObject)A)); 1458 1459 /* Grid needs to be shared between multiple Mats on the same communicator, implement by attribute caching on the MPI_Comm */ 1460 if (Petsc_Elemental_keyval == MPI_KEYVAL_INVALID) { 1461 ierr = MPI_Keyval_create(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Elemental_keyval,(void*)0); 1462 /* ierr = MPI_Comm_create_Keyval(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Elemental_keyval,(void*)0); -- new version? */ 1463 } 1464 ierr = PetscCommDuplicate(cxxcomm.comm,&icomm,NULL);CHKERRQ(ierr); 1465 ierr = MPI_Attr_get(icomm,Petsc_Elemental_keyval,(void**)&commgrid,(int*)&flg);CHKERRQ(ierr); 1466 if (!flg) { 1467 ierr = PetscNewLog(A,&commgrid);CHKERRQ(ierr); 1468 1469 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"Elemental Options","Mat");CHKERRQ(ierr); 1470 /* displayed default grid sizes (CommSize,1) are set by us arbitrarily until El::Grid() is called */ 1471 ierr = PetscOptionsInt("-mat_elemental_grid_height","Grid Height","None",El::mpi::Size(cxxcomm),&optv1,&flg1);CHKERRQ(ierr); 1472 if (flg1) { 1473 if (El::mpi::Size(cxxcomm) % optv1 != 0) { 1474 SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Grid Height %D must evenly divide CommSize %D",optv1,(PetscInt)El::mpi::Size(cxxcomm)); 1475 } 1476 commgrid->grid = new El::Grid(cxxcomm,optv1); /* use user-provided grid height */ 1477 } else { 1478 commgrid->grid = new El::Grid(cxxcomm); /* use Elemental default grid sizes */ 1479 /* printf("new commgrid->grid = %p\n",commgrid->grid); -- memory leak revealed by valgrind? */ 1480 } 1481 commgrid->grid_refct = 1; 1482 ierr = MPI_Attr_put(icomm,Petsc_Elemental_keyval,(void*)commgrid);CHKERRQ(ierr); 1483 PetscOptionsEnd(); 1484 } else { 1485 commgrid->grid_refct++; 1486 } 1487 ierr = PetscCommDestroy(&icomm);CHKERRQ(ierr); 1488 a->grid = commgrid->grid; 1489 a->emat = new El::DistMatrix<PetscElemScalar>(*a->grid); 1490 a->pivot = new El::DistMatrix<PetscInt,El::VC,El::STAR>(*a->grid); 1491 a->roworiented = PETSC_TRUE; 1492 1493 ierr = PetscObjectComposeFunction((PetscObject)A,"MatGetOwnershipIS_C",MatGetOwnershipIS_Elemental);CHKERRQ(ierr); 1494 ierr = PetscObjectComposeFunction((PetscObject)A,"MatElementalHermitianGenDefEig_C",MatElementalHermitianGenDefEig_Elemental);CHKERRQ(ierr); 1495 1496 ierr = PetscObjectChangeTypeName((PetscObject)A,MATELEMENTAL);CHKERRQ(ierr); 1497 PetscFunctionReturn(0); 1498 } 1499