1 /* 2 GAMG geometric-algebric multigrid PC - Mark Adams 2011 3 */ 4 #include "private/matimpl.h" 5 #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 6 #include <private/kspimpl.h> 7 8 #if defined PETSC_USE_LOG 9 PetscLogEvent gamg_setup_events[NUM_SET]; 10 #endif 11 #define GAMG_MAXLEVELS 30 12 13 /* #define GAMG_STAGES */ 14 #if (defined PETSC_USE_LOG && defined GAMG_STAGES) 15 static PetscLogStage gamg_stages[GAMG_MAXLEVELS]; 16 #endif 17 18 static PetscFList GAMGList = 0; 19 20 /* ----------------------------------------------------------------------------- */ 21 #undef __FUNCT__ 22 #define __FUNCT__ "PCReset_GAMG" 23 PetscErrorCode PCReset_GAMG(PC pc) 24 { 25 PetscErrorCode ierr; 26 PC_MG *mg = (PC_MG*)pc->data; 27 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 28 29 PetscFunctionBegin; 30 if( pc_gamg->data != 0 ) { /* this should not happen, cleaned up in SetUp */ 31 ierr = PetscFree(pc_gamg->data); CHKERRQ(ierr); 32 } 33 pc_gamg->data = 0; pc_gamg->data_sz = 0; 34 PetscFunctionReturn(0); 35 } 36 37 /* -------------------------------------------------------------------------- */ 38 /* 39 createLevel: create coarse op with RAP. repartition and/or reduce number 40 of active processors. 41 42 Input Parameter: 43 . pc - parameters 44 . Amat_fine - matrix on this fine (k) level 45 . cbs - coarse block size 46 . isLast - 47 In/Output Parameter: 48 . a_P_inout - prolongation operator to the next level (k-1) 49 . a_nactive_proc - number of active procs 50 Output Parameter: 51 . a_Amat_crs - coarse matrix that is created (k-1) 52 */ 53 54 #undef __FUNCT__ 55 #define __FUNCT__ "createLevel" 56 PetscErrorCode createLevel( const PC pc, 57 const Mat Amat_fine, 58 const PetscInt cbs, 59 const PetscBool isLast, 60 Mat *a_P_inout, 61 PetscMPIInt *a_nactive_proc, 62 Mat *a_Amat_crs 63 ) 64 { 65 PC_MG *mg = (PC_MG*)pc->data; 66 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 67 const PetscBool repart = pc_gamg->repart; 68 const PetscInt min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit; 69 PetscErrorCode ierr; 70 Mat Cmat,Pnew,Pold=*a_P_inout; 71 IS new_indices,isnum; 72 MPI_Comm wcomm = ((PetscObject)Amat_fine)->comm; 73 PetscMPIInt mype,npe,new_npe,nactive = *a_nactive_proc; 74 PetscInt neq,NN,Istart,Iend,Istart0,Iend0,ncrs_new,ncrs0,rfactor; 75 76 PetscFunctionBegin; 77 ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr); 78 ierr = MPI_Comm_size( wcomm, &npe ); CHKERRQ(ierr); 79 80 /* RAP */ 81 #ifdef USE_R 82 /* make R wih brute force for now */ 83 ierr = MatTranspose( Pold, Pnew ); 84 ierr = MatDestroy( &Pold ); CHKERRQ(ierr); 85 ierr = MatRARt( Amat_fine, Pnew, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr); 86 Pold = Pnew; 87 #else 88 ierr = MatPtAP( Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr); 89 #endif 90 ierr = MatSetBlockSize( Cmat, cbs ); CHKERRQ(ierr); 91 ierr = MatGetOwnershipRange( Cmat, &Istart0, &Iend0 ); CHKERRQ(ierr); 92 ncrs0 = (Iend0-Istart0)/cbs; assert((Iend0-Istart0)%cbs == 0); 93 94 /* get number of PEs to make active, reduce */ 95 ierr = MatGetSize( Cmat, &neq, &NN ); CHKERRQ(ierr); 96 new_npe = (PetscMPIInt)((float)neq/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 97 if( new_npe == 0 || neq < coarse_max ) new_npe = 1; 98 else if (new_npe >= nactive ) new_npe = nactive; /* no change, rare */ 99 if( isLast ) new_npe = 1; 100 101 if( !repart && new_npe==nactive ) { 102 *a_Amat_crs = Cmat; /* output - no repartitioning or reduction */ 103 } 104 else { 105 /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */ 106 const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,data_sz=ndata_rows*ndata_cols; 107 const PetscInt stride0=ncrs0*pc_gamg->data_cell_rows; 108 PetscInt *counts,is_sz,*newproc_idx,ii,jj,kk,strideNew,*tidx,inpe,targetPE; 109 IS isnewproc; 110 VecScatter vecscat; 111 PetscScalar *array; 112 Vec src_crd, dest_crd; 113 IS isscat; 114 115 ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr); 116 #if defined PETSC_USE_LOG 117 ierr = PetscLogEventBegin(gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 118 #endif 119 if( repart ) { 120 /* create sub communicator */ 121 Mat adj; 122 123 if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,neq); 124 125 /* MatPartitioningApply call MatConvert, which is collective */ 126 if( cbs == 1 ) { 127 ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 128 } 129 else{ 130 /* make a scalar matrix to partition */ 131 Mat tMat; 132 PetscInt ncols,jj,Ii; 133 const PetscScalar *vals; 134 const PetscInt *idx; 135 PetscInt *d_nnz, *o_nnz; 136 static PetscInt llev = 0; 137 138 ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr); 139 ierr = PetscMalloc( ncrs0*sizeof(PetscInt), &o_nnz ); CHKERRQ(ierr); 140 for ( Ii = Istart0, jj = 0 ; Ii < Iend0 ; Ii += cbs, jj++ ) { 141 ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr); 142 d_nnz[jj] = ncols/cbs; 143 o_nnz[jj] = ncols/cbs; 144 ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr); 145 if( d_nnz[jj] > ncrs0 ) d_nnz[jj] = ncrs0; 146 if( o_nnz[jj] > (neq/cbs-ncrs0) ) o_nnz[jj] = neq/cbs-ncrs0; 147 } 148 149 ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0, 150 PETSC_DETERMINE, PETSC_DETERMINE, 151 0, d_nnz, 0, o_nnz, 152 &tMat ); 153 CHKERRQ(ierr); 154 ierr = PetscFree( d_nnz ); CHKERRQ(ierr); 155 ierr = PetscFree( o_nnz ); CHKERRQ(ierr); 156 157 for ( ii = Istart0; ii < Iend0; ii++ ) { 158 PetscInt dest_row = ii/cbs; 159 ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr); 160 for( jj = 0 ; jj < ncols ; jj++ ){ 161 PetscInt dest_col = idx[jj]/cbs; 162 PetscScalar v = 1.0; 163 ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr); 164 } 165 ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr); 166 } 167 ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168 ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 169 170 if( llev++ == -1 ) { 171 PetscViewer viewer; char fname[32]; 172 sprintf(fname,"part_mat_%d.mat",llev); 173 PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer); 174 ierr = MatView( tMat, viewer ); CHKERRQ(ierr); 175 ierr = PetscViewerDestroy( &viewer ); 176 } 177 178 ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 179 180 ierr = MatDestroy( &tMat ); CHKERRQ(ierr); 181 } 182 183 { /* partition: get newproc_idx, set is_sz */ 184 char prefix[256]; 185 const char *pcpre; 186 const PetscInt *is_idx; 187 MatPartitioning mpart; 188 IS proc_is; 189 190 ierr = MatPartitioningCreate( wcomm, &mpart ); CHKERRQ(ierr); 191 ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr); 192 ierr = PCGetOptionsPrefix( pc, &pcpre );CHKERRQ(ierr); 193 ierr = PetscSNPrintf(prefix,sizeof prefix,"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr); 194 ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 195 ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr); 196 ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr); 197 ierr = MatPartitioningApply( mpart, &proc_is ); CHKERRQ(ierr); 198 ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr); 199 200 /* collect IS info */ 201 ierr = ISGetLocalSize( proc_is, &is_sz ); CHKERRQ(ierr); 202 ierr = PetscMalloc( cbs*is_sz*sizeof(PetscInt), &newproc_idx ); CHKERRQ(ierr); 203 ierr = ISGetIndices( proc_is, &is_idx ); CHKERRQ(ierr); 204 NN = 1; /* bring to "front" of machine */ 205 /*NN = npe/new_npe;*/ /* spread partitioning across machine */ 206 for( kk = jj = 0 ; kk < is_sz ; kk++ ){ 207 for( ii = 0 ; ii < cbs ; ii++, jj++ ){ 208 newproc_idx[jj] = is_idx[kk] * NN; /* distribution */ 209 } 210 } 211 ierr = ISRestoreIndices( proc_is, &is_idx ); CHKERRQ(ierr); 212 ierr = ISDestroy( &proc_is ); CHKERRQ(ierr); 213 214 is_sz *= cbs; 215 } 216 ierr = MatDestroy( &adj ); CHKERRQ(ierr); 217 218 ierr = ISCreateGeneral( wcomm, is_sz, newproc_idx, PETSC_COPY_VALUES, &isnewproc ); 219 CHKERRQ(ierr); 220 if( newproc_idx != 0 ) { 221 ierr = PetscFree( newproc_idx ); CHKERRQ(ierr); 222 } 223 } 224 else { /* simple aggreagtion of parts */ 225 /* find factor */ 226 if( new_npe == 1 ) rfactor = npe; /* easy */ 227 else { 228 PetscReal best_fact = 0.; 229 jj = -1; 230 for( kk = 1 ; kk <= npe ; kk++ ){ 231 if( npe%kk==0 ) { /* a candidate */ 232 PetscReal nactpe = (PetscReal)npe/(PetscReal)kk, fact = nactpe/(PetscReal)new_npe; 233 if(fact > 1.0) fact = 1./fact; /* keep fact < 1 */ 234 if( fact > best_fact ) { 235 best_fact = fact; jj = kk; 236 } 237 } 238 } 239 if(jj!=-1) rfactor = jj; 240 else rfactor = 1; /* prime? */ 241 } 242 new_npe = npe/rfactor; 243 244 if( new_npe==nactive ) { 245 *a_Amat_crs = Cmat; /* output - no repartitioning or reduction */ 246 ierr = PetscFree( counts ); CHKERRQ(ierr); 247 *a_nactive_proc = new_npe; /* output */ 248 if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s aggregate processors noop: new_npe=%d, neq=%d\n",mype,__FUNCT__,new_npe,neq); 249 PetscFunctionReturn(0); 250 } 251 252 /* reduce - isnewproc */ 253 targetPE = mype/rfactor; 254 255 if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s aggregate processors: npe: %d --> %d, neq=%d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,neq); 256 is_sz = ncrs0*cbs; 257 ierr = ISCreateStride( wcomm, is_sz, targetPE, 0, &isnewproc ); 258 CHKERRQ(ierr); 259 } 260 261 /* 262 Create an index set from the isnewproc index set to indicate the mapping TO 263 */ 264 ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr); 265 /* 266 Determine how many elements are assigned to each processor 267 */ 268 inpe = npe; 269 ierr = ISPartitioningCount( isnewproc, inpe, counts ); CHKERRQ(ierr); 270 ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); 271 ncrs_new = counts[mype]/cbs; 272 strideNew = ncrs_new*ndata_rows; 273 #if defined PETSC_USE_LOG 274 ierr = PetscLogEventEnd(gamg_setup_events[SET12],0,0,0,0); CHKERRQ(ierr); 275 #endif 276 /* Create a vector to contain the newly ordered element information */ 277 ierr = VecCreate( wcomm, &dest_crd ); 278 ierr = VecSetSizes( dest_crd, data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr); 279 ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */ 280 /* 281 There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 282 a block size of ...). Note, ISs are expanded into equation space by 'cbs'. 283 */ 284 ierr = PetscMalloc( (ncrs0*data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr); 285 ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr); 286 for(ii=0,jj=0; ii<ncrs0 ; ii++) { 287 PetscInt id = idx[ii*cbs]/cbs; /* get node back */ 288 for( kk=0; kk<data_sz ; kk++, jj++) tidx[jj] = id*data_sz + kk; 289 } 290 ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr); 291 ierr = ISCreateGeneral( wcomm, data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat ); 292 CHKERRQ(ierr); 293 ierr = PetscFree( tidx ); CHKERRQ(ierr); 294 /* 295 Create a vector to contain the original vertex information for each element 296 */ 297 ierr = VecCreateSeq( PETSC_COMM_SELF, data_sz*ncrs0, &src_crd ); CHKERRQ(ierr); 298 for( jj=0; jj<ndata_cols ; jj++ ) { 299 for( ii=0 ; ii<ncrs0 ; ii++) { 300 for( kk=0; kk<ndata_rows ; kk++ ) { 301 PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*data_sz + kk*ndata_cols + jj; 302 PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 303 ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES ); CHKERRQ(ierr); 304 } 305 } 306 } 307 ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr); 308 ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr); 309 /* 310 Scatter the element vertex information (still in the original vertex ordering) 311 to the correct processor 312 */ 313 ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat); 314 CHKERRQ(ierr); 315 ierr = ISDestroy( &isscat ); CHKERRQ(ierr); 316 ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 317 ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 318 ierr = VecScatterDestroy( &vecscat ); CHKERRQ(ierr); 319 ierr = VecDestroy( &src_crd ); CHKERRQ(ierr); 320 /* 321 Put the element vertex data into a new allocation of the gdata->ele 322 */ 323 ierr = PetscFree( pc_gamg->data ); CHKERRQ(ierr); 324 ierr = PetscMalloc( data_sz*ncrs_new*sizeof(PetscReal), &pc_gamg->data ); CHKERRQ(ierr); 325 pc_gamg->data_sz = data_sz*ncrs_new; 326 ierr = VecGetArray( dest_crd, &array ); CHKERRQ(ierr); 327 for( jj=0; jj<ndata_cols ; jj++ ) { 328 for( ii=0 ; ii<ncrs_new ; ii++) { 329 for( kk=0; kk<ndata_rows ; kk++ ) { 330 PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*data_sz + kk*ndata_cols + jj; 331 pc_gamg->data[ix] = PetscRealPart(array[jx]); 332 array[jx] = 1.e300; 333 } 334 } 335 } 336 ierr = VecRestoreArray( dest_crd, &array ); CHKERRQ(ierr); 337 ierr = VecDestroy( &dest_crd ); CHKERRQ(ierr); 338 #if defined PETSC_USE_LOG 339 ierr = PetscLogEventBegin(gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 340 #endif 341 /* 342 Invert for MatGetSubMatrix 343 */ 344 ierr = ISInvertPermutation( isnum, ncrs_new*cbs, &new_indices ); CHKERRQ(ierr); 345 ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */ 346 ierr = ISDestroy( &isnum ); CHKERRQ(ierr); 347 #if defined PETSC_USE_LOG 348 ierr = PetscLogEventEnd(gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 349 ierr = PetscLogEventBegin(gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 350 #endif 351 /* A_crs output */ 352 ierr = MatGetSubMatrix( Cmat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs ); 353 CHKERRQ(ierr); 354 355 ierr = MatDestroy( &Cmat ); CHKERRQ(ierr); 356 Cmat = *a_Amat_crs; /* output */ 357 ierr = MatSetBlockSize( Cmat, cbs ); CHKERRQ(ierr); 358 #if defined PETSC_USE_LOG 359 ierr = PetscLogEventEnd(gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 360 #endif 361 /* prolongator */ 362 ierr = MatGetOwnershipRange( Pold, &Istart, &Iend ); CHKERRQ(ierr); 363 { 364 IS findices; 365 #if defined PETSC_USE_LOG 366 ierr = PetscLogEventBegin(gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 367 #endif 368 ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices); CHKERRQ(ierr); 369 #ifdef USE_R 370 ierr = MatGetSubMatrix( Pold, new_indices, findices, MAT_INITIAL_MATRIX, &Pnew ); 371 #else 372 ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew ); 373 #endif 374 CHKERRQ(ierr); 375 ierr = ISDestroy( &findices ); CHKERRQ(ierr); 376 #if defined PETSC_USE_LOG 377 ierr = PetscLogEventEnd(gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 378 #endif 379 } 380 ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr); 381 *a_P_inout = Pnew; /* output - repartitioned */ 382 ierr = ISDestroy( &new_indices ); CHKERRQ(ierr); 383 ierr = PetscFree( counts ); CHKERRQ(ierr); 384 } 385 386 *a_nactive_proc = new_npe; /* output */ 387 388 if( !PETSC_TRUE ) { 389 PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs; 390 if(llev==0) { 391 sprintf(fname,"Cmat_%d.m",llev++); 392 PetscViewerASCIIOpen(wcomm,fname,&viewer); 393 ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 394 ierr = MatView(Amat_fine, viewer ); CHKERRQ(ierr); 395 ierr = PetscViewerDestroy( &viewer ); 396 } 397 sprintf(fname,"Cmat_%d.m",llev++); 398 PetscViewerASCIIOpen(wcomm,fname,&viewer); 399 ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 400 ierr = MatView(Cmat, viewer ); CHKERRQ(ierr); 401 ierr = PetscViewerDestroy( &viewer ); 402 } 403 404 PetscFunctionReturn(0); 405 } 406 407 /* -------------------------------------------------------------------------- */ 408 /* 409 PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 410 by setting data structures and options. 411 412 Input Parameter: 413 . pc - the preconditioner context 414 415 Application Interface Routine: PCSetUp() 416 417 Notes: 418 The interface routine PCSetUp() is not usually called directly by 419 the user, but instead is called by PCApply() if necessary. 420 */ 421 #undef __FUNCT__ 422 #define __FUNCT__ "PCSetUp_GAMG" 423 PetscErrorCode PCSetUp_GAMG( PC pc ) 424 { 425 PetscErrorCode ierr; 426 PC_MG *mg = (PC_MG*)pc->data; 427 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 428 PC_MG_Levels **mglevels = mg->levels; 429 Mat Pmat = pc->pmat; 430 PetscInt fine_level,level,level1,M,N,bs,nloc,lidx,Istart,Iend; 431 MPI_Comm wcomm = ((PetscObject)pc)->comm; 432 PetscMPIInt mype,npe,nactivepe; 433 Mat Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS]; 434 PetscReal emaxs[GAMG_MAXLEVELS]; 435 PetscLogDouble nnz0,nnztot; 436 MatInfo info; 437 438 PetscFunctionBegin; 439 pc_gamg->setup_count++; 440 if( pc->setupcalled > 0 ) { 441 /* just do Galerkin grids */ 442 Mat B,dA,dB; 443 444 /* PCSetUp_MG seems to insists on setting this to GMRES */ 445 ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr); 446 447 if( pc_gamg->Nlevels > 1 ) { 448 /* currently only handle case where mat and pmat are the same on coarser levels */ 449 ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr); 450 /* (re)set to get dirty flag */ 451 ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 452 ierr = KSPSetUp( mglevels[pc_gamg->Nlevels-1]->smoothd ); CHKERRQ(ierr); 453 454 for (level=pc_gamg->Nlevels-2; level>-1; level--) { 455 ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 456 /* the first time through the matrix structure has changed from repartitioning */ 457 if( pc_gamg->setup_count == 2 ) { 458 ierr = MatDestroy( &B ); CHKERRQ(ierr); 459 ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 460 mglevels[level]->A = B; 461 } 462 else { 463 ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 464 } 465 ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr); 466 dB = B; 467 /* setup KSP/PC */ 468 ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr); 469 } 470 } 471 472 pc->setupcalled = 2; 473 474 PetscFunctionReturn(0); 475 } 476 477 ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 478 ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 479 /* GAMG requires input of fine-grid matrix. It determines nlevels. */ 480 ierr = MatGetBlockSize( Pmat, &bs ); CHKERRQ(ierr); 481 ierr = MatGetSize( Pmat, &M, &N );CHKERRQ(ierr); 482 ierr = MatGetOwnershipRange( Pmat, &Istart, &Iend ); CHKERRQ(ierr); 483 nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0); 484 485 if( pc_gamg->data == 0 && nloc > 0 ) { 486 if(!pc_gamg->createdefaultdata){ 487 SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_LIB,"'createdefaultdata' not set?!?! need to support NULL data!!!"); 488 } 489 ierr = pc_gamg->createdefaultdata( pc ); CHKERRQ(ierr); 490 } 491 492 /* Get A_i and R_i */ 493 if (pc_gamg->verbose) { 494 if(pc_gamg->verbose==1) ierr = MatGetInfo(Pmat,MAT_LOCAL,&info); 495 else ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info); 496 CHKERRQ(ierr); 497 nnz0 = info.nz_used; 498 nnztot = info.nz_used; 499 PetscPrintf(wcomm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n", 500 mype,__FUNCT__,0,N,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols, 501 (int)(nnz0/(PetscReal)N),npe); 502 } 503 for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */ 504 level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */ 505 level++ ){ 506 level1 = level + 1; 507 #if (defined PETSC_USE_LOG && defined GAMG_STAGES) 508 ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr ); 509 #endif 510 #if defined PETSC_USE_LOG 511 ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0); CHKERRQ(ierr); 512 #endif 513 { /* construct prolongator */ 514 Mat Gmat; assert(pc_gamg->graph); 515 IS selected, llist; assert(pc_gamg->coarsen); 516 517 ierr = pc_gamg->graph( pc, Aarr[level], &Gmat ); CHKERRQ(ierr); 518 ierr = pc_gamg->coarsen( pc, Gmat, &selected, &llist ); CHKERRQ(ierr); 519 ierr = pc_gamg->prolongator( pc, Aarr[level], Gmat, selected, llist, &Parr[level1] ); 520 CHKERRQ(ierr); 521 522 if( Parr[level1] ){ 523 /* get new block size of coarse matrices */ 524 if( pc_gamg->col_bs_id != -1 ){ 525 PetscBool flag; 526 ierr = PetscObjectComposedDataGetInt( (PetscObject)Parr[level1], pc_gamg->col_bs_id, bs, flag ); 527 CHKERRQ( ierr ); assert(flag); 528 } 529 } 530 531 if( pc_gamg->optprol && Parr[level1] ){ 532 /* smooth */ 533 ierr = pc_gamg->optprol( pc, Aarr[level], &Parr[level1] ); CHKERRQ(ierr); 534 } 535 536 ierr = MatDestroy( &Gmat ); CHKERRQ(ierr); 537 ierr = ISDestroy( &selected ); CHKERRQ(ierr); 538 ierr = ISDestroy( &llist ); CHKERRQ(ierr); 539 } 540 #if defined PETSC_USE_LOG 541 ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 542 #endif 543 /* cache eigen estimate */ 544 if( pc_gamg->emax_id != -1 ){ 545 PetscBool flag; 546 ierr = PetscObjectComposedDataGetReal( (PetscObject)Aarr[level], pc_gamg->emax_id, emaxs[level], flag ); 547 CHKERRQ( ierr ); 548 if( !flag ) emaxs[level] = -1.; 549 } 550 else emaxs[level] = -1.; 551 if(level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 552 if( !Parr[level1] ) { 553 if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level); 554 break; 555 } 556 #if defined PETSC_USE_LOG 557 ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 558 #endif 559 ierr = createLevel( pc, Aarr[level], bs, 560 (PetscBool)(level==pc_gamg->Nlevels-2), 561 &Parr[level1], &nactivepe, &Aarr[level1] ); 562 CHKERRQ(ierr); 563 #if defined PETSC_USE_LOG 564 ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 565 #endif 566 ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr); 567 ierr = MatGetInfo( Aarr[level1], MAT_GLOBAL_SUM, &info ); CHKERRQ(ierr); 568 if (pc_gamg->verbose){ 569 PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n", 570 mype,__FUNCT__,(int)level1,N,pc_gamg->data_cell_cols, 571 (int)(info.nz_used/(PetscReal)N),nactivepe); 572 } 573 /* stop if one node */ 574 if( M/pc_gamg->data_cell_cols < 2 ) { 575 level++; 576 break; 577 } 578 /* Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id; */ 579 /* v = 1.e-10; /\* LU factor has hard wired numbers for small diags so this needs to match (yuk) *\/ */ 580 /* ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr); */ 581 /* nloceq = Iend-Istart; */ 582 /* ierr = MatGetVecs( Aarr[level1], &diag, 0 ); CHKERRQ(ierr); */ 583 /* ierr = MatGetDiagonal( Aarr[level1], diag ); CHKERRQ(ierr); */ 584 /* ierr = VecGetArray( diag, &data_arr ); CHKERRQ(ierr); */ 585 /* for(kk=0;kk<nloceq;kk++){ */ 586 /* if(data_arr[kk]==0.0) { */ 587 /* id = kk + Istart; */ 588 /* ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES); */ 589 /* CHKERRQ(ierr); */ 590 /* PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1); */ 591 /* } */ 592 /* } */ 593 /* ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr); */ 594 /* ierr = VecDestroy( &diag ); CHKERRQ(ierr); */ 595 ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 596 ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 597 598 #if (defined PETSC_USE_LOG && defined GAMG_STAGES) 599 ierr = PetscLogStagePop(); CHKERRQ( ierr ); 600 #endif 601 if(pc_gamg->verbose){ 602 if(pc_gamg->verbose==1) ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info); 603 else ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); 604 CHKERRQ(ierr); 605 nnztot += info.nz_used; 606 } 607 } /* levels */ 608 609 if( pc_gamg->data ) { 610 ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr ); 611 pc_gamg->data = 0; 612 } 613 614 if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid compexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0); 615 pc_gamg->Nlevels = level + 1; 616 fine_level = level; 617 ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr); 618 619 if( pc_gamg->Nlevels > 1 ) { /* don't setup MG if */ 620 /* set default smoothers */ 621 for ( lidx = 1, level = pc_gamg->Nlevels-2; 622 lidx <= fine_level; 623 lidx++, level--) { 624 PetscReal emax, emin; 625 KSP smoother; PC subpc; 626 PetscBool isCheb; 627 /* set defaults */ 628 ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr); 629 ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr); 630 ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 631 /* ierr = KSPSetTolerances(smoother,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,2); CHKERRQ(ierr); */ 632 ierr = PCSetType( subpc, PCJACOBI ); CHKERRQ(ierr); 633 ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 634 /* overide defaults with input parameters */ 635 ierr = KSPSetFromOptions( smoother ); CHKERRQ(ierr); 636 637 ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 638 /* do my own cheby */ 639 ierr = PetscTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &isCheb ); CHKERRQ(ierr); 640 641 if( isCheb ) { 642 ierr = PetscTypeCompare( (PetscObject)subpc, PCJACOBI, &isCheb ); CHKERRQ(ierr); 643 if( isCheb && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */ 644 else{ /* eigen estimate 'emax' */ 645 KSP eksp; Mat Lmat = Aarr[level]; 646 Vec bb, xx; 647 648 ierr = MatGetVecs( Lmat, &bb, 0 ); CHKERRQ(ierr); 649 ierr = MatGetVecs( Lmat, &xx, 0 ); CHKERRQ(ierr); 650 { 651 PetscRandom rctx; 652 ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr); 653 ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 654 ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 655 ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr); 656 } 657 ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr); 658 ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 ); 659 CHKERRQ(ierr); 660 ierr = KSPSetNormType( eksp, KSP_NORM_NONE ); CHKERRQ(ierr); 661 662 ierr = KSPAppendOptionsPrefix( eksp, "est_"); CHKERRQ(ierr); 663 ierr = KSPSetFromOptions( eksp ); CHKERRQ(ierr); 664 665 ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr); 666 ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr ); 667 ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr); 668 669 { /* set PC type to be same as smoother - does not get all parameters!!! */ 670 const PCType type; PC pc; 671 ierr = PCGetType( subpc, &type ); CHKERRQ(ierr); 672 ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr ); 673 ierr = PCSetType( pc, type ); CHKERRQ(ierr); 674 } 675 676 /* solve - keep stuff out of logging */ 677 ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 678 ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 679 ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr); 680 ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 681 ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 682 683 ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr); 684 685 ierr = VecDestroy( &xx ); CHKERRQ(ierr); 686 ierr = VecDestroy( &bb ); CHKERRQ(ierr); 687 ierr = KSPDestroy( &eksp ); CHKERRQ(ierr); 688 689 if (pc_gamg->verbose) { 690 PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e\n",__FUNCT__,emax,emin); 691 } 692 } 693 { 694 PetscInt N1, N0, tt; 695 ierr = MatGetSize( Aarr[level], &N1, &tt ); CHKERRQ(ierr); 696 ierr = MatGetSize( Aarr[level+1], &N0, &tt ); CHKERRQ(ierr); 697 /* heuristic - is this crap? */ 698 emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); 699 emax *= 1.05; 700 } 701 ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr); 702 } 703 } 704 { 705 /* coarse grid */ 706 KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 707 Mat Lmat = Aarr[pc_gamg->Nlevels-1]; 708 ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr); 709 ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 710 ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 711 ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 712 ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr); 713 ierr = PCSetUp( subpc ); CHKERRQ(ierr); 714 ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); 715 assert(ii==1); 716 ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 717 ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr); 718 } 719 720 /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 721 ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 722 ierr = PCSetFromOptions_MG(pc); CHKERRQ(ierr); 723 ierr = PetscOptionsEnd();CHKERRQ(ierr); 724 if (mg->galerkin != 2) SETERRQ(wcomm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used."); 725 726 /* set interpolation between the levels, clean up */ 727 for (lidx=0,level=pc_gamg->Nlevels-1; 728 lidx<fine_level; 729 lidx++, level--){ 730 ierr = PCMGSetInterpolation( pc, lidx+1, Parr[level] );CHKERRQ(ierr); 731 ierr = MatDestroy( &Parr[level] ); CHKERRQ(ierr); 732 ierr = MatDestroy( &Aarr[level] ); CHKERRQ(ierr); 733 } 734 735 /* setupcalled is set to 0 so that MG is setup from scratch */ 736 pc->setupcalled = 0; 737 ierr = PCSetUp_MG( pc );CHKERRQ( ierr ); 738 pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */ 739 740 { 741 KSP smoother; /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */ 742 ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr); 743 ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr); 744 ierr = KSPSetUp( smoother ); CHKERRQ(ierr); 745 } 746 } 747 else { 748 KSP smoother; 749 if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",mype,__FUNCT__); 750 ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr); 751 ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 752 ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr); 753 /* setupcalled is set to 0 so that MG is setup from scratch */ 754 pc->setupcalled = 0; 755 ierr = PCSetUp_MG( pc );CHKERRQ( ierr ); 756 pc->setupcalled = 1; /* use 1 as signal that this has not been re-setup */ 757 } 758 PetscFunctionReturn(0); 759 } 760 761 /* ------------------------------------------------------------------------- */ 762 /* 763 PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 764 that was created with PCCreate_GAMG(). 765 766 Input Parameter: 767 . pc - the preconditioner context 768 769 Application Interface Routine: PCDestroy() 770 */ 771 #undef __FUNCT__ 772 #define __FUNCT__ "PCDestroy_GAMG" 773 PetscErrorCode PCDestroy_GAMG( PC pc ) 774 { 775 PetscErrorCode ierr; 776 PC_MG *mg = (PC_MG*)pc->data; 777 PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 778 779 PetscFunctionBegin; 780 ierr = PCReset_GAMG( pc );CHKERRQ(ierr); 781 ierr = PetscFree( pc_gamg );CHKERRQ(ierr); 782 ierr = PCDestroy_MG( pc );CHKERRQ(ierr); 783 PetscFunctionReturn(0); 784 } 785 786 787 #undef __FUNCT__ 788 #define __FUNCT__ "PCGAMGSetProcEqLim" 789 /*@ 790 PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via 791 processor reduction. 792 793 Not Collective on PC 794 795 Input Parameters: 796 . pc - the preconditioner context 797 798 799 Options Database Key: 800 . -pc_gamg_process_eq_limit 801 802 Level: intermediate 803 804 Concepts: Unstructured multrigrid preconditioner 805 806 .seealso: () 807 @*/ 808 PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 809 { 810 PetscErrorCode ierr; 811 812 PetscFunctionBegin; 813 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 814 ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 815 PetscFunctionReturn(0); 816 } 817 818 EXTERN_C_BEGIN 819 #undef __FUNCT__ 820 #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG" 821 PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 822 { 823 PC_MG *mg = (PC_MG*)pc->data; 824 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 825 826 PetscFunctionBegin; 827 if(n>0) pc_gamg->min_eq_proc = n; 828 PetscFunctionReturn(0); 829 } 830 EXTERN_C_END 831 832 #undef __FUNCT__ 833 #define __FUNCT__ "PCGAMGSetCoarseEqLim" 834 /*@ 835 PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids. 836 837 Collective on PC 838 839 Input Parameters: 840 . pc - the preconditioner context 841 842 843 Options Database Key: 844 . -pc_gamg_coarse_eq_limit 845 846 Level: intermediate 847 848 Concepts: Unstructured multrigrid preconditioner 849 850 .seealso: () 851 @*/ 852 PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 853 { 854 PetscErrorCode ierr; 855 856 PetscFunctionBegin; 857 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 858 ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 859 PetscFunctionReturn(0); 860 } 861 862 EXTERN_C_BEGIN 863 #undef __FUNCT__ 864 #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG" 865 PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 866 { 867 PC_MG *mg = (PC_MG*)pc->data; 868 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 869 870 PetscFunctionBegin; 871 if(n>0) pc_gamg->coarse_eq_limit = n; 872 PetscFunctionReturn(0); 873 } 874 EXTERN_C_END 875 876 #undef __FUNCT__ 877 #define __FUNCT__ "PCGAMGSetRepartitioning" 878 /*@ 879 PCGAMGSetRepartitioning - Repartition the coarse grids 880 881 Collective on PC 882 883 Input Parameters: 884 . pc - the preconditioner context 885 886 887 Options Database Key: 888 . -pc_gamg_repartition 889 890 Level: intermediate 891 892 Concepts: Unstructured multrigrid preconditioner 893 894 .seealso: () 895 @*/ 896 PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n) 897 { 898 PetscErrorCode ierr; 899 900 PetscFunctionBegin; 901 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 902 ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 903 PetscFunctionReturn(0); 904 } 905 906 EXTERN_C_BEGIN 907 #undef __FUNCT__ 908 #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG" 909 PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n) 910 { 911 PC_MG *mg = (PC_MG*)pc->data; 912 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 913 914 PetscFunctionBegin; 915 pc_gamg->repart = n; 916 PetscFunctionReturn(0); 917 } 918 EXTERN_C_END 919 920 #undef __FUNCT__ 921 #define __FUNCT__ "PCGAMGSetNlevels" 922 /*@ 923 PCGAMGSetNlevels - 924 925 Not collective on PC 926 927 Input Parameters: 928 . pc - the preconditioner context 929 930 931 Options Database Key: 932 . -pc_mg_levels 933 934 Level: intermediate 935 936 Concepts: Unstructured multrigrid preconditioner 937 938 .seealso: () 939 @*/ 940 PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 941 { 942 PetscErrorCode ierr; 943 944 PetscFunctionBegin; 945 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 946 ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 947 PetscFunctionReturn(0); 948 } 949 950 EXTERN_C_BEGIN 951 #undef __FUNCT__ 952 #define __FUNCT__ "PCGAMGSetNlevels_GAMG" 953 PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 954 { 955 PC_MG *mg = (PC_MG*)pc->data; 956 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 957 958 PetscFunctionBegin; 959 pc_gamg->Nlevels = n; 960 PetscFunctionReturn(0); 961 } 962 EXTERN_C_END 963 964 #undef __FUNCT__ 965 #define __FUNCT__ "PCGAMGSetThreshold" 966 /*@ 967 PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 968 969 Not collective on PC 970 971 Input Parameters: 972 . pc - the preconditioner context 973 974 975 Options Database Key: 976 . -pc_gamg_threshold 977 978 Level: intermediate 979 980 Concepts: Unstructured multrigrid preconditioner 981 982 .seealso: () 983 @*/ 984 PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n) 985 { 986 PetscErrorCode ierr; 987 988 PetscFunctionBegin; 989 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 990 ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr); 991 PetscFunctionReturn(0); 992 } 993 994 EXTERN_C_BEGIN 995 #undef __FUNCT__ 996 #define __FUNCT__ "PCGAMGSetThreshold_GAMG" 997 PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n) 998 { 999 PC_MG *mg = (PC_MG*)pc->data; 1000 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1001 1002 PetscFunctionBegin; 1003 pc_gamg->threshold = n; 1004 PetscFunctionReturn(0); 1005 } 1006 EXTERN_C_END 1007 1008 #undef __FUNCT__ 1009 #define __FUNCT__ "PCGAMGSetType" 1010 /*@ 1011 PCGAMGSetType - Set solution method - calls sub create method 1012 1013 Collective on PC 1014 1015 Input Parameters: 1016 . pc - the preconditioner context 1017 1018 1019 Options Database Key: 1020 . -pc_gamg_type 1021 1022 Level: intermediate 1023 1024 Concepts: Unstructured multrigrid preconditioner 1025 1026 .seealso: () 1027 @*/ 1028 PetscErrorCode PCGAMGSetType( PC pc, const PCGAMGType type ) 1029 { 1030 PetscErrorCode ierr; 1031 1032 PetscFunctionBegin; 1033 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1034 ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,const PCGAMGType),(pc,type)); 1035 CHKERRQ(ierr); 1036 PetscFunctionReturn(0); 1037 } 1038 1039 EXTERN_C_BEGIN 1040 #undef __FUNCT__ 1041 #define __FUNCT__ "PCGAMGSetType_GAMG" 1042 PetscErrorCode PCGAMGSetType_GAMG( PC pc, const PCGAMGType type ) 1043 { 1044 PetscErrorCode ierr,(*r)(PC); 1045 1046 PetscFunctionBegin; 1047 ierr = PetscFListFind(GAMGList,((PetscObject)pc)->comm,type,PETSC_FALSE,(PetscVoidStarFunction)&r); 1048 CHKERRQ(ierr); 1049 1050 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 1051 1052 /* call sub create method */ 1053 ierr = (*r)(pc); CHKERRQ(ierr); 1054 1055 PetscFunctionReturn(0); 1056 } 1057 EXTERN_C_END 1058 1059 #undef __FUNCT__ 1060 #define __FUNCT__ "PCSetFromOptions_GAMG" 1061 PetscErrorCode PCSetFromOptions_GAMG( PC pc ) 1062 { 1063 PetscErrorCode ierr; 1064 PC_MG *mg = (PC_MG*)pc->data; 1065 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1066 PetscBool flag; 1067 MPI_Comm wcomm = ((PetscObject)pc)->comm; 1068 1069 PetscFunctionBegin; 1070 ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr); 1071 { 1072 /* -pc_gamg_verbose */ 1073 ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG", 1074 "none", pc_gamg->verbose, 1075 &pc_gamg->verbose, PETSC_NULL ); 1076 CHKERRQ(ierr); 1077 1078 /* -pc_gamg_repartition */ 1079 ierr = PetscOptionsBool("-pc_gamg_repartition", 1080 "Repartion coarse grids (false)", 1081 "PCGAMGRepartitioning", 1082 pc_gamg->repart, 1083 &pc_gamg->repart, 1084 &flag); 1085 CHKERRQ(ierr); 1086 1087 /* -pc_gamg_process_eq_limit */ 1088 ierr = PetscOptionsInt("-pc_gamg_process_eq_limit", 1089 "Limit (goal) on number of equations per process on coarse grids", 1090 "PCGAMGSetProcEqLim", 1091 pc_gamg->min_eq_proc, 1092 &pc_gamg->min_eq_proc, 1093 &flag ); 1094 CHKERRQ(ierr); 1095 1096 /* -pc_gamg_coarse_eq_limit */ 1097 ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit", 1098 "Limit on number of equations for the coarse grid", 1099 "PCGAMGSetCoarseEqLim", 1100 pc_gamg->coarse_eq_limit, 1101 &pc_gamg->coarse_eq_limit, 1102 &flag ); 1103 CHKERRQ(ierr); 1104 1105 /* -pc_gamg_threshold */ 1106 ierr = PetscOptionsReal("-pc_gamg_threshold", 1107 "Relative threshold to use for dropping edges in aggregation graph", 1108 "PCGAMGSetThreshold", 1109 pc_gamg->threshold, 1110 &pc_gamg->threshold, 1111 &flag ); 1112 CHKERRQ(ierr); 1113 if(flag && pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold); 1114 1115 ierr = PetscOptionsInt("-pc_mg_levels", 1116 "Set number of MG levels", 1117 "PCGAMGSetNlevels", 1118 pc_gamg->Nlevels, 1119 &pc_gamg->Nlevels, 1120 &flag ); 1121 } 1122 ierr = PetscOptionsTail();CHKERRQ(ierr); 1123 1124 PetscFunctionReturn(0); 1125 } 1126 1127 /* -------------------------------------------------------------------------- */ 1128 /*MC 1129 PCCreate_GAMG - Geometric algebraic multigrid (AMG) preconditioning framework. 1130 - This is the entry point to GAMG, registered in pcregis.c 1131 1132 Options Database Keys: 1133 Multigrid options(inherited) 1134 + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType) 1135 . -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp) 1136 . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown) 1137 - -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 1138 1139 Level: intermediate 1140 1141 Concepts: multigrid 1142 1143 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 1144 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 1145 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 1146 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 1147 PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 1148 M*/ 1149 EXTERN_C_BEGIN 1150 #undef __FUNCT__ 1151 #define __FUNCT__ "PCCreate_GAMG" 1152 PetscErrorCode PCCreate_GAMG( PC pc ) 1153 { 1154 PetscErrorCode ierr; 1155 PC_GAMG *pc_gamg; 1156 PC_MG *mg; 1157 PetscClassId cookie; 1158 1159 #if defined PETSC_USE_LOG 1160 static long count = 0; 1161 #endif 1162 1163 PetscFunctionBegin; 1164 /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 1165 ierr = PCSetType( pc, PCMG ); CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 1166 ierr = PetscObjectChangeTypeName( (PetscObject)pc, PCGAMG ); CHKERRQ(ierr); 1167 1168 /* create a supporting struct and attach it to pc */ 1169 ierr = PetscNewLog( pc, PC_GAMG, &pc_gamg ); CHKERRQ(ierr); 1170 mg = (PC_MG*)pc->data; 1171 mg->galerkin = 2; /* Use Galerkin, but it is computed externally */ 1172 mg->innerctx = pc_gamg; 1173 1174 pc_gamg->setup_count = 0; 1175 /* these should be in subctx but repartitioning needs simple arrays */ 1176 pc_gamg->data_sz = 0; 1177 pc_gamg->data = 0; 1178 1179 /* register AMG type */ 1180 if( !GAMGList ){ 1181 ierr = PetscFListAdd(&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr); 1182 ierr = PetscFListAdd(&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr); 1183 } 1184 1185 /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 1186 pc->ops->setfromoptions = PCSetFromOptions_GAMG; 1187 pc->ops->setup = PCSetUp_GAMG; 1188 pc->ops->reset = PCReset_GAMG; 1189 pc->ops->destroy = PCDestroy_GAMG; 1190 1191 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1192 "PCGAMGSetProcEqLim_C", 1193 "PCGAMGSetProcEqLim_GAMG", 1194 PCGAMGSetProcEqLim_GAMG); 1195 CHKERRQ(ierr); 1196 1197 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1198 "PCGAMGSetCoarseEqLim_C", 1199 "PCGAMGSetCoarseEqLim_GAMG", 1200 PCGAMGSetCoarseEqLim_GAMG); 1201 CHKERRQ(ierr); 1202 1203 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1204 "PCGAMGSetRepartitioning_C", 1205 "PCGAMGSetRepartitioning_GAMG", 1206 PCGAMGSetRepartitioning_GAMG); 1207 CHKERRQ(ierr); 1208 1209 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1210 "PCGAMGSetThreshold_C", 1211 "PCGAMGSetThreshold_GAMG", 1212 PCGAMGSetThreshold_GAMG); 1213 CHKERRQ(ierr); 1214 1215 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1216 "PCGAMGSetType_C", 1217 "PCGAMGSetType_GAMG", 1218 PCGAMGSetType_GAMG); 1219 CHKERRQ(ierr); 1220 1221 pc_gamg->repart = PETSC_FALSE; 1222 pc_gamg->min_eq_proc = 100; 1223 pc_gamg->coarse_eq_limit = 800; 1224 pc_gamg->threshold = 0.05; 1225 pc_gamg->Nlevels = GAMG_MAXLEVELS; 1226 pc_gamg->verbose = 0; 1227 pc_gamg->emax_id = -1; 1228 pc_gamg->col_bs_id = -1; 1229 1230 /* instantiate derived type */ 1231 ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr); 1232 { 1233 char tname[256] = GAMGAGG; 1234 ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType", 1235 GAMGList, tname, tname, sizeof(tname), PETSC_NULL ); 1236 CHKERRQ(ierr); 1237 ierr = PCGAMGSetType( pc, tname ); CHKERRQ(ierr); 1238 } 1239 ierr = PetscOptionsTail();CHKERRQ(ierr); 1240 1241 #if defined PETSC_USE_LOG 1242 if( count++ == 0 ) { 1243 PetscClassIdRegister("GAMG Setup",&cookie); 1244 PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]); 1245 PetscLogEventRegister(" Graph", cookie, &gamg_setup_events[GRAPH]); 1246 /* PetscLogEventRegister(" G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]); */ 1247 /* PetscLogEventRegister(" G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]); */ 1248 /* PetscLogEventRegister(" G.Square", cookie, &gamg_setup_events[GRAPH_SQR]); */ 1249 PetscLogEventRegister(" MIS/Agg", cookie, &gamg_setup_events[SET4]); 1250 PetscLogEventRegister(" geo: growSupp", cookie, &gamg_setup_events[SET5]); 1251 PetscLogEventRegister(" geo: triangle", cookie, &gamg_setup_events[SET6]); 1252 PetscLogEventRegister(" search&set", cookie, &gamg_setup_events[FIND_V]); 1253 PetscLogEventRegister(" SA: colect data", cookie, &gamg_setup_events[SET7]); 1254 PetscLogEventRegister(" SA: frmProl0", cookie, &gamg_setup_events[SET8]); 1255 PetscLogEventRegister(" SA: smooth", cookie, &gamg_setup_events[SET9]); 1256 PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]); 1257 PetscLogEventRegister(" repartition", cookie, &gamg_setup_events[SET12]); 1258 PetscLogEventRegister(" Invert-Sort", cookie, &gamg_setup_events[SET13]); 1259 PetscLogEventRegister(" Move A", cookie, &gamg_setup_events[SET14]); 1260 PetscLogEventRegister(" Move P", cookie, &gamg_setup_events[SET15]); 1261 1262 /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */ 1263 /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */ 1264 /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */ 1265 /* create timer stages */ 1266 #if defined GAMG_STAGES 1267 { 1268 char str[32]; 1269 sprintf(str,"MG Level %d (finest)",0); 1270 PetscLogStageRegister(str, &gamg_stages[0]); 1271 PetscInt lidx; 1272 for (lidx=1;lidx<9;lidx++){ 1273 sprintf(str,"MG Level %d",lidx); 1274 PetscLogStageRegister(str, &gamg_stages[lidx]); 1275 } 1276 } 1277 #endif 1278 } 1279 #endif 1280 PetscFunctionReturn(0); 1281 } 1282 EXTERN_C_END 1283