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