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