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