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