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