1 /* 2 GAMG geometric-algebric multiogrid PC - Mark Adams 2011 3 */ 4 #include <../src/ksp/pc/impls/gamg/gamg.h> 5 6 PetscLogEvent gamg_setup_stages[NUM_SET]; 7 8 /* Private context for the GAMG preconditioner */ 9 typedef struct gamg_TAG { 10 PetscInt m_dim; 11 PetscInt m_Nlevels; 12 PetscInt m_data_sz; 13 PetscReal *m_data; /* blocked vector of vertex data on fine grid (coordinates) */ 14 } PC_GAMG; 15 16 #define TOP_GRID_LIM 100 17 18 /* -----------------------------------------------------------------------------*/ 19 #undef __FUNCT__ 20 #define __FUNCT__ "PCReset_GAMG" 21 PetscErrorCode PCReset_GAMG(PC pc) 22 { 23 PetscErrorCode ierr; 24 PC_MG *mg = (PC_MG*)pc->data; 25 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 26 27 PetscFunctionBegin; 28 ierr = PetscFree(pc_gamg->m_data);CHKERRQ(ierr); 29 pc_gamg->m_data = 0; pc_gamg->m_data_sz = 0; 30 PetscFunctionReturn(0); 31 } 32 33 /* -------------------------------------------------------------------------- */ 34 /* 35 PCSetCoordinates_GAMG 36 37 Input Parameter: 38 . pc - the preconditioner context 39 */ 40 EXTERN_C_BEGIN 41 #undef __FUNCT__ 42 #define __FUNCT__ "PCSetCoordinates_GAMG" 43 PetscErrorCode PCSetCoordinates_GAMG( PC a_pc, PetscInt a_ndm, PetscReal *a_coords ) 44 { 45 PC_MG *mg = (PC_MG*)a_pc->data; 46 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 47 PetscErrorCode ierr; 48 PetscInt arrsz, bs, my0, tt, ii, nloc, sz, kk; 49 Mat mat = a_pc->pmat; 50 PetscBool useSA = PETSC_FALSE, flag; 51 char str[16]; 52 53 PetscFunctionBegin; 54 ierr = PetscOptionsGetString(PETSC_NULL,"-pc_gamg_type",str,16,&flag); CHKERRQ( ierr ); 55 useSA = (PetscBool)(flag && strcmp(str,"sa") == 0); 56 ierr = MatGetBlockSize( mat, &bs ); CHKERRQ( ierr ); 57 ierr = MatGetOwnershipRange( mat, &my0, &tt ); CHKERRQ(ierr); 58 nloc = (tt-my0)/bs; sz = (useSA ? (a_coords==0 ? a_ndm*a_ndm: (a_ndm==2 ? 3*a_ndm : 6*a_ndm)) : a_ndm ); 59 arrsz = nloc*sz; 60 61 // put coordinates 62 if (!pc_gamg->m_data || (pc_gamg->m_data_sz != arrsz)) { 63 ierr = PetscFree( pc_gamg->m_data ); CHKERRQ(ierr); 64 ierr = PetscMalloc((arrsz+1)*sizeof(double), &pc_gamg->m_data ); CHKERRQ(ierr); 65 } 66 for(tt=0;tt<arrsz;tt++)pc_gamg->m_data[tt] = 0.; 67 pc_gamg->m_data[arrsz] = -99.; 68 /* copy data in */ 69 if( useSA ) { 70 for(tt=0,kk=0;tt<nloc;tt++){ 71 PetscReal *data = &pc_gamg->m_data[tt*sz]; 72 for(ii=0;ii<a_ndm;ii++) data[a_ndm*ii + ii] = 1.0; /* translational mode */ 73 if( a_coords != 0 ) { 74 if( a_ndm == 2 ){ 75 data[4] = a_coords[2*tt+1]; data[5] = -a_coords[2*tt]; 76 } 77 else { 78 data[10] = -a_coords[2*tt+2]; data[12] = a_coords[2*tt+2]; data[15] = a_coords[2*tt+1]; 79 data[11] = -a_coords[2*tt+1]; data[14] = -a_coords[2*tt]; data[16] = a_coords[2*tt]; 80 } 81 } 82 } 83 } 84 else { 85 for(tt=0;tt<arrsz;tt++){ 86 pc_gamg->m_data[tt] = a_coords[tt]; 87 } 88 } 89 assert(pc_gamg->m_data[arrsz] == -99.); 90 pc_gamg->m_data_sz = arrsz; 91 pc_gamg->m_dim = a_ndm; 92 93 PetscFunctionReturn(0); 94 } 95 EXTERN_C_END 96 97 /* -------------------------------------------------------------------------- */ 98 /* 99 partitionLevel 100 101 Input Parameter: 102 . a_Amat_fine - matrix on this fine (k) level 103 . a_data_sz - size of data to move 104 In/Output Parameter: 105 . a_P_inout - prolongation operator to the next level (k-1) 106 . a_coarse_data - data that need to be moved 107 . a_active_proc - number of active procs 108 Output Parameter: 109 . a_Amat_crs - coarse matrix that is created (k-1) 110 */ 111 112 #undef __FUNCT__ 113 #define __FUNCT__ "partitionLevel" 114 PetscErrorCode partitionLevel( Mat a_Amat_fine, 115 PetscInt a_data_sz, 116 Mat *a_P_inout, 117 PetscReal **a_coarse_data, 118 PetscMPIInt *a_active_proc, 119 Mat *a_Amat_crs 120 ) 121 { 122 PetscErrorCode ierr; 123 Mat Amat, Pnew, Pold = *a_P_inout; 124 IS new_indices,isnum; 125 MPI_Comm wcomm = ((PetscObject)a_Amat_fine)->comm; 126 PetscMPIInt nactive_procs,mype,npe; 127 PetscInt Istart,Iend,Istart0,Iend0,ncrs0,ncrs_new,bs,bs2; 128 129 PetscFunctionBegin; 130 ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 131 ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 132 ierr = MatGetBlockSize( a_Amat_fine, &bs ); CHKERRQ(ierr); 133 /* RAP */ 134 ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Amat ); CHKERRQ(ierr); 135 ierr = MatSetBlockSize( Amat, bs ); CHKERRQ(ierr); 136 ierr = MatGetOwnershipRange( Amat, &Istart0, &Iend0 ); CHKERRQ(ierr); 137 ncrs0 = (Iend0-Istart0)/bs; assert( (Iend0-Istart0)%bs == 0 ); 138 139 /* Repartition Amat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */ 140 { 141 PetscInt neq,N,counts[npe]; 142 IS isnewproc; 143 PetscMPIInt new_npe,targ_npe; 144 145 ierr = MatGetSize( Amat, &neq, &N );CHKERRQ(ierr); 146 #define MIN_EQ_PROC 100 147 nactive_procs = *a_active_proc; 148 targ_npe = neq/MIN_EQ_PROC; /* hardwire min. number of eq/proc */ 149 if( targ_npe == 0 || neq < TOP_GRID_LIM ) new_npe = 1; /* chop coarsest grid */ 150 else if (targ_npe >= nactive_procs ) new_npe = nactive_procs; /* no change */ 151 else { 152 PetscMPIInt factstart,fact; 153 new_npe = -9999; 154 factstart = nactive_procs; 155 for(fact=factstart;fact>0;fact--){ /* try to find a better number of procs */ 156 if( nactive_procs%fact==0 && neq/(nactive_procs/fact) > MIN_EQ_PROC ) { 157 new_npe = nactive_procs/fact; 158 } 159 } 160 assert(new_npe != -9999); 161 } 162 *a_active_proc = new_npe; /* output for next time */ 163 { /* partition: get 'isnewproc' */ 164 MatPartitioning mpart; 165 Mat adj; 166 const PetscInt *is_idx; 167 PetscInt is_sz,kk,jj,ii,old_fact=(npe/nactive_procs),new_fact=(npe/new_npe),*isnewproc_idx; 168 /* create sub communicator */ 169 MPI_Comm cm,new_comm; 170 int membershipKey = mype % old_fact; 171 ierr = MPI_Comm_split(wcomm, membershipKey, mype, &cm); CHKERRQ(ierr); 172 ierr = PetscCommDuplicate( cm, &new_comm, PETSC_NULL ); CHKERRQ(ierr); 173 ierr = MPI_Comm_free( &cm ); CHKERRQ(ierr); 174 175 /* MatPartitioningApply call MatConvert, which is collective */ 176 if( bs==1) { 177 ierr = MatConvert( Amat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 178 } 179 else { 180 Mat tMat; 181 PetscInt Ii,ncols; const PetscScalar *vals; const PetscInt *idx; 182 ierr = MatCreateMPIAIJ( wcomm, ncrs0, ncrs0, 183 PETSC_DETERMINE, PETSC_DETERMINE, 184 25, PETSC_NULL, 10, PETSC_NULL, 185 &tMat ); 186 187 for ( Ii = Istart0; Ii < Iend0; Ii++ ) { 188 PetscInt dest_row = Ii/bs; 189 ierr = MatGetRow(Amat,Ii,&ncols,&idx,&vals); CHKERRQ(ierr); 190 for( jj = 0 ; jj < ncols ; jj++ ){ 191 PetscInt dest_col = idx[jj]/bs; 192 PetscScalar v = 1.0; 193 ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES); CHKERRQ(ierr); 194 } 195 ierr = MatRestoreRow(Amat,Ii,&ncols,&idx,&vals); CHKERRQ(ierr); 196 } 197 ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 198 ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 199 200 ierr = MatConvert( tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj ); CHKERRQ(ierr); 201 202 ierr = MatDestroy( &tMat ); CHKERRQ(ierr); 203 } 204 if( membershipKey == 0 ) { 205 /* hack to fix global data that pmetis.c uses in 'adj' */ 206 for( kk=0 , jj=0 ; kk<=npe ; jj++, kk += old_fact ) { 207 adj->rmap->range[jj] = adj->rmap->range[kk]; 208 } 209 ierr = MatPartitioningCreate( new_comm, &mpart ); CHKERRQ(ierr); 210 ierr = MatPartitioningSetAdjacency( mpart, adj ); CHKERRQ(ierr); 211 ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr); 212 ierr = MatPartitioningSetNParts( mpart, new_npe ); CHKERRQ(ierr); 213 ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr); 214 ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr); 215 /* collect IS info */ 216 ierr = ISGetLocalSize( isnewproc, &is_sz ); CHKERRQ(ierr); 217 ierr = PetscMalloc( bs*is_sz*sizeof(PetscInt), &isnewproc_idx ); CHKERRQ(ierr); 218 ierr = ISGetIndices( isnewproc, &is_idx ); CHKERRQ(ierr); 219 /* spread partitioning across machine - probably the right thing to do but machine spec. */ 220 for(kk=0,jj=0;kk<is_sz;kk++){ 221 for(ii=0 ; ii<bs ; ii++, jj++ ) { /* expand for equation level by 'bs' */ 222 isnewproc_idx[jj] = is_idx[kk] * new_fact; 223 } 224 } 225 ierr = ISRestoreIndices( isnewproc, &is_idx ); CHKERRQ(ierr); 226 ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); 227 is_sz *= bs; 228 } 229 else { 230 isnewproc_idx = 0; 231 is_sz = 0; 232 } 233 ierr = MatDestroy( &adj ); CHKERRQ(ierr); 234 ierr = MPI_Comm_free( &new_comm ); CHKERRQ(ierr); 235 ierr = ISCreateGeneral( wcomm, is_sz, isnewproc_idx, PETSC_COPY_VALUES, &isnewproc ); 236 if( membershipKey == 0 ) { 237 ierr = PetscFree( isnewproc_idx ); CHKERRQ(ierr); 238 } 239 } 240 /* 241 Create an index set from the isnewproc index set to indicate the mapping TO 242 */ 243 ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr); 244 /* 245 Determine how many elements are assigned to each processor 246 */ 247 ierr = ISPartitioningCount( isnewproc, npe, counts ); CHKERRQ(ierr); 248 ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); 249 ncrs_new = counts[mype]/bs; 250 } 251 252 { /* Create a vector to contain the newly ordered element information */ 253 const PetscInt *idx; 254 PetscInt i,j,k; 255 IS isscat; 256 PetscScalar *array; 257 Vec src_crd, dest_crd; 258 PetscReal *data = *a_coarse_data; 259 VecScatter vecscat; 260 PetscInt tidx[ncrs0*a_data_sz]; 261 262 ierr = VecCreate( wcomm, &dest_crd ); 263 ierr = VecSetSizes( dest_crd, a_data_sz*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr); 264 ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /*funny vector-get global options?*/ 265 /* 266 There are 'a_data_sz' data items per node, (one can think of the vectors of having a 267 block size of 'a_data_sz'). Note, ISs are expanded into equation space by 'bs'. 268 */ 269 ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr); 270 for(i=0,j=0; i<ncrs0 ; i++) { 271 PetscInt lid = idx[i*bs]/bs; assert(idx[i*bs]%bs==0); 272 for( k=0; k<a_data_sz ; k++, j++) tidx[j] = lid*a_data_sz + k; 273 } 274 ierr = ISCreateGeneral( wcomm, a_data_sz*ncrs0, tidx, PETSC_COPY_VALUES, &isscat ); 275 CHKERRQ(ierr); 276 ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr); 277 /* 278 Create a vector to contain the original vertex information for each element 279 */ 280 ierr = VecCreateSeq( PETSC_COMM_SELF, a_data_sz*ncrs0, &src_crd ); CHKERRQ(ierr); 281 for (i=0; i<a_data_sz*ncrs0; i++) { 282 ierr = VecSetValues(src_crd, 1, &i, &data[i], INSERT_VALUES ); CHKERRQ(ierr); 283 } 284 ierr = VecAssemblyBegin(src_crd); CHKERRQ(ierr); 285 ierr = VecAssemblyEnd(src_crd); CHKERRQ(ierr); 286 /* 287 Scatter the element vertex information (still in the original vertex ordering) 288 to the correct processor 289 */ 290 ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat); 291 CHKERRQ(ierr); 292 ierr = ISDestroy( &isscat ); CHKERRQ(ierr); 293 ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 294 ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 295 ierr = VecScatterDestroy( &vecscat ); CHKERRQ(ierr); 296 ierr = VecDestroy( &src_crd ); CHKERRQ(ierr); 297 /* 298 Put the element vertex data into a new allocation of the gdata->ele 299 */ 300 ierr = PetscFree( *a_coarse_data ); CHKERRQ(ierr); 301 ierr = PetscMalloc( a_data_sz*ncrs_new*sizeof(PetscReal), a_coarse_data ); CHKERRQ(ierr); 302 VecGetLocalSize( dest_crd, &i ); assert(i==a_data_sz*ncrs_new); 303 304 ierr = VecGetArray( dest_crd, &array ); CHKERRQ(ierr); 305 data = *a_coarse_data; 306 for (i=0; i<a_data_sz*ncrs_new; i++) data[i] = PetscRealPart(array[i]); 307 ierr = VecRestoreArray( dest_crd, &array ); CHKERRQ(ierr); 308 ierr = VecDestroy( &dest_crd ); CHKERRQ(ierr); 309 } 310 /* 311 Invert for MatGetSubMatrix 312 */ 313 ierr = ISInvertPermutation( isnum, ncrs_new*bs, &new_indices ); CHKERRQ(ierr); 314 ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */ 315 ierr = ISDestroy( &isnum ); CHKERRQ(ierr); 316 /* A_crs output */ 317 ierr = MatGetSubMatrix( Amat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs ); 318 CHKERRQ(ierr); 319 320 ierr = MatDestroy( &Amat ); CHKERRQ(ierr); 321 Amat = *a_Amat_crs; 322 ierr = MatSetBlockSize( Amat, bs ); CHKERRQ(ierr); 323 324 /* prolongator */ 325 ierr = MatGetOwnershipRange( Pold, &Istart, &Iend ); CHKERRQ(ierr); 326 ierr = MatGetBlockSize( Pold, &bs2 ); CHKERRQ(ierr); assert(bs==bs2); 327 { 328 IS findices; 329 ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices); CHKERRQ(ierr); 330 ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew ); 331 CHKERRQ(ierr); 332 ierr = ISDestroy( &findices ); CHKERRQ(ierr); 333 ierr = MatSetBlockSize( Pnew, bs ); CHKERRQ(ierr); 334 } 335 ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr); 336 *a_P_inout = Pnew; /* output */ 337 ierr = MatGetBlockSize( Pnew, &bs2 ); CHKERRQ(ierr); assert(bs==bs2); 338 ierr = ISDestroy( &new_indices ); CHKERRQ(ierr); 339 340 PetscFunctionReturn(0); 341 } 342 343 #define GAMG_MAXLEVELS 30 344 #if defined(PETSC_USE_LOG) 345 PetscLogStage gamg_stages[20]; 346 #endif 347 /* -------------------------------------------------------------------------- */ 348 /* 349 PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 350 by setting data structures and options. 351 352 Input Parameter: 353 . pc - the preconditioner context 354 355 Application Interface Routine: PCSetUp() 356 357 Notes: 358 The interface routine PCSetUp() is not usually called directly by 359 the user, but instead is called by PCApply() if necessary. 360 */ 361 #undef __FUNCT__ 362 #define __FUNCT__ "PCSetUp_GAMG" 363 PetscErrorCode PCSetUp_GAMG( PC a_pc ) 364 { 365 PetscErrorCode ierr; 366 PC_MG *mg = (PC_MG*)a_pc->data; 367 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 368 Mat Amat = a_pc->mat, Pmat = a_pc->pmat; 369 PetscBool isSeq, isMPI; 370 PetscInt fine_level, level, level1, M, N, bs, nloc, lidx, data_sz, Istart, Iend; 371 MPI_Comm wcomm = ((PetscObject)a_pc)->comm; 372 PetscMPIInt mype,npe,nactivepe; 373 PetscBool isOK, useSA = PETSC_FALSE, flag; 374 Mat Aarr[GAMG_MAXLEVELS], Parr[GAMG_MAXLEVELS]; 375 PetscReal *coarse_data = 0, *data, emaxs[GAMG_MAXLEVELS]; 376 char str[16]; 377 378 PetscFunctionBegin; 379 if( a_pc->setupcalled ) { 380 /* no state data in GAMG to destroy */ 381 ierr = PCReset_MG( a_pc ); CHKERRQ(ierr); 382 } 383 ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 384 ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 385 /* setup special features of PCGAMG */ 386 ierr = PetscTypeCompare((PetscObject)Amat, MATSEQAIJ, &isSeq);CHKERRQ(ierr); 387 ierr = PetscTypeCompare((PetscObject)Amat, MATMPIAIJ, &isMPI);CHKERRQ(ierr); 388 if (isMPI) { 389 } else if (isSeq) { 390 } else SETERRQ1(wcomm,PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with GAMG. GAMG can only handle AIJ matrices.",((PetscObject)Amat)->type_name); 391 392 /* GAMG requires input of fine-grid matrix. It determines nlevels. */ 393 ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr); 394 ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr); 395 ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr); 396 nloc = (Iend-Istart)/bs; assert((Iend-Istart)%bs == 0); 397 398 ierr = PetscOptionsGetString(PETSC_NULL,"-pc_gamg_type",str,16,&flag); CHKERRQ( ierr ); 399 useSA = (PetscBool)(flag && strcmp(str,"sa") == 0); 400 if( pc_gamg->m_data == 0 ) { 401 useSA = PETSC_TRUE; /* use SA if no data */ 402 ierr = PCSetCoordinates_GAMG( a_pc, 1, 0 ); CHKERRQ( ierr ); 403 } 404 data = pc_gamg->m_data; 405 /* Get A_i and R_i */ 406 data_sz = pc_gamg->m_data_sz/nloc; 407 PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s level %d N=%d, data size=%d m_data_sz=%d nloc=%d\n",0,__FUNCT__,0,N,data_sz,pc_gamg->m_data_sz,nloc); 408 for ( level=0, Aarr[0] = Pmat, nactivepe = npe; /* hard wired stopping logic */ 409 level < GAMG_MAXLEVELS-1 && (level==0 || M/bs>TOP_GRID_LIM) && (npe==1 || nactivepe>1); 410 level++ ) { 411 level1 = level + 1; 412 ierr = PetscLogEventBegin(gamg_setup_stages[SET1],0,0,0,0);CHKERRQ(ierr); 413 ierr = createProlongation( Aarr[level], data, pc_gamg->m_dim, &data_sz, 414 &Parr[level1], &coarse_data, &isOK, &emaxs[level] ); 415 CHKERRQ(ierr); 416 ierr = PetscLogEventEnd(gamg_setup_stages[SET1],0,0,0,0);CHKERRQ(ierr); 417 418 ierr = PetscFree( data ); CHKERRQ( ierr ); 419 if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */ 420 if( isOK ) { 421 ierr = PetscLogEventBegin(gamg_setup_stages[SET2],0,0,0,0);CHKERRQ(ierr); 422 ierr = partitionLevel( Aarr[level], pc_gamg->m_data_sz/nloc, 423 &Parr[level1], &coarse_data, &nactivepe, &Aarr[level1] ); 424 CHKERRQ(ierr); 425 ierr = PetscLogEventEnd(gamg_setup_stages[SET2],0,0,0,0);CHKERRQ(ierr); 426 ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr); 427 PetscPrintf(PETSC_COMM_WORLD,"\t\t[%d]%s %d) N=%d, %d data size, %d active pes\n",0,__FUNCT__,level1,N,data_sz,nactivepe); 428 } 429 else{ 430 break; 431 } 432 data = coarse_data; 433 } 434 ierr = PetscFree( coarse_data ); CHKERRQ( ierr ); 435 PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1); 436 pc_gamg->m_data = 0; /* destroyed coordinate data */ 437 pc_gamg->m_Nlevels = level + 1; 438 fine_level = level; 439 ierr = PCMGSetLevels(a_pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr); 440 441 /* set default smoothers */ 442 for ( lidx=1, level = pc_gamg->m_Nlevels-2; 443 lidx <= fine_level; 444 lidx++, level--) { 445 PetscReal emax, emin; 446 KSP smoother; PC subpc; 447 ierr = PCMGGetSmoother( a_pc, lidx, &smoother ); CHKERRQ(ierr); 448 ierr = KSPSetType( smoother, KSPCHEBYCHEV );CHKERRQ(ierr); 449 if( emaxs[level] > 0.0 ) { 450 emax = 1.05*emaxs[level]; 451 } 452 else{ /* eigen estimate 'emax' */ 453 KSP eksp; Mat Lmat = Aarr[level]; 454 Vec bb, xx; PC pc; 455 PetscInt N1, N0, tt; 456 ierr = MatGetVecs( Lmat, &bb, 0 ); CHKERRQ(ierr); 457 ierr = MatGetVecs( Lmat, &xx, 0 ); CHKERRQ(ierr); 458 { 459 PetscRandom rctx; 460 ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr); 461 ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 462 ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 463 ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr); 464 } 465 ierr = KSPCreate(wcomm,&eksp);CHKERRQ(ierr); 466 ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr); 467 ierr = KSPSetOperators( eksp, Lmat, Lmat, DIFFERENT_NONZERO_PATTERN ); CHKERRQ( ierr ); 468 ierr = KSPGetPC( eksp, &pc );CHKERRQ( ierr ); 469 ierr = PCSetType( pc, PCJACOBI ); CHKERRQ(ierr); /* should be same as above */ 470 ierr = KSPSetTolerances( eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 5 ); 471 CHKERRQ(ierr); 472 ierr = KSPSetConvergenceTest( eksp, KSPSkipConverged, 0, 0 ); CHKERRQ(ierr); 473 ierr = KSPSetNormType( eksp, KSP_NORM_NONE ); CHKERRQ(ierr); 474 475 ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr); 476 ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr); 477 ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr); 478 ierr = MatGetSize( Lmat, &N1, &tt ); CHKERRQ(ierr); 479 ierr = MatGetSize( Aarr[level+1], &N0, &tt );CHKERRQ(ierr); 480 PetscPrintf(PETSC_COMM_WORLD,"\t\t\t%s max eigen = %e (N=%d)\n",__FUNCT__,emax,N1/bs); 481 emax *= 1.05; 482 483 ierr = VecDestroy( &xx ); CHKERRQ(ierr); 484 ierr = VecDestroy( &bb ); CHKERRQ(ierr); 485 ierr = KSPDestroy( &eksp ); CHKERRQ(ierr); 486 487 emin = emax/((PetscReal)N1/(PetscReal)N0); /* this should be about the coarsening rate */ 488 } 489 ierr = KSPSetOperators( smoother, Aarr[level], Aarr[level], DIFFERENT_NONZERO_PATTERN ); 490 ierr = KSPChebychevSetEigenvalues( smoother, emax, emin );CHKERRQ(ierr); 491 ierr = KSPGetPC( smoother, &subpc ); CHKERRQ(ierr); 492 ierr = PCSetType( subpc, PCJACOBI ); CHKERRQ(ierr); 493 ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 494 } 495 { 496 KSP smoother; /* coarse grid */ 497 Mat Lmat = Aarr[pc_gamg->m_Nlevels-1]; 498 ierr = PCMGGetSmoother( a_pc, 0, &smoother ); CHKERRQ(ierr); 499 ierr = KSPSetOperators( smoother, Lmat, Lmat, DIFFERENT_NONZERO_PATTERN ); 500 CHKERRQ(ierr); 501 ierr = KSPSetNormType( smoother, KSP_NORM_NONE ); CHKERRQ(ierr); 502 } 503 /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 504 ierr = PCSetFromOptions_MG(a_pc); CHKERRQ(ierr); 505 { 506 PetscBool galerkin; 507 ierr = PCMGGetGalerkin( a_pc, &galerkin); CHKERRQ(ierr); 508 if(galerkin){ 509 SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG."); 510 } 511 } 512 513 /* set interpolation between the levels, create timer stages, clean up */ 514 if( PETSC_FALSE ) { 515 char str[32]; 516 sprintf(str,"MG Level %d (%d)",0,pc_gamg->m_Nlevels-1); 517 PetscLogStageRegister(str, &gamg_stages[fine_level]); 518 } 519 for (lidx=0,level=pc_gamg->m_Nlevels-1; 520 lidx<fine_level; 521 lidx++, level--){ 522 ierr = PCMGSetInterpolation( a_pc, lidx+1, Parr[level] );CHKERRQ(ierr); 523 if( !PETSC_TRUE ) { 524 PetscViewer viewer; char fname[32]; 525 sprintf(fname,"Amat_%d.m",level); 526 ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); 527 ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 528 ierr = MatView( Aarr[level], viewer ); CHKERRQ(ierr); 529 ierr = PetscViewerDestroy( &viewer ); 530 } 531 ierr = MatDestroy( &Parr[level] ); CHKERRQ(ierr); 532 ierr = MatDestroy( &Aarr[level] ); CHKERRQ(ierr); 533 if( PETSC_FALSE ) { 534 char str[32]; 535 sprintf(str,"MG Level %d (%d)",lidx+1,level-1); 536 PetscLogStageRegister(str, &gamg_stages[level-1]); 537 } 538 } 539 540 /* setupcalled is set to 0 so that MG is setup from scratch */ 541 a_pc->setupcalled = 0; 542 ierr = PCSetUp_MG(a_pc);CHKERRQ(ierr); 543 PetscFunctionReturn(0); 544 } 545 546 /* ------------------------------------------------------------------------- */ 547 /* 548 PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 549 that was created with PCCreate_GAMG(). 550 551 Input Parameter: 552 . pc - the preconditioner context 553 554 Application Interface Routine: PCDestroy() 555 */ 556 #undef __FUNCT__ 557 #define __FUNCT__ "PCDestroy_GAMG" 558 PetscErrorCode PCDestroy_GAMG(PC pc) 559 { 560 PetscErrorCode ierr; 561 PC_MG *mg = (PC_MG*)pc->data; 562 PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 563 564 PetscFunctionBegin; 565 ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 566 ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 567 ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 568 PetscFunctionReturn(0); 569 } 570 571 #undef __FUNCT__ 572 #define __FUNCT__ "PCSetFromOptions_GAMG" 573 PetscErrorCode PCSetFromOptions_GAMG(PC pc) 574 { 575 /* PetscErrorCode ierr; */ 576 /* PC_MG *mg = (PC_MG*)pc->data; */ 577 /* PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; */ 578 /* MPI_Comm comm = ((PetscObject)pc)->comm; */ 579 580 PetscFunctionBegin; 581 PetscFunctionReturn(0); 582 } 583 584 /* -------------------------------------------------------------------------- */ 585 /* 586 PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG 587 588 Input Parameter: 589 . pc - the preconditioner context 590 591 Application Interface Routine: PCCreate() 592 593 */ 594 /* MC 595 PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide 596 fine grid discretization matrix and coordinates on the fine grid. 597 598 Options Database Key: 599 Multigrid options(inherited) 600 + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) 601 . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) 602 . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) 603 -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 604 GAMG options: 605 606 Level: intermediate 607 Concepts: multigrid 608 609 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 610 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), 611 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 612 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 613 PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 614 M */ 615 616 EXTERN_C_BEGIN 617 #undef __FUNCT__ 618 #define __FUNCT__ "PCCreate_GAMG" 619 PetscErrorCode PCCreate_GAMG(PC pc) 620 { 621 PetscErrorCode ierr; 622 PC_GAMG *pc_gamg; 623 PC_MG *mg; 624 PetscClassId cookie; 625 626 PetscFunctionBegin; 627 /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 628 ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 629 ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr); 630 631 /* create a supporting struct and attach it to pc */ 632 ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr); 633 mg = (PC_MG*)pc->data; 634 mg->innerctx = pc_gamg; 635 636 pc_gamg->m_Nlevels = -1; 637 638 /* overwrite the pointers of PCMG by the functions of PCGAMG */ 639 pc->ops->setfromoptions = PCSetFromOptions_GAMG; 640 pc->ops->setup = PCSetUp_GAMG; 641 pc->ops->reset = PCReset_GAMG; 642 pc->ops->destroy = PCDestroy_GAMG; 643 644 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 645 "PCSetCoordinates_C", 646 "PCSetCoordinates_GAMG", 647 PCSetCoordinates_GAMG);CHKERRQ(ierr); 648 649 PetscClassIdRegister("GAMG Setup",&cookie); 650 PetscLogEventRegister("GAMG-createProl", cookie, &gamg_setup_stages[SET1]); 651 PetscLogEventRegister("GAMG-partLevel", cookie, &gamg_setup_stages[SET2]); 652 PetscLogEventRegister("GAMG-MIS Graph", cookie, &gamg_setup_stages[SET3]); 653 PetscLogEventRegister("GAMG-MIS-Agg", cookie, &gamg_setup_stages[SET4]); 654 PetscLogEventRegister("GAMG-growSupp", cookie, &gamg_setup_stages[SET5]); 655 PetscLogEventRegister("GAMG-tri-Prol", cookie, &gamg_setup_stages[SET6]); 656 PetscLogEventRegister("GAMG-find-prol", cookie, &gamg_setup_stages[FIND_V]); 657 658 PetscFunctionReturn(0); 659 } 660 EXTERN_C_END 661