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