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