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