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