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 Mat Pmat = pc->pmat; 429 PetscInt fine_level,level,level1,M,N,bs,nloc,lidx,Istart,Iend; 430 MPI_Comm wcomm = ((PetscObject)pc)->comm; 431 PetscMPIInt mype,npe,nactivepe; 432 Mat Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS]; 433 PetscReal emaxs[GAMG_MAXLEVELS]; 434 PetscLogDouble nnz0,nnztot; 435 MatInfo info; 436 437 PetscFunctionBegin; 438 ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 439 ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 440 if( pc_gamg->setup_count++ > 0 ) { 441 PC_MG_Levels **mglevels = mg->levels; 442 /* just do Galerkin grids */ 443 Mat B,dA,dB; 444 assert(pc->setupcalled); 445 446 if( pc_gamg->Nlevels > 1 ) { 447 /* currently only handle case where mat and pmat are the same on coarser levels */ 448 ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr); 449 /* (re)set to get dirty flag */ 450 ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 451 ierr = KSPSetUp( mglevels[pc_gamg->Nlevels-1]->smoothd ); CHKERRQ(ierr); 452 453 for (level=pc_gamg->Nlevels-2; level>-1; level--) { 454 ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 455 /* the first time through the matrix structure has changed from repartitioning */ 456 if( pc_gamg->setup_count==2 /*&& (pc_gamg->repart || level==0)*/) { 457 ierr = MatDestroy( &B ); CHKERRQ(ierr); 458 ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 459 mglevels[level]->A = B; 460 } 461 else { 462 ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 463 } 464 ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr); 465 dB = B; 466 /* setup KSP/PC */ 467 ierr = KSPSetUp( mglevels[level]->smoothd ); CHKERRQ(ierr); 468 } 469 } 470 471 ierr = PCSetUp_MG( pc );CHKERRQ( ierr ); 472 473 /* PCSetUp_MG seems to insists on setting this to GMRES */ 474 ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr); 475 476 PetscFunctionReturn(0); 477 } 478 assert(pc->setupcalled == 0); 479 480 /* GAMG requires input of fine-grid matrix. It determines nlevels. */ 481 ierr = MatGetBlockSize( Pmat, &bs ); CHKERRQ(ierr); 482 ierr = MatGetSize( Pmat, &M, &N );CHKERRQ(ierr); 483 ierr = MatGetOwnershipRange( Pmat, &Istart, &Iend ); CHKERRQ(ierr); 484 nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0); 485 486 if( pc_gamg->data == 0 && nloc > 0 ) { 487 if(!pc_gamg->createdefaultdata){ 488 SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_LIB,"'createdefaultdata' not set?!?! need to support NULL data!!!"); 489 } 490 ierr = pc_gamg->createdefaultdata( pc ); CHKERRQ(ierr); 491 } 492 493 /* Get A_i and R_i */ 494 if (pc_gamg->verbose) { 495 if(pc_gamg->verbose==1) ierr = MatGetInfo(Pmat,MAT_LOCAL,&info); 496 else ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info); 497 CHKERRQ(ierr); 498 nnz0 = info.nz_used; 499 nnztot = info.nz_used; 500 PetscPrintf(wcomm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n", 501 mype,__FUNCT__,0,N,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols, 502 (int)(nnz0/(PetscReal)N),npe); 503 } 504 505 for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */ 506 level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */ 507 level++ ){ 508 level1 = level + 1; 509 #if (defined PETSC_USE_LOG && defined GAMG_STAGES) 510 ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr ); 511 #endif 512 #if defined PETSC_USE_LOG 513 ierr = PetscLogEventBegin(gamg_setup_events[SET1],0,0,0,0); CHKERRQ(ierr); 514 #endif 515 { /* construct prolongator */ 516 Mat Gmat; 517 IS selected, llist; 518 assert(pc_gamg->graph); 519 assert(pc_gamg->coarsen); 520 521 ierr = pc_gamg->graph( pc, Aarr[level], &Gmat ); CHKERRQ(ierr); 522 ierr = pc_gamg->coarsen( pc, Gmat, &selected, &llist ); CHKERRQ(ierr); 523 ierr = pc_gamg->prolongator( pc, Aarr[level], Gmat, selected, llist, &Parr[level1] ); 524 CHKERRQ(ierr); 525 526 if( Parr[level1] ){ 527 /* get new block size of coarse matrices */ 528 if( pc_gamg->col_bs_id != -1 ){ 529 PetscBool flag; 530 ierr = PetscObjectComposedDataGetInt( (PetscObject)Parr[level1], pc_gamg->col_bs_id, bs, flag ); 531 CHKERRQ( ierr ); assert(flag); 532 } 533 } 534 535 if( pc_gamg->optprol && Parr[level1] ){ 536 /* smooth */ 537 ierr = pc_gamg->optprol( pc, Aarr[level], &Parr[level1] ); CHKERRQ(ierr); 538 } 539 540 ierr = MatDestroy( &Gmat ); CHKERRQ(ierr); 541 ierr = ISDestroy( &selected ); CHKERRQ(ierr); 542 ierr = ISDestroy( &llist ); CHKERRQ(ierr); 543 } 544 #if defined PETSC_USE_LOG 545 ierr = PetscLogEventEnd(gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 546 #endif 547 /* cache eigen estimate */ 548 if( pc_gamg->emax_id != -1 ){ 549 PetscBool flag; 550 ierr = PetscObjectComposedDataGetReal( (PetscObject)Aarr[level], pc_gamg->emax_id, emaxs[level], flag ); 551 CHKERRQ( ierr ); 552 if( !flag ) emaxs[level] = -1.; 553 } 554 else emaxs[level] = -1.; 555 if(level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 556 if( !Parr[level1] ) { 557 if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level); 558 break; 559 } 560 #if defined PETSC_USE_LOG 561 ierr = PetscLogEventBegin(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 562 #endif 563 ierr = createLevel( pc, Aarr[level], bs, 564 (PetscBool)(level==pc_gamg->Nlevels-2), 565 &Parr[level1], &nactivepe, &Aarr[level1] ); 566 CHKERRQ(ierr); 567 #if defined PETSC_USE_LOG 568 ierr = PetscLogEventEnd(gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 569 #endif 570 ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr); 571 if (pc_gamg->verbose){ 572 ierr = MatGetInfo( Aarr[level1], MAT_GLOBAL_SUM, &info ); CHKERRQ(ierr); 573 PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n", 574 mype,__FUNCT__,(int)level1,N,pc_gamg->data_cell_cols, 575 (int)(info.nz_used/(PetscReal)N),nactivepe); 576 } 577 /* stop if one node */ 578 if( M/pc_gamg->data_cell_cols < 2 ) { 579 level++; 580 break; 581 } 582 /* Vec diag; PetscScalar *data_arr,v; PetscInt Istart,Iend,kk,nloceq,id; */ 583 /* v = 1.e-10; /\* LU factor has hard wired numbers for small diags so this needs to match (yuk) *\/ */ 584 /* ierr = MatGetOwnershipRange(Aarr[level1], &Istart, &Iend); CHKERRQ(ierr); */ 585 /* nloceq = Iend-Istart; */ 586 /* ierr = MatGetVecs( Aarr[level1], &diag, 0 ); CHKERRQ(ierr); */ 587 /* ierr = MatGetDiagonal( Aarr[level1], diag ); CHKERRQ(ierr); */ 588 /* ierr = VecGetArray( diag, &data_arr ); CHKERRQ(ierr); */ 589 /* for(kk=0;kk<nloceq;kk++){ */ 590 /* if(data_arr[kk]==0.0) { */ 591 /* id = kk + Istart; */ 592 /* ierr = MatSetValues(Aarr[level1],1,&id,1,&id,&v,INSERT_VALUES); */ 593 /* CHKERRQ(ierr); */ 594 /* PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s warning: added zero to diag (%d) on level %d \n",mype,__FUNCT__,id,level1); */ 595 /* } */ 596 /* } */ 597 /* ierr = VecRestoreArray( diag, &data_arr ); CHKERRQ(ierr); */ 598 /* ierr = VecDestroy( &diag ); CHKERRQ(ierr); */ 599 /* ierr = MatAssemblyBegin(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); */ 600 /* ierr = MatAssemblyEnd(Aarr[level1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); */ 601 602 #if (defined PETSC_USE_LOG && defined GAMG_STAGES) 603 ierr = PetscLogStagePop(); CHKERRQ( ierr ); 604 #endif 605 if(pc_gamg->verbose){ 606 if(pc_gamg->verbose==1) ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info); 607 else ierr = MatGetInfo(Aarr[level1],MAT_GLOBAL_SUM,&info); 608 CHKERRQ(ierr); 609 nnztot += info.nz_used; 610 } 611 } /* levels */ 612 613 if( pc_gamg->data ) { 614 ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr ); 615 pc_gamg->data = 0; 616 } 617 618 if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid compexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0); 619 pc_gamg->Nlevels = level + 1; 620 fine_level = level; 621 ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr); 622 623 /* simple setup */ 624 if( !PETSC_TRUE ){ 625 PC_MG_Levels **mglevels = mg->levels; 626 for (lidx=0,level=pc_gamg->Nlevels-1; 627 lidx<fine_level; 628 lidx++, level--){ 629 ierr = PCMGSetInterpolation( pc, lidx+1, Parr[level] );CHKERRQ(ierr); 630 ierr = KSPSetOperators( mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 631 ierr = MatDestroy( &Parr[level] ); CHKERRQ(ierr); 632 ierr = MatDestroy( &Aarr[level] ); CHKERRQ(ierr); 633 } 634 ierr = KSPSetOperators( mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 635 636 ierr = PCSetUp_MG( pc ); CHKERRQ( ierr ); 637 } 638 else if( pc_gamg->Nlevels > 1 ) { /* don't setup MG if one level */ 639 /* set default smoothers & set operators */ 640 for ( lidx = 1, level = pc_gamg->Nlevels-2; 641 lidx <= fine_level; 642 lidx++, level--) { 643 KSP smoother; PC subpc; 644 /* set defaults */ 645 ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr); 646 ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr); 647 ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 648 ierr = PCSetType( subpc, PCJACOBI ); CHKERRQ(ierr); 649 ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 650 /* set ops */ 651 ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 652 ierr = PCMGSetInterpolation( pc, lidx, Parr[level+1] );CHKERRQ(ierr); 653 } 654 { 655 /* coarse grid */ 656 KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 657 Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0; 658 ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr); 659 ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 660 ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 661 ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 662 ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr); 663 ierr = PCSetUp( subpc ); CHKERRQ(ierr); 664 ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); assert(ii==1); 665 ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 666 ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr); 667 /* set ops */ 668 ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 669 } 670 671 /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 672 ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 673 ierr = PCSetFromOptions_MG(pc); CHKERRQ(ierr); 674 ierr = PetscOptionsEnd();CHKERRQ(ierr); 675 if (mg->galerkin != 2) SETERRQ(wcomm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used."); 676 677 /* create cheby smoothers */ 678 for ( lidx = 1, level = pc_gamg->Nlevels-2; 679 lidx <= fine_level; 680 lidx++, level--) { 681 KSP smoother; 682 PetscBool isCheb; 683 684 ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr); 685 /* do my own cheby */ 686 ierr = PetscTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &isCheb ); CHKERRQ(ierr); 687 if( isCheb ) { 688 PetscReal emax, emin; 689 PC subpc; 690 691 ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 692 ierr = PetscTypeCompare( (PetscObject)subpc, PCJACOBI, &isCheb ); CHKERRQ(ierr); 693 if( isCheb && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */ 694 else{ /* eigen estimate 'emax' */ 695 KSP eksp; Mat Lmat = Aarr[level]; 696 Vec bb, xx; 697 698 ierr = MatGetVecs( Lmat, &bb, 0 ); CHKERRQ(ierr); 699 ierr = MatGetVecs( Lmat, &xx, 0 ); CHKERRQ(ierr); 700 { 701 PetscRandom rctx; 702 ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr); 703 ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 704 ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 705 ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr); 706 } 707 ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr); 708 ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 ); 709 CHKERRQ(ierr); 710 ierr = KSPSetNormType( eksp, KSP_NORM_NONE ); CHKERRQ(ierr); 711 712 ierr = KSPAppendOptionsPrefix( eksp, "est_"); CHKERRQ(ierr); 713 ierr = KSPSetFromOptions( eksp ); CHKERRQ(ierr); 714 715 ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr); 716 ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr ); 717 ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr); 718 719 { /* set PC type to be same as smoother - does not get all parameters!!! */ 720 const PCType type; PC pc; 721 ierr = PCGetType( subpc, &type ); CHKERRQ(ierr); 722 ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr ); 723 ierr = PCSetType( pc, type ); CHKERRQ(ierr); 724 } 725 726 /* solve - keep stuff out of logging */ 727 ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 728 ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 729 ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr); 730 ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 731 ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 732 733 ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr); 734 735 ierr = VecDestroy( &xx ); CHKERRQ(ierr); 736 ierr = VecDestroy( &bb ); CHKERRQ(ierr); 737 ierr = KSPDestroy( &eksp ); CHKERRQ(ierr); 738 739 if (pc_gamg->verbose) { 740 PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e\n",__FUNCT__,emax,emin); 741 } 742 } 743 { 744 PetscInt N1, N0, tt; 745 ierr = MatGetSize( Aarr[level], &N1, &tt ); CHKERRQ(ierr); 746 ierr = MatGetSize( Aarr[level+1], &N0, &tt ); CHKERRQ(ierr); 747 /* heuristic - is this crap? */ 748 emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); 749 emax *= 1.05; 750 } 751 ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr); 752 } 753 } 754 755 /* clean up */ 756 for(level=1;level<pc_gamg->Nlevels;level++){ 757 ierr = MatDestroy( &Parr[level] ); CHKERRQ(ierr); 758 ierr = MatDestroy( &Aarr[level] ); CHKERRQ(ierr); 759 } 760 761 ierr = PCSetUp_MG( pc );CHKERRQ( ierr ); 762 763 { 764 KSP smoother; /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */ 765 ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr); 766 ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr); 767 //ierr = KSPSetUp( smoother ); CHKERRQ(ierr); 768 } 769 } 770 else { 771 KSP smoother; 772 if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",mype,__FUNCT__); 773 ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr); 774 ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 775 ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr); 776 ierr = PCSetUp_MG( pc );CHKERRQ( ierr ); 777 } 778 779 PetscFunctionReturn(0); 780 } 781 782 /* ------------------------------------------------------------------------- */ 783 /* 784 PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 785 that was created with PCCreate_GAMG(). 786 787 Input Parameter: 788 . pc - the preconditioner context 789 790 Application Interface Routine: PCDestroy() 791 */ 792 #undef __FUNCT__ 793 #define __FUNCT__ "PCDestroy_GAMG" 794 PetscErrorCode PCDestroy_GAMG( PC pc ) 795 { 796 PetscErrorCode ierr; 797 PC_MG *mg = (PC_MG*)pc->data; 798 PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 799 800 PetscFunctionBegin; 801 ierr = PCReset_GAMG( pc );CHKERRQ(ierr); 802 ierr = PetscFree( pc_gamg );CHKERRQ(ierr); 803 ierr = PCDestroy_MG( pc );CHKERRQ(ierr); 804 PetscFunctionReturn(0); 805 } 806 807 808 #undef __FUNCT__ 809 #define __FUNCT__ "PCGAMGSetProcEqLim" 810 /*@ 811 PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via 812 processor reduction. 813 814 Not Collective on PC 815 816 Input Parameters: 817 . pc - the preconditioner context 818 819 820 Options Database Key: 821 . -pc_gamg_process_eq_limit 822 823 Level: intermediate 824 825 Concepts: Unstructured multrigrid preconditioner 826 827 .seealso: () 828 @*/ 829 PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 830 { 831 PetscErrorCode ierr; 832 833 PetscFunctionBegin; 834 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 835 ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 836 PetscFunctionReturn(0); 837 } 838 839 EXTERN_C_BEGIN 840 #undef __FUNCT__ 841 #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG" 842 PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 843 { 844 PC_MG *mg = (PC_MG*)pc->data; 845 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 846 847 PetscFunctionBegin; 848 if(n>0) pc_gamg->min_eq_proc = n; 849 PetscFunctionReturn(0); 850 } 851 EXTERN_C_END 852 853 #undef __FUNCT__ 854 #define __FUNCT__ "PCGAMGSetCoarseEqLim" 855 /*@ 856 PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids. 857 858 Collective on PC 859 860 Input Parameters: 861 . pc - the preconditioner context 862 863 864 Options Database Key: 865 . -pc_gamg_coarse_eq_limit 866 867 Level: intermediate 868 869 Concepts: Unstructured multrigrid preconditioner 870 871 .seealso: () 872 @*/ 873 PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 874 { 875 PetscErrorCode ierr; 876 877 PetscFunctionBegin; 878 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 879 ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 880 PetscFunctionReturn(0); 881 } 882 883 EXTERN_C_BEGIN 884 #undef __FUNCT__ 885 #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG" 886 PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 887 { 888 PC_MG *mg = (PC_MG*)pc->data; 889 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 890 891 PetscFunctionBegin; 892 if(n>0) pc_gamg->coarse_eq_limit = n; 893 PetscFunctionReturn(0); 894 } 895 EXTERN_C_END 896 897 #undef __FUNCT__ 898 #define __FUNCT__ "PCGAMGSetRepartitioning" 899 /*@ 900 PCGAMGSetRepartitioning - Repartition the coarse grids 901 902 Collective on PC 903 904 Input Parameters: 905 . pc - the preconditioner context 906 907 908 Options Database Key: 909 . -pc_gamg_repartition 910 911 Level: intermediate 912 913 Concepts: Unstructured multrigrid preconditioner 914 915 .seealso: () 916 @*/ 917 PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n) 918 { 919 PetscErrorCode ierr; 920 921 PetscFunctionBegin; 922 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 923 ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 924 PetscFunctionReturn(0); 925 } 926 927 EXTERN_C_BEGIN 928 #undef __FUNCT__ 929 #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG" 930 PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n) 931 { 932 PC_MG *mg = (PC_MG*)pc->data; 933 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 934 935 PetscFunctionBegin; 936 pc_gamg->repart = n; 937 PetscFunctionReturn(0); 938 } 939 EXTERN_C_END 940 941 #undef __FUNCT__ 942 #define __FUNCT__ "PCGAMGSetNlevels" 943 /*@ 944 PCGAMGSetNlevels - 945 946 Not collective on PC 947 948 Input Parameters: 949 . pc - the preconditioner context 950 951 952 Options Database Key: 953 . -pc_mg_levels 954 955 Level: intermediate 956 957 Concepts: Unstructured multrigrid preconditioner 958 959 .seealso: () 960 @*/ 961 PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 962 { 963 PetscErrorCode ierr; 964 965 PetscFunctionBegin; 966 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 967 ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 968 PetscFunctionReturn(0); 969 } 970 971 EXTERN_C_BEGIN 972 #undef __FUNCT__ 973 #define __FUNCT__ "PCGAMGSetNlevels_GAMG" 974 PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 975 { 976 PC_MG *mg = (PC_MG*)pc->data; 977 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 978 979 PetscFunctionBegin; 980 pc_gamg->Nlevels = n; 981 PetscFunctionReturn(0); 982 } 983 EXTERN_C_END 984 985 #undef __FUNCT__ 986 #define __FUNCT__ "PCGAMGSetThreshold" 987 /*@ 988 PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 989 990 Not collective on PC 991 992 Input Parameters: 993 . pc - the preconditioner context 994 995 996 Options Database Key: 997 . -pc_gamg_threshold 998 999 Level: intermediate 1000 1001 Concepts: Unstructured multrigrid preconditioner 1002 1003 .seealso: () 1004 @*/ 1005 PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n) 1006 { 1007 PetscErrorCode ierr; 1008 1009 PetscFunctionBegin; 1010 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1011 ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr); 1012 PetscFunctionReturn(0); 1013 } 1014 1015 EXTERN_C_BEGIN 1016 #undef __FUNCT__ 1017 #define __FUNCT__ "PCGAMGSetThreshold_GAMG" 1018 PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n) 1019 { 1020 PC_MG *mg = (PC_MG*)pc->data; 1021 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1022 1023 PetscFunctionBegin; 1024 pc_gamg->threshold = n; 1025 PetscFunctionReturn(0); 1026 } 1027 EXTERN_C_END 1028 1029 #undef __FUNCT__ 1030 #define __FUNCT__ "PCGAMGSetType" 1031 /*@ 1032 PCGAMGSetType - Set solution method - calls sub create method 1033 1034 Collective on PC 1035 1036 Input Parameters: 1037 . pc - the preconditioner context 1038 1039 1040 Options Database Key: 1041 . -pc_gamg_type 1042 1043 Level: intermediate 1044 1045 Concepts: Unstructured multrigrid preconditioner 1046 1047 .seealso: () 1048 @*/ 1049 PetscErrorCode PCGAMGSetType( PC pc, const PCGAMGType type ) 1050 { 1051 PetscErrorCode ierr; 1052 1053 PetscFunctionBegin; 1054 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1055 ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,const PCGAMGType),(pc,type)); 1056 CHKERRQ(ierr); 1057 PetscFunctionReturn(0); 1058 } 1059 1060 EXTERN_C_BEGIN 1061 #undef __FUNCT__ 1062 #define __FUNCT__ "PCGAMGSetType_GAMG" 1063 PetscErrorCode PCGAMGSetType_GAMG( PC pc, const PCGAMGType type ) 1064 { 1065 PetscErrorCode ierr,(*r)(PC); 1066 1067 PetscFunctionBegin; 1068 ierr = PetscFListFind(GAMGList,((PetscObject)pc)->comm,type,PETSC_FALSE,(PetscVoidStarFunction)&r); 1069 CHKERRQ(ierr); 1070 1071 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 1072 1073 /* call sub create method */ 1074 ierr = (*r)(pc); CHKERRQ(ierr); 1075 1076 PetscFunctionReturn(0); 1077 } 1078 EXTERN_C_END 1079 1080 #undef __FUNCT__ 1081 #define __FUNCT__ "PCSetFromOptions_GAMG" 1082 PetscErrorCode PCSetFromOptions_GAMG( PC pc ) 1083 { 1084 PetscErrorCode ierr; 1085 PC_MG *mg = (PC_MG*)pc->data; 1086 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1087 PetscBool flag; 1088 MPI_Comm wcomm = ((PetscObject)pc)->comm; 1089 1090 PetscFunctionBegin; 1091 ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr); 1092 { 1093 /* -pc_gamg_verbose */ 1094 ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG", 1095 "none", pc_gamg->verbose, 1096 &pc_gamg->verbose, PETSC_NULL ); 1097 CHKERRQ(ierr); 1098 1099 /* -pc_gamg_repartition */ 1100 ierr = PetscOptionsBool("-pc_gamg_repartition", 1101 "Repartion coarse grids (false)", 1102 "PCGAMGRepartitioning", 1103 pc_gamg->repart, 1104 &pc_gamg->repart, 1105 &flag); 1106 CHKERRQ(ierr); 1107 1108 /* -pc_gamg_process_eq_limit */ 1109 ierr = PetscOptionsInt("-pc_gamg_process_eq_limit", 1110 "Limit (goal) on number of equations per process on coarse grids", 1111 "PCGAMGSetProcEqLim", 1112 pc_gamg->min_eq_proc, 1113 &pc_gamg->min_eq_proc, 1114 &flag ); 1115 CHKERRQ(ierr); 1116 1117 /* -pc_gamg_coarse_eq_limit */ 1118 ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit", 1119 "Limit on number of equations for the coarse grid", 1120 "PCGAMGSetCoarseEqLim", 1121 pc_gamg->coarse_eq_limit, 1122 &pc_gamg->coarse_eq_limit, 1123 &flag ); 1124 CHKERRQ(ierr); 1125 1126 /* -pc_gamg_threshold */ 1127 ierr = PetscOptionsReal("-pc_gamg_threshold", 1128 "Relative threshold to use for dropping edges in aggregation graph", 1129 "PCGAMGSetThreshold", 1130 pc_gamg->threshold, 1131 &pc_gamg->threshold, 1132 &flag ); 1133 CHKERRQ(ierr); 1134 if(flag && pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold); 1135 1136 ierr = PetscOptionsInt("-pc_mg_levels", 1137 "Set number of MG levels", 1138 "PCGAMGSetNlevels", 1139 pc_gamg->Nlevels, 1140 &pc_gamg->Nlevels, 1141 &flag ); 1142 } 1143 ierr = PetscOptionsTail();CHKERRQ(ierr); 1144 1145 PetscFunctionReturn(0); 1146 } 1147 1148 /* -------------------------------------------------------------------------- */ 1149 /*MC 1150 PCCreate_GAMG - Geometric algebraic multigrid (AMG) preconditioning framework. 1151 - This is the entry point to GAMG, registered in pcregis.c 1152 1153 Options Database Keys: 1154 Multigrid options(inherited) 1155 + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType) 1156 . -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp) 1157 . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown) 1158 - -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 1159 1160 Level: intermediate 1161 1162 Concepts: multigrid 1163 1164 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 1165 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 1166 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 1167 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 1168 PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 1169 M*/ 1170 EXTERN_C_BEGIN 1171 #undef __FUNCT__ 1172 #define __FUNCT__ "PCCreate_GAMG" 1173 PetscErrorCode PCCreate_GAMG( PC pc ) 1174 { 1175 PetscErrorCode ierr; 1176 PC_GAMG *pc_gamg; 1177 PC_MG *mg; 1178 PetscClassId cookie; 1179 1180 #if defined PETSC_USE_LOG 1181 static long count = 0; 1182 #endif 1183 1184 PetscFunctionBegin; 1185 /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 1186 ierr = PCSetType( pc, PCMG ); CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 1187 ierr = PetscObjectChangeTypeName( (PetscObject)pc, PCGAMG ); CHKERRQ(ierr); 1188 1189 /* create a supporting struct and attach it to pc */ 1190 ierr = PetscNewLog( pc, PC_GAMG, &pc_gamg ); CHKERRQ(ierr); 1191 mg = (PC_MG*)pc->data; 1192 mg->galerkin = 2; /* Use Galerkin, but it is computed externally */ 1193 mg->innerctx = pc_gamg; 1194 1195 pc_gamg->setup_count = 0; 1196 /* these should be in subctx but repartitioning needs simple arrays */ 1197 pc_gamg->data_sz = 0; 1198 pc_gamg->data = 0; 1199 1200 /* register AMG type */ 1201 if( !GAMGList ){ 1202 ierr = PetscFListAdd(&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr); 1203 ierr = PetscFListAdd(&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr); 1204 } 1205 1206 /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 1207 pc->ops->setfromoptions = PCSetFromOptions_GAMG; 1208 pc->ops->setup = PCSetUp_GAMG; 1209 pc->ops->reset = PCReset_GAMG; 1210 pc->ops->destroy = PCDestroy_GAMG; 1211 1212 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1213 "PCGAMGSetProcEqLim_C", 1214 "PCGAMGSetProcEqLim_GAMG", 1215 PCGAMGSetProcEqLim_GAMG); 1216 CHKERRQ(ierr); 1217 1218 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1219 "PCGAMGSetCoarseEqLim_C", 1220 "PCGAMGSetCoarseEqLim_GAMG", 1221 PCGAMGSetCoarseEqLim_GAMG); 1222 CHKERRQ(ierr); 1223 1224 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1225 "PCGAMGSetRepartitioning_C", 1226 "PCGAMGSetRepartitioning_GAMG", 1227 PCGAMGSetRepartitioning_GAMG); 1228 CHKERRQ(ierr); 1229 1230 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1231 "PCGAMGSetThreshold_C", 1232 "PCGAMGSetThreshold_GAMG", 1233 PCGAMGSetThreshold_GAMG); 1234 CHKERRQ(ierr); 1235 1236 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1237 "PCGAMGSetType_C", 1238 "PCGAMGSetType_GAMG", 1239 PCGAMGSetType_GAMG); 1240 CHKERRQ(ierr); 1241 1242 pc_gamg->repart = PETSC_FALSE; 1243 pc_gamg->min_eq_proc = 100; 1244 pc_gamg->coarse_eq_limit = 800; 1245 pc_gamg->threshold = 0.05; 1246 pc_gamg->Nlevels = GAMG_MAXLEVELS; 1247 pc_gamg->verbose = 0; 1248 pc_gamg->emax_id = -1; 1249 pc_gamg->col_bs_id = -1; 1250 1251 /* instantiate derived type */ 1252 ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr); 1253 { 1254 char tname[256] = GAMGAGG; 1255 ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType", 1256 GAMGList, tname, tname, sizeof(tname), PETSC_NULL ); 1257 CHKERRQ(ierr); 1258 ierr = PCGAMGSetType( pc, tname ); CHKERRQ(ierr); 1259 } 1260 ierr = PetscOptionsTail();CHKERRQ(ierr); 1261 1262 #if defined PETSC_USE_LOG 1263 if( count++ == 0 ) { 1264 PetscClassIdRegister("GAMG Setup",&cookie); 1265 PetscLogEventRegister("GAMG: createProl", cookie, &gamg_setup_events[SET1]); 1266 PetscLogEventRegister(" Graph", cookie, &gamg_setup_events[GRAPH]); 1267 /* PetscLogEventRegister(" G.Mat", cookie, &gamg_setup_events[GRAPH_MAT]); */ 1268 /* PetscLogEventRegister(" G.Filter", cookie, &gamg_setup_events[GRAPH_FILTER]); */ 1269 /* PetscLogEventRegister(" G.Square", cookie, &gamg_setup_events[GRAPH_SQR]); */ 1270 PetscLogEventRegister(" MIS/Agg", cookie, &gamg_setup_events[SET4]); 1271 PetscLogEventRegister(" geo: growSupp", cookie, &gamg_setup_events[SET5]); 1272 PetscLogEventRegister(" geo: triangle", cookie, &gamg_setup_events[SET6]); 1273 PetscLogEventRegister(" search&set", cookie, &gamg_setup_events[FIND_V]); 1274 PetscLogEventRegister(" SA: colect data", cookie, &gamg_setup_events[SET7]); 1275 PetscLogEventRegister(" SA: frmProl0", cookie, &gamg_setup_events[SET8]); 1276 PetscLogEventRegister(" SA: smooth", cookie, &gamg_setup_events[SET9]); 1277 PetscLogEventRegister("GAMG: partLevel", cookie, &gamg_setup_events[SET2]); 1278 PetscLogEventRegister(" repartition", cookie, &gamg_setup_events[SET12]); 1279 PetscLogEventRegister(" Invert-Sort", cookie, &gamg_setup_events[SET13]); 1280 PetscLogEventRegister(" Move A", cookie, &gamg_setup_events[SET14]); 1281 PetscLogEventRegister(" Move P", cookie, &gamg_setup_events[SET15]); 1282 1283 /* PetscLogEventRegister(" PL move data", cookie, &gamg_setup_events[SET13]); */ 1284 /* PetscLogEventRegister("GAMG: fix", cookie, &gamg_setup_events[SET10]); */ 1285 /* PetscLogEventRegister("GAMG: set levels", cookie, &gamg_setup_events[SET11]); */ 1286 /* create timer stages */ 1287 #if defined GAMG_STAGES 1288 { 1289 char str[32]; 1290 sprintf(str,"MG Level %d (finest)",0); 1291 PetscLogStageRegister(str, &gamg_stages[0]); 1292 PetscInt lidx; 1293 for (lidx=1;lidx<9;lidx++){ 1294 sprintf(str,"MG Level %d",lidx); 1295 PetscLogStageRegister(str, &gamg_stages[lidx]); 1296 } 1297 } 1298 #endif 1299 } 1300 #endif 1301 PetscFunctionReturn(0); 1302 } 1303 EXTERN_C_END 1304