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