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