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