1 /* 2 GAMG geometric-algebric multigrid PC - Mark Adams 2011 3 */ 4 #include "petsc-private/matimpl.h" 5 #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 6 #include <petsc-private/kspimpl.h> 7 8 #if defined PETSC_GAMG_USE_LOG 9 PetscLogEvent petsc_gamg_setup_events[NUM_SET]; 10 #endif 11 12 #if defined PETSC_USE_LOG 13 PetscLogEvent PC_GAMGGgraph_AGG; 14 PetscLogEvent PC_GAMGGgraph_GEO; 15 PetscLogEvent PC_GAMGCoarsen_AGG; 16 PetscLogEvent PC_GAMGCoarsen_GEO; 17 PetscLogEvent PC_GAMGProlongator_AGG; 18 PetscLogEvent PC_GAMGProlongator_GEO; 19 PetscLogEvent PC_GAMGOptprol_AGG; 20 PetscLogEvent PC_GAMGKKTProl_AGG; 21 #endif 22 23 #define GAMG_MAXLEVELS 30 24 25 /* #define GAMG_STAGES */ 26 #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 27 static PetscLogStage gamg_stages[GAMG_MAXLEVELS]; 28 #endif 29 30 static PetscFList GAMGList = 0; 31 32 /* ----------------------------------------------------------------------------- */ 33 #undef __FUNCT__ 34 #define __FUNCT__ "PCReset_GAMG" 35 PetscErrorCode PCReset_GAMG(PC pc) 36 { 37 PetscErrorCode ierr; 38 PC_MG *mg = (PC_MG*)pc->data; 39 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 40 41 PetscFunctionBegin; 42 if( pc_gamg->data ) { /* this should not happen, cleaned up in SetUp */ 43 ierr = PetscFree( pc_gamg->data ); CHKERRQ(ierr); 44 } 45 pc_gamg->data = PETSC_NULL; pc_gamg->data_sz = 0; 46 PetscFunctionReturn(0); 47 } 48 49 /* private 2x2 Mat Nest for Stokes */ 50 typedef struct{ 51 Mat A11,A21,A12,Amat; 52 IS prim_is,constr_is; 53 }GAMGKKTMat; 54 55 #undef __FUNCT__ 56 #define __FUNCT__ "GAMGKKTMatCreate" 57 static PetscErrorCode GAMGKKTMatCreate( Mat A, PetscBool iskkt, GAMGKKTMat *out ) 58 { 59 PetscFunctionBegin; 60 out->Amat = A; 61 if( iskkt ){ 62 PetscErrorCode ierr; 63 IS is_constraint, is_prime; 64 PetscInt nmin,nmax; 65 66 ierr = MatGetOwnershipRange( A, &nmin, &nmax ); CHKERRQ(ierr); 67 ierr = MatFindZeroDiagonals( A, &is_constraint ); CHKERRQ(ierr); 68 ierr = ISComplement( is_constraint, nmin, nmax, &is_prime ); CHKERRQ(ierr); 69 out->prim_is = is_prime; 70 out->constr_is = is_constraint; 71 72 ierr = MatGetSubMatrix( A, is_prime, is_prime, MAT_INITIAL_MATRIX, &out->A11); CHKERRQ(ierr); 73 ierr = MatGetSubMatrix( A, is_prime, is_constraint, MAT_INITIAL_MATRIX, &out->A12); CHKERRQ(ierr); 74 ierr = MatGetSubMatrix( A, is_constraint, is_prime, MAT_INITIAL_MATRIX, &out->A21); CHKERRQ(ierr); 75 } 76 else { 77 out->A11 = A; 78 out->A21 = PETSC_NULL; 79 out->A12 = PETSC_NULL; 80 out->prim_is = PETSC_NULL; 81 out->constr_is = PETSC_NULL; 82 } 83 PetscFunctionReturn(0); 84 } 85 86 #undef __FUNCT__ 87 #define __FUNCT__ "GAMGKKTMatDestroy" 88 static PetscErrorCode GAMGKKTMatDestroy( GAMGKKTMat *mat ) 89 { 90 PetscErrorCode ierr; 91 92 PetscFunctionBegin; 93 if( mat->A11 && mat->A11 != mat->Amat ) { 94 ierr = MatDestroy( &mat->A11 ); CHKERRQ(ierr); 95 } 96 ierr = MatDestroy( &mat->A21 ); CHKERRQ(ierr); 97 ierr = MatDestroy( &mat->A12 ); CHKERRQ(ierr); 98 99 ierr = ISDestroy( &mat->prim_is ); CHKERRQ(ierr); 100 ierr = ISDestroy( &mat->constr_is ); CHKERRQ(ierr); 101 102 PetscFunctionReturn(0); 103 } 104 105 /* -------------------------------------------------------------------------- */ 106 /* 107 createLevel: create coarse op with RAP. repartition and/or reduce number 108 of active processors. 109 110 Input Parameter: 111 . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 112 'pc_gamg->data_sz' are changed via repartitioning/reduction. 113 . Amat_fine - matrix on this fine (k) level 114 . cbs - coarse block size 115 . isLast - 116 . stokes - 117 In/Output Parameter: 118 . a_P_inout - prolongation operator to the next level (k-->k-1) 119 . a_nactive_proc - number of active procs 120 Output Parameter: 121 . a_Amat_crs - coarse matrix that is created (k-1) 122 */ 123 124 #undef __FUNCT__ 125 #define __FUNCT__ "createLevel" 126 static PetscErrorCode createLevel( const PC pc, 127 const Mat Amat_fine, 128 const PetscInt cbs, 129 const PetscBool isLast, 130 const PetscBool stokes, 131 Mat *a_P_inout, 132 Mat *a_Amat_crs, 133 PetscMPIInt *a_nactive_proc 134 ) 135 { 136 PetscErrorCode ierr; 137 PC_MG *mg = (PC_MG*)pc->data; 138 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 139 const PetscBool repart = pc_gamg->repart; 140 const PetscInt min_eq_proc = pc_gamg->min_eq_proc, coarse_max = pc_gamg->coarse_eq_limit; 141 Mat Cmat,Pold=*a_P_inout; 142 MPI_Comm wcomm = ((PetscObject)Amat_fine)->comm; 143 PetscMPIInt mype,npe,new_npe,nactive=*a_nactive_proc; 144 PetscInt ncrs_eq,ncrs_prim; 145 146 PetscFunctionBegin; 147 ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr); 148 ierr = MPI_Comm_size( wcomm, &npe ); CHKERRQ(ierr); 149 /* RAP */ 150 ierr = MatPtAP( Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat ); CHKERRQ(ierr); 151 152 /* set 'ncrs_prim' (nodes), 'ncrs_eq' (equations)*/ 153 ncrs_prim = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows; 154 assert(pc_gamg->data_sz%(pc_gamg->data_cell_cols*pc_gamg->data_cell_rows)==0); 155 ierr = MatGetLocalSize( Cmat, &ncrs_eq, PETSC_NULL ); CHKERRQ(ierr); 156 157 /* get number of PEs to make active 'new_npe', reduce, can be any integer 1-P */ 158 { 159 PetscInt ncrs_eq_glob,ncrs_eq_ave; 160 ierr = MatGetSize( Cmat, &ncrs_eq_glob, PETSC_NULL ); CHKERRQ(ierr); 161 ncrs_eq_ave = ncrs_eq_glob/npe; 162 new_npe = (PetscMPIInt)((float)ncrs_eq_ave/(float)min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 163 if( new_npe == 0 || ncrs_eq_ave < coarse_max ) new_npe = 1; 164 else if ( new_npe >= nactive ) new_npe = nactive; /* no change, rare */ 165 if( isLast ) new_npe = 1; 166 } 167 168 if( !repart && new_npe==nactive ) { 169 *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 170 } 171 else { 172 const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols; 173 PetscInt *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_prim_new,ncrs_eq_new,nloc_old; 174 IS is_eq_newproc,is_eq_newproc_prim,is_eq_num,is_eq_num_prim,isscat,new_eq_indices; 175 VecScatter vecscat; 176 PetscScalar *array; 177 Vec src_crd, dest_crd; 178 179 nloc_old = ncrs_eq/cbs; assert(ncrs_eq%cbs==0); 180 #if defined PETSC_GAMG_USE_LOG 181 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 182 #endif 183 /* make 'is_eq_newproc' */ 184 ierr = PetscMalloc( npe*sizeof(PetscInt), &counts ); CHKERRQ(ierr); 185 if( repart && !stokes ) { 186 /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */ 187 Mat adj; 188 189 if (pc_gamg->verbose>0) { 190 if (pc_gamg->verbose==1) PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,ncrs_eq); 191 else { 192 PetscInt n; 193 ierr = MPI_Allreduce( &ncrs_eq, &n, 1, MPIU_INT, MPI_SUM, wcomm );CHKERRQ(ierr); 194 PetscPrintf(wcomm,"\t[%d]%s repartition: npe (active): %d --> %d, neq = %d\n",mype,__FUNCT__,*a_nactive_proc,new_npe,n); 195 } 196 } 197 198 /* get 'adj' */ 199 if( cbs == 1 ) { 200 ierr = MatConvert( Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 201 } 202 else{ 203 /* make a scalar matrix to partition (no Stokes here) */ 204 Mat tMat; 205 PetscInt Istart_crs,Iend_crs,ncols,jj,Ii; 206 const PetscScalar *vals; 207 const PetscInt *idx; 208 PetscInt *d_nnz, *o_nnz, M, N; 209 static PetscInt llev = 0; 210 211 ierr = PetscMalloc( ncrs_prim*sizeof(PetscInt), &d_nnz ); CHKERRQ(ierr); 212 ierr = PetscMalloc( ncrs_prim*sizeof(PetscInt), &o_nnz ); CHKERRQ(ierr); 213 ierr = MatGetOwnershipRange( Cmat, &Istart_crs, &Iend_crs ); CHKERRQ(ierr); 214 ierr = MatGetSize( Cmat, &M, &N );CHKERRQ(ierr); 215 for ( Ii = Istart_crs, jj = 0 ; Ii < Iend_crs ; Ii += cbs, jj++ ) { 216 ierr = MatGetRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr); 217 d_nnz[jj] = ncols/cbs; 218 o_nnz[jj] = ncols/cbs; 219 ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0); CHKERRQ(ierr); 220 if( d_nnz[jj] > ncrs_prim ) d_nnz[jj] = ncrs_prim; 221 if( o_nnz[jj] > (M/cbs-ncrs_prim) ) o_nnz[jj] = M/cbs-ncrs_prim; 222 } 223 224 ierr = MatCreate( wcomm, &tMat ); CHKERRQ(ierr); 225 ierr = MatSetSizes( tMat, ncrs_prim, ncrs_prim, 226 PETSC_DETERMINE, PETSC_DETERMINE ); 227 CHKERRQ(ierr); 228 ierr = MatSetType(tMat,MATAIJ); CHKERRQ(ierr); 229 ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr); 230 ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 231 ierr = PetscFree( d_nnz ); CHKERRQ(ierr); 232 ierr = PetscFree( o_nnz ); CHKERRQ(ierr); 233 234 for ( ii = Istart_crs; ii < Iend_crs; ii++ ) { 235 PetscInt dest_row = ii/cbs; 236 ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr); 237 for( jj = 0 ; jj < ncols ; jj++ ){ 238 PetscInt dest_col = idx[jj]/cbs; 239 PetscScalar v = 1.0; 240 ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr); 241 } 242 ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals); CHKERRQ(ierr); 243 } 244 ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 245 ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 246 247 if( llev++ == -1 ) { 248 PetscViewer viewer; char fname[32]; 249 ierr = PetscSNPrintf(fname,sizeof fname,"part_mat_%D.mat",llev);CHKERRQ(ierr); 250 PetscViewerBinaryOpen(wcomm,fname,FILE_MODE_WRITE,&viewer); 251 ierr = MatView( tMat, viewer ); CHKERRQ(ierr); 252 ierr = PetscViewerDestroy( &viewer ); 253 } 254 255 ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 256 257 ierr = MatDestroy( &tMat ); CHKERRQ(ierr); 258 } /* create 'adj' */ 259 260 { /* partition: get newproc_idx */ 261 char prefix[256]; 262 const char *pcpre; 263 const PetscInt *is_idx; 264 MatPartitioning mpart; 265 IS proc_is; 266 PetscInt targetPE; 267 268 ierr = MatPartitioningCreate( wcomm, &mpart ); CHKERRQ(ierr); 269 ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr); 270 ierr = PCGetOptionsPrefix( pc, &pcpre );CHKERRQ(ierr); 271 ierr = PetscSNPrintf(prefix,sizeof prefix,"%spc_gamg_",pcpre?pcpre:"");CHKERRQ(ierr); 272 ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 273 ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr); 274 ierr = MatPartitioningSetNParts( mpart, new_npe );CHKERRQ(ierr); 275 ierr = MatPartitioningApply( mpart, &proc_is ); CHKERRQ(ierr); 276 ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr); 277 278 /* collect IS info */ 279 ierr = PetscMalloc( ncrs_eq*sizeof(PetscInt), &newproc_idx ); CHKERRQ(ierr); 280 ierr = ISGetIndices( proc_is, &is_idx ); CHKERRQ(ierr); 281 targetPE = 1; /* bring to "front" of machine */ 282 /*targetPE = npe/new_npe;*/ /* spread partitioning across machine */ 283 for( kk = jj = 0 ; kk < nloc_old ; kk++ ){ 284 for( ii = 0 ; ii < cbs ; ii++, jj++ ){ 285 newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */ 286 } 287 } 288 ierr = ISRestoreIndices( proc_is, &is_idx ); CHKERRQ(ierr); 289 ierr = ISDestroy( &proc_is ); CHKERRQ(ierr); 290 } 291 ierr = MatDestroy( &adj ); CHKERRQ(ierr); 292 293 ierr = ISCreateGeneral( wcomm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc ); 294 CHKERRQ(ierr); 295 if( newproc_idx != 0 ) { 296 ierr = PetscFree( newproc_idx ); CHKERRQ(ierr); 297 } 298 } /* repartitioning */ 299 else { /* simple aggreagtion of parts -- 'is_eq_newproc' */ 300 301 PetscInt rfactor,targetPE; 302 /* find factor */ 303 if( new_npe == 1 ) rfactor = npe; /* easy */ 304 else { 305 PetscReal best_fact = 0.; 306 jj = -1; 307 for( kk = 1 ; kk <= npe ; kk++ ){ 308 if( npe%kk==0 ) { /* a candidate */ 309 PetscReal nactpe = (PetscReal)npe/(PetscReal)kk, fact = nactpe/(PetscReal)new_npe; 310 if(fact > 1.0) fact = 1./fact; /* keep fact < 1 */ 311 if( fact > best_fact ) { 312 best_fact = fact; jj = kk; 313 } 314 } 315 } 316 if( jj != -1 ) rfactor = jj; 317 else rfactor = 1; /* does this happen .. a prime */ 318 } 319 new_npe = npe/rfactor; 320 321 if( new_npe==nactive ) { 322 *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */ 323 ierr = PetscFree( counts ); CHKERRQ(ierr); 324 if (pc_gamg->verbose>0){ 325 PetscPrintf(wcomm,"\t[%d]%s aggregate processors noop: new_npe=%d, neq(loc)=%d\n",mype,__FUNCT__,new_npe,ncrs_eq); 326 } 327 PetscFunctionReturn(0); 328 } 329 330 if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s number of equations (loc) %d with simple aggregation\n",mype,__FUNCT__,ncrs_eq); 331 targetPE = mype/rfactor; 332 ierr = ISCreateStride( wcomm, ncrs_eq, targetPE, 0, &is_eq_newproc ); CHKERRQ(ierr); 333 334 if( stokes ) { 335 ierr = ISCreateStride( wcomm, ncrs_prim*cbs, targetPE, 0, &is_eq_newproc_prim ); CHKERRQ(ierr); 336 } 337 } /* end simple 'is_eq_newproc' */ 338 339 /* 340 Create an index set from the is_eq_newproc index set to indicate the mapping TO 341 */ 342 ierr = ISPartitioningToNumbering( is_eq_newproc, &is_eq_num ); CHKERRQ(ierr); 343 if( stokes ) { 344 ierr = ISPartitioningToNumbering( is_eq_newproc_prim, &is_eq_num_prim ); CHKERRQ(ierr); 345 } 346 else is_eq_num_prim = is_eq_num; 347 /* 348 Determine how many equations/vertices are assigned to each processor 349 */ 350 ierr = ISPartitioningCount( is_eq_newproc, npe, counts ); CHKERRQ(ierr); 351 ncrs_eq_new = counts[mype]; 352 ierr = ISDestroy( &is_eq_newproc ); CHKERRQ(ierr); 353 if( stokes ) { 354 ierr = ISPartitioningCount( is_eq_newproc_prim, npe, counts ); CHKERRQ(ierr); 355 ierr = ISDestroy( &is_eq_newproc_prim ); CHKERRQ(ierr); 356 ncrs_prim_new = counts[mype]/cbs; /* nodes */ 357 } 358 else ncrs_prim_new = ncrs_eq_new/cbs; /* eqs */ 359 360 ierr = PetscFree( counts ); CHKERRQ(ierr); 361 #if defined PETSC_GAMG_USE_LOG 362 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0); CHKERRQ(ierr); 363 #endif 364 365 /* move data (for primal equations only) */ 366 /* Create a vector to contain the newly ordered element information */ 367 ierr = VecCreate( wcomm, &dest_crd ); 368 ierr = VecSetSizes( dest_crd, node_data_sz*ncrs_prim_new, PETSC_DECIDE ); CHKERRQ(ierr); 369 ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /* this is needed! */ 370 /* 371 There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 372 a block size of ...). Note, ISs are expanded into equation space by 'cbs'. 373 */ 374 ierr = PetscMalloc( (ncrs_prim*node_data_sz)*sizeof(PetscInt), &tidx ); CHKERRQ(ierr); 375 ierr = ISGetIndices( is_eq_num_prim, &idx ); CHKERRQ(ierr); 376 for(ii=0,jj=0; ii<ncrs_prim ; ii++) { 377 PetscInt id = idx[ii*cbs]/cbs; /* get node back */ 378 for( kk=0; kk<node_data_sz ; kk++, jj++) tidx[jj] = id*node_data_sz + kk; 379 } 380 ierr = ISRestoreIndices( is_eq_num_prim, &idx ); CHKERRQ(ierr); 381 ierr = ISCreateGeneral( wcomm, node_data_sz*ncrs_prim, tidx, PETSC_COPY_VALUES, &isscat ); 382 CHKERRQ(ierr); 383 ierr = PetscFree( tidx ); CHKERRQ(ierr); 384 /* 385 Create a vector to contain the original vertex information for each element 386 */ 387 ierr = VecCreateSeq( PETSC_COMM_SELF, node_data_sz*ncrs_prim, &src_crd ); CHKERRQ(ierr); 388 for( jj=0; jj<ndata_cols ; jj++ ) { 389 const PetscInt stride0=ncrs_prim*pc_gamg->data_cell_rows; 390 for( ii=0 ; ii<ncrs_prim ; ii++) { 391 for( kk=0; kk<ndata_rows ; kk++ ) { 392 PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj; 393 PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 394 ierr = VecSetValues( src_crd, 1, &jx, &tt, INSERT_VALUES ); CHKERRQ(ierr); 395 } 396 } 397 } 398 ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr); 399 ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr); 400 /* 401 Scatter the element vertex information (still in the original vertex ordering) 402 to the correct processor 403 */ 404 ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat); 405 CHKERRQ(ierr); 406 ierr = ISDestroy( &isscat ); CHKERRQ(ierr); 407 ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 408 ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 409 ierr = VecScatterDestroy( &vecscat ); CHKERRQ(ierr); 410 ierr = VecDestroy( &src_crd ); CHKERRQ(ierr); 411 /* 412 Put the element vertex data into a new allocation of the gdata->ele 413 */ 414 ierr = PetscFree( pc_gamg->data ); CHKERRQ(ierr); 415 ierr = PetscMalloc( node_data_sz*ncrs_prim_new*sizeof(PetscReal), &pc_gamg->data ); CHKERRQ(ierr); 416 pc_gamg->data_sz = node_data_sz*ncrs_prim_new; 417 strideNew = ncrs_prim_new*ndata_rows; 418 ierr = VecGetArray( dest_crd, &array ); CHKERRQ(ierr); 419 for( jj=0; jj<ndata_cols ; jj++ ) { 420 for( ii=0 ; ii<ncrs_prim_new ; ii++) { 421 for( kk=0; kk<ndata_rows ; kk++ ) { 422 PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj; 423 pc_gamg->data[ix] = PetscRealPart(array[jx]); 424 } 425 } 426 } 427 ierr = VecRestoreArray( dest_crd, &array ); CHKERRQ(ierr); 428 ierr = VecDestroy( &dest_crd ); CHKERRQ(ierr); 429 430 /* move A and P (columns) with new layout */ 431 #if defined PETSC_GAMG_USE_LOG 432 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 433 #endif 434 435 /* 436 Invert for MatGetSubMatrix 437 */ 438 ierr = ISInvertPermutation( is_eq_num, ncrs_eq_new, &new_eq_indices ); CHKERRQ(ierr); 439 ierr = ISSort( new_eq_indices ); CHKERRQ(ierr); /* is this needed? */ 440 if(is_eq_num != is_eq_num_prim) { 441 ierr = ISDestroy( &is_eq_num_prim ); CHKERRQ(ierr); /* could be same as 'is_eq_num' */ 442 } 443 ierr = ISDestroy( &is_eq_num ); CHKERRQ(ierr); 444 #if defined PETSC_GAMG_USE_LOG 445 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 446 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 447 #endif 448 /* 'a_Amat_crs' output */ 449 { 450 Mat mat; 451 ierr = MatGetSubMatrix( Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat ); 452 CHKERRQ(ierr); 453 *a_Amat_crs = mat; 454 } 455 ierr = MatDestroy( &Cmat ); CHKERRQ(ierr); 456 457 #if defined PETSC_GAMG_USE_LOG 458 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 459 #endif 460 /* prolongator */ 461 { 462 IS findices; 463 PetscInt Istart,Iend; 464 Mat Pnew; 465 ierr = MatGetOwnershipRange( Pold, &Istart, &Iend ); CHKERRQ(ierr); 466 #if defined PETSC_GAMG_USE_LOG 467 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 468 #endif 469 ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices); CHKERRQ(ierr); 470 ierr = MatGetSubMatrix( Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew ); 471 CHKERRQ(ierr); 472 ierr = ISDestroy( &findices ); CHKERRQ(ierr); 473 #if defined PETSC_GAMG_USE_LOG 474 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 475 #endif 476 ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr); 477 478 /* output - repartitioned */ 479 *a_P_inout = Pnew; 480 } 481 ierr = ISDestroy( &new_eq_indices ); CHKERRQ(ierr); 482 483 *a_nactive_proc = new_npe; /* output */ 484 } 485 486 /* outout matrix data */ 487 if( !PETSC_TRUE ) { 488 PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs; 489 if(llev==0) { 490 sprintf(fname,"Cmat_%d.m",llev++); 491 PetscViewerASCIIOpen(wcomm,fname,&viewer); 492 ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 493 ierr = MatView(Amat_fine, viewer ); CHKERRQ(ierr); 494 ierr = PetscViewerDestroy( &viewer ); 495 } 496 sprintf(fname,"Cmat_%d.m",llev++); 497 PetscViewerASCIIOpen(wcomm,fname,&viewer); 498 ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 499 ierr = MatView(Cmat, viewer ); CHKERRQ(ierr); 500 ierr = PetscViewerDestroy( &viewer ); 501 } 502 503 PetscFunctionReturn(0); 504 } 505 506 /* -------------------------------------------------------------------------- */ 507 /* 508 PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 509 by setting data structures and options. 510 511 Input Parameter: 512 . pc - the preconditioner context 513 514 Application Interface Routine: PCSetUp() 515 516 Notes: 517 The interface routine PCSetUp() is not usually called directly by 518 the user, but instead is called by PCApply() if necessary. 519 */ 520 #undef __FUNCT__ 521 #define __FUNCT__ "PCSetUp_GAMG" 522 PetscErrorCode PCSetUp_GAMG( PC pc ) 523 { 524 PetscErrorCode ierr; 525 PC_MG *mg = (PC_MG*)pc->data; 526 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 527 Mat Pmat = pc->pmat; 528 PetscInt fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS]; 529 MPI_Comm wcomm = ((PetscObject)pc)->comm; 530 PetscMPIInt mype,npe,nactivepe; 531 Mat Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS]; 532 PetscReal emaxs[GAMG_MAXLEVELS]; 533 IS *ASMLocalIDsArr[GAMG_MAXLEVELS],removedEqs[GAMG_MAXLEVELS]; 534 PetscInt level_bs[GAMG_MAXLEVELS]; 535 GAMGKKTMat kktMatsArr[GAMG_MAXLEVELS]; 536 PetscLogDouble nnz0=0.,nnztot=0.; 537 MatInfo info; 538 PetscBool stokes = PETSC_FALSE; 539 540 PetscFunctionBegin; 541 ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 542 ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 543 if (pc_gamg->verbose>2) PetscPrintf(wcomm,"[%d]%s pc_gamg->setup_count=%d pc->setupcalled=%d\n",mype,__FUNCT__,pc_gamg->setup_count,pc->setupcalled); 544 if( pc_gamg->setup_count++ > 0 ) { 545 PC_MG_Levels **mglevels = mg->levels; 546 /* just do Galerkin grids */ 547 Mat B,dA,dB; 548 assert(pc->setupcalled); 549 550 if( pc_gamg->Nlevels > 1 ) { 551 /* currently only handle case where mat and pmat are the same on coarser levels */ 552 ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,PETSC_NULL);CHKERRQ(ierr); 553 /* (re)set to get dirty flag */ 554 ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 555 556 for (level=pc_gamg->Nlevels-2; level>-1; level--) { 557 /* the first time through the matrix structure has changed from repartitioning */ 558 if( pc_gamg->setup_count==2 /*&& (pc_gamg->repart || level==0)*/) { 559 ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 560 ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr); 561 mglevels[level]->A = B; 562 } 563 else { 564 ierr = KSPGetOperators(mglevels[level]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 565 ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 566 } 567 ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN); CHKERRQ(ierr); 568 dB = B; 569 } 570 } 571 572 ierr = PCSetUp_MG( pc );CHKERRQ( ierr ); 573 574 /* PCSetUp_MG seems to insists on setting this to GMRES */ 575 ierr = KSPSetType( mglevels[0]->smoothd, KSPPREONLY ); CHKERRQ(ierr); 576 577 PetscFunctionReturn(0); 578 } 579 assert(pc->setupcalled == 0); 580 581 ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 582 583 if( pc_gamg->data == 0 ) { 584 if( !pc_gamg->createdefaultdata ){ 585 SETERRQ(wcomm,PETSC_ERR_LIB,"'createdefaultdata' not set(?) need to support NULL data"); 586 } 587 if( stokes ) { 588 SETERRQ(wcomm,PETSC_ERR_LIB,"Need data (eg, PCSetCoordinates) for Stokes problems"); 589 } 590 ierr = pc_gamg->createdefaultdata( pc, Pmat ); CHKERRQ(ierr); 591 } 592 593 /* get basic dims */ 594 if( stokes ) { 595 bs = pc_gamg->data_cell_rows; /* this is agg-mg specific */ 596 /* nloc = pc_gamg->data_sz/pc_gamg->data_cell_cols/bs; assert(pc_gamg->data_sz%(pc_gamg->data_cell_cols*bs)==0); */ 597 } 598 else { 599 /* PetscInt Istart,Iend,M,N; */ 600 ierr = MatGetBlockSize( Pmat, &bs ); CHKERRQ(ierr); 601 /* ierr = MatGetSize( Pmat, &M, &N );CHKERRQ(ierr); */ 602 /* ierr = MatGetOwnershipRange( Pmat, &Istart, &Iend ); CHKERRQ(ierr); */ 603 /* nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0); */ 604 } 605 606 ierr = MatGetSize( Pmat, &M, &qq );CHKERRQ(ierr); 607 if (pc_gamg->verbose) { 608 if(pc_gamg->verbose==1) ierr = MatGetInfo(Pmat,MAT_LOCAL,&info); 609 else ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info); 610 CHKERRQ(ierr); 611 nnz0 = info.nz_used; 612 nnztot = info.nz_used; 613 PetscPrintf(wcomm,"\t[%d]%s level %d N=%d, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n", 614 mype,__FUNCT__,0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols, 615 (int)(nnz0/(PetscReal)M),npe); 616 } 617 618 /* Get A_i and R_i */ 619 for ( level=0, Aarr[0]=Pmat, nactivepe = npe; /* hard wired stopping logic */ 620 level < (pc_gamg->Nlevels-1) && (level==0 || M>pc_gamg->coarse_eq_limit); /* && (npe==1 || nactivepe>1); */ 621 level++ ){ 622 level1 = level + 1; 623 #if defined PETSC_GAMG_USE_LOG 624 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0); CHKERRQ(ierr); 625 #if (defined GAMG_STAGES) 626 ierr = PetscLogStagePush(gamg_stages[level]); CHKERRQ( ierr ); 627 #endif 628 #endif 629 /* deal with Stokes, get sub matrices */ 630 ierr = GAMGKKTMatCreate( Aarr[level], stokes, &kktMatsArr[level] ); CHKERRQ(ierr); 631 632 { /* construct prolongator */ 633 Mat Gmat; 634 PetscCoarsenData *agg_lists; 635 Mat Prol11,Prol22; 636 637 level_bs[level] = bs; 638 ierr = pc_gamg->graph( pc,kktMatsArr[level].A11, &Gmat ); CHKERRQ(ierr); 639 ierr = pc_gamg->coarsen( pc, &Gmat, &agg_lists ); CHKERRQ(ierr); 640 ierr = pc_gamg->prolongator( pc, kktMatsArr[level].A11, Gmat, agg_lists, &Prol11 ); CHKERRQ(ierr); 641 642 /* could have failed to create new level */ 643 if( Prol11 ){ 644 /* get new block size of coarse matrices */ 645 if( pc_gamg->col_bs_id != -1 ){ 646 PetscBool flag; 647 ierr = PetscObjectComposedDataGetInt( (PetscObject)Prol11, pc_gamg->col_bs_id, bs, flag ); 648 CHKERRQ( ierr ); assert(flag); 649 } 650 651 if( stokes ) { 652 if(!pc_gamg->formkktprol) SETERRQ(wcomm,PETSC_ERR_USER,"Stokes not supportd by AMG method."); 653 /* R A12 == (T = A21 P)'; G = T' T; coarsen G; form plain agg with G */ 654 ierr = pc_gamg->formkktprol( pc, Prol11, kktMatsArr[level].A21, &Prol22 ); CHKERRQ(ierr); 655 } 656 657 if( pc_gamg->optprol ){ 658 /* smooth */ 659 ierr = pc_gamg->optprol( pc, kktMatsArr[level].A11, &Prol11 ); CHKERRQ(ierr); 660 } 661 662 if( stokes ) { 663 IS is_row[2] = {kktMatsArr[level].prim_is,kktMatsArr[level].constr_is}; 664 Mat a[4] = {Prol11, PETSC_NULL, PETSC_NULL, Prol22 }; 665 ierr = MatCreateNest(wcomm,2,is_row, 2, is_row, a, &Parr[level1] ); CHKERRQ(ierr); 666 } 667 else { 668 Parr[level1] = Prol11; 669 } 670 } 671 else Parr[level1] = PETSC_NULL; 672 673 if ( pc_gamg->use_aggs_in_gasm ) { 674 ierr = PetscCDGetASMBlocks(agg_lists, level_bs[level], &nASMBlocksArr[level], &ASMLocalIDsArr[level] ); CHKERRQ(ierr); 675 } 676 677 ierr = PetscCDGetRemovedIS( agg_lists, &removedEqs[level] ); CHKERRQ(ierr); 678 679 ierr = MatDestroy( &Gmat ); CHKERRQ(ierr); 680 ierr = PetscCDDestroy( agg_lists ); CHKERRQ(ierr); 681 } /* construct prolongator scope */ 682 #if defined PETSC_GAMG_USE_LOG 683 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 684 #endif 685 /* cache eigen estimate */ 686 if( pc_gamg->emax_id != -1 ){ 687 PetscBool flag; 688 ierr = PetscObjectComposedDataGetReal( (PetscObject)kktMatsArr[level].A11, pc_gamg->emax_id, emaxs[level], flag ); 689 CHKERRQ( ierr ); 690 if( !flag ) emaxs[level] = -1.; 691 } 692 else emaxs[level] = -1.; 693 if(level==0) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 694 if( !Parr[level1] ) { 695 if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s stop gridding, level %d\n",mype,__FUNCT__,level); 696 break; 697 } 698 #if defined PETSC_GAMG_USE_LOG 699 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 700 #endif 701 702 ierr = createLevel( pc, Aarr[level], bs, (PetscBool)(level==pc_gamg->Nlevels-2), 703 stokes, &Parr[level1], &Aarr[level1], &nactivepe ); 704 CHKERRQ(ierr); 705 706 #if defined PETSC_GAMG_USE_LOG 707 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 708 #endif 709 ierr = MatGetSize( Aarr[level1], &M, &qq );CHKERRQ(ierr); 710 711 if (pc_gamg->verbose > 0){ 712 PetscInt NN = M; 713 if(pc_gamg->verbose==1) { 714 ierr = MatGetInfo(Aarr[level1],MAT_LOCAL,&info); CHKERRQ(ierr); 715 ierr = MatGetLocalSize( Aarr[level1], &NN, &qq ); 716 } 717 else ierr = MatGetInfo( Aarr[level1], MAT_GLOBAL_SUM, &info ); 718 719 CHKERRQ(ierr); 720 nnztot += info.nz_used; 721 PetscPrintf(wcomm,"\t\t[%d]%s %d) N=%d, n data cols=%d, nnz/row (ave)=%d, %d active pes\n", 722 mype,__FUNCT__,(int)level1,M,pc_gamg->data_cell_cols, 723 (int)(info.nz_used/(PetscReal)NN), nactivepe ); 724 CHKERRQ(ierr); 725 } 726 727 /* stop if one node -- could pull back for singular problems */ 728 if( M/pc_gamg->data_cell_cols < 2 ) { 729 level++; 730 break; 731 } 732 #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 733 ierr = PetscLogStagePop(); CHKERRQ( ierr ); 734 #endif 735 } /* levels */ 736 737 if( pc_gamg->data ) { 738 ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr ); 739 pc_gamg->data = PETSC_NULL; 740 } 741 742 if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s %d levels, grid complexity = %g\n",0,__FUNCT__,level+1,nnztot/nnz0); 743 pc_gamg->Nlevels = level + 1; 744 fine_level = level; 745 ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,PETSC_NULL);CHKERRQ(ierr); 746 747 /* simple setup */ 748 if( !PETSC_TRUE ){ 749 PC_MG_Levels **mglevels = mg->levels; 750 for (lidx=0,level=pc_gamg->Nlevels-1; 751 lidx<fine_level; 752 lidx++, level--){ 753 ierr = PCMGSetInterpolation( pc, lidx+1, Parr[level] );CHKERRQ(ierr); 754 ierr = KSPSetOperators( mglevels[lidx]->smoothd, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN );CHKERRQ(ierr); 755 ierr = MatDestroy( &Parr[level] ); CHKERRQ(ierr); 756 ierr = MatDestroy( &Aarr[level] ); CHKERRQ(ierr); 757 } 758 ierr = KSPSetOperators( mglevels[fine_level]->smoothd, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 759 760 ierr = PCSetUp_MG( pc ); CHKERRQ( ierr ); 761 } 762 else if( pc_gamg->Nlevels > 1 ) { /* don't setup MG if one level */ 763 /* set default smoothers & set operators */ 764 for ( lidx = 1, level = pc_gamg->Nlevels-2; 765 lidx <= fine_level; 766 lidx++, level--) { 767 KSP smoother; 768 PC subpc; 769 770 ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr); 771 ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 772 773 ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 774 /* set ops */ 775 ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 776 ierr = PCMGSetInterpolation( pc, lidx, Parr[level+1] );CHKERRQ(ierr); 777 778 /* create field split PC, get subsmoother */ 779 if( stokes ) { 780 KSP *ksps; 781 PetscInt nn; 782 ierr = PCFieldSplitSetIS(subpc,"0",kktMatsArr[level].prim_is); CHKERRQ(ierr); 783 ierr = PCFieldSplitSetIS(subpc,"1",kktMatsArr[level].constr_is); CHKERRQ(ierr); 784 ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps); CHKERRQ(ierr); 785 smoother = ksps[0]; 786 ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 787 ierr = PetscFree( ksps ); CHKERRQ(ierr); 788 } 789 ierr = GAMGKKTMatDestroy( &kktMatsArr[level] ); CHKERRQ(ierr); 790 791 /* set defaults */ 792 ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr); 793 794 /* override defaults and command line args (!) */ 795 if ( pc_gamg->use_aggs_in_gasm ) { 796 PetscInt sz; 797 IS *is; 798 799 sz = nASMBlocksArr[level]; 800 is = ASMLocalIDsArr[level]; 801 ierr = PCSetType( subpc, PCGASM ); CHKERRQ(ierr); 802 if(sz==0){ 803 IS is; 804 PetscInt my0,kk; 805 ierr = MatGetOwnershipRange( Aarr[level], &my0, &kk ); CHKERRQ(ierr); 806 ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is ); CHKERRQ(ierr); 807 ierr = PCGASMSetLocalSubdomains( subpc, 1, &is, PETSC_NULL ); CHKERRQ(ierr); 808 ierr = ISDestroy( &is ); CHKERRQ(ierr); 809 } 810 else { 811 PetscInt kk; 812 ierr = PCGASMSetLocalSubdomains( subpc, sz, is, PETSC_NULL ); CHKERRQ(ierr); 813 for(kk=0;kk<sz;kk++){ 814 ierr = ISDestroy( &is[kk] ); CHKERRQ(ierr); 815 } 816 ierr = PetscFree( is ); CHKERRQ(ierr); 817 } 818 ierr = PCGASMSetOverlap( subpc, 0 ); CHKERRQ(ierr); 819 820 ASMLocalIDsArr[level] = PETSC_NULL; 821 nASMBlocksArr[level] = 0; 822 ierr = PCGASMSetType( subpc, PC_GASM_BASIC ); CHKERRQ(ierr); 823 } 824 else { 825 ierr = PCSetType( subpc, PCJACOBI ); CHKERRQ(ierr); 826 } 827 } 828 { 829 /* coarse grid */ 830 KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 831 Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0; 832 ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr); 833 ierr = KSPSetOperators( smoother, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 834 ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 835 ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 836 ierr = PCSetType( subpc, PCBJACOBI ); CHKERRQ(ierr); 837 ierr = PCSetUp( subpc ); CHKERRQ(ierr); 838 ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); assert(ii==1); 839 ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 840 ierr = PCSetType( pc2, PCLU ); CHKERRQ(ierr); 841 } 842 843 /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 844 ierr = PetscObjectOptionsBegin( (PetscObject)pc );CHKERRQ(ierr); 845 ierr = PCSetFromOptions_MG( pc ); CHKERRQ(ierr); 846 ierr = PetscOptionsEnd(); CHKERRQ(ierr); 847 if (mg->galerkin != 2) SETERRQ(wcomm,PETSC_ERR_USER,"GAMG does Galerkin manually so the -pc_mg_galerkin option must not be used."); 848 849 /* create cheby smoothers */ 850 for ( lidx = 1, level = pc_gamg->Nlevels-2; 851 lidx <= fine_level; 852 lidx++, level--) { 853 KSP smoother; 854 PetscBool flag; 855 PC subpc; 856 857 ierr = PCMGGetSmoother( pc, lidx, &smoother ); CHKERRQ(ierr); 858 ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 859 860 /* create field split PC, get subsmoother */ 861 if( stokes ) { 862 KSP *ksps; 863 PetscInt nn; 864 ierr = PCFieldSplitGetSubKSP(subpc,&nn,&ksps); CHKERRQ(ierr); 865 smoother = ksps[0]; 866 ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 867 ierr = PetscFree( ksps ); CHKERRQ(ierr); 868 } 869 870 /* do my own cheby */ 871 ierr = PetscObjectTypeCompare( (PetscObject)smoother, KSPCHEBYCHEV, &flag ); CHKERRQ(ierr); 872 if( flag ) { 873 PetscReal emax, emin; 874 ierr = PetscObjectTypeCompare( (PetscObject)subpc, PCJACOBI, &flag ); CHKERRQ(ierr); 875 if( flag && emaxs[level] > 0.0 ) emax=emaxs[level]; /* eigen estimate only for diagnal PC */ 876 else{ /* eigen estimate 'emax' */ 877 KSP eksp; Mat Lmat = Aarr[level]; 878 Vec bb, xx; 879 880 ierr = MatGetVecs( Lmat, &bb, 0 ); CHKERRQ(ierr); 881 ierr = MatGetVecs( Lmat, &xx, 0 ); CHKERRQ(ierr); 882 { 883 PetscRandom rctx; 884 ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr); 885 ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 886 ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 887 ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr); 888 } 889 890 if( removedEqs[level] ) { 891 /* being very careful - zeroing out BC rows (this is not done in agg.c estimates) */ 892 PetscScalar *zeros; 893 PetscInt ii,jj, *idx_bs, sz, bs=level_bs[level]; 894 const PetscInt *idx; 895 ierr = ISGetLocalSize( removedEqs[level], &sz ); CHKERRQ(ierr); 896 ierr = PetscMalloc( bs*sz*sizeof(PetscScalar), &zeros ); CHKERRQ(ierr); 897 for(ii=0;ii<bs*sz;ii++) zeros[ii] = 0.; 898 ierr = PetscMalloc( bs*sz*sizeof(PetscInt), &idx_bs ); CHKERRQ(ierr); 899 ierr = ISGetIndices( removedEqs[level], &idx); CHKERRQ(ierr); 900 for(ii=0;ii<sz;ii++) { 901 for(jj=0;jj<bs;jj++) { 902 idx_bs[ii] = bs*idx[ii]+jj; 903 } 904 } 905 ierr = ISRestoreIndices( removedEqs[level], &idx ); CHKERRQ(ierr); 906 if( sz > 0 ) { 907 ierr = VecSetValues( bb, sz, idx_bs, zeros, INSERT_VALUES ); CHKERRQ(ierr); 908 } 909 ierr = PetscFree( idx_bs ); CHKERRQ(ierr); 910 ierr = PetscFree( zeros ); CHKERRQ(ierr); 911 ierr = VecAssemblyBegin(bb); CHKERRQ(ierr); 912 ierr = VecAssemblyEnd(bb); CHKERRQ(ierr); 913 } 914 ierr = KSPCreate( wcomm, &eksp );CHKERRQ(ierr); 915 ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10 ); 916 CHKERRQ(ierr); 917 ierr = KSPSetNormType( eksp, KSP_NORM_NONE ); CHKERRQ(ierr); 918 ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr); 919 ierr = KSPAppendOptionsPrefix( eksp, "gamg_est_"); CHKERRQ(ierr); 920 ierr = KSPSetFromOptions( eksp ); CHKERRQ(ierr); 921 922 ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr); 923 ierr = KSPSetOperators( eksp, Lmat, Lmat, SAME_NONZERO_PATTERN ); CHKERRQ( ierr ); 924 ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr); 925 926 /* set PC type to be same as smoother */ 927 ierr = KSPSetPC( eksp, subpc ); CHKERRQ( ierr ); 928 929 /* solve - keep stuff out of logging */ 930 ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 931 ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 932 ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr); 933 ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 934 ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 935 936 ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr); 937 938 ierr = VecDestroy( &xx ); CHKERRQ(ierr); 939 ierr = VecDestroy( &bb ); CHKERRQ(ierr); 940 ierr = KSPDestroy( &eksp ); CHKERRQ(ierr); 941 942 if( pc_gamg->verbose > 0 ) { 943 PetscInt N1, tt; 944 ierr = MatGetSize( Aarr[level], &N1, &tt ); CHKERRQ(ierr); 945 PetscPrintf(wcomm,"\t\t\t%s PC setup max eigen=%e min=%e on level %d (N=%d)\n",__FUNCT__,emax,emin,lidx,N1); 946 } 947 } 948 { 949 PetscInt N1, N0, tt; 950 ierr = MatGetSize( Aarr[level], &N1, &tt ); CHKERRQ(ierr); 951 ierr = MatGetSize( Aarr[level+1], &N0, &tt ); CHKERRQ(ierr); 952 /* heuristic - is this crap? */ 953 emin = 1.*emax/((PetscReal)N1/(PetscReal)N0); 954 emax *= 1.05; 955 } 956 ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr); 957 } /* setup checby flag */ 958 959 if( removedEqs[level] ) { 960 ierr = ISDestroy( &removedEqs[level] ); CHKERRQ(ierr); 961 } 962 } /* non-coarse levels */ 963 964 /* clean up */ 965 for(level=1;level<pc_gamg->Nlevels;level++){ 966 ierr = MatDestroy( &Parr[level] ); CHKERRQ(ierr); 967 ierr = MatDestroy( &Aarr[level] ); CHKERRQ(ierr); 968 } 969 970 ierr = PCSetUp_MG( pc );CHKERRQ( ierr ); 971 972 if( PETSC_FALSE ){ 973 KSP smoother; /* PCSetUp_MG seems to insists on setting this to GMRES on coarse grid */ 974 ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr); 975 ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr); 976 } 977 } 978 else { 979 KSP smoother; 980 if (pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s one level solver used (system is seen as DD). Using default solver.\n",mype,__FUNCT__); 981 ierr = PCMGGetSmoother( pc, 0, &smoother ); CHKERRQ(ierr); 982 ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], SAME_NONZERO_PATTERN ); CHKERRQ(ierr); 983 ierr = KSPSetType( smoother, KSPPREONLY ); CHKERRQ(ierr); 984 ierr = PCSetUp_MG( pc );CHKERRQ( ierr ); 985 } 986 987 PetscFunctionReturn(0); 988 } 989 990 /* ------------------------------------------------------------------------- */ 991 /* 992 PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 993 that was created with PCCreate_GAMG(). 994 995 Input Parameter: 996 . pc - the preconditioner context 997 998 Application Interface Routine: PCDestroy() 999 */ 1000 #undef __FUNCT__ 1001 #define __FUNCT__ "PCDestroy_GAMG" 1002 PetscErrorCode PCDestroy_GAMG( PC pc ) 1003 { 1004 PetscErrorCode ierr; 1005 PC_MG *mg = (PC_MG*)pc->data; 1006 PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 1007 1008 PetscFunctionBegin; 1009 ierr = PCReset_GAMG( pc );CHKERRQ(ierr); 1010 ierr = PetscFree( pc_gamg );CHKERRQ(ierr); 1011 ierr = PCDestroy_MG( pc );CHKERRQ(ierr); 1012 PetscFunctionReturn(0); 1013 } 1014 1015 1016 #undef __FUNCT__ 1017 #define __FUNCT__ "PCGAMGSetProcEqLim" 1018 /*@ 1019 PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via 1020 processor reduction. 1021 1022 Not Collective on PC 1023 1024 Input Parameters: 1025 . pc - the preconditioner context 1026 1027 1028 Options Database Key: 1029 . -pc_gamg_process_eq_limit 1030 1031 Level: intermediate 1032 1033 Concepts: Unstructured multrigrid preconditioner 1034 1035 .seealso: () 1036 @*/ 1037 PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 1038 { 1039 PetscErrorCode ierr; 1040 1041 PetscFunctionBegin; 1042 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1043 ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 1044 PetscFunctionReturn(0); 1045 } 1046 1047 EXTERN_C_BEGIN 1048 #undef __FUNCT__ 1049 #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG" 1050 PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 1051 { 1052 PC_MG *mg = (PC_MG*)pc->data; 1053 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1054 1055 PetscFunctionBegin; 1056 if(n>0) pc_gamg->min_eq_proc = n; 1057 PetscFunctionReturn(0); 1058 } 1059 EXTERN_C_END 1060 1061 #undef __FUNCT__ 1062 #define __FUNCT__ "PCGAMGSetCoarseEqLim" 1063 /*@ 1064 PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids. 1065 1066 Collective on PC 1067 1068 Input Parameters: 1069 . pc - the preconditioner context 1070 1071 1072 Options Database Key: 1073 . -pc_gamg_coarse_eq_limit 1074 1075 Level: intermediate 1076 1077 Concepts: Unstructured multrigrid preconditioner 1078 1079 .seealso: () 1080 @*/ 1081 PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 1082 { 1083 PetscErrorCode ierr; 1084 1085 PetscFunctionBegin; 1086 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1087 ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 1088 PetscFunctionReturn(0); 1089 } 1090 1091 EXTERN_C_BEGIN 1092 #undef __FUNCT__ 1093 #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG" 1094 PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 1095 { 1096 PC_MG *mg = (PC_MG*)pc->data; 1097 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1098 1099 PetscFunctionBegin; 1100 if(n>0) pc_gamg->coarse_eq_limit = n; 1101 PetscFunctionReturn(0); 1102 } 1103 EXTERN_C_END 1104 1105 #undef __FUNCT__ 1106 #define __FUNCT__ "PCGAMGSetRepartitioning" 1107 /*@ 1108 PCGAMGSetRepartitioning - Repartition the coarse grids 1109 1110 Collective on PC 1111 1112 Input Parameters: 1113 . pc - the preconditioner context 1114 1115 1116 Options Database Key: 1117 . -pc_gamg_repartition 1118 1119 Level: intermediate 1120 1121 Concepts: Unstructured multrigrid preconditioner 1122 1123 .seealso: () 1124 @*/ 1125 PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n) 1126 { 1127 PetscErrorCode ierr; 1128 1129 PetscFunctionBegin; 1130 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1131 ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1132 PetscFunctionReturn(0); 1133 } 1134 1135 EXTERN_C_BEGIN 1136 #undef __FUNCT__ 1137 #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG" 1138 PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n) 1139 { 1140 PC_MG *mg = (PC_MG*)pc->data; 1141 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1142 1143 PetscFunctionBegin; 1144 pc_gamg->repart = n; 1145 PetscFunctionReturn(0); 1146 } 1147 EXTERN_C_END 1148 1149 #undef __FUNCT__ 1150 #define __FUNCT__ "PCGAMGSetUseASMAggs" 1151 /*@ 1152 PCGAMGSetUseASMAggs - 1153 1154 Collective on PC 1155 1156 Input Parameters: 1157 . pc - the preconditioner context 1158 1159 1160 Options Database Key: 1161 . -pc_gamg_use_agg_gasm 1162 1163 Level: intermediate 1164 1165 Concepts: Unstructured multrigrid preconditioner 1166 1167 .seealso: () 1168 @*/ 1169 PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n) 1170 { 1171 PetscErrorCode ierr; 1172 1173 PetscFunctionBegin; 1174 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1175 ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1176 PetscFunctionReturn(0); 1177 } 1178 1179 EXTERN_C_BEGIN 1180 #undef __FUNCT__ 1181 #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG" 1182 PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n) 1183 { 1184 PC_MG *mg = (PC_MG*)pc->data; 1185 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1186 1187 PetscFunctionBegin; 1188 pc_gamg->use_aggs_in_gasm = n; 1189 PetscFunctionReturn(0); 1190 } 1191 EXTERN_C_END 1192 1193 #undef __FUNCT__ 1194 #define __FUNCT__ "PCGAMGSetNlevels" 1195 /*@ 1196 PCGAMGSetNlevels - 1197 1198 Not collective on PC 1199 1200 Input Parameters: 1201 . pc - the preconditioner context 1202 1203 1204 Options Database Key: 1205 . -pc_mg_levels 1206 1207 Level: intermediate 1208 1209 Concepts: Unstructured multrigrid preconditioner 1210 1211 .seealso: () 1212 @*/ 1213 PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 1214 { 1215 PetscErrorCode ierr; 1216 1217 PetscFunctionBegin; 1218 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1219 ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 1220 PetscFunctionReturn(0); 1221 } 1222 1223 EXTERN_C_BEGIN 1224 #undef __FUNCT__ 1225 #define __FUNCT__ "PCGAMGSetNlevels_GAMG" 1226 PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 1227 { 1228 PC_MG *mg = (PC_MG*)pc->data; 1229 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1230 1231 PetscFunctionBegin; 1232 pc_gamg->Nlevels = n; 1233 PetscFunctionReturn(0); 1234 } 1235 EXTERN_C_END 1236 1237 #undef __FUNCT__ 1238 #define __FUNCT__ "PCGAMGSetThreshold" 1239 /*@ 1240 PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 1241 1242 Not collective on PC 1243 1244 Input Parameters: 1245 . pc - the preconditioner context 1246 1247 1248 Options Database Key: 1249 . -pc_gamg_threshold 1250 1251 Level: intermediate 1252 1253 Concepts: Unstructured multrigrid preconditioner 1254 1255 .seealso: () 1256 @*/ 1257 PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n) 1258 { 1259 PetscErrorCode ierr; 1260 1261 PetscFunctionBegin; 1262 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1263 ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr); 1264 PetscFunctionReturn(0); 1265 } 1266 1267 EXTERN_C_BEGIN 1268 #undef __FUNCT__ 1269 #define __FUNCT__ "PCGAMGSetThreshold_GAMG" 1270 PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n) 1271 { 1272 PC_MG *mg = (PC_MG*)pc->data; 1273 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1274 1275 PetscFunctionBegin; 1276 pc_gamg->threshold = n; 1277 PetscFunctionReturn(0); 1278 } 1279 EXTERN_C_END 1280 1281 #undef __FUNCT__ 1282 #define __FUNCT__ "PCGAMGSetType" 1283 /*@ 1284 PCGAMGSetType - Set solution method - calls sub create method 1285 1286 Collective on PC 1287 1288 Input Parameters: 1289 . pc - the preconditioner context 1290 1291 1292 Options Database Key: 1293 . -pc_gamg_type 1294 1295 Level: intermediate 1296 1297 Concepts: Unstructured multrigrid preconditioner 1298 1299 .seealso: () 1300 @*/ 1301 PetscErrorCode PCGAMGSetType( PC pc, const PCGAMGType type ) 1302 { 1303 PetscErrorCode ierr; 1304 1305 PetscFunctionBegin; 1306 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1307 ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,const PCGAMGType),(pc,type)); 1308 CHKERRQ(ierr); 1309 PetscFunctionReturn(0); 1310 } 1311 1312 EXTERN_C_BEGIN 1313 #undef __FUNCT__ 1314 #define __FUNCT__ "PCGAMGSetType_GAMG" 1315 PetscErrorCode PCGAMGSetType_GAMG( PC pc, const PCGAMGType type ) 1316 { 1317 PetscErrorCode ierr,(*r)(PC); 1318 1319 PetscFunctionBegin; 1320 ierr = PetscFListFind(GAMGList,((PetscObject)pc)->comm,type,PETSC_FALSE,(PetscVoidStarFunction)&r); 1321 CHKERRQ(ierr); 1322 1323 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 1324 1325 /* call sub create method */ 1326 ierr = (*r)(pc); CHKERRQ(ierr); 1327 1328 PetscFunctionReturn(0); 1329 } 1330 EXTERN_C_END 1331 1332 #undef __FUNCT__ 1333 #define __FUNCT__ "PCSetFromOptions_GAMG" 1334 PetscErrorCode PCSetFromOptions_GAMG( PC pc ) 1335 { 1336 PetscErrorCode ierr; 1337 PC_MG *mg = (PC_MG*)pc->data; 1338 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1339 PetscBool flag; 1340 MPI_Comm wcomm = ((PetscObject)pc)->comm; 1341 1342 PetscFunctionBegin; 1343 ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr); 1344 { 1345 /* -pc_gamg_verbose */ 1346 ierr = PetscOptionsInt("-pc_gamg_verbose","Verbose (debugging) output for PCGAMG", 1347 "none", pc_gamg->verbose, 1348 &pc_gamg->verbose, PETSC_NULL ); 1349 CHKERRQ(ierr); 1350 1351 /* -pc_gamg_repartition */ 1352 ierr = PetscOptionsBool("-pc_gamg_repartition", 1353 "Repartion coarse grids (false)", 1354 "PCGAMGRepartitioning", 1355 pc_gamg->repart, 1356 &pc_gamg->repart, 1357 &flag); 1358 CHKERRQ(ierr); 1359 1360 /* -pc_gamg_use_agg_gasm */ 1361 ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm", 1362 "Use aggregation agragates for GASM smoother (false)", 1363 "PCGAMGUseASMAggs", 1364 pc_gamg->use_aggs_in_gasm, 1365 &pc_gamg->use_aggs_in_gasm, 1366 &flag); 1367 CHKERRQ(ierr); 1368 1369 /* -pc_gamg_process_eq_limit */ 1370 ierr = PetscOptionsInt("-pc_gamg_process_eq_limit", 1371 "Limit (goal) on number of equations per process on coarse grids", 1372 "PCGAMGSetProcEqLim", 1373 pc_gamg->min_eq_proc, 1374 &pc_gamg->min_eq_proc, 1375 &flag ); 1376 CHKERRQ(ierr); 1377 1378 /* -pc_gamg_coarse_eq_limit */ 1379 ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit", 1380 "Limit on number of equations for the coarse grid", 1381 "PCGAMGSetCoarseEqLim", 1382 pc_gamg->coarse_eq_limit, 1383 &pc_gamg->coarse_eq_limit, 1384 &flag ); 1385 CHKERRQ(ierr); 1386 1387 /* -pc_gamg_threshold */ 1388 ierr = PetscOptionsReal("-pc_gamg_threshold", 1389 "Relative threshold to use for dropping edges in aggregation graph", 1390 "PCGAMGSetThreshold", 1391 pc_gamg->threshold, 1392 &pc_gamg->threshold, 1393 &flag ); 1394 CHKERRQ(ierr); 1395 if(flag && pc_gamg->verbose) PetscPrintf(wcomm,"\t[%d]%s threshold set %e\n",0,__FUNCT__,pc_gamg->threshold); 1396 1397 ierr = PetscOptionsInt("-pc_mg_levels", 1398 "Set number of MG levels", 1399 "PCGAMGSetNlevels", 1400 pc_gamg->Nlevels, 1401 &pc_gamg->Nlevels, 1402 &flag ); 1403 } 1404 ierr = PetscOptionsTail();CHKERRQ(ierr); 1405 1406 PetscFunctionReturn(0); 1407 } 1408 1409 /* -------------------------------------------------------------------------- */ 1410 /*MC 1411 PCGAMG - Geometric algebraic multigrid (AMG) preconditioning framework. 1412 - This is the entry point to GAMG, registered in pcregis.c 1413 1414 Options Database Keys: 1415 Multigrid options(inherited) 1416 + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (PCMGSetCycleType) 1417 . -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp) 1418 . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown) 1419 - -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 1420 1421 Level: intermediate 1422 1423 Concepts: multigrid 1424 1425 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 1426 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 1427 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 1428 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 1429 PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 1430 M*/ 1431 EXTERN_C_BEGIN 1432 #undef __FUNCT__ 1433 #define __FUNCT__ "PCCreate_GAMG" 1434 PetscErrorCode PCCreate_GAMG( PC pc ) 1435 { 1436 PetscErrorCode ierr; 1437 PC_GAMG *pc_gamg; 1438 PC_MG *mg; 1439 #if defined PETSC_GAMG_USE_LOG 1440 static long count = 0; 1441 #endif 1442 1443 PetscFunctionBegin; 1444 1445 /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 1446 ierr = PCSetType( pc, PCMG ); CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 1447 ierr = PetscObjectChangeTypeName( (PetscObject)pc, PCGAMG ); CHKERRQ(ierr); 1448 1449 /* create a supporting struct and attach it to pc */ 1450 ierr = PetscNewLog( pc, PC_GAMG, &pc_gamg ); CHKERRQ(ierr); 1451 mg = (PC_MG*)pc->data; 1452 mg->galerkin = 2; /* Use Galerkin, but it is computed externally */ 1453 mg->innerctx = pc_gamg; 1454 1455 pc_gamg->setup_count = 0; 1456 /* these should be in subctx but repartitioning needs simple arrays */ 1457 pc_gamg->data_sz = 0; 1458 pc_gamg->data = 0; 1459 1460 /* register AMG type */ 1461 if( !GAMGList ){ 1462 ierr = PetscFListAdd(&GAMGList,GAMGGEO,"PCCreateGAMG_GEO",(void(*)(void))PCCreateGAMG_GEO);CHKERRQ(ierr); 1463 ierr = PetscFListAdd(&GAMGList,GAMGAGG,"PCCreateGAMG_AGG",(void(*)(void))PCCreateGAMG_AGG);CHKERRQ(ierr); 1464 } 1465 1466 /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 1467 pc->ops->setfromoptions = PCSetFromOptions_GAMG; 1468 pc->ops->setup = PCSetUp_GAMG; 1469 pc->ops->reset = PCReset_GAMG; 1470 pc->ops->destroy = PCDestroy_GAMG; 1471 1472 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1473 "PCGAMGSetProcEqLim_C", 1474 "PCGAMGSetProcEqLim_GAMG", 1475 PCGAMGSetProcEqLim_GAMG); 1476 CHKERRQ(ierr); 1477 1478 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1479 "PCGAMGSetCoarseEqLim_C", 1480 "PCGAMGSetCoarseEqLim_GAMG", 1481 PCGAMGSetCoarseEqLim_GAMG); 1482 CHKERRQ(ierr); 1483 1484 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1485 "PCGAMGSetRepartitioning_C", 1486 "PCGAMGSetRepartitioning_GAMG", 1487 PCGAMGSetRepartitioning_GAMG); 1488 CHKERRQ(ierr); 1489 1490 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1491 "PCGAMGSetUseASMAggs_C", 1492 "PCGAMGSetUseASMAggs_GAMG", 1493 PCGAMGSetUseASMAggs_GAMG); 1494 CHKERRQ(ierr); 1495 1496 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1497 "PCGAMGSetThreshold_C", 1498 "PCGAMGSetThreshold_GAMG", 1499 PCGAMGSetThreshold_GAMG); 1500 CHKERRQ(ierr); 1501 1502 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1503 "PCGAMGSetType_C", 1504 "PCGAMGSetType_GAMG", 1505 PCGAMGSetType_GAMG); 1506 CHKERRQ(ierr); 1507 1508 pc_gamg->repart = PETSC_FALSE; 1509 pc_gamg->use_aggs_in_gasm = PETSC_FALSE; 1510 pc_gamg->min_eq_proc = 100; 1511 pc_gamg->coarse_eq_limit = 800; 1512 pc_gamg->threshold = 0.001; 1513 pc_gamg->Nlevels = GAMG_MAXLEVELS; 1514 pc_gamg->verbose = 0; 1515 pc_gamg->emax_id = -1; 1516 pc_gamg->col_bs_id = -1; 1517 1518 /* private events */ 1519 #if defined PETSC_GAMG_USE_LOG 1520 if( count++ == 0 ) { 1521 PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]); 1522 PetscLogEventRegister(" Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]); 1523 /* PetscLogEventRegister(" G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */ 1524 /* PetscLogEventRegister(" G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */ 1525 /* PetscLogEventRegister(" G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */ 1526 PetscLogEventRegister(" MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]); 1527 PetscLogEventRegister(" geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]); 1528 PetscLogEventRegister(" geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]); 1529 PetscLogEventRegister(" search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]); 1530 PetscLogEventRegister(" SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]); 1531 PetscLogEventRegister(" SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]); 1532 PetscLogEventRegister(" SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]); 1533 PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]); 1534 PetscLogEventRegister(" repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]); 1535 PetscLogEventRegister(" Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]); 1536 PetscLogEventRegister(" Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]); 1537 PetscLogEventRegister(" Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]); 1538 1539 /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */ 1540 /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */ 1541 /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */ 1542 /* create timer stages */ 1543 #if defined GAMG_STAGES 1544 { 1545 char str[32]; 1546 sprintf(str,"MG Level %d (finest)",0); 1547 PetscLogStageRegister(str, &gamg_stages[0]); 1548 PetscInt lidx; 1549 for (lidx=1;lidx<9;lidx++){ 1550 sprintf(str,"MG Level %d",lidx); 1551 PetscLogStageRegister(str, &gamg_stages[lidx]); 1552 } 1553 } 1554 #endif 1555 } 1556 #endif 1557 /* general events */ 1558 #if defined PETSC_USE_LOG 1559 PetscLogEventRegister("PCGAMGgraph_AGG", 0, &PC_GAMGGgraph_AGG); 1560 PetscLogEventRegister("PCGAMGgraph_GEO", PC_CLASSID, &PC_GAMGGgraph_GEO); 1561 PetscLogEventRegister("PCGAMGcoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG); 1562 PetscLogEventRegister("PCGAMGcoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO); 1563 PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG); 1564 PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO); 1565 PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptprol_AGG); 1566 PetscLogEventRegister("GAMGKKTProl_AGG", PC_CLASSID, &PC_GAMGKKTProl_AGG); 1567 #endif 1568 1569 /* instantiate derived type */ 1570 ierr = PetscOptionsHead("GAMG options"); CHKERRQ(ierr); 1571 { 1572 char tname[256] = GAMGAGG; 1573 ierr = PetscOptionsList("-pc_gamg_type","Type of GAMG method","PCGAMGSetType", 1574 GAMGList, tname, tname, sizeof(tname), PETSC_NULL ); 1575 CHKERRQ(ierr); 1576 ierr = PCGAMGSetType( pc, tname ); CHKERRQ(ierr); 1577 } 1578 ierr = PetscOptionsTail(); CHKERRQ(ierr); 1579 1580 PetscFunctionReturn(0); 1581 } 1582 EXTERN_C_END 1583