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