1 /* 2 GAMG geometric-algebric multiogrid PC - Mark Adams 2011 3 */ 4 5 #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 6 #include <petsc-private/kspimpl.h> 7 8 #include <petscblaslapack.h> 9 10 typedef struct { 11 PetscInt nsmooths; 12 PetscBool sym_graph; 13 PetscBool square_graph; 14 } PC_GAMG_AGG; 15 16 #undef __FUNCT__ 17 #define __FUNCT__ "PCGAMGSetNSmooths" 18 /*@ 19 PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) 20 21 Not Collective on PC 22 23 Input Parameters: 24 . pc - the preconditioner context 25 26 Options Database Key: 27 . -pc_gamg_agg_nsmooths 28 29 Level: intermediate 30 31 Concepts: Aggregation AMG preconditioner 32 33 .seealso: () 34 @*/ 35 PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n) 36 { 37 PetscErrorCode ierr; 38 39 PetscFunctionBegin; 40 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 41 ierr = PetscTryMethod(pc,"PCGAMGSetNSmooths_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 42 PetscFunctionReturn(0); 43 } 44 45 #undef __FUNCT__ 46 #define __FUNCT__ "PCGAMGSetNSmooths_GAMG" 47 PetscErrorCode PCGAMGSetNSmooths_GAMG(PC pc, PetscInt n) 48 { 49 PC_MG *mg = (PC_MG*)pc->data; 50 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 51 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 52 53 PetscFunctionBegin; 54 pc_gamg_agg->nsmooths = n; 55 PetscFunctionReturn(0); 56 } 57 58 #undef __FUNCT__ 59 #define __FUNCT__ "PCGAMGSetSymGraph" 60 /*@ 61 PCGAMGSetSymGraph - 62 63 Not Collective on PC 64 65 Input Parameters: 66 . pc - the preconditioner context 67 68 Options Database Key: 69 . -pc_gamg_sym_graph 70 71 Level: intermediate 72 73 Concepts: Aggregation AMG preconditioner 74 75 .seealso: () 76 @*/ 77 PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n) 78 { 79 PetscErrorCode ierr; 80 81 PetscFunctionBegin; 82 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 83 ierr = PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 84 PetscFunctionReturn(0); 85 } 86 87 #undef __FUNCT__ 88 #define __FUNCT__ "PCGAMGSetSymGraph_GAMG" 89 PetscErrorCode PCGAMGSetSymGraph_GAMG(PC pc, PetscBool n) 90 { 91 PC_MG *mg = (PC_MG*)pc->data; 92 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 93 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 94 95 PetscFunctionBegin; 96 pc_gamg_agg->sym_graph = n; 97 PetscFunctionReturn(0); 98 } 99 100 #undef __FUNCT__ 101 #define __FUNCT__ "PCGAMGSetSquareGraph" 102 /*@ 103 PCGAMGSetSquareGraph - 104 105 Not Collective on PC 106 107 Input Parameters: 108 . pc - the preconditioner context 109 110 Options Database Key: 111 . -pc_gamg_square_graph 112 113 Level: intermediate 114 115 Concepts: Aggregation AMG preconditioner 116 117 .seealso: () 118 @*/ 119 PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscBool n) 120 { 121 PetscErrorCode ierr; 122 123 PetscFunctionBegin; 124 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 125 ierr = PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 126 PetscFunctionReturn(0); 127 } 128 129 #undef __FUNCT__ 130 #define __FUNCT__ "PCGAMGSetSquareGraph_GAMG" 131 PetscErrorCode PCGAMGSetSquareGraph_GAMG(PC pc, PetscBool n) 132 { 133 PC_MG *mg = (PC_MG*)pc->data; 134 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 135 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 136 137 PetscFunctionBegin; 138 pc_gamg_agg->square_graph = n; 139 PetscFunctionReturn(0); 140 } 141 142 /* -------------------------------------------------------------------------- */ 143 /* 144 PCSetFromOptions_GAMG_AGG 145 146 Input Parameter: 147 . pc - 148 */ 149 #undef __FUNCT__ 150 #define __FUNCT__ "PCSetFromOptions_GAMG_AGG" 151 PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc) 152 { 153 PetscErrorCode ierr; 154 PC_MG *mg = (PC_MG*)pc->data; 155 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 156 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 157 PetscBool flag; 158 159 PetscFunctionBegin; 160 161 ierr = PetscOptionsHead("GAMG-AGG options");CHKERRQ(ierr); 162 { 163 /* -pc_gamg_agg_nsmooths */ 164 pc_gamg_agg->nsmooths = 1; 165 166 ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", 167 "smoothing steps for smoothed aggregation, usually 1 (0)", 168 "PCGAMGSetNSmooths", 169 pc_gamg_agg->nsmooths, 170 &pc_gamg_agg->nsmooths, 171 &flag);CHKERRQ(ierr); 172 173 /* -pc_gamg_sym_graph */ 174 pc_gamg_agg->sym_graph = PETSC_FALSE; 175 176 ierr = PetscOptionsBool("-pc_gamg_sym_graph", 177 "Set for asymmetric matrices", 178 "PCGAMGSetSymGraph", 179 pc_gamg_agg->sym_graph, 180 &pc_gamg_agg->sym_graph, 181 &flag);CHKERRQ(ierr); 182 183 /* -pc_gamg_square_graph */ 184 pc_gamg_agg->square_graph = PETSC_TRUE; 185 186 ierr = PetscOptionsBool("-pc_gamg_square_graph", 187 "For faster coarsening and lower coarse grid complexity", 188 "PCGAMGSetSquareGraph", 189 pc_gamg_agg->square_graph, 190 &pc_gamg_agg->square_graph, 191 &flag);CHKERRQ(ierr); 192 } 193 ierr = PetscOptionsTail();CHKERRQ(ierr); 194 PetscFunctionReturn(0); 195 } 196 197 /* -------------------------------------------------------------------------- */ 198 /* 199 PCDestroy_AGG 200 201 Input Parameter: 202 . pc - 203 */ 204 #undef __FUNCT__ 205 #define __FUNCT__ "PCDestroy_GAMG_AGG" 206 PetscErrorCode PCDestroy_GAMG_AGG(PC pc) 207 { 208 PetscErrorCode ierr; 209 PC_MG *mg = (PC_MG*)pc->data; 210 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 211 212 PetscFunctionBegin; 213 ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr); 214 PetscFunctionReturn(0); 215 } 216 217 /* -------------------------------------------------------------------------- */ 218 /* 219 PCSetCoordinates_AGG 220 - collective 221 222 Input Parameter: 223 . pc - the preconditioner context 224 . ndm - dimesion of data (used for dof/vertex for Stokes) 225 . a_nloc - number of vertices local 226 . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...} 227 */ 228 229 #undef __FUNCT__ 230 #define __FUNCT__ "PCSetCoordinates_AGG" 231 static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 232 { 233 PC_MG *mg = (PC_MG*)pc->data; 234 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 235 PetscErrorCode ierr; 236 PetscInt arrsz,kk,ii,jj,nloc,ndatarows,ndf; 237 Mat mat = pc->pmat; 238 239 PetscFunctionBegin; 240 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 241 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 242 nloc = a_nloc; 243 244 /* SA: null space vectors */ 245 ierr = MatGetBlockSize(mat, &ndf);CHKERRQ(ierr); /* this does not work for Stokes */ 246 if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */ 247 else if (coords) { 248 if (ndm > ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %d > block size %d",ndm,ndf); 249 pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */ 250 if (ndm != ndf) { 251 if (pc_gamg->data_cell_cols != ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Don't know how to create null space for ndm=%d, ndf=%d. Use MatSetNearNullSpace.",ndm,ndf); 252 } 253 } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */ 254 pc_gamg->data_cell_rows = ndatarows = ndf; 255 if (pc_gamg->data_cell_cols <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->data_cell_cols %D <= 0",pc_gamg->data_cell_cols); 256 arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols; 257 258 /* create data - syntactic sugar that should be refactored at some point */ 259 if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) { 260 ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 261 ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->data);CHKERRQ(ierr); 262 } 263 /* copy data in - column oriented */ 264 for (kk=0; kk<nloc; kk++) { 265 const PetscInt M = nloc*pc_gamg->data_cell_rows; /* stride into data */ 266 PetscReal *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */ 267 if (pc_gamg->data_cell_cols==1) *data = 1.0; 268 else { 269 /* translational modes */ 270 for (ii=0;ii<ndatarows;ii++) { 271 for (jj=0;jj<ndatarows;jj++) { 272 if (ii==jj)data[ii*M + jj] = 1.0; 273 else data[ii*M + jj] = 0.0; 274 } 275 } 276 277 /* rotational modes */ 278 if (coords) { 279 if (ndm == 2) { 280 data += 2*M; 281 data[0] = -coords[2*kk+1]; 282 data[1] = coords[2*kk]; 283 } else { 284 data += 3*M; 285 data[0] = 0.0; data[M+0] = coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1]; 286 data[1] = -coords[3*kk+2]; data[M+1] = 0.0; data[2*M+1] = coords[3*kk]; 287 data[2] = coords[3*kk+1]; data[M+2] = -coords[3*kk]; data[2*M+2] = 0.0; 288 } 289 } 290 } 291 } 292 293 pc_gamg->data_sz = arrsz; 294 PetscFunctionReturn(0); 295 } 296 297 typedef PetscInt NState; 298 static const NState NOT_DONE=-2; 299 static const NState DELETED =-1; 300 static const NState REMOVED =-3; 301 #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED) 302 303 /* -------------------------------------------------------------------------- */ 304 /* 305 smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific 306 - AGG-MG specific: clears singletons out of 'selected_2' 307 308 Input Parameter: 309 . Gmat_2 - glabal matrix of graph (data not defined) 310 . Gmat_1 - base graph to grab with 311 Input/Output Parameter: 312 . aggs_2 - linked list of aggs with gids) 313 */ 314 #undef __FUNCT__ 315 #define __FUNCT__ "smoothAggs" 316 static PetscErrorCode smoothAggs(const Mat Gmat_2, /* base (squared) graph */ 317 const Mat Gmat_1, /* base graph */ 318 /* const IS selected_2, [nselected local] selected vertices */ 319 PetscCoarsenData *aggs_2) /* [nselected local] global ID of aggregate */ 320 { 321 PetscErrorCode ierr; 322 PetscBool isMPI; 323 Mat_SeqAIJ *matA_1, *matB_1=0, *matA_2, *matB_2=0; 324 MPI_Comm comm; 325 PetscMPIInt rank,size; 326 PetscInt lid,*ii,*idx,ix,Iend,my0,kk,n,j; 327 Mat_MPIAIJ *mpimat_2 = 0, *mpimat_1=0; 328 const PetscInt nloc = Gmat_2->rmap->n; 329 PetscScalar *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid; 330 PetscInt *lid_cprowID_1; 331 NState *lid_state; 332 Vec ghost_par_orig2; 333 334 PetscFunctionBegin; 335 ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr); 336 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 337 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 338 ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr); 339 340 if (PETSC_FALSE) { 341 PetscViewer viewer; char fname[32]; static int llev=0; 342 sprintf(fname,"Gmat2_%d.m",llev++); 343 PetscViewerASCIIOpen(comm,fname,&viewer); 344 ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 345 ierr = MatView(Gmat_2, viewer);CHKERRQ(ierr); 346 ierr = PetscViewerDestroy(&viewer); 347 } 348 349 /* get submatrices */ 350 ierr = PetscObjectTypeCompare((PetscObject)Gmat_1, MATMPIAIJ, &isMPI);CHKERRQ(ierr); 351 if (isMPI) { 352 /* grab matrix objects */ 353 mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data; 354 mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data; 355 matA_1 = (Mat_SeqAIJ*)mpimat_1->A->data; 356 matB_1 = (Mat_SeqAIJ*)mpimat_1->B->data; 357 matA_2 = (Mat_SeqAIJ*)mpimat_2->A->data; 358 matB_2 = (Mat_SeqAIJ*)mpimat_2->B->data; 359 360 /* force compressed row storage for B matrix in AuxMat */ 361 ierr = MatCheckCompressedRow(mpimat_1->B,matB_1->nonzerorowcnt,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr); 362 363 ierr = PetscMalloc(nloc*sizeof(PetscInt), &lid_cprowID_1);CHKERRQ(ierr); 364 for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1; 365 for (ix=0; ix<matB_1->compressedrow.nrows; ix++) { 366 PetscInt lid = matB_1->compressedrow.rindex[ix]; 367 lid_cprowID_1[lid] = ix; 368 } 369 } else { 370 matA_1 = (Mat_SeqAIJ*)Gmat_1->data; 371 matA_2 = (Mat_SeqAIJ*)Gmat_2->data; 372 lid_cprowID_1 = NULL; 373 } 374 if (!(matA_1 && !matA_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_1 && !matA_1->compressedrow.use)"); 375 if (!(matB_1==0 || matB_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_1==0 || matB_1->compressedrow.use)"); 376 if (!(matA_2 && !matA_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_2 && !matA_2->compressedrow.use)"); 377 if (!(matB_2==0 || matB_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_2==0 || matB_2->compressedrow.use)"); 378 379 /* get state of locals and selected gid for deleted */ 380 ierr = PetscMalloc(nloc*sizeof(NState), &lid_state);CHKERRQ(ierr); 381 ierr = PetscMalloc(nloc*sizeof(PetscScalar), &lid_parent_gid);CHKERRQ(ierr); 382 for (lid = 0; lid < nloc; lid++) { 383 lid_parent_gid[lid] = -1.0; 384 lid_state[lid] = DELETED; 385 } 386 387 /* set lid_state */ 388 for (lid = 0; lid < nloc; lid++) { 389 PetscCDPos pos; 390 ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 391 if (pos) { 392 PetscInt gid1; 393 394 ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr); 395 if (gid1 != lid+my0) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %D != lid %D + my0 %D",gid1,lid,my0); 396 lid_state[lid] = gid1; 397 } 398 } 399 400 /* map local to selected local, DELETED means a ghost owns it */ 401 for (lid=kk=0; lid<nloc; lid++) { 402 NState state = lid_state[lid]; 403 if (IS_SELECTED(state)) { 404 PetscCDPos pos; 405 ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 406 while (pos) { 407 PetscInt gid1; 408 ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr); 409 ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr); 410 411 if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0); 412 } 413 } 414 } 415 /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */ 416 if (isMPI) { 417 Vec tempVec; 418 /* get 'cpcol_1_state' */ 419 ierr = MatGetVecs(Gmat_1, &tempVec, 0);CHKERRQ(ierr); 420 for (kk=0,j=my0; kk<nloc; kk++,j++) { 421 PetscScalar v = (PetscScalar)lid_state[kk]; 422 ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 423 } 424 ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 425 ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 426 ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 427 ierr = VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 428 ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr); 429 /* get 'cpcol_2_state' */ 430 ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 431 ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 432 ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr); 433 /* get 'cpcol_2_par_orig' */ 434 for (kk=0,j=my0; kk<nloc; kk++,j++) { 435 PetscScalar v = (PetscScalar)lid_parent_gid[kk]; 436 ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 437 } 438 ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 439 ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 440 ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr); 441 ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 442 ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 443 ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr); 444 445 ierr = VecDestroy(&tempVec);CHKERRQ(ierr); 446 } /* ismpi */ 447 448 /* doit */ 449 for (lid=0; lid<nloc; lid++) { 450 NState state = lid_state[lid]; 451 if (IS_SELECTED(state)) { 452 /* steal locals */ 453 ii = matA_1->i; n = ii[lid+1] - ii[lid]; 454 idx = matA_1->j + ii[lid]; 455 for (j=0; j<n; j++) { 456 PetscInt lidj = idx[j], sgid; 457 NState statej = lid_state[lidj]; 458 if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */ 459 lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */ 460 if (sgid >= my0 && sgid < Iend) { /* I'm stealing this local from a local sgid */ 461 PetscInt hav=0,slid=sgid-my0,gidj=lidj+my0; 462 PetscCDPos pos,last=NULL; 463 /* looking for local from local so id_llist_2 works */ 464 ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr); 465 while (pos) { 466 PetscInt gid; 467 ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr); 468 if (gid == gidj) { 469 if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null"); 470 ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr); 471 ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr); 472 hav = 1; 473 break; 474 } else last = pos; 475 476 ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr); 477 } 478 if (hav!=1) { 479 if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 480 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav); 481 } 482 } else { /* I'm stealing this local, owned by a ghost */ 483 if (sgid != -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Have un-symmetric graph (apparently). Use '-pc_gamg_sym_graph true' to symetrize the graph or '-pc_gamg_threshold 0.0' if the matrix is structurally symmetric."); 484 ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr); 485 } 486 } 487 } /* local neighbors */ 488 } else if (state == DELETED && lid_cprowID_1) { 489 PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]); 490 /* see if I have a selected ghost neighbor that will steal me */ 491 if ((ix=lid_cprowID_1[lid]) != -1) { 492 ii = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix]; 493 idx = matB_1->j + ii[ix]; 494 for (j=0; j<n; j++) { 495 PetscInt cpid = idx[j]; 496 NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]); 497 if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */ 498 lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */ 499 if (sgidold>=my0 && sgidold<Iend) { /* this was mine */ 500 PetscInt hav=0,oldslidj=sgidold-my0; 501 PetscCDPos pos,last=NULL; 502 /* remove from 'oldslidj' list */ 503 ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr); 504 while (pos) { 505 PetscInt gid; 506 ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr); 507 if (lid+my0 == gid) { 508 /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove lid from oldslidj list *\/ */ 509 if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null"); 510 ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr); 511 /* ghost (PetscScalar)statej will add this later */ 512 hav = 1; 513 break; 514 } else last = pos; 515 516 ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr); 517 } 518 if (hav!=1) { 519 if (hav==0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 520 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav); 521 } 522 } else { 523 /* ghosts remove this later */ 524 } 525 } 526 } 527 } 528 } /* selected/deleted */ 529 } /* node loop */ 530 531 if (isMPI) { 532 PetscScalar *cpcol_2_parent,*cpcol_2_gid; 533 Vec tempVec,ghostgids2,ghostparents2; 534 PetscInt cpid,nghost_2; 535 GAMGHashTable gid_cpid; 536 537 ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr); 538 ierr = MatGetVecs(Gmat_2, &tempVec, 0);CHKERRQ(ierr); 539 540 /* get 'cpcol_2_parent' */ 541 for (kk=0,j=my0; kk<nloc; kk++,j++) { 542 ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr); 543 } 544 ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 545 ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 546 ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr); 547 ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 548 ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 549 ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr); 550 551 /* get 'cpcol_2_gid' */ 552 for (kk=0,j=my0; kk<nloc; kk++,j++) { 553 PetscScalar v = (PetscScalar)j; 554 ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 555 } 556 ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 557 ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 558 ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr); 559 ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 560 ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 561 ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr); 562 563 ierr = VecDestroy(&tempVec);CHKERRQ(ierr); 564 565 /* look for deleted ghosts and add to table */ 566 ierr = GAMGTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr); 567 for (cpid = 0; cpid < nghost_2; cpid++) { 568 NState state = (NState)PetscRealPart(cpcol_2_state[cpid]); 569 if (state==DELETED) { 570 PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 571 PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]); 572 if (sgid_old == -1 && sgid_new != -1) { 573 PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 574 ierr = GAMGTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr); 575 } 576 } 577 } 578 579 /* look for deleted ghosts and see if they moved - remove it */ 580 for (lid=0; lid<nloc; lid++) { 581 NState state = lid_state[lid]; 582 if (IS_SELECTED(state)) { 583 PetscCDPos pos,last=NULL; 584 /* look for deleted ghosts and see if they moved */ 585 ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 586 while (pos) { 587 PetscInt gid; 588 ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr); 589 590 if (gid < my0 || gid >= Iend) { 591 ierr = GAMGTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr); 592 if (cpid != -1) { 593 /* a moved ghost - */ 594 /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove 'flid' from list *\/ */ 595 ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr); 596 } else last = pos; 597 } else last = pos; 598 599 ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr); 600 } /* loop over list of deleted */ 601 } /* selected */ 602 } 603 ierr = GAMGTableDestroy(&gid_cpid);CHKERRQ(ierr); 604 605 /* look at ghosts, see if they changed - and it */ 606 for (cpid = 0; cpid < nghost_2; cpid++) { 607 PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 608 if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */ 609 PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 610 PetscInt slid_new=sgid_new-my0,hav=0; 611 PetscCDPos pos; 612 /* search for this gid to see if I have it */ 613 ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr); 614 while (pos) { 615 PetscInt gidj; 616 ierr = PetscLLNGetID(pos, &gidj);CHKERRQ(ierr); 617 ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr); 618 619 if (gidj == gid) { hav = 1; break; } 620 } 621 if (hav != 1) { 622 /* insert 'flidj' into head of llist */ 623 ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr); 624 } 625 } 626 } 627 628 ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr); 629 ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr); 630 ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr); 631 ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr); 632 ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr); 633 ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr); 634 ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr); 635 ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr); 636 } 637 638 ierr = PetscFree(lid_parent_gid);CHKERRQ(ierr); 639 ierr = PetscFree(lid_state);CHKERRQ(ierr); 640 PetscFunctionReturn(0); 641 } 642 643 /* -------------------------------------------------------------------------- */ 644 /* 645 PCSetData_AGG - called if data is not set with PCSetCoordinates. 646 Looks in Mat for near null space. 647 Does not work for Stokes 648 649 Input Parameter: 650 . pc - 651 . a_A - matrix to get (near) null space out of. 652 */ 653 #undef __FUNCT__ 654 #define __FUNCT__ "PCSetData_AGG" 655 PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) 656 { 657 PetscErrorCode ierr; 658 PC_MG *mg = (PC_MG*)pc->data; 659 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 660 MatNullSpace mnull; 661 662 PetscFunctionBegin; 663 ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr); 664 if (!mnull) { 665 PetscInt bs,NN,MM; 666 ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); 667 ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr); 668 if (MM % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %D must be divisible by bs %D",MM,bs); 669 ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, NULL);CHKERRQ(ierr); 670 } else { 671 PetscReal *nullvec; 672 PetscBool has_const; 673 PetscInt i,j,mlocal,nvec,bs; 674 const Vec *vecs; const PetscScalar *v; 675 ierr = MatGetLocalSize(a_A,&mlocal,NULL);CHKERRQ(ierr); 676 ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr); 677 pc_gamg->data_sz = (nvec+!!has_const)*mlocal; 678 ierr = PetscMalloc((nvec+!!has_const)*mlocal*sizeof(*nullvec),&nullvec);CHKERRQ(ierr); 679 if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0; 680 for (i=0; i<nvec; i++) { 681 ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr); 682 for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]); 683 ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr); 684 } 685 pc_gamg->data = nullvec; 686 pc_gamg->data_cell_cols = (nvec+!!has_const); 687 688 ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); 689 690 pc_gamg->data_cell_rows = bs; 691 } 692 PetscFunctionReturn(0); 693 } 694 695 /* -------------------------------------------------------------------------- */ 696 /* 697 formProl0 698 699 Input Parameter: 700 . agg_llists - list of arrays with aggregates 701 . bs - block size 702 . nSAvec - column bs of new P 703 . my0crs - global index of start of locals 704 . data_stride - bs*(nloc nodes + ghost nodes) 705 . data_in[data_stride*nSAvec] - local data on fine grid 706 . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 707 Output Parameter: 708 . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 709 . a_Prol - prolongation operator 710 */ 711 #undef __FUNCT__ 712 #define __FUNCT__ "formProl0" 713 static PetscErrorCode formProl0(const PetscCoarsenData *agg_llists, /* list from selected vertices of aggregate unselected vertices */ 714 const PetscInt bs, /* (row) block size */ 715 const PetscInt nSAvec, /* column bs */ 716 const PetscInt my0crs, /* global index of start of locals */ 717 const PetscInt data_stride, /* (nloc+nghost)*bs */ 718 PetscReal data_in[], /* [data_stride][nSAvec] */ 719 const PetscInt flid_fgid[], /* [data_stride/bs] */ 720 PetscReal **a_data_out, 721 Mat a_Prol) /* prolongation operator (output)*/ 722 { 723 PetscErrorCode ierr; 724 PetscInt Istart,my0,Iend,nloc,clid,flid,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride; 725 MPI_Comm comm; 726 PetscMPIInt rank, size; 727 PetscReal *out_data; 728 PetscCDPos pos; 729 GAMGHashTable fgid_flid; 730 731 /* #define OUT_AGGS */ 732 #if defined(OUT_AGGS) 733 static PetscInt llev = 0; char fname[32]; FILE *file = NULL; PetscInt pM; 734 #endif 735 736 PetscFunctionBegin; 737 ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr); 738 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 739 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 740 ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr); 741 nloc = (Iend-Istart)/bs; my0 = Istart/bs; 742 if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Iend %D - Istart %d must be divisible by bs %D",Iend,Istart,bs); 743 Iend /= bs; 744 nghosts = data_stride/bs - nloc; 745 746 ierr = GAMGTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr); 747 for (kk=0; kk<nghosts; kk++) { 748 ierr = GAMGTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr); 749 } 750 751 #if defined(OUT_AGGS) 752 sprintf(fname,"aggs_%d_%d.m",llev++,rank); 753 if (llev==1) file = fopen(fname,"w"); 754 MatGetSize(a_Prol, &pM, &jj); 755 #endif 756 757 /* count selected -- same as number of cols of P */ 758 for (nSelected=mm=0; mm<nloc; mm++) { 759 PetscBool ise; 760 ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr); 761 if (!ise) nSelected++; 762 } 763 ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr); 764 if ((ii/nSAvec) != my0crs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D != my0crs %D",ii,nSAvec,my0crs); 765 if (nSelected != (jj-ii)/nSAvec) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nSelected %D != (jj %D - ii %D)/nSAvec %D",nSelected,jj,ii,nSAvec); 766 767 /* aloc space for coarse point data (output) */ 768 out_data_stride = nSelected*nSAvec; 769 770 ierr = PetscMalloc(out_data_stride*nSAvec*sizeof(PetscReal), &out_data);CHKERRQ(ierr); 771 for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=PETSC_MAX_REAL; 772 *a_data_out = out_data; /* output - stride nSelected*nSAvec */ 773 774 /* find points and set prolongation */ 775 minsz = 100; 776 ndone = 0; 777 for (mm = clid = 0; mm < nloc; mm++) { 778 ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr); 779 if (jj > 0) { 780 const PetscInt lid = mm, cgid = my0crs + clid; 781 PetscInt cids[100]; /* max bs */ 782 PetscBLASInt asz =jj,M=asz*bs,N=nSAvec,INFO; 783 PetscBLASInt Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs; 784 PetscScalar *qqc,*qqr,*TAU,*WORK; 785 PetscInt *fids; 786 PetscReal *data; 787 /* count agg */ 788 if (asz<minsz) minsz = asz; 789 790 /* get block */ 791 ierr = PetscMalloc((Mdata*N)*sizeof(PetscScalar), &qqc);CHKERRQ(ierr); 792 ierr = PetscMalloc((M*N)*sizeof(PetscScalar), &qqr);CHKERRQ(ierr); 793 ierr = PetscMalloc(N*sizeof(PetscScalar), &TAU);CHKERRQ(ierr); 794 ierr = PetscMalloc(LWORK*sizeof(PetscScalar), &WORK);CHKERRQ(ierr); 795 ierr = PetscMalloc(M*sizeof(PetscInt), &fids);CHKERRQ(ierr); 796 797 aggID = 0; 798 ierr = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr); 799 while (pos) { 800 PetscInt gid1; 801 ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr); 802 ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr); 803 804 if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0; 805 else { 806 ierr = GAMGTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr); 807 if (flid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table"); 808 } 809 /* copy in B_i matrix - column oriented */ 810 data = &data_in[flid*bs]; 811 for (kk = ii = 0; ii < bs; ii++) { 812 for (jj = 0; jj < N; jj++) { 813 PetscReal d = data[jj*data_stride + ii]; 814 qqc[jj*Mdata + aggID*bs + ii] = d; 815 } 816 } 817 #if defined(OUT_AGGS) 818 if (llev==1) { 819 char str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^"; 820 PetscInt MM,pi,pj; 821 str[12] = col[(clid+17*rank)%6]; str[13] = sim[(clid+17*rank)%11]; 822 MM = (PetscInt)(PetscSqrtReal((PetscReal)pM)); 823 pj = gid1/MM; pi = gid1%MM; 824 fprintf(file,str,(double)pi,(double)pj); 825 /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */ 826 } 827 #endif 828 /* set fine IDs */ 829 for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk; 830 831 aggID++; 832 } 833 834 /* pad with zeros */ 835 for (ii = asz*bs; ii < Mdata; ii++) { 836 for (jj = 0; jj < N; jj++, kk++) { 837 qqc[jj*Mdata + ii] = .0; 838 } 839 } 840 841 ndone += aggID; 842 /* QR */ 843 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 844 PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 845 ierr = PetscFPTrapPop();CHKERRQ(ierr); 846 if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error"); 847 /* get R - column oriented - output B_{i+1} */ 848 { 849 PetscReal *data = &out_data[clid*nSAvec]; 850 for (jj = 0; jj < nSAvec; jj++) { 851 for (ii = 0; ii < nSAvec; ii++) { 852 if (data[jj*out_data_stride + ii] != PETSC_MAX_REAL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"data[jj*out_data_stride + ii] != %e",PETSC_MAX_REAL); 853 if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]); 854 else data[jj*out_data_stride + ii] = 0.; 855 } 856 } 857 } 858 859 /* get Q - row oriented */ 860 PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 861 if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO); 862 863 for (ii = 0; ii < M; ii++) { 864 for (jj = 0; jj < N; jj++) { 865 qqr[N*ii + jj] = qqc[jj*Mdata + ii]; 866 } 867 } 868 869 /* add diagonal block of P0 */ 870 for (kk=0; kk<N; kk++) { 871 cids[kk] = N*cgid + kk; /* global col IDs in P0 */ 872 } 873 ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr); 874 875 ierr = PetscFree(qqc);CHKERRQ(ierr); 876 ierr = PetscFree(qqr);CHKERRQ(ierr); 877 ierr = PetscFree(TAU);CHKERRQ(ierr); 878 ierr = PetscFree(WORK);CHKERRQ(ierr); 879 ierr = PetscFree(fids);CHKERRQ(ierr); 880 clid++; 881 } /* coarse agg */ 882 } /* for all fine nodes */ 883 ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 884 ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 885 886 /* ierr = MPI_Allreduce(&ndone, &ii, 1, MPIU_INT, MPIU_SUM, comm); */ 887 /* MatGetSize(a_Prol, &kk, &jj); */ 888 /* ierr = MPI_Allreduce(&minsz, &jj, 1, MPIU_INT, MPIU_MIN, comm); */ 889 /* PetscPrintf(comm," **** [%d]%s %d total done, %d nodes (%d local done), min agg. size = %d\n",rank,__FUNCT__,ii,kk/bs,ndone,jj); */ 890 891 #if defined(OUT_AGGS) 892 if (llev==1) fclose(file); 893 #endif 894 ierr = GAMGTableDestroy(&fgid_flid);CHKERRQ(ierr); 895 PetscFunctionReturn(0); 896 } 897 898 /* -------------------------------------------------------------------------- */ 899 /* 900 PCGAMGgraph_AGG 901 902 Input Parameter: 903 . pc - this 904 . Amat - matrix on this fine level 905 Output Parameter: 906 . a_Gmat - 907 */ 908 #undef __FUNCT__ 909 #define __FUNCT__ "PCGAMGgraph_AGG" 910 PetscErrorCode PCGAMGgraph_AGG(PC pc,const Mat Amat,Mat *a_Gmat) 911 { 912 PetscErrorCode ierr; 913 PC_MG *mg = (PC_MG*)pc->data; 914 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 915 const PetscInt verbose = pc_gamg->verbose; 916 const PetscReal vfilter = pc_gamg->threshold; 917 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 918 PetscMPIInt rank,size; 919 Mat Gmat; 920 MPI_Comm comm; 921 PetscBool /* set,flg , */ symm; 922 923 PetscFunctionBegin; 924 ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 925 #if defined PETSC_USE_LOG 926 ierr = PetscLogEventBegin(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr); 927 #endif 928 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 929 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 930 931 /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */ 932 symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */ 933 934 ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr); 935 ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm, verbose);CHKERRQ(ierr); 936 937 *a_Gmat = Gmat; 938 939 #if defined PETSC_USE_LOG 940 ierr = PetscLogEventEnd(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr); 941 #endif 942 PetscFunctionReturn(0); 943 } 944 945 /* -------------------------------------------------------------------------- */ 946 /* 947 PCGAMGCoarsen_AGG 948 949 Input Parameter: 950 . a_pc - this 951 Input/Output Parameter: 952 . a_Gmat1 - graph on this fine level - coarsening can change this (squares it) 953 Output Parameter: 954 . agg_lists - list of aggregates 955 */ 956 #undef __FUNCT__ 957 #define __FUNCT__ "PCGAMGCoarsen_AGG" 958 PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists) 959 { 960 PetscErrorCode ierr; 961 PC_MG *mg = (PC_MG*)a_pc->data; 962 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 963 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 964 Mat mat,Gmat2, Gmat1 = *a_Gmat1; /* squared graph */ 965 const PetscInt verbose = pc_gamg->verbose; 966 IS perm; 967 PetscInt Ii,nloc,bs,n,m; 968 PetscInt *permute; 969 PetscBool *bIndexSet; 970 MatCoarsen crs; 971 MPI_Comm comm; 972 PetscMPIInt rank,size; 973 974 PetscFunctionBegin; 975 #if defined PETSC_USE_LOG 976 ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr); 977 #endif 978 ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr); 979 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 980 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 981 ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr); 982 ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr); 983 if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs); 984 nloc = n/bs; 985 986 if (pc_gamg_agg->square_graph) { 987 if (verbose > 1) PetscPrintf(comm,"[%d]%s square graph\n",rank,__FUNCT__); 988 /* ierr = MatMatTransposeMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2); */ 989 ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr); 990 if (verbose > 2) { 991 ierr = PetscPrintf(comm,"[%d]%s square graph done\n",rank,__FUNCT__);CHKERRQ(ierr); 992 /* check for symetry */ 993 if (verbose > 4) { 994 995 } 996 } 997 } else Gmat2 = Gmat1; 998 999 /* get MIS aggs */ 1000 /* randomize */ 1001 ierr = PetscMalloc(nloc*sizeof(PetscInt), &permute);CHKERRQ(ierr); 1002 ierr = PetscMalloc(nloc*sizeof(PetscBool), &bIndexSet);CHKERRQ(ierr); 1003 for (Ii = 0; Ii < nloc; Ii++) { 1004 bIndexSet[Ii] = PETSC_FALSE; 1005 permute[Ii] = Ii; 1006 } 1007 srand(1); /* make deterministic */ 1008 for (Ii = 0; Ii < nloc; Ii++) { 1009 PetscInt iSwapIndex = rand()%nloc; 1010 if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 1011 PetscInt iTemp = permute[iSwapIndex]; 1012 permute[iSwapIndex] = permute[Ii]; 1013 permute[Ii] = iTemp; 1014 bIndexSet[iSwapIndex] = PETSC_TRUE; 1015 } 1016 } 1017 ierr = PetscFree(bIndexSet);CHKERRQ(ierr); 1018 1019 if (verbose > 1) PetscPrintf(comm,"[%d]%s coarsen graph\n",rank,__FUNCT__); 1020 1021 ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr); 1022 #if defined PETSC_GAMG_USE_LOG 1023 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 1024 #endif 1025 ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr); 1026 /* ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr); */ 1027 ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr); 1028 ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr); 1029 ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr); 1030 ierr = MatCoarsenSetVerbose(crs, pc_gamg->verbose);CHKERRQ(ierr); 1031 ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr); 1032 ierr = MatCoarsenApply(crs);CHKERRQ(ierr); 1033 ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */ 1034 ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr); 1035 1036 ierr = ISDestroy(&perm);CHKERRQ(ierr); 1037 ierr = PetscFree(permute);CHKERRQ(ierr); 1038 #if defined PETSC_GAMG_USE_LOG 1039 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 1040 #endif 1041 if (verbose > 2) PetscPrintf(comm,"[%d]%s coarsen graph done\n",rank,__FUNCT__); 1042 1043 /* smooth aggs */ 1044 if (Gmat2 != Gmat1) { 1045 const PetscCoarsenData *llist = *agg_lists; 1046 ierr = smoothAggs(Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr); 1047 ierr = MatDestroy(&Gmat1);CHKERRQ(ierr); 1048 *a_Gmat1 = Gmat2; /* output */ 1049 ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr); 1050 if (mat) SETERRQ(comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????"); 1051 } else { 1052 const PetscCoarsenData *llist = *agg_lists; 1053 /* see if we have a matrix that takes precedence (returned from MatCoarsenAppply) */ 1054 ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr); 1055 if (mat) { 1056 ierr = MatDestroy(&Gmat1);CHKERRQ(ierr); 1057 *a_Gmat1 = mat; /* output */ 1058 } 1059 } 1060 #if defined PETSC_USE_LOG 1061 ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr); 1062 #endif 1063 if (verbose > 2) PetscPrintf(comm,"[%d]%s PCGAMGCoarsen_AGG done\n",rank,__FUNCT__); 1064 PetscFunctionReturn(0); 1065 } 1066 1067 /* -------------------------------------------------------------------------- */ 1068 /* 1069 PCGAMGProlongator_AGG 1070 1071 Input Parameter: 1072 . pc - this 1073 . Amat - matrix on this fine level 1074 . Graph - used to get ghost data for nodes in 1075 . agg_lists - list of aggregates 1076 Output Parameter: 1077 . a_P_out - prolongation operator to the next level 1078 */ 1079 #undef __FUNCT__ 1080 #define __FUNCT__ "PCGAMGProlongator_AGG" 1081 PetscErrorCode PCGAMGProlongator_AGG(PC pc,const Mat Amat,const Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out) 1082 { 1083 PC_MG *mg = (PC_MG*)pc->data; 1084 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1085 const PetscInt verbose = pc_gamg->verbose; 1086 const PetscInt data_cols = pc_gamg->data_cell_cols; 1087 PetscErrorCode ierr; 1088 PetscInt Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs; 1089 Mat Prol; 1090 PetscMPIInt rank, size; 1091 MPI_Comm comm; 1092 const PetscInt col_bs = data_cols; 1093 PetscReal *data_w_ghost; 1094 PetscInt myCrs0, nbnodes=0, *flid_fgid; 1095 1096 PetscFunctionBegin; 1097 ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 1098 #if defined PETSC_USE_LOG 1099 ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1100 #endif 1101 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1102 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 1103 ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr); 1104 ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr); 1105 nloc = (Iend-Istart)/bs; my0 = Istart/bs; 1106 if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs); 1107 1108 /* get 'nLocalSelected' */ 1109 for (ii=0, nLocalSelected = 0; ii < nloc; ii++) { 1110 PetscBool ise; 1111 /* filter out singletons 0 or 1? */ 1112 ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr); 1113 if (!ise) nLocalSelected++; 1114 } 1115 1116 /* create prolongator, create P matrix */ 1117 ierr = MatCreate(comm, &Prol);CHKERRQ(ierr); 1118 ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1119 ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr); 1120 ierr = MatSetType(Prol, MATAIJ);CHKERRQ(ierr); 1121 ierr = MatSeqAIJSetPreallocation(Prol, data_cols, NULL);CHKERRQ(ierr); 1122 ierr = MatMPIAIJSetPreallocation(Prol,data_cols, NULL,data_cols, NULL);CHKERRQ(ierr); 1123 /* nloc*bs, nLocalSelected*col_bs, */ 1124 /* PETSC_DETERMINE, PETSC_DETERMINE, */ 1125 /* data_cols, NULL, data_cols, NULL, */ 1126 /* &Prol); */ 1127 1128 /* can get all points "removed" */ 1129 ierr = MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr); 1130 if (ii==0) { 1131 if (verbose > 0) PetscPrintf(comm,"[%d]%s no selected points on coarse grid\n",rank,__FUNCT__); 1132 ierr = MatDestroy(&Prol);CHKERRQ(ierr); 1133 *a_P_out = NULL; /* out */ 1134 PetscFunctionReturn(0); 1135 } 1136 if (verbose > 0) PetscPrintf(comm,"\t\t[%d]%s New grid %d nodes\n",rank,__FUNCT__,ii/col_bs); 1137 ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr); 1138 1139 if ((kk-myCrs0) % col_bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(kk %D -myCrs0 %D) not divisible by col_bs %D",kk,myCrs0,col_bs); 1140 myCrs0 = myCrs0/col_bs; 1141 if ((kk/col_bs-myCrs0) != nLocalSelected) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(kk %D/col_bs %D - myCrs0 %D) != nLocalSelected %D)",kk,col_bs,myCrs0,nLocalSelected); 1142 1143 /* create global vector of data in 'data_w_ghost' */ 1144 #if defined PETSC_GAMG_USE_LOG 1145 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 1146 #endif 1147 if (size > 1) { /* */ 1148 PetscReal *tmp_gdata,*tmp_ldata,*tp2; 1149 ierr = PetscMalloc(nloc*sizeof(PetscReal), &tmp_ldata);CHKERRQ(ierr); 1150 for (jj = 0; jj < data_cols; jj++) { 1151 for (kk = 0; kk < bs; kk++) { 1152 PetscInt ii,stride; 1153 const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk; 1154 for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp; 1155 1156 ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr); 1157 1158 if (jj==0 && kk==0) { /* now I know how many todal nodes - allocate */ 1159 ierr = PetscMalloc(stride*bs*data_cols*sizeof(PetscReal), &data_w_ghost);CHKERRQ(ierr); 1160 nbnodes = bs*stride; 1161 } 1162 tp2 = data_w_ghost + jj*bs*stride + kk; 1163 for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii]; 1164 ierr = PetscFree(tmp_gdata);CHKERRQ(ierr); 1165 } 1166 } 1167 ierr = PetscFree(tmp_ldata);CHKERRQ(ierr); 1168 } else { 1169 nbnodes = bs*nloc; 1170 data_w_ghost = (PetscReal*)pc_gamg->data; 1171 } 1172 1173 /* get P0 */ 1174 if (size > 1) { 1175 PetscReal *fid_glid_loc,*fiddata; 1176 PetscInt stride; 1177 1178 ierr = PetscMalloc(nloc*sizeof(PetscReal), &fid_glid_loc);CHKERRQ(ierr); 1179 for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk); 1180 ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr); 1181 ierr = PetscMalloc(stride*sizeof(PetscInt), &flid_fgid);CHKERRQ(ierr); 1182 for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 1183 ierr = PetscFree(fiddata);CHKERRQ(ierr); 1184 1185 if (stride != nbnodes/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs); 1186 ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr); 1187 } else { 1188 ierr = PetscMalloc(nloc*sizeof(PetscInt), &flid_fgid);CHKERRQ(ierr); 1189 for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk; 1190 } 1191 #if defined PETSC_GAMG_USE_LOG 1192 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 1193 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 1194 #endif 1195 { 1196 PetscReal *data_out = NULL; 1197 ierr = formProl0(agg_lists, bs, data_cols, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr); 1198 ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 1199 1200 pc_gamg->data = data_out; 1201 pc_gamg->data_cell_rows = data_cols; 1202 pc_gamg->data_sz = data_cols*data_cols*nLocalSelected; 1203 } 1204 #if defined PETSC_GAMG_USE_LOG 1205 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 1206 #endif 1207 if (size > 1) ierr = PetscFree(data_w_ghost);CHKERRQ(ierr); 1208 ierr = PetscFree(flid_fgid);CHKERRQ(ierr); 1209 1210 *a_P_out = Prol; /* out */ 1211 #if defined PETSC_USE_LOG 1212 ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1213 #endif 1214 PetscFunctionReturn(0); 1215 } 1216 1217 /* -------------------------------------------------------------------------- */ 1218 /* 1219 PCGAMGOptprol_AGG 1220 1221 Input Parameter: 1222 . pc - this 1223 . Amat - matrix on this fine level 1224 In/Output Parameter: 1225 . a_P_out - prolongation operator to the next level 1226 */ 1227 #undef __FUNCT__ 1228 #define __FUNCT__ "PCGAMGOptprol_AGG" 1229 PetscErrorCode PCGAMGOptprol_AGG(PC pc,const Mat Amat,Mat *a_P) 1230 { 1231 PetscErrorCode ierr; 1232 PC_MG *mg = (PC_MG*)pc->data; 1233 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1234 const PetscInt verbose = pc_gamg->verbose; 1235 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 1236 PetscInt jj; 1237 PetscMPIInt rank,size; 1238 Mat Prol = *a_P; 1239 MPI_Comm comm; 1240 1241 PetscFunctionBegin; 1242 ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 1243 #if defined PETSC_USE_LOG 1244 ierr = PetscLogEventBegin(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr); 1245 #endif 1246 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1247 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 1248 1249 /* smooth P0 */ 1250 for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) { 1251 Mat tMat; 1252 Vec diag; 1253 PetscReal alpha, emax, emin; 1254 #if defined PETSC_GAMG_USE_LOG 1255 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 1256 #endif 1257 if (jj == 0) { 1258 KSP eksp; 1259 Vec bb, xx; 1260 PC epc; 1261 ierr = MatGetVecs(Amat, &bb, 0);CHKERRQ(ierr); 1262 ierr = MatGetVecs(Amat, &xx, 0);CHKERRQ(ierr); 1263 { 1264 PetscRandom rctx; 1265 ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr); 1266 ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 1267 ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 1268 ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 1269 } 1270 1271 /* zeroing out BC rows -- needed for crazy matrices */ 1272 { 1273 PetscInt Istart,Iend,ncols,jj,Ii; 1274 PetscScalar zero = 0.0; 1275 ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr); 1276 for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) { 1277 ierr = MatGetRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr); 1278 if (ncols <= 1) { 1279 ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr); 1280 } 1281 ierr = MatRestoreRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr); 1282 } 1283 ierr = VecAssemblyBegin(bb);CHKERRQ(ierr); 1284 ierr = VecAssemblyEnd(bb);CHKERRQ(ierr); 1285 } 1286 1287 ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr); 1288 ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr); 1289 ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr); 1290 ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr); 1291 ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr); 1292 ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr); 1293 1294 ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr); 1295 ierr = KSPSetOperators(eksp, Amat, Amat, SAME_NONZERO_PATTERN);CHKERRQ(ierr); 1296 ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr); 1297 1298 ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr); 1299 ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr); /* smoother in smoothed agg. */ 1300 1301 /* solve - keep stuff out of logging */ 1302 ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 1303 ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 1304 ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr); 1305 ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 1306 ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 1307 1308 ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr); 1309 if (verbose > 0) { 1310 PetscPrintf(comm,"\t\t\t%s smooth P0: max eigen=%e min=%e PC=%s\n", 1311 __FUNCT__,emax,emin,PCJACOBI); 1312 } 1313 ierr = VecDestroy(&xx);CHKERRQ(ierr); 1314 ierr = VecDestroy(&bb);CHKERRQ(ierr); 1315 ierr = KSPDestroy(&eksp);CHKERRQ(ierr); 1316 1317 if (pc_gamg->emax_id == -1) { 1318 ierr = PetscObjectComposedDataRegister(&pc_gamg->emax_id);CHKERRQ(ierr); 1319 if (pc_gamg->emax_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->emax_id == -1"); 1320 } 1321 ierr = PetscObjectComposedDataSetScalar((PetscObject)Amat, pc_gamg->emax_id, emax);CHKERRQ(ierr); 1322 } 1323 1324 /* smooth P1 := (I - omega/lam D^{-1}A)P0 */ 1325 ierr = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr); 1326 ierr = MatGetVecs(Amat, &diag, 0);CHKERRQ(ierr); 1327 ierr = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */ 1328 ierr = VecReciprocal(diag);CHKERRQ(ierr); 1329 ierr = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr); 1330 ierr = VecDestroy(&diag);CHKERRQ(ierr); 1331 alpha = -1.4/emax; 1332 ierr = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 1333 ierr = MatDestroy(&Prol);CHKERRQ(ierr); 1334 Prol = tMat; 1335 #if defined PETSC_GAMG_USE_LOG 1336 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 1337 #endif 1338 } 1339 #if defined PETSC_USE_LOG 1340 ierr = PetscLogEventEnd(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr); 1341 #endif 1342 *a_P = Prol; 1343 PetscFunctionReturn(0); 1344 } 1345 1346 /* -------------------------------------------------------------------------- */ 1347 /* 1348 PCCreateGAMG_AGG 1349 1350 Input Parameter: 1351 . pc - 1352 */ 1353 #undef __FUNCT__ 1354 #define __FUNCT__ "PCCreateGAMG_AGG" 1355 PetscErrorCode PCCreateGAMG_AGG(PC pc) 1356 { 1357 PetscErrorCode ierr; 1358 PC_MG *mg = (PC_MG*)pc->data; 1359 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1360 PC_GAMG_AGG *pc_gamg_agg; 1361 1362 PetscFunctionBegin; 1363 /* create sub context for SA */ 1364 ierr = PetscNewLog(pc, PC_GAMG_AGG, &pc_gamg_agg);CHKERRQ(ierr); 1365 pc_gamg->subctx = pc_gamg_agg; 1366 1367 pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 1368 pc_gamg->ops->destroy = PCDestroy_GAMG_AGG; 1369 /* reset does not do anything; setup not virtual */ 1370 1371 /* set internal function pointers */ 1372 pc_gamg->ops->graph = PCGAMGgraph_AGG; 1373 pc_gamg->ops->coarsen = PCGAMGCoarsen_AGG; 1374 pc_gamg->ops->prolongator = PCGAMGProlongator_AGG; 1375 pc_gamg->ops->optprol = PCGAMGOptprol_AGG; 1376 1377 pc_gamg->ops->createdefaultdata = PCSetData_AGG; 1378 1379 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG);CHKERRQ(ierr); 1380 PetscFunctionReturn(0); 1381 } 1382