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