15b89ad90SMark F. Adams /* 20cd22d39SHong Zhang GAMG geometric-algebric multigrid PC - Mark Adams 2011 35b89ad90SMark F. Adams */ 4af0996ceSBarry Smith #include <petsc/private/matimpl.h> 5389730f3SMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 65b42dca8SJed Brown #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /* Hack to access same_local_solves */ 7*18c3aa7eSMark #include <../src/ksp/ksp/impls/cheby/chebyshevimpl.h> /*I "petscksp.h" I*/ 8f96513f1SMatthew G Knepley 90cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 100cbbd2e1SMark F. Adams PetscLogEvent petsc_gamg_setup_events[NUM_SET]; 11b4fbaa2aSMark F. Adams #endif 120cbbd2e1SMark F. Adams 130cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 14fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_AGG; 15fd1112cbSBarry Smith PetscLogEvent PC_GAMGGraph_GEO; 160cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_AGG; 170cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGCoarsen_GEO; 180cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_AGG; 190cbbd2e1SMark F. Adams PetscLogEvent PC_GAMGProlongator_GEO; 20fd1112cbSBarry Smith PetscLogEvent PC_GAMGOptProlongator_AGG; 210cbbd2e1SMark F. Adams #endif 220cbbd2e1SMark F. Adams 23b8fd24d8SMark F. Adams /* #define GAMG_STAGES */ 240cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 25*18c3aa7eSMark static PetscLogStage gamg_stages[PETSC_MG_MAXLEVELS]; 26b4fbaa2aSMark F. Adams #endif 27f96513f1SMatthew G Knepley 28140e18c1SBarry Smith static PetscFunctionList GAMGList = 0; 293e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized; 309d5b6da9SMark F. Adams 31d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */ 32d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 33d3d6bff4SMark F. Adams { 34*18c3aa7eSMark PetscErrorCode ierr, level; 35d3d6bff4SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 36d3d6bff4SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 37d3d6bff4SMark F. Adams 38d3d6bff4SMark F. Adams PetscFunctionBegin; 3922a233eaSStefano Zampini ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 401c1aac46SBarry Smith pc_gamg->data_sz = 0; 41878e152fSMark F. Adams ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr); 42*18c3aa7eSMark for (level = 0; level < PETSC_MG_MAXLEVELS ; level++) { 43*18c3aa7eSMark mg->min_eigen_DinvA[level] = 0; 44*18c3aa7eSMark mg->max_eigen_DinvA[level] = 0; 45*18c3aa7eSMark } 46*18c3aa7eSMark pc_gamg->emin = 0; 47*18c3aa7eSMark pc_gamg->emax = 0; 48a2f3521dSMark F. Adams PetscFunctionReturn(0); 49a2f3521dSMark F. Adams } 50a2f3521dSMark F. Adams 515b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 525b89ad90SMark F. Adams /* 53c238b0ebSToby Isaac PCGAMGCreateLevel_GAMG: create coarse op with RAP. repartition and/or reduce number 54a147abb0SMark F. Adams of active processors. 555b89ad90SMark F. Adams 565b89ad90SMark F. Adams Input Parameter: 57a2f3521dSMark F. Adams . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 58a2f3521dSMark F. Adams 'pc_gamg->data_sz' are changed via repartitioning/reduction. 599d5b6da9SMark F. Adams . Amat_fine - matrix on this fine (k) level 60c5bfad50SMark F. Adams . cr_bs - coarse block size 613530afc2SMark F. Adams In/Output Parameter: 62a2f3521dSMark F. Adams . a_P_inout - prolongation operator to the next level (k-->k-1) 63afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 6411e60469SMark F. Adams Output Parameter: 653530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 665b89ad90SMark F. Adams */ 675cb416c2SMark F. Adams 68171cca9aSMark Adams static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc,Mat Amat_fine,PetscInt cr_bs,Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc,IS * Pcolumnperm, PetscBool is_last) 695b89ad90SMark F. Adams { 70a2f3521dSMark F. Adams PetscErrorCode ierr; 719d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 72486a8d0bSJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 73a2f3521dSMark F. Adams Mat Cmat,Pold=*a_P_inout; 743b4367a7SBarry Smith MPI_Comm comm; 75c5df96a5SBarry Smith PetscMPIInt rank,size,new_size,nactive=*a_nactive_proc; 763ae0bb68SMark Adams PetscInt ncrs_eq,ncrs,f_bs; 775b89ad90SMark F. Adams 785b89ad90SMark F. Adams PetscFunctionBegin; 793b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr); 803b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 813b4367a7SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 82c5bfad50SMark F. Adams ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr); 839d5b6da9SMark F. Adams ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr); 84038e3b61SMark F. Adams 85ce7c7f2fSMark Adams if (Pcolumnperm) *Pcolumnperm = NULL; 86ce7c7f2fSMark Adams 873ae0bb68SMark Adams /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/ 880298fd71SBarry Smith ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr); 893ae0bb68SMark Adams if (pc_gamg->data_cell_rows>0) { 903ae0bb68SMark Adams ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows; 9173911c69SBarry Smith } else { 923ae0bb68SMark Adams PetscInt bs; 933ae0bb68SMark Adams ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr); 943ae0bb68SMark Adams ncrs = ncrs_eq/bs; 953ae0bb68SMark Adams } 96c5df96a5SBarry Smith /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */ 97171cca9aSMark Adams if (is_last && !pc_gamg->use_parallel_coarse_grid_solver) new_size = 1; 98171cca9aSMark Adams else { 99472110cdSMark F. Adams PetscInt ncrs_eq_glob; 1000298fd71SBarry Smith ierr = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr); 101a90e85d9SMark Adams new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 1027f66b68fSMark Adams if (!new_size) new_size = 1; /* not likely, posible? */ 103c5df96a5SBarry Smith else if (new_size >= nactive) new_size = nactive; /* no change, rare */ 104a2f3521dSMark F. Adams } 105f852f58cSMark F. Adams 1062e3501ffSMark Adams if (new_size==nactive) { 107ef3f0257SMark Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 108ce7c7f2fSMark Adams if (new_size < size) { 109ce7c7f2fSMark Adams /* odd case where multiple coarse grids are on one processor or no coarsening ... */ 11031cb4603SMark Adams ierr = PetscInfo1(pc,"reduced grid using same number of processors (%d) as last grid (use larger coarse grid)\n",nactive);CHKERRQ(ierr); 111ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 112ce7c7f2fSMark Adams ierr = MatPinToCPU(*a_Amat_crs,PETSC_TRUE);CHKERRQ(ierr); 113ce7c7f2fSMark Adams ierr = MatPinToCPU(*a_P_inout,PETSC_TRUE);CHKERRQ(ierr); 114ce7c7f2fSMark Adams } 115ce7c7f2fSMark Adams } 116ef3f0257SMark Adams /* we know that the grid structure can be reused in MatPtAP */ 1172e3501ffSMark Adams } else { /* reduce active processors - we know that the grid structure can NOT be reused in MatPtAP */ 118192c0e8bSMark Adams PetscInt *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old,expand_factor=1,rfactor=1; 119885364a3SMark Adams IS is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices; 12071959b99SBarry Smith nloc_old = ncrs_eq/cr_bs; 12171959b99SBarry Smith if (ncrs_eq % cr_bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ncrs_eq %D not divisible by cr_bs %D",ncrs_eq,cr_bs); 122ce7c7f2fSMark Adams /* get new_size and rfactor */ 123ce7c7f2fSMark Adams if (pc_gamg->layout_type==PCGAMG_LAYOUT_SPREAD || !pc_gamg->repart) { 124ce7c7f2fSMark Adams /* find factor */ 125ce7c7f2fSMark Adams if (new_size == 1) rfactor = size; /* don't modify */ 126ce7c7f2fSMark Adams else { 127ce7c7f2fSMark Adams PetscReal best_fact = 0.; 128ce7c7f2fSMark Adams jj = -1; 129ce7c7f2fSMark Adams for (kk = 1 ; kk <= size ; kk++) { 130ce7c7f2fSMark Adams if (!(size%kk)) { /* a candidate */ 131ce7c7f2fSMark Adams PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size; 132ce7c7f2fSMark Adams if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */ 133ce7c7f2fSMark Adams if (fact > best_fact) { 134ce7c7f2fSMark Adams best_fact = fact; jj = kk; 135ce7c7f2fSMark Adams } 136ce7c7f2fSMark Adams } 137ce7c7f2fSMark Adams } 138ce7c7f2fSMark Adams if (jj != -1) rfactor = jj; 139ce7c7f2fSMark Adams else rfactor = 1; /* a prime */ 140ce7c7f2fSMark Adams if (pc_gamg->layout_type == PCGAMG_LAYOUT_COMPACT) expand_factor = 1; 141ce7c7f2fSMark Adams else expand_factor = rfactor; 142ce7c7f2fSMark Adams } 143ce7c7f2fSMark Adams new_size = size/rfactor; /* make new size one that is factor */ 1444cdfd227SMark if (new_size==nactive) { /* no repartitioning or reduction, bail out because nested here (rare) */ 1454cdfd227SMark *a_Amat_crs = Cmat; 14631cb4603SMark Adams ierr = PetscInfo2(pc,"Finding factorable processor set stopped reduction: new_size=%d, neq(loc)=%D\n",new_size,ncrs_eq);CHKERRQ(ierr); 147ce7c7f2fSMark Adams PetscFunctionReturn(0); 148ce7c7f2fSMark Adams } 149ce7c7f2fSMark Adams } 1504cdfd227SMark #if defined PETSC_GAMG_USE_LOG 1514cdfd227SMark ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 1524cdfd227SMark #endif 153a2f3521dSMark F. Adams /* make 'is_eq_newproc' */ 154785e854fSJed Brown ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr); 1552e3501ffSMark Adams if (pc_gamg->repart) { 156a2f3521dSMark F. Adams /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */ 1575a9b9e01SMark F. Adams Mat adj; 15831cb4603SMark Adams ierr = PetscInfo4(pc,"Repartition: size (active): %d --> %d, %D local equations, using %s process layout\n",*a_nactive_proc, new_size, ncrs_eq, (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) ? "compact" : "spread");CHKERRQ(ierr); 159a2f3521dSMark F. Adams /* get 'adj' */ 160c5bfad50SMark F. Adams if (cr_bs == 1) { 161038e3b61SMark F. Adams ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 162806fa848SBarry Smith } else { 163a2f3521dSMark F. Adams /* make a scalar matrix to partition (no Stokes here) */ 164eb07cef2SMark F. Adams Mat tMat; 165a2f3521dSMark F. Adams PetscInt Istart_crs,Iend_crs,ncols,jj,Ii; 166b4fbaa2aSMark F. Adams const PetscScalar *vals; 167b4fbaa2aSMark F. Adams const PetscInt *idx; 168a2f3521dSMark F. Adams PetscInt *d_nnz, *o_nnz, M, N; 16939d09545SMark Adams static PetscInt llev = 0; /* ugly but just used for debugging */ 170d9558ea9SBarry Smith MatType mtype; 171b4fbaa2aSMark F. Adams 172e632b94dSBarry Smith ierr = PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);CHKERRQ(ierr); 173a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr); 174a2f3521dSMark F. Adams ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr); 175c5bfad50SMark F. Adams for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) { 17658471d46SMark F. Adams ierr = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr); 177c5bfad50SMark F. Adams d_nnz[jj] = ncols/cr_bs; 178c5bfad50SMark F. Adams o_nnz[jj] = ncols/cr_bs; 17958471d46SMark F. Adams ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr); 1803ae0bb68SMark Adams if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs; 1813ae0bb68SMark Adams if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs; 18258471d46SMark F. Adams } 1836876a03eSMark F. Adams 184d9558ea9SBarry Smith ierr = MatGetType(Amat_fine,&mtype);CHKERRQ(ierr); 1853b4367a7SBarry Smith ierr = MatCreate(comm, &tMat);CHKERRQ(ierr); 1863ae0bb68SMark Adams ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 187d9558ea9SBarry Smith ierr = MatSetType(tMat,mtype);CHKERRQ(ierr); 188a2f3521dSMark F. Adams ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr); 189a2f3521dSMark F. Adams ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 190e632b94dSBarry Smith ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr); 191eb07cef2SMark F. Adams 192a2f3521dSMark F. Adams for (ii = Istart_crs; ii < Iend_crs; ii++) { 193c5bfad50SMark F. Adams PetscInt dest_row = ii/cr_bs; 19422063be5SMark F. Adams ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 195eb07cef2SMark F. Adams for (jj = 0; jj < ncols; jj++) { 196c5bfad50SMark F. Adams PetscInt dest_col = idx[jj]/cr_bs; 197eb07cef2SMark F. Adams PetscScalar v = 1.0; 198eb07cef2SMark F. Adams ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr); 199eb07cef2SMark F. Adams } 20022063be5SMark F. Adams ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 201eb07cef2SMark F. Adams } 202eb07cef2SMark F. Adams ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 203eb07cef2SMark F. Adams ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 204eb07cef2SMark F. Adams 205b4fbaa2aSMark F. Adams if (llev++ == -1) { 206b4fbaa2aSMark F. Adams PetscViewer viewer; char fname[32]; 2078caf3d72SBarry Smith ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr); 2083b4367a7SBarry Smith PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer); 209b4fbaa2aSMark F. Adams ierr = MatView(tMat, viewer);CHKERRQ(ierr); 2103bf036e2SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 211b4fbaa2aSMark F. Adams } 212eb07cef2SMark F. Adams ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 213eb07cef2SMark F. Adams ierr = MatDestroy(&tMat);CHKERRQ(ierr); 214a2f3521dSMark F. Adams } /* create 'adj' */ 215f150b916SMark F. Adams 216a2f3521dSMark F. Adams { /* partition: get newproc_idx */ 2175a9b9e01SMark F. Adams char prefix[256]; 2185a9b9e01SMark F. Adams const char *pcpre; 219b4fbaa2aSMark F. Adams const PetscInt *is_idx; 220b4fbaa2aSMark F. Adams MatPartitioning mpart; 221a4b7d37bSMark F. Adams IS proc_is; 2222f03bc48SMark F. Adams 2233b4367a7SBarry Smith ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr); 2245ef31b24SMark F. Adams ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr); 2259d5b6da9SMark F. Adams ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr); 2268caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr); 22759a0be82SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 22811e60469SMark F. Adams ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 229c5df96a5SBarry Smith ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr); 230a4b7d37bSMark F. Adams ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr); 23111e60469SMark F. Adams ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 2325a9b9e01SMark F. Adams 2335ef31b24SMark F. Adams /* collect IS info */ 234785e854fSJed Brown ierr = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr); 235a4b7d37bSMark F. Adams ierr = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr); 236a2f3521dSMark F. Adams for (kk = jj = 0 ; kk < nloc_old ; kk++) { 237c5bfad50SMark F. Adams for (ii = 0 ; ii < cr_bs ; ii++, jj++) { 238ce7c7f2fSMark Adams newproc_idx[jj] = is_idx[kk] * expand_factor; /* distribution */ 239eb07cef2SMark F. Adams } 2405ef31b24SMark F. Adams } 241a4b7d37bSMark F. Adams ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr); 242a4b7d37bSMark F. Adams ierr = ISDestroy(&proc_is);CHKERRQ(ierr); 2435ef31b24SMark F. Adams } 2445ef31b24SMark F. Adams ierr = MatDestroy(&adj);CHKERRQ(ierr); 2455a9b9e01SMark F. Adams 2463b4367a7SBarry Smith ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr); 2478263b398SMark F. Adams ierr = PetscFree(newproc_idx);CHKERRQ(ierr); 24831cb4603SMark Adams } else { /* simple aggregation of parts -- 'is_eq_newproc' */ 249ce7c7f2fSMark Adams PetscInt targetPE; 2504cdfd227SMark if (new_size==nactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"new_size==nactive. Should not happen"); 251302440fdSBarry Smith ierr = PetscInfo1(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);CHKERRQ(ierr); 252ce7c7f2fSMark Adams targetPE = (rank/rfactor)*expand_factor; 2533b4367a7SBarry Smith ierr = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr); 254a2f3521dSMark F. Adams } /* end simple 'is_eq_newproc' */ 255e33ef3b1SMark F. Adams 25611e60469SMark F. Adams /* 257a2f3521dSMark F. Adams Create an index set from the is_eq_newproc index set to indicate the mapping TO 25811e60469SMark F. Adams */ 259a2f3521dSMark F. Adams ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr); 2607700e67bSMark Adams is_eq_num_prim = is_eq_num; 26111e60469SMark F. Adams /* 262a2f3521dSMark F. Adams Determine how many equations/vertices are assigned to each processor 26311e60469SMark F. Adams */ 264c5df96a5SBarry Smith ierr = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr); 265c5df96a5SBarry Smith ncrs_eq_new = counts[rank]; 266a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_newproc);CHKERRQ(ierr); 267ce7c7f2fSMark Adams ncrs_new = ncrs_eq_new/cr_bs; 268a2f3521dSMark F. Adams 269a2f3521dSMark F. Adams ierr = PetscFree(counts);CHKERRQ(ierr); 270885364a3SMark Adams /* data movement scope -- this could be moved to subclasses so that we don't try to cram all auxilary data into some complex abstracted thing */ 271885364a3SMark Adams { 272885364a3SMark Adams Vec src_crd, dest_crd; 273885364a3SMark Adams const PetscInt *idx,ndata_rows=pc_gamg->data_cell_rows,ndata_cols=pc_gamg->data_cell_cols,node_data_sz=ndata_rows*ndata_cols; 274885364a3SMark Adams VecScatter vecscat; 275885364a3SMark Adams PetscScalar *array; 276885364a3SMark Adams IS isscat; 277a2f3521dSMark F. Adams /* move data (for primal equations only) */ 27822063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 2793b4367a7SBarry Smith ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr); 2803ae0bb68SMark Adams ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr); 281c0dedaeaSBarry Smith ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */ 28211e60469SMark F. Adams /* 2839d5b6da9SMark F. Adams There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 284c5bfad50SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'cr_bs'. 28511e60469SMark F. Adams */ 286854ce69bSBarry Smith ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr); 287a2f3521dSMark F. Adams ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 2883ae0bb68SMark Adams for (ii=0,jj=0; ii<ncrs; ii++) { 289c5bfad50SMark F. Adams PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */ 290a2f3521dSMark F. Adams for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk; 29111e60469SMark F. Adams } 292a2f3521dSMark F. Adams ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 2933ae0bb68SMark Adams ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr); 29492a756f0SMark F. Adams ierr = PetscFree(tidx);CHKERRQ(ierr); 29511e60469SMark F. Adams /* 29611e60469SMark F. Adams Create a vector to contain the original vertex information for each element 29711e60469SMark F. Adams */ 2983ae0bb68SMark Adams ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr); 2999d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols; jj++) { 3003ae0bb68SMark Adams const PetscInt stride0=ncrs*pc_gamg->data_cell_rows; 3013ae0bb68SMark Adams for (ii=0; ii<ncrs; ii++) { 3029d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows; kk++) { 303a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj; 304c8b0795cSMark F. Adams PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 305676e1743SMark F. Adams ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr); 306d3d6bff4SMark F. Adams } 307038e3b61SMark F. Adams } 308eb07cef2SMark F. Adams } 309eb07cef2SMark F. Adams ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr); 310eb07cef2SMark F. Adams ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr); 31111e60469SMark F. Adams /* 31211e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 31311e60469SMark F. Adams to the correct processor 31411e60469SMark F. Adams */ 3159448b7f1SJunchao Zhang ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr); 31611e60469SMark F. Adams ierr = ISDestroy(&isscat);CHKERRQ(ierr); 31711e60469SMark F. Adams ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 31811e60469SMark F. Adams ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 31911e60469SMark F. Adams ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 32011e60469SMark F. Adams ierr = VecDestroy(&src_crd);CHKERRQ(ierr); 32111e60469SMark F. Adams /* 32211e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 32311e60469SMark F. Adams */ 324c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 325578f55a3SPeter Brune ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr); 3262fa5cd67SKarl Rupp 3273ae0bb68SMark Adams pc_gamg->data_sz = node_data_sz*ncrs_new; 3283ae0bb68SMark Adams strideNew = ncrs_new*ndata_rows; 3292fa5cd67SKarl Rupp 33011e60469SMark F. Adams ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr); 3319d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols; jj++) { 3323ae0bb68SMark Adams for (ii=0; ii<ncrs_new; ii++) { 3339d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows; kk++) { 334a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj; 335c8b0795cSMark F. Adams pc_gamg->data[ix] = PetscRealPart(array[jx]); 336d3d6bff4SMark F. Adams } 337038e3b61SMark F. Adams } 338038e3b61SMark F. Adams } 33911e60469SMark F. Adams ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr); 34011e60469SMark F. Adams ierr = VecDestroy(&dest_crd);CHKERRQ(ierr); 341885364a3SMark Adams } 342a2f3521dSMark F. Adams /* move A and P (columns) with new layout */ 3430cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3440cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 345ed3f9983SMark F. Adams #endif 34611e60469SMark F. Adams /* 3477dae84e0SHong Zhang Invert for MatCreateSubMatrix 34811e60469SMark F. Adams */ 349a2f3521dSMark F. Adams ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr); 350a2f3521dSMark F. Adams ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */ 351c5bfad50SMark F. Adams ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr); 352a2f3521dSMark F. Adams if (is_eq_num != is_eq_num_prim) { 353a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */ 354a2f3521dSMark F. Adams } 3553cb8563fSToby Isaac if (Pcolumnperm) { 3563cb8563fSToby Isaac ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr); 3573cb8563fSToby Isaac *Pcolumnperm = new_eq_indices; 3583cb8563fSToby Isaac } 359a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr); 3600cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3610cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 3620cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 363ed3f9983SMark F. Adams #endif 364a2f3521dSMark F. Adams /* 'a_Amat_crs' output */ 365a2f3521dSMark F. Adams { 366a2f3521dSMark F. Adams Mat mat; 3677dae84e0SHong Zhang ierr = MatCreateSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr); 368a2f3521dSMark F. Adams *a_Amat_crs = mat; 369a2f3521dSMark F. Adams } 370038e3b61SMark F. Adams ierr = MatDestroy(&Cmat);CHKERRQ(ierr); 371a2f3521dSMark F. Adams 3720cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3730cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 374ed3f9983SMark F. Adams #endif 37511e60469SMark F. Adams /* prolongator */ 37611e60469SMark F. Adams { 37711e60469SMark F. Adams IS findices; 378a2f3521dSMark F. Adams PetscInt Istart,Iend; 379a2f3521dSMark F. Adams Mat Pnew; 38062294041SBarry Smith 381a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr); 3820cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3830cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 384ed3f9983SMark F. Adams #endif 3853b4367a7SBarry Smith ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr); 386c5bfad50SMark F. Adams ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr); 3877dae84e0SHong Zhang ierr = MatCreateSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr); 38811e60469SMark F. Adams ierr = ISDestroy(&findices);CHKERRQ(ierr); 389c5bfad50SMark F. Adams 3900cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3910cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 392ed3f9983SMark F. Adams #endif 3933530afc2SMark F. Adams ierr = MatDestroy(a_P_inout);CHKERRQ(ierr); 394a2f3521dSMark F. Adams 395a2f3521dSMark F. Adams /* output - repartitioned */ 396a2f3521dSMark F. Adams *a_P_inout = Pnew; 397e33ef3b1SMark F. Adams } 398a2f3521dSMark F. Adams ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr); 3995b89ad90SMark F. Adams 400c5df96a5SBarry Smith *a_nactive_proc = new_size; /* output */ 401ce7c7f2fSMark Adams 402ce7c7f2fSMark Adams /* pinning on reduced grids, not a bad heuristic and optimization gets folded into process reduction optimization */ 403ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 404ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 4058bca76a6SMark Adams static PetscInt llev = 2; 40639d09545SMark Adams ierr = PetscInfo1(pc,"Pinning level %D to the CPU\n",llev++);CHKERRQ(ierr); 407ce7c7f2fSMark Adams #endif 408ce7c7f2fSMark Adams ierr = MatPinToCPU(*a_Amat_crs,PETSC_TRUE);CHKERRQ(ierr); 409ce7c7f2fSMark Adams ierr = MatPinToCPU(*a_P_inout,PETSC_TRUE);CHKERRQ(ierr); 410ce7c7f2fSMark Adams if (1) { /* lvec is created, need to pin it, this is done in MatSetUpMultiply_MPIAIJ. Hack */ 411ce7c7f2fSMark Adams Mat A = *a_Amat_crs, P = *a_P_inout; 412ce7c7f2fSMark Adams PetscMPIInt size; 413ce7c7f2fSMark Adams ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 414ce7c7f2fSMark Adams if (size > 1) { 415ce7c7f2fSMark Adams Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data, *p = (Mat_MPIAIJ*)P->data; 416ce7c7f2fSMark Adams ierr = VecPinToCPU(a->lvec,PETSC_TRUE);CHKERRQ(ierr); 417ce7c7f2fSMark Adams ierr = VecPinToCPU(p->lvec,PETSC_TRUE);CHKERRQ(ierr); 418ce7c7f2fSMark Adams } 419ce7c7f2fSMark Adams } 420ce7c7f2fSMark Adams } 4214cdfd227SMark #if defined PETSC_GAMG_USE_LOG 4224cdfd227SMark ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 4234cdfd227SMark #endif 424a2f3521dSMark F. Adams } 4255b89ad90SMark F. Adams PetscFunctionReturn(0); 4265b89ad90SMark F. Adams } 4275b89ad90SMark F. Adams 4285b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 4295b89ad90SMark F. Adams /* 4305b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 4315b89ad90SMark F. Adams by setting data structures and options. 4325b89ad90SMark F. Adams 4335b89ad90SMark F. Adams Input Parameter: 4345b89ad90SMark F. Adams . pc - the preconditioner context 4355b89ad90SMark F. Adams 4365b89ad90SMark F. Adams */ 4379d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc) 4385b89ad90SMark F. Adams { 4395b89ad90SMark F. Adams PetscErrorCode ierr; 4409d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 4415b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 4422adcac29SMark F. Adams Mat Pmat = pc->pmat; 443*18c3aa7eSMark PetscInt fine_level,level,level1,bs,M,N,qq,lidx,nASMBlocksArr[PETSC_MG_MAXLEVELS]; 4443b4367a7SBarry Smith MPI_Comm comm; 445c5df96a5SBarry Smith PetscMPIInt rank,size,nactivepe; 446*18c3aa7eSMark Mat Aarr[PETSC_MG_MAXLEVELS],Parr[PETSC_MG_MAXLEVELS]; 447*18c3aa7eSMark IS *ASMLocalIDsArr[PETSC_MG_MAXLEVELS]; 448a2f3521dSMark F. Adams PetscLogDouble nnz0=0.,nnztot=0.; 449569f4572SMark Adams MatInfo info; 450171cca9aSMark Adams PetscBool is_last = PETSC_FALSE; 4515ef31b24SMark F. Adams 4525b89ad90SMark F. Adams PetscFunctionBegin; 4533b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 4543b4367a7SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 4553b4367a7SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 456dfd5c07aSMark F. Adams 45784d3f75bSMark F. Adams if (pc_gamg->setup_count++ > 0) { 4581c1aac46SBarry Smith if ((PetscBool)(!pc_gamg->reuse_prol)) { 459878e152fSMark F. Adams /* reset everything */ 460878e152fSMark F. Adams ierr = PCReset_MG(pc);CHKERRQ(ierr); 461878e152fSMark F. Adams pc->setupcalled = 0; 462806fa848SBarry Smith } else { 46384d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 46403a628feSMark F. Adams /* just do Galerkin grids */ 46558471d46SMark F. Adams Mat B,dA,dB; 46658471d46SMark F. Adams 46771959b99SBarry Smith if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet"); 4689d5b6da9SMark F. Adams if (pc_gamg->Nlevels > 1) { 46958471d46SMark F. Adams /* currently only handle case where mat and pmat are the same on coarser levels */ 47023ee1639SBarry Smith ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 47158471d46SMark F. Adams /* (re)set to get dirty flag */ 47223ee1639SBarry Smith ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr); 47358471d46SMark F. Adams 4742fb0b348SMark F. Adams for (level=pc_gamg->Nlevels-2; level>=0; level--) { 475ef3f0257SMark Adams /* 2nd solve, matrix structure can change from repartitioning or process reduction but don't know if we have process reduction here. Should fix */ 476ef3f0257SMark Adams if (pc_gamg->setup_count==2 /* && pc_gamg->repart||reduction */) { 477ef3f0257SMark Adams ierr = PetscInfo2(pc,"new RAP after first solve level %D, %D setup\n",level,pc_gamg->setup_count);CHKERRQ(ierr); 478ef3f0257SMark Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,2.0,&B);CHKERRQ(ierr); 479084a8fe3SJed Brown ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr); 48003a628feSMark F. Adams mglevels[level]->A = B; 481806fa848SBarry Smith } else { 482ef3f0257SMark Adams ierr = PetscInfo2(pc,"RAP after first solve reusing matrix level %D, %D setup\n",level,pc_gamg->setup_count);CHKERRQ(ierr); 48323ee1639SBarry Smith ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr); 48458471d46SMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 48503a628feSMark F. Adams } 48623ee1639SBarry Smith ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr); 48758471d46SMark F. Adams dB = B; 48858471d46SMark F. Adams } 4895f8cf99dSMark F. Adams } 490d5280255SMark F. Adams 491d5280255SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 49258471d46SMark F. Adams PetscFunctionReturn(0); 493eb07cef2SMark F. Adams } 494878e152fSMark F. Adams } 495f6536408SMark F. Adams 496878e152fSMark F. Adams if (!pc_gamg->data) { 497878e152fSMark F. Adams if (pc_gamg->orig_data) { 498878e152fSMark F. Adams ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 4990298fd71SBarry Smith ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr); 5002fa5cd67SKarl Rupp 501878e152fSMark F. Adams pc_gamg->data_sz = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols; 502878e152fSMark F. Adams pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 503878e152fSMark F. Adams pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 5042fa5cd67SKarl Rupp 505785e854fSJed Brown ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr); 506878e152fSMark F. Adams for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 507806fa848SBarry Smith } else { 5081ab5ffc9SJed Brown if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data"); 5097700e67bSMark Adams ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr); 5109d5b6da9SMark F. Adams } 511878e152fSMark F. Adams } 512878e152fSMark F. Adams 513878e152fSMark F. Adams /* cache original data for reuse */ 5141c1aac46SBarry Smith if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) { 515785e854fSJed Brown ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr); 516878e152fSMark F. Adams for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 517878e152fSMark F. Adams pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 518878e152fSMark F. Adams pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 519878e152fSMark F. Adams } 520038e3b61SMark F. Adams 521302f38e8SMark F. Adams /* get basic dims */ 522302f38e8SMark F. Adams ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 523171cca9aSMark Adams ierr = MatGetSize(Pmat, &M, &N);CHKERRQ(ierr); 52484d3f75bSMark F. Adams 525569f4572SMark Adams ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */ 526569f4572SMark Adams nnz0 = info.nz_used; 527569f4572SMark Adams nnztot = info.nz_used; 52862294041SBarry Smith ierr = PetscInfo6(pc,"level %d) N=%D, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n",0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols,(int)(nnz0/(PetscReal)M+0.5),size);CHKERRQ(ierr); 529569f4572SMark Adams 530a2f3521dSMark F. Adams /* Get A_i and R_i */ 53162294041SBarry Smith for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) { 5329ab59c8bSMark Adams pc_gamg->current_level = level; 533*18c3aa7eSMark if (level >= PETSC_MG_MAXLEVELS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Too many levels %D",level); 5345b89ad90SMark F. Adams level1 = level + 1; 5350cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 5360cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 537a2f3521dSMark F. Adams #if (defined GAMG_STAGES) 538a2f3521dSMark F. Adams ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr); 539b4fbaa2aSMark F. Adams #endif 540a2f3521dSMark F. Adams #endif 541c8b0795cSMark F. Adams { /* construct prolongator */ 542725b86d8SJed Brown Mat Gmat; 5430cbbd2e1SMark F. Adams PetscCoarsenData *agg_lists; 5447700e67bSMark Adams Mat Prol11; 545c8b0795cSMark F. Adams 5467700e67bSMark Adams ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr); 5471ab5ffc9SJed Brown ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr); 5487700e67bSMark Adams ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr); 549c8b0795cSMark F. Adams 550a2f3521dSMark F. Adams /* could have failed to create new level */ 551a2f3521dSMark F. Adams if (Prol11) { 5529d5b6da9SMark F. Adams /* get new block size of coarse matrices */ 5530298fd71SBarry Smith ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr); 554a2f3521dSMark F. Adams 555fd1112cbSBarry Smith if (pc_gamg->ops->optprolongator) { 556c8b0795cSMark F. Adams /* smooth */ 557fd1112cbSBarry Smith ierr = pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);CHKERRQ(ierr); 558c8b0795cSMark F. Adams } 559c8b0795cSMark F. Adams 5600c3bc534SBarry Smith if (pc_gamg->use_aggs_in_asm) { 5611b18a24aSMark Adams PetscInt bs; 5621b18a24aSMark Adams ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr); 5630a3c815dSMark Adams ierr = PetscCDGetASMBlocks(agg_lists, bs, Gmat, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr); 564ffc955d6SMark F. Adams } 565ffc955d6SMark F. Adams 5664bde40a0SMark Adams Parr[level1] = Prol11; 5674bde40a0SMark Adams } else Parr[level1] = NULL; /* failed to coarsen */ 5684bde40a0SMark Adams 569a2f3521dSMark F. Adams ierr = MatDestroy(&Gmat);CHKERRQ(ierr); 57041b27cdeSMark F. Adams ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr); 571a2f3521dSMark F. Adams } /* construct prolongator scope */ 5720cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 5730cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 574c8b0795cSMark F. Adams #endif 5757f66b68fSMark Adams if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 576171cca9aSMark Adams if (!Parr[level1]) { /* failed to coarsen */ 577569f4572SMark Adams ierr = PetscInfo1(pc,"Stop gridding, level %D\n",level);CHKERRQ(ierr); 57862294041SBarry Smith #if defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES 579a90e85d9SMark Adams ierr = PetscLogStagePop();CHKERRQ(ierr); 580a90e85d9SMark Adams #endif 581c8b0795cSMark F. Adams break; 582c8b0795cSMark F. Adams } 5830cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 5840cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 585b4fbaa2aSMark F. Adams #endif 586171cca9aSMark Adams ierr = MatGetSize(Parr[level1], &M, &N);CHKERRQ(ierr); /* N is next M, a loop test variables */ 587171cca9aSMark Adams if (is_last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Is last ????????"); 588171cca9aSMark Adams if (N <= pc_gamg->coarse_eq_limit) is_last = PETSC_TRUE; 5890e2909e1SMark Adams if (level1 == pc_gamg->Nlevels-1) is_last = PETSC_TRUE; 590171cca9aSMark Adams ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs, &Parr[level1], &Aarr[level1], &nactivepe, NULL, is_last);CHKERRQ(ierr); 591a2f3521dSMark F. Adams 5920cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 5930cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 594b4fbaa2aSMark F. Adams #endif 595171cca9aSMark Adams ierr = MatGetSize(Aarr[level1], &M, &N);CHKERRQ(ierr); /* M is loop test variables */ 596569f4572SMark Adams ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr); 597569f4572SMark Adams nnztot += info.nz_used; 5981d5b2942SMark Adams ierr = PetscInfo5(pc,"%d) N=%D, n data cols=%d, nnz/row (ave)=%d, %d active pes\n",level1,M,pc_gamg->data_cell_cols,(int)(info.nz_used/(PetscReal)M),nactivepe);CHKERRQ(ierr); 599569f4572SMark Adams 6000cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 601b4fbaa2aSMark F. Adams ierr = PetscLogStagePop();CHKERRQ(ierr); 602b4fbaa2aSMark F. Adams #endif 603a90e85d9SMark Adams /* stop if one node or one proc -- could pull back for singular problems */ 6049ab59c8bSMark Adams if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2) ) { 6059ab59c8bSMark Adams ierr = PetscInfo2(pc,"HARD stop of coarsening on level %D. Grid too small: %D block nodes\n",level,M/bs);CHKERRQ(ierr); 606a90e85d9SMark Adams level++; 607a90e85d9SMark Adams break; 608a90e85d9SMark Adams } 609c8b0795cSMark F. Adams } /* levels */ 610c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 611c8b0795cSMark F. Adams 612569f4572SMark Adams ierr = PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);CHKERRQ(ierr); 6139d5b6da9SMark F. Adams pc_gamg->Nlevels = level + 1; 6145b89ad90SMark F. Adams fine_level = level; 6150298fd71SBarry Smith ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr); 6165b89ad90SMark F. Adams 61762294041SBarry Smith if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */ 618d5280255SMark F. Adams /* set default smoothers & set operators */ 61962294041SBarry Smith for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) { 620ffc955d6SMark F. Adams KSP smoother; 621ffc955d6SMark F. Adams PC subpc; 622a2f3521dSMark F. Adams 6239d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 624f6536408SMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 625ffc955d6SMark F. Adams 626a2f3521dSMark F. Adams ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 627a2f3521dSMark F. Adams /* set ops */ 62823ee1639SBarry Smith ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr); 629a2f3521dSMark F. Adams ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr); 630a2f3521dSMark F. Adams 631a2f3521dSMark F. Adams /* set defaults */ 6326c9de887SHong Zhang ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr); 633a2f3521dSMark F. Adams 6340c3bc534SBarry Smith /* set blocks for ASM smoother that uses the 'aggregates' */ 6350c3bc534SBarry Smith if (pc_gamg->use_aggs_in_asm) { 6362d3561bbSSatish Balay PetscInt sz; 6377a28f3e5SMark Adams IS *iss; 638a2f3521dSMark F. Adams 6392d3561bbSSatish Balay sz = nASMBlocksArr[level]; 6407a28f3e5SMark Adams iss = ASMLocalIDsArr[level]; 6410c3bc534SBarry Smith ierr = PCSetType(subpc, PCASM);CHKERRQ(ierr); 6420a3c815dSMark Adams ierr = PCASMSetOverlap(subpc, 0);CHKERRQ(ierr); 6430c3bc534SBarry Smith ierr = PCASMSetType(subpc,PC_ASM_BASIC);CHKERRQ(ierr); 6447f66b68fSMark Adams if (!sz) { 645ffc955d6SMark F. Adams IS is; 6460a3c815dSMark Adams ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);CHKERRQ(ierr); 6477a28f3e5SMark Adams ierr = PCASMSetLocalSubdomains(subpc, 1, NULL, &is);CHKERRQ(ierr); 648a94c3b12SMark F. Adams ierr = ISDestroy(&is);CHKERRQ(ierr); 649806fa848SBarry Smith } else { 650a94c3b12SMark F. Adams PetscInt kk; 6517a28f3e5SMark Adams ierr = PCASMSetLocalSubdomains(subpc, sz, NULL, iss);CHKERRQ(ierr); 652a94c3b12SMark F. Adams for (kk=0; kk<sz; kk++) { 6537a28f3e5SMark Adams ierr = ISDestroy(&iss[kk]);CHKERRQ(ierr); 654a94c3b12SMark F. Adams } 6557a28f3e5SMark Adams ierr = PetscFree(iss);CHKERRQ(ierr); 656ffc955d6SMark F. Adams } 6570298fd71SBarry Smith ASMLocalIDsArr[level] = NULL; 658ffc955d6SMark F. Adams nASMBlocksArr[level] = 0; 659806fa848SBarry Smith } else { 660890ffe84SMark Adams ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr); 661ffc955d6SMark F. Adams } 662d5280255SMark F. Adams } 663d5280255SMark F. Adams { 664d5280255SMark F. Adams /* coarse grid */ 665d5280255SMark F. Adams KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 666d5280255SMark F. Adams Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0; 667d5280255SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 66823ee1639SBarry Smith ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr); 669cf8ae1d3SMark Adams if (!pc_gamg->use_parallel_coarse_grid_solver) { 670d5280255SMark F. Adams ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 671d5280255SMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 672d5280255SMark F. Adams ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr); 673d5280255SMark F. Adams ierr = PCSetUp(subpc);CHKERRQ(ierr); 67471959b99SBarry Smith ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); 67571959b99SBarry Smith if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii); 676d5280255SMark F. Adams ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 677d5280255SMark F. Adams ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr); 6789dbfc187SHong Zhang ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 6792fb0b348SMark F. Adams ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr); 68008e36f19SMark Adams ierr = KSPSetType(k2[0], KSPPREONLY);CHKERRQ(ierr); 6815b42dca8SJed Brown /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in 6825b42dca8SJed Brown * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to 6835b42dca8SJed Brown * view every subdomain as though they were different. */ 6845b42dca8SJed Brown ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE; 685d5280255SMark F. Adams } 686cf8ae1d3SMark Adams } 687d5280255SMark F. Adams 688d5280255SMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 689d5280255SMark F. Adams ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 690e55864a3SBarry Smith ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr); 691d5280255SMark F. Adams ierr = PetscOptionsEnd();CHKERRQ(ierr); 69269aca0b8SBarry Smith ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr); 693d5280255SMark F. Adams 694*18c3aa7eSMark ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 695*18c3aa7eSMark 696*18c3aa7eSMark /* setup cheby eigen estimates from SA */ 697*18c3aa7eSMark for (lidx = 1, level = pc_gamg->Nlevels-2; level >= 0 ; lidx++, level--) { 698*18c3aa7eSMark KSP smoother; 699*18c3aa7eSMark PetscBool ischeb; 700*18c3aa7eSMark ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 701*18c3aa7eSMark ierr = PetscObjectTypeCompare((PetscObject)smoother,KSPCHEBYSHEV,&ischeb);CHKERRQ(ierr); 702*18c3aa7eSMark if (ischeb) { 703*18c3aa7eSMark KSP_Chebyshev *cheb = (KSP_Chebyshev*)smoother->data; 704*18c3aa7eSMark if (mg->max_eigen_DinvA[level] > 0 && cheb->emax == 0.) { /* let command line emax override using SA's eigenvalues */ 705*18c3aa7eSMark PC subpc; 706*18c3aa7eSMark PetscBool isjac; 707*18c3aa7eSMark ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 708*18c3aa7eSMark ierr = PetscObjectTypeCompare((PetscObject)subpc,PCJACOBI,&isjac);CHKERRQ(ierr); 709*18c3aa7eSMark if ( (isjac && pc_gamg->use_sa_esteig==-1) || pc_gamg->use_sa_esteig==1) { 710*18c3aa7eSMark PetscReal emax,emin; 711*18c3aa7eSMark emin = mg->min_eigen_DinvA[level]; 712*18c3aa7eSMark emax = mg->max_eigen_DinvA[level]; 713*18c3aa7eSMark ierr = PetscInfo4(pc,"PCSetUp_GAMG: call KSPChebyshevSetEigenvalues on level %D (N=%D) with emax = %g emin = %g\n",level,Aarr[level]->rmap->N,(double)emax,(double)emin);CHKERRQ(ierr); 714*18c3aa7eSMark cheb->emin_computed = emin; 715*18c3aa7eSMark cheb->emax_computed = emax; 716*18c3aa7eSMark ierr = KSPChebyshevSetEigenvalues(smoother, cheb->tform[2]*emin + cheb->tform[3]*emax, cheb->tform[0]*emin + cheb->tform[1]*emax);CHKERRQ(ierr); 717*18c3aa7eSMark } 718*18c3aa7eSMark } 719*18c3aa7eSMark } 720*18c3aa7eSMark } 721*18c3aa7eSMark 722d5280255SMark F. Adams /* clean up */ 723d5280255SMark F. Adams for (level=1; level<pc_gamg->Nlevels; level++) { 724587fa25dSMark F. Adams ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr); 725587fa25dSMark F. Adams ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr); 7265b89ad90SMark F. Adams } 727*18c3aa7eSMark 728806fa848SBarry Smith } else { 7295f8cf99dSMark F. Adams KSP smoother; 730302440fdSBarry Smith ierr = PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");CHKERRQ(ierr); 7319d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr); 73223ee1639SBarry Smith ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr); 7335f8cf99dSMark F. Adams ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr); 7349d5b6da9SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 7355f8cf99dSMark F. Adams } 7365b89ad90SMark F. Adams PetscFunctionReturn(0); 7375b89ad90SMark F. Adams } 7385b89ad90SMark F. Adams 739eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */ 7405b89ad90SMark F. Adams /* 7415b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 7425b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 7435b89ad90SMark F. Adams 7445b89ad90SMark F. Adams Input Parameter: 7455b89ad90SMark F. Adams . pc - the preconditioner context 7465b89ad90SMark F. Adams 7475b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 7485b89ad90SMark F. Adams */ 7495b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc) 7505b89ad90SMark F. Adams { 7515b89ad90SMark F. Adams PetscErrorCode ierr; 7525b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 7535b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 7545b89ad90SMark F. Adams 7555b89ad90SMark F. Adams PetscFunctionBegin; 7565b89ad90SMark F. Adams ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 7579b8ffb57SJed Brown if (pc_gamg->ops->destroy) { 7589b8ffb57SJed Brown ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr); 7599b8ffb57SJed Brown } 7601ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr); 7611ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr); 7625b89ad90SMark F. Adams ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 7635b89ad90SMark F. Adams ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 7645b89ad90SMark F. Adams PetscFunctionReturn(0); 7655b89ad90SMark F. Adams } 7665b89ad90SMark F. Adams 767676e1743SMark F. Adams /*@ 768cab9ed1eSBarry Smith PCGAMGSetProcEqLim - Set number of equations to aim for per process on the coarse grids via processor reduction. 769676e1743SMark F. Adams 7701cc46a46SBarry Smith Logically Collective on PC 771676e1743SMark F. Adams 772676e1743SMark F. Adams Input Parameters: 7731cc46a46SBarry Smith + pc - the preconditioner context 7741cc46a46SBarry Smith - n - the number of equations 775676e1743SMark F. Adams 776676e1743SMark F. Adams 777676e1743SMark F. Adams Options Database Key: 7781cc46a46SBarry Smith . -pc_gamg_process_eq_limit <limit> 779676e1743SMark F. Adams 78095452b02SPatrick Sanan Notes: 78195452b02SPatrick Sanan GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> equations on each process 782cab9ed1eSBarry Smith that has degrees of freedom 783cab9ed1eSBarry Smith 784676e1743SMark F. Adams Level: intermediate 785676e1743SMark F. Adams 7861c1aac46SBarry Smith .seealso: PCGAMGSetCoarseEqLim() 787676e1743SMark F. Adams @*/ 788676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 789676e1743SMark F. Adams { 790676e1743SMark F. Adams PetscErrorCode ierr; 791676e1743SMark F. Adams 792676e1743SMark F. Adams PetscFunctionBegin; 793676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 794676e1743SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 795676e1743SMark F. Adams PetscFunctionReturn(0); 796676e1743SMark F. Adams } 797676e1743SMark F. Adams 7981e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 799676e1743SMark F. Adams { 800c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 801c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 802676e1743SMark F. Adams 803676e1743SMark F. Adams PetscFunctionBegin; 8049d5b6da9SMark F. Adams if (n>0) pc_gamg->min_eq_proc = n; 805676e1743SMark F. Adams PetscFunctionReturn(0); 806676e1743SMark F. Adams } 807676e1743SMark F. Adams 808389730f3SMark F. Adams /*@ 809cab9ed1eSBarry Smith PCGAMGSetCoarseEqLim - Set maximum number of equations on coarsest grid. 810389730f3SMark F. Adams 811389730f3SMark F. Adams Collective on PC 812389730f3SMark F. Adams 813389730f3SMark F. Adams Input Parameters: 8141cc46a46SBarry Smith + pc - the preconditioner context 8151cc46a46SBarry Smith - n - maximum number of equations to aim for 816389730f3SMark F. Adams 817389730f3SMark F. Adams Options Database Key: 8181cc46a46SBarry Smith . -pc_gamg_coarse_eq_limit <limit> 819389730f3SMark F. Adams 82074329af1SBarry Smith Notes: For example -pc_gamg_coarse_eq_limit 1000 will stop coarsening once the coarse grid 82174329af1SBarry Smith has less than 1000 unknowns. 82274329af1SBarry Smith 823389730f3SMark F. Adams Level: intermediate 824389730f3SMark F. Adams 8251c1aac46SBarry Smith .seealso: PCGAMGSetProcEqLim() 826389730f3SMark F. Adams @*/ 827389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 828389730f3SMark F. Adams { 829389730f3SMark F. Adams PetscErrorCode ierr; 830389730f3SMark F. Adams 831389730f3SMark F. Adams PetscFunctionBegin; 832389730f3SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 833389730f3SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 834389730f3SMark F. Adams PetscFunctionReturn(0); 835389730f3SMark F. Adams } 836389730f3SMark F. Adams 8371e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 838389730f3SMark F. Adams { 839389730f3SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 840389730f3SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 841389730f3SMark F. Adams 842389730f3SMark F. Adams PetscFunctionBegin; 8439d5b6da9SMark F. Adams if (n>0) pc_gamg->coarse_eq_limit = n; 844389730f3SMark F. Adams PetscFunctionReturn(0); 845389730f3SMark F. Adams } 846389730f3SMark F. Adams 847676e1743SMark F. Adams /*@ 848cab9ed1eSBarry Smith PCGAMGSetRepartition - Repartition the degrees of freedom across the processors on the coarser grids 849676e1743SMark F. Adams 850676e1743SMark F. Adams Collective on PC 851676e1743SMark F. Adams 852676e1743SMark F. Adams Input Parameters: 8531cc46a46SBarry Smith + pc - the preconditioner context 8541cc46a46SBarry Smith - n - PETSC_TRUE or PETSC_FALSE 855676e1743SMark F. Adams 856676e1743SMark F. Adams Options Database Key: 8571cc46a46SBarry Smith . -pc_gamg_repartition <true,false> 858676e1743SMark F. Adams 85995452b02SPatrick Sanan Notes: 86095452b02SPatrick Sanan this will generally improve the loading balancing of the work on each level 861cab9ed1eSBarry Smith 862676e1743SMark F. Adams Level: intermediate 863676e1743SMark F. Adams 864676e1743SMark F. Adams .seealso: () 865676e1743SMark F. Adams @*/ 866cab9ed1eSBarry Smith PetscErrorCode PCGAMGSetRepartition(PC pc, PetscBool n) 867676e1743SMark F. Adams { 868676e1743SMark F. Adams PetscErrorCode ierr; 869676e1743SMark F. Adams 870676e1743SMark F. Adams PetscFunctionBegin; 871676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 872cab9ed1eSBarry Smith ierr = PetscTryMethod(pc,"PCGAMGSetRepartition_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 873676e1743SMark F. Adams PetscFunctionReturn(0); 874676e1743SMark F. Adams } 875676e1743SMark F. Adams 876cab9ed1eSBarry Smith static PetscErrorCode PCGAMGSetRepartition_GAMG(PC pc, PetscBool n) 877676e1743SMark F. Adams { 878c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 879c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 880676e1743SMark F. Adams 881676e1743SMark F. Adams PetscFunctionBegin; 8829d5b6da9SMark F. Adams pc_gamg->repart = n; 883676e1743SMark F. Adams PetscFunctionReturn(0); 884676e1743SMark F. Adams } 885676e1743SMark F. Adams 886dfd5c07aSMark F. Adams /*@ 887*18c3aa7eSMark PCGAMGSetEstEigKSPMaxIt - Set number of KSP iterations in eigen estimator (for Cheby) 888*18c3aa7eSMark 889*18c3aa7eSMark Collective on PC 890*18c3aa7eSMark 891*18c3aa7eSMark Input Parameters: 892*18c3aa7eSMark + pc - the preconditioner context 893*18c3aa7eSMark - n - number of its 894*18c3aa7eSMark 895*18c3aa7eSMark Options Database Key: 896*18c3aa7eSMark . -pc_gamg_esteig_ksp_max_it <its> 897*18c3aa7eSMark 898*18c3aa7eSMark Notes: 899*18c3aa7eSMark 900*18c3aa7eSMark Level: intermediate 901*18c3aa7eSMark 902*18c3aa7eSMark .seealso: () 903*18c3aa7eSMark @*/ 904*18c3aa7eSMark PetscErrorCode PCGAMGSetEstEigKSPMaxIt(PC pc, PetscInt n) 905*18c3aa7eSMark { 906*18c3aa7eSMark PetscErrorCode ierr; 907*18c3aa7eSMark 908*18c3aa7eSMark PetscFunctionBegin; 909*18c3aa7eSMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 910*18c3aa7eSMark ierr = PetscTryMethod(pc,"PCGAMGSetEstEigKSPMaxIt_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 911*18c3aa7eSMark PetscFunctionReturn(0); 912*18c3aa7eSMark } 913*18c3aa7eSMark 914*18c3aa7eSMark static PetscErrorCode PCGAMGSetEstEigKSPMaxIt_GAMG(PC pc, PetscInt n) 915*18c3aa7eSMark { 916*18c3aa7eSMark PC_MG *mg = (PC_MG*)pc->data; 917*18c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 918*18c3aa7eSMark 919*18c3aa7eSMark PetscFunctionBegin; 920*18c3aa7eSMark pc_gamg->esteig_max_it = n; 921*18c3aa7eSMark PetscFunctionReturn(0); 922*18c3aa7eSMark } 923*18c3aa7eSMark 924*18c3aa7eSMark /*@ 925*18c3aa7eSMark PCGAMGSetUseSAEstEig - Use eigen estimate from smoothed aggregation for Cheby smoother 926*18c3aa7eSMark 927*18c3aa7eSMark Collective on PC 928*18c3aa7eSMark 929*18c3aa7eSMark Input Parameters: 930*18c3aa7eSMark + pc - the preconditioner context 931*18c3aa7eSMark - n - number of its 932*18c3aa7eSMark 933*18c3aa7eSMark Options Database Key: 934*18c3aa7eSMark . -pc_gamg_use_sa_esteig <true,false> 935*18c3aa7eSMark 936*18c3aa7eSMark Notes: 937*18c3aa7eSMark 938*18c3aa7eSMark Level: intermediate 939*18c3aa7eSMark 940*18c3aa7eSMark .seealso: () 941*18c3aa7eSMark @*/ 942*18c3aa7eSMark PetscErrorCode PCGAMGSetUseSAEstEig(PC pc, PetscBool n) 943*18c3aa7eSMark { 944*18c3aa7eSMark PetscErrorCode ierr; 945*18c3aa7eSMark 946*18c3aa7eSMark PetscFunctionBegin; 947*18c3aa7eSMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 948*18c3aa7eSMark ierr = PetscTryMethod(pc,"PCGAMGSetUseSAEstEig_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 949*18c3aa7eSMark PetscFunctionReturn(0); 950*18c3aa7eSMark } 951*18c3aa7eSMark 952*18c3aa7eSMark static PetscErrorCode PCGAMGSetUseSAEstEig_GAMG(PC pc, PetscInt n) 953*18c3aa7eSMark { 954*18c3aa7eSMark PC_MG *mg = (PC_MG*)pc->data; 955*18c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 956*18c3aa7eSMark 957*18c3aa7eSMark PetscFunctionBegin; 958*18c3aa7eSMark pc_gamg->use_sa_esteig = n ? 1 : 0; 959*18c3aa7eSMark PetscFunctionReturn(0); 960*18c3aa7eSMark } 961*18c3aa7eSMark 962*18c3aa7eSMark /*@C 963*18c3aa7eSMark PCGAMGSetEstEigKSPType - Set type of KSP in eigen estimator (for Cheby) 964*18c3aa7eSMark 965*18c3aa7eSMark Collective on PC 966*18c3aa7eSMark 967*18c3aa7eSMark Input Parameters: 968*18c3aa7eSMark + pc - the preconditioner context 969*18c3aa7eSMark - t - ksp type 970*18c3aa7eSMark 971*18c3aa7eSMark Options Database Key: 972*18c3aa7eSMark . -pc_gamg_esteig_ksp_type <type> 973*18c3aa7eSMark 974*18c3aa7eSMark Notes: 975*18c3aa7eSMark 976*18c3aa7eSMark Level: intermediate 977*18c3aa7eSMark 978*18c3aa7eSMark .seealso: () 979*18c3aa7eSMark @*/ 980*18c3aa7eSMark PetscErrorCode PCGAMGSetEstEigKSPType(PC pc, char t[]) 981*18c3aa7eSMark { 982*18c3aa7eSMark PetscErrorCode ierr; 983*18c3aa7eSMark 984*18c3aa7eSMark PetscFunctionBegin; 985*18c3aa7eSMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 986*18c3aa7eSMark ierr = PetscTryMethod(pc,"PCGAMGSetEstEigKSPType_C",(PC,char[]),(pc,t));CHKERRQ(ierr); 987*18c3aa7eSMark PetscFunctionReturn(0); 988*18c3aa7eSMark } 989*18c3aa7eSMark 990*18c3aa7eSMark static PetscErrorCode PCGAMGSetEstEigKSPType_GAMG(PC pc, char t[]) 991*18c3aa7eSMark { 992*18c3aa7eSMark PetscErrorCode ierr; 993*18c3aa7eSMark PC_MG *mg = (PC_MG*)pc->data; 994*18c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 995*18c3aa7eSMark 996*18c3aa7eSMark PetscFunctionBegin; 997*18c3aa7eSMark ierr = PetscStrcpy(pc_gamg->esteig_type,t);CHKERRQ(ierr); 998*18c3aa7eSMark PetscFunctionReturn(0); 999*18c3aa7eSMark } 1000*18c3aa7eSMark 1001*18c3aa7eSMark /*@ 1002*18c3aa7eSMark PCGAMGSetEigenvalues - Set eigenvalues 1003*18c3aa7eSMark 1004*18c3aa7eSMark Collective on PC 1005*18c3aa7eSMark 1006*18c3aa7eSMark Input Parameters: 1007*18c3aa7eSMark + pc - the preconditioner context 1008*18c3aa7eSMark - emax - max eigenvalue 1009*18c3aa7eSMark . emin - min eigenvalue 1010*18c3aa7eSMark 1011*18c3aa7eSMark Options Database Key: 1012*18c3aa7eSMark . -pc_gamg_eigenvalues 1013*18c3aa7eSMark 1014*18c3aa7eSMark Level: intermediate 1015*18c3aa7eSMark 1016*18c3aa7eSMark .seealso: PCGAMGSetEstEigKSPMaxIt(), PCGAMGSetUseSAEstEig(), PCGAMGSetEstEigKSPType() 1017*18c3aa7eSMark @*/ 1018*18c3aa7eSMark PetscErrorCode PCGAMGSetEigenvalues(PC pc, PetscReal emax,PetscReal emin) 1019*18c3aa7eSMark { 1020*18c3aa7eSMark PetscErrorCode ierr; 1021*18c3aa7eSMark 1022*18c3aa7eSMark PetscFunctionBegin; 1023*18c3aa7eSMark PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1024*18c3aa7eSMark ierr = PetscTryMethod(pc,"PCGAMGSetEigenvalues_C",(PC,PetscReal,PetscReal),(pc,emax,emin));CHKERRQ(ierr); 1025*18c3aa7eSMark PetscFunctionReturn(0); 1026*18c3aa7eSMark } 1027*18c3aa7eSMark static PetscErrorCode PCGAMGSetEigenvalues_GAMG(PC pc,PetscReal emax,PetscReal emin) 1028*18c3aa7eSMark { 1029*18c3aa7eSMark PC_MG *mg = (PC_MG*)pc->data; 1030*18c3aa7eSMark PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1031*18c3aa7eSMark 1032*18c3aa7eSMark PetscFunctionBegin; 1033*18c3aa7eSMark if (emax <= emin) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Maximum eigenvalue must be larger than minimum: max %g min %g",(double)emax,(double)emin); 1034*18c3aa7eSMark if (emax*emin <= 0.0) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Both eigenvalues must be of the same sign: max %g min %g",(double)emax,(double)emin); 1035*18c3aa7eSMark pc_gamg->emax = emax; 1036*18c3aa7eSMark pc_gamg->emin = emin; 1037*18c3aa7eSMark 1038*18c3aa7eSMark PetscFunctionReturn(0); 1039*18c3aa7eSMark } 1040*18c3aa7eSMark 1041*18c3aa7eSMark /*@ 1042cab9ed1eSBarry Smith PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding algebraic multigrid preconditioner 1043dfd5c07aSMark F. Adams 1044dfd5c07aSMark F. Adams Collective on PC 1045dfd5c07aSMark F. Adams 1046dfd5c07aSMark F. Adams Input Parameters: 10471cc46a46SBarry Smith + pc - the preconditioner context 10481cc46a46SBarry Smith - n - PETSC_TRUE or PETSC_FALSE 1049dfd5c07aSMark F. Adams 1050dfd5c07aSMark F. Adams Options Database Key: 10511cc46a46SBarry Smith . -pc_gamg_reuse_interpolation <true,false> 1052dfd5c07aSMark F. Adams 1053dfd5c07aSMark F. Adams Level: intermediate 1054dfd5c07aSMark F. Adams 105595452b02SPatrick Sanan Notes: 105695452b02SPatrick Sanan this may negatively affect the convergence rate of the method on new matrices if the matrix entries change a great deal, but allows 1057cab9ed1eSBarry Smith rebuilding the preconditioner quicker. 1058cab9ed1eSBarry Smith 1059dfd5c07aSMark F. Adams .seealso: () 1060dfd5c07aSMark F. Adams @*/ 10611cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n) 1062dfd5c07aSMark F. Adams { 1063dfd5c07aSMark F. Adams PetscErrorCode ierr; 1064dfd5c07aSMark F. Adams 1065dfd5c07aSMark F. Adams PetscFunctionBegin; 1066dfd5c07aSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10671cc46a46SBarry Smith ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 1068dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1069dfd5c07aSMark F. Adams } 1070dfd5c07aSMark F. Adams 10711cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n) 1072dfd5c07aSMark F. Adams { 1073dfd5c07aSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1074dfd5c07aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1075dfd5c07aSMark F. Adams 1076dfd5c07aSMark F. Adams PetscFunctionBegin; 1077dfd5c07aSMark F. Adams pc_gamg->reuse_prol = n; 1078dfd5c07aSMark F. Adams PetscFunctionReturn(0); 1079dfd5c07aSMark F. Adams } 1080dfd5c07aSMark F. Adams 1081ffc955d6SMark F. Adams /*@ 1082cab9ed1eSBarry Smith PCGAMGASMSetUseAggs - Have the PCGAMG smoother on each level use the aggregates defined by the coarsening process as the subdomains for the additive Schwarz preconditioner. 1083ffc955d6SMark F. Adams 1084ffc955d6SMark F. Adams Collective on PC 1085ffc955d6SMark F. Adams 1086ffc955d6SMark F. Adams Input Parameters: 1087cab9ed1eSBarry Smith + pc - the preconditioner context 1088cab9ed1eSBarry Smith - flg - PETSC_TRUE to use aggregates, PETSC_FALSE to not 1089ffc955d6SMark F. Adams 1090ffc955d6SMark F. Adams Options Database Key: 1091cab9ed1eSBarry Smith . -pc_gamg_asm_use_agg 1092ffc955d6SMark F. Adams 1093ffc955d6SMark F. Adams Level: intermediate 1094ffc955d6SMark F. Adams 1095ffc955d6SMark F. Adams .seealso: () 1096ffc955d6SMark F. Adams @*/ 1097cab9ed1eSBarry Smith PetscErrorCode PCGAMGASMSetUseAggs(PC pc, PetscBool flg) 1098ffc955d6SMark F. Adams { 1099ffc955d6SMark F. Adams PetscErrorCode ierr; 1100ffc955d6SMark F. Adams 1101ffc955d6SMark F. Adams PetscFunctionBegin; 1102ffc955d6SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1103cab9ed1eSBarry Smith ierr = PetscTryMethod(pc,"PCGAMGASMSetUseAggs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 1104ffc955d6SMark F. Adams PetscFunctionReturn(0); 1105ffc955d6SMark F. Adams } 1106ffc955d6SMark F. Adams 1107cab9ed1eSBarry Smith static PetscErrorCode PCGAMGASMSetUseAggs_GAMG(PC pc, PetscBool flg) 1108ffc955d6SMark F. Adams { 1109ffc955d6SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1110ffc955d6SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1111ffc955d6SMark F. Adams 1112ffc955d6SMark F. Adams PetscFunctionBegin; 1113cab9ed1eSBarry Smith pc_gamg->use_aggs_in_asm = flg; 1114ffc955d6SMark F. Adams PetscFunctionReturn(0); 1115ffc955d6SMark F. Adams } 1116ffc955d6SMark F. Adams 1117171cca9aSMark Adams /*@ 1118cf8ae1d3SMark Adams PCGAMGSetUseParallelCoarseGridSolve - allow a parallel coarse grid solver 1119171cca9aSMark Adams 1120171cca9aSMark Adams Collective on PC 1121171cca9aSMark Adams 1122171cca9aSMark Adams Input Parameters: 1123171cca9aSMark Adams + pc - the preconditioner context 1124cf8ae1d3SMark Adams - flg - PETSC_TRUE to not force coarse grid onto one processor 1125171cca9aSMark Adams 1126171cca9aSMark Adams Options Database Key: 1127cf8ae1d3SMark Adams . -pc_gamg_use_parallel_coarse_grid_solver 1128171cca9aSMark Adams 1129171cca9aSMark Adams Level: intermediate 1130171cca9aSMark Adams 113139d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetCpuPinCoarseGrids() 1132171cca9aSMark Adams @*/ 1133171cca9aSMark Adams PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC pc, PetscBool flg) 1134171cca9aSMark Adams { 1135171cca9aSMark Adams PetscErrorCode ierr; 1136171cca9aSMark Adams 1137171cca9aSMark Adams PetscFunctionBegin; 1138171cca9aSMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1139171cca9aSMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetUseParallelCoarseGridSolve_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 1140171cca9aSMark Adams PetscFunctionReturn(0); 1141171cca9aSMark Adams } 1142171cca9aSMark Adams 1143171cca9aSMark Adams static PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve_GAMG(PC pc, PetscBool flg) 1144171cca9aSMark Adams { 1145171cca9aSMark Adams PC_MG *mg = (PC_MG*)pc->data; 1146171cca9aSMark Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1147171cca9aSMark Adams 1148171cca9aSMark Adams PetscFunctionBegin; 1149171cca9aSMark Adams pc_gamg->use_parallel_coarse_grid_solver = flg; 1150ffc955d6SMark F. Adams PetscFunctionReturn(0); 1151ffc955d6SMark F. Adams } 1152ffc955d6SMark F. Adams 11534ef23d27SMark F. Adams /*@ 1154ce7c7f2fSMark Adams PCGAMGSetCpuPinCoarseGrids - pin reduced grids to CPU 1155ce7c7f2fSMark Adams 1156ce7c7f2fSMark Adams Collective on PC 1157ce7c7f2fSMark Adams 1158ce7c7f2fSMark Adams Input Parameters: 1159ce7c7f2fSMark Adams + pc - the preconditioner context 1160ce7c7f2fSMark Adams - flg - PETSC_TRUE to pin coarse grids to CPU 1161ce7c7f2fSMark Adams 1162ce7c7f2fSMark Adams Options Database Key: 1163ce7c7f2fSMark Adams . -pc_gamg_cpu_pin_coarse_grids 1164ce7c7f2fSMark Adams 1165ce7c7f2fSMark Adams Level: intermediate 1166ce7c7f2fSMark Adams 116739d09545SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType(), PCGAMGSetUseParallelCoarseGridSolve() 1168ce7c7f2fSMark Adams @*/ 1169ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC pc, PetscBool flg) 1170ce7c7f2fSMark Adams { 1171ce7c7f2fSMark Adams PetscErrorCode ierr; 1172ce7c7f2fSMark Adams 1173ce7c7f2fSMark Adams PetscFunctionBegin; 1174ce7c7f2fSMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1175ce7c7f2fSMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetCpuPinCoarseGrids_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 1176ce7c7f2fSMark Adams PetscFunctionReturn(0); 1177ce7c7f2fSMark Adams } 1178ce7c7f2fSMark Adams 1179ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCpuPinCoarseGrids_GAMG(PC pc, PetscBool flg) 1180ce7c7f2fSMark Adams { 1181ce7c7f2fSMark Adams PC_MG *mg = (PC_MG*)pc->data; 1182ce7c7f2fSMark Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1183ce7c7f2fSMark Adams 1184ce7c7f2fSMark Adams PetscFunctionBegin; 1185ce7c7f2fSMark Adams pc_gamg->cpu_pin_coarse_grids = flg; 1186ce7c7f2fSMark Adams PetscFunctionReturn(0); 1187ce7c7f2fSMark Adams } 1188ce7c7f2fSMark Adams 1189ce7c7f2fSMark Adams /*@ 1190ce7c7f2fSMark Adams PCGAMGSetCoarseGridLayoutType - place reduce grids on processors with natural order (compact type) 1191ce7c7f2fSMark Adams 1192ce7c7f2fSMark Adams Collective on PC 1193ce7c7f2fSMark Adams 1194ce7c7f2fSMark Adams Input Parameters: 1195ce7c7f2fSMark Adams + pc - the preconditioner context 1196ce7c7f2fSMark Adams - flg - Layout type 1197ce7c7f2fSMark Adams 1198ce7c7f2fSMark Adams Options Database Key: 1199ce7c7f2fSMark Adams . -pc_gamg_coarse_grid_layout_type 1200ce7c7f2fSMark Adams 1201ce7c7f2fSMark Adams Level: intermediate 1202ce7c7f2fSMark Adams 120339d09545SMark Adams .seealso: PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetCpuPinCoarseGrids() 1204ce7c7f2fSMark Adams @*/ 1205ce7c7f2fSMark Adams PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC pc, PCGAMGLayoutType flg) 1206ce7c7f2fSMark Adams { 1207ce7c7f2fSMark Adams PetscErrorCode ierr; 1208ce7c7f2fSMark Adams 1209ce7c7f2fSMark Adams PetscFunctionBegin; 1210ce7c7f2fSMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1211ce7c7f2fSMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetCoarseGridLayoutType_C",(PC,PCGAMGLayoutType),(pc,flg));CHKERRQ(ierr); 1212ce7c7f2fSMark Adams PetscFunctionReturn(0); 1213ce7c7f2fSMark Adams } 1214ce7c7f2fSMark Adams 1215ce7c7f2fSMark Adams static PetscErrorCode PCGAMGSetCoarseGridLayoutType_GAMG(PC pc, PCGAMGLayoutType flg) 1216ce7c7f2fSMark Adams { 1217ce7c7f2fSMark Adams PC_MG *mg = (PC_MG*)pc->data; 1218ce7c7f2fSMark Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1219ce7c7f2fSMark Adams 1220ce7c7f2fSMark Adams PetscFunctionBegin; 1221ce7c7f2fSMark Adams pc_gamg->layout_type = flg; 1222ce7c7f2fSMark Adams PetscFunctionReturn(0); 1223ce7c7f2fSMark Adams } 1224ce7c7f2fSMark Adams 1225ce7c7f2fSMark Adams /*@ 12261cc46a46SBarry Smith PCGAMGSetNlevels - Sets the maximum number of levels PCGAMG will use 12274ef23d27SMark F. Adams 12284ef23d27SMark F. Adams Not collective on PC 12294ef23d27SMark F. Adams 12304ef23d27SMark F. Adams Input Parameters: 12311cc46a46SBarry Smith + pc - the preconditioner 12321cc46a46SBarry Smith - n - the maximum number of levels to use 12334ef23d27SMark F. Adams 12344ef23d27SMark F. Adams Options Database Key: 12354ef23d27SMark F. Adams . -pc_mg_levels 12364ef23d27SMark F. Adams 12374ef23d27SMark F. Adams Level: intermediate 12384ef23d27SMark F. Adams 12394ef23d27SMark F. Adams .seealso: () 12404ef23d27SMark F. Adams @*/ 12414ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 12424ef23d27SMark F. Adams { 12434ef23d27SMark F. Adams PetscErrorCode ierr; 12444ef23d27SMark F. Adams 12454ef23d27SMark F. Adams PetscFunctionBegin; 12464ef23d27SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12474ef23d27SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 12484ef23d27SMark F. Adams PetscFunctionReturn(0); 12494ef23d27SMark F. Adams } 12504ef23d27SMark F. Adams 12511e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 12524ef23d27SMark F. Adams { 12534ef23d27SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 12544ef23d27SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 12554ef23d27SMark F. Adams 12564ef23d27SMark F. Adams PetscFunctionBegin; 12579d5b6da9SMark F. Adams pc_gamg->Nlevels = n; 12584ef23d27SMark F. Adams PetscFunctionReturn(0); 12594ef23d27SMark F. Adams } 12604ef23d27SMark F. Adams 12613542efc5SMark F. Adams /*@ 12623542efc5SMark F. Adams PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 12633542efc5SMark F. Adams 12643542efc5SMark F. Adams Not collective on PC 12653542efc5SMark F. Adams 12663542efc5SMark F. Adams Input Parameters: 12671cc46a46SBarry Smith + pc - the preconditioner context 1268b001cb0fSBarry Smith - threshold - the threshold value, 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph 12693542efc5SMark F. Adams 12703542efc5SMark F. Adams Options Database Key: 12711cc46a46SBarry Smith . -pc_gamg_threshold <threshold> 12723542efc5SMark F. Adams 127395452b02SPatrick Sanan Notes: 1274af3c827dSMark Adams Increasing the threshold decreases the rate of coarsening. Conversely reducing the threshold increases the rate of coarsening (aggressive coarsening) and thereby reduces the complexity of the coarse grids, and generally results in slower solver converge rates. Reducing coarse grid complexity reduced the complexity of Galerkin coarse grid construction considerably. 1275af3c827dSMark Adams Before coarsening or aggregating the graph, GAMG removes small values from the graph with this threshold, and thus reducing the coupling in the graph and a different (perhaps better) coarser set of points. 1276cab9ed1eSBarry Smith 12773542efc5SMark F. Adams Level: intermediate 12783542efc5SMark F. Adams 1279af3c827dSMark Adams .seealso: PCGAMGFilterGraph(), PCGAMGSetSquareGraph() 12803542efc5SMark F. Adams @*/ 1281c1eae691SMark Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal v[], PetscInt n) 12823542efc5SMark F. Adams { 12833542efc5SMark F. Adams PetscErrorCode ierr; 12843542efc5SMark F. Adams 12853542efc5SMark F. Adams PetscFunctionBegin; 12863542efc5SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1287c1eae691SMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal[],PetscInt),(pc,v,n));CHKERRQ(ierr); 12883542efc5SMark F. Adams PetscFunctionReturn(0); 12893542efc5SMark F. Adams } 12903542efc5SMark F. Adams 1291c1eae691SMark Adams static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal v[], PetscInt n) 12923542efc5SMark F. Adams { 1293c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1294c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1295c1eae691SMark Adams PetscInt i; 1296c1eae691SMark Adams PetscFunctionBegin; 1297c1eae691SMark Adams for (i=0;i<n;i++) pc_gamg->threshold[i] = v[i]; 1298*18c3aa7eSMark do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS); 1299c1eae691SMark Adams PetscFunctionReturn(0); 1300c1eae691SMark Adams } 1301c1eae691SMark Adams 1302c1eae691SMark Adams /*@ 1303c1eae691SMark Adams PCGAMGSetThresholdScale - Relative threshold reduction at each level 1304c1eae691SMark Adams 1305c1eae691SMark Adams Not collective on PC 1306c1eae691SMark Adams 1307c1eae691SMark Adams Input Parameters: 1308c1eae691SMark Adams + pc - the preconditioner context 1309c1eae691SMark Adams - scale - the threshold value reduction, ussually < 1.0 1310c1eae691SMark Adams 1311c1eae691SMark Adams Options Database Key: 1312c1eae691SMark Adams . -pc_gamg_threshold_scale <v> 1313c1eae691SMark Adams 1314c1eae691SMark Adams Level: advanced 1315c1eae691SMark Adams 1316c1eae691SMark Adams .seealso: () 1317c1eae691SMark Adams @*/ 1318c1eae691SMark Adams PetscErrorCode PCGAMGSetThresholdScale(PC pc, PetscReal v) 1319c1eae691SMark Adams { 1320c1eae691SMark Adams PetscErrorCode ierr; 13213542efc5SMark F. Adams 13223542efc5SMark F. Adams PetscFunctionBegin; 1323c1eae691SMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1324c1eae691SMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetThresholdScale_C",(PC,PetscReal),(pc,v));CHKERRQ(ierr); 1325c1eae691SMark Adams PetscFunctionReturn(0); 1326c1eae691SMark Adams } 1327c1eae691SMark Adams 1328c1eae691SMark Adams static PetscErrorCode PCGAMGSetThresholdScale_GAMG(PC pc, PetscReal v) 1329c1eae691SMark Adams { 1330c1eae691SMark Adams PC_MG *mg = (PC_MG*)pc->data; 1331c1eae691SMark Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1332c1eae691SMark Adams PetscFunctionBegin; 1333c1eae691SMark Adams pc_gamg->threshold_scale = v; 13343542efc5SMark F. Adams PetscFunctionReturn(0); 13353542efc5SMark F. Adams } 13363542efc5SMark F. Adams 1337e20c40e8SBarry Smith /*@C 1338c60c7ad4SBarry Smith PCGAMGSetType - Set solution method 1339676e1743SMark F. Adams 1340676e1743SMark F. Adams Collective on PC 1341676e1743SMark F. Adams 1342676e1743SMark F. Adams Input Parameters: 1343c60c7ad4SBarry Smith + pc - the preconditioner context 1344c60c7ad4SBarry Smith - type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL 1345676e1743SMark F. Adams 1346676e1743SMark F. Adams Options Database Key: 1347cab9ed1eSBarry Smith . -pc_gamg_type <agg,geo,classical> - type of algebraic multigrid to apply 1348676e1743SMark F. Adams 1349676e1743SMark F. Adams Level: intermediate 1350676e1743SMark F. Adams 1351cab9ed1eSBarry Smith .seealso: PCGAMGGetType(), PCGAMG, PCGAMGType 1352676e1743SMark F. Adams @*/ 135319fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) 1354676e1743SMark F. Adams { 1355676e1743SMark F. Adams PetscErrorCode ierr; 1356676e1743SMark F. Adams 1357676e1743SMark F. Adams PetscFunctionBegin; 1358676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1359806fa848SBarry Smith ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr); 1360676e1743SMark F. Adams PetscFunctionReturn(0); 1361676e1743SMark F. Adams } 1362676e1743SMark F. Adams 1363e20c40e8SBarry Smith /*@C 1364c60c7ad4SBarry Smith PCGAMGGetType - Get solution method 1365c60c7ad4SBarry Smith 1366c60c7ad4SBarry Smith Collective on PC 1367c60c7ad4SBarry Smith 1368c60c7ad4SBarry Smith Input Parameter: 1369c60c7ad4SBarry Smith . pc - the preconditioner context 1370c60c7ad4SBarry Smith 1371c60c7ad4SBarry Smith Output Parameter: 1372c60c7ad4SBarry Smith . type - the type of algorithm used 1373c60c7ad4SBarry Smith 1374c60c7ad4SBarry Smith Level: intermediate 1375c60c7ad4SBarry Smith 13761c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType 1377c60c7ad4SBarry Smith @*/ 1378c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type) 1379c60c7ad4SBarry Smith { 1380c60c7ad4SBarry Smith PetscErrorCode ierr; 1381c60c7ad4SBarry Smith 1382c60c7ad4SBarry Smith PetscFunctionBegin; 1383c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1384c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr); 1385c60c7ad4SBarry Smith PetscFunctionReturn(0); 1386c60c7ad4SBarry Smith } 1387c60c7ad4SBarry Smith 1388c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type) 1389c60c7ad4SBarry Smith { 1390c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1391c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1392c60c7ad4SBarry Smith 1393c60c7ad4SBarry Smith PetscFunctionBegin; 1394c60c7ad4SBarry Smith *type = pc_gamg->type; 1395c60c7ad4SBarry Smith PetscFunctionReturn(0); 1396c60c7ad4SBarry Smith } 1397c60c7ad4SBarry Smith 13981e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) 1399676e1743SMark F. Adams { 14009d5b6da9SMark F. Adams PetscErrorCode ierr,(*r)(PC); 14011ab5ffc9SJed Brown PC_MG *mg = (PC_MG*)pc->data; 14021ab5ffc9SJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1403676e1743SMark F. Adams 1404676e1743SMark F. Adams PetscFunctionBegin; 1405c60c7ad4SBarry Smith pc_gamg->type = type; 14061c9cd337SJed Brown ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr); 14079d5b6da9SMark F. Adams if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 14081ab5ffc9SJed Brown if (pc_gamg->ops->destroy) { 14091ab5ffc9SJed Brown ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr); 14101ab5ffc9SJed Brown ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr); 1411e616c208SToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 14123ae0bb68SMark Adams /* cleaning up common data in pc_gamg - this should disapear someday */ 14133ae0bb68SMark Adams pc_gamg->data_cell_cols = 0; 14143ae0bb68SMark Adams pc_gamg->data_cell_rows = 0; 14153ae0bb68SMark Adams pc_gamg->orig_data_cell_cols = 0; 14163ae0bb68SMark Adams pc_gamg->orig_data_cell_rows = 0; 14173ae0bb68SMark Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 14183ae0bb68SMark Adams pc_gamg->data_sz = 0; 14191ab5ffc9SJed Brown } 14201ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr); 14211ab5ffc9SJed Brown ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr); 14229d5b6da9SMark F. Adams ierr = (*r)(pc);CHKERRQ(ierr); 1423676e1743SMark F. Adams PetscFunctionReturn(0); 1424676e1743SMark F. Adams } 1425676e1743SMark F. Adams 14269d766c59SMark Adams /* -------------------------------------------------------------------------- */ 14279d766c59SMark Adams /* 14289d766c59SMark Adams PCMGGetGridComplexity - compute coarse grid complexity of MG hierarchy 14299d766c59SMark Adams 14309d766c59SMark Adams Input Parameter: 14319d766c59SMark Adams . pc - the preconditioner context 14329d766c59SMark Adams 14339d766c59SMark Adams Output Parameter: 14349d766c59SMark Adams . gc - grid complexity = sum_i(nnz_i) / nnz_0 14359d766c59SMark Adams 14369d766c59SMark Adams Level: advanced 14379d766c59SMark Adams */ 14389d766c59SMark Adams static PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc) 14399d766c59SMark Adams { 14409d766c59SMark Adams PetscErrorCode ierr; 14419d766c59SMark Adams PC_MG *mg = (PC_MG*)pc->data; 14429d766c59SMark Adams PC_MG_Levels **mglevels = mg->levels; 14433966268fSBarry Smith PetscInt lev; 14443966268fSBarry Smith PetscLogDouble nnz0 = 0, sgc = 0; 14459d766c59SMark Adams MatInfo info; 14463966268fSBarry Smith 14479d766c59SMark Adams PetscFunctionBegin; 14489d766c59SMark Adams if (!mg->nlevels) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MG has no levels"); 14493966268fSBarry Smith for (lev=0; lev<mg->nlevels; lev++) { 145062a6e064SMark Adams Mat dB; 145162a6e064SMark Adams ierr = KSPGetOperators(mglevels[lev]->smoothd,NULL,&dB);CHKERRQ(ierr); 145262a6e064SMark Adams ierr = MatGetInfo(dB,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */ 14533966268fSBarry Smith sgc += info.nz_used; 14549d766c59SMark Adams if (lev==mg->nlevels-1) nnz0 = info.nz_used; 14559d766c59SMark Adams } 14563966268fSBarry Smith if (nnz0 > 0) *gc = (PetscReal)(sgc/nnz0); 14573966268fSBarry Smith else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number for grid points on finest level is not available"); 14589d766c59SMark Adams PetscFunctionReturn(0); 14599d766c59SMark Adams } 14609d766c59SMark Adams 14615adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer) 14625adeb434SBarry Smith { 1463c1eae691SMark Adams PetscErrorCode ierr,i; 14645adeb434SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 14655adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 146623b2d91dSMark Adams PetscReal gc=0; 14675adeb434SBarry Smith PetscFunctionBegin; 14685adeb434SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," GAMG specific options\n");CHKERRQ(ierr); 1469459726d8SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," Threshold for dropping small values in graph on each level =");CHKERRQ(ierr); 1470c1eae691SMark Adams for (i=0;i<pc_gamg->current_level;i++) { 1471c1eae691SMark Adams ierr = PetscViewerASCIIPrintf(viewer," %g",(double)pc_gamg->threshold[i]);CHKERRQ(ierr); 1472c1eae691SMark Adams } 1473459726d8SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1474459726d8SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," Threshold scaling factor for each level not specified = %g\n",(double)pc_gamg->threshold_scale);CHKERRQ(ierr); 1475cab9ed1eSBarry Smith if (pc_gamg->use_aggs_in_asm) { 1476cab9ed1eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using aggregates from coarsening process to define subdomains for PCASM\n");CHKERRQ(ierr); 1477cab9ed1eSBarry Smith } 1478171cca9aSMark Adams if (pc_gamg->use_parallel_coarse_grid_solver) { 1479171cca9aSMark Adams ierr = PetscViewerASCIIPrintf(viewer," Using parallel coarse grid solver (all coarse grid equations not put on one process)\n");CHKERRQ(ierr); 1480171cca9aSMark Adams } 1481ce7c7f2fSMark Adams #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 1482ce7c7f2fSMark Adams if (pc_gamg->cpu_pin_coarse_grids) { 1483ce7c7f2fSMark Adams /* ierr = PetscViewerASCIIPrintf(viewer," Pinning coarse grids to the CPU)\n");CHKERRQ(ierr); */ 1484ce7c7f2fSMark Adams } 1485ce7c7f2fSMark Adams #endif 1486ce7c7f2fSMark Adams /* if (pc_gamg->layout_type==PCGAMG_LAYOUT_COMPACT) { */ 1487ce7c7f2fSMark Adams /* ierr = PetscViewerASCIIPrintf(viewer," Put reduced grids on processes in natural order (ie, 0,1,2...)\n");CHKERRQ(ierr); */ 1488ce7c7f2fSMark Adams /* } else { */ 1489ce7c7f2fSMark Adams /* ierr = PetscViewerASCIIPrintf(viewer," Put reduced grids on whole machine (ie, 0,1*f,2*f...,np-f)\n");CHKERRQ(ierr); */ 1490ce7c7f2fSMark Adams /* } */ 14915adeb434SBarry Smith if (pc_gamg->ops->view) { 14925adeb434SBarry Smith ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr); 14935adeb434SBarry Smith } 14949d766c59SMark Adams ierr = PCMGGetGridComplexity(pc,&gc);CHKERRQ(ierr); 14959d766c59SMark Adams ierr = PetscViewerASCIIPrintf(viewer," Complexity: grid = %g\n",gc);CHKERRQ(ierr); 14965adeb434SBarry Smith PetscFunctionReturn(0); 14975adeb434SBarry Smith } 14985adeb434SBarry Smith 14994416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc) 15005b89ad90SMark F. Adams { 1501676e1743SMark F. Adams PetscErrorCode ierr; 1502676e1743SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1503676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1504*18c3aa7eSMark PetscBool flag,f2; 15053b4367a7SBarry Smith MPI_Comm comm; 1506*18c3aa7eSMark char prefix[256],tname[32]; 1507c1eae691SMark Adams PetscInt i,n; 150814a9496bSBarry Smith const char *pcpre; 1509d605200eSMark Adams static const char *LayoutTypes[] = {"compact","spread","PCGAMGLayoutType","PC_GAMG_LAYOUT",0}; 15105b89ad90SMark F. Adams PetscFunctionBegin; 15113b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1512e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr); 15131a1c1e04SBarry Smith ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr); 1514bd94a7aaSJed Brown if (flag) { 1515bd94a7aaSJed Brown ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr); 15161ab5ffc9SJed Brown } 1517*18c3aa7eSMark ierr = PetscOptionsFList("-pc_gamg_esteig_ksp_type","Krylov method for eigen estimator","PCGAMGSetEstEigKSPType",KSPList,pc_gamg->esteig_type,tname,sizeof(tname),&flag);CHKERRQ(ierr); 1518*18c3aa7eSMark if (flag) { 1519*18c3aa7eSMark ierr = PCGAMGSetEstEigKSPType(pc,tname);CHKERRQ(ierr); 1520*18c3aa7eSMark } 1521cab9ed1eSBarry Smith ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGSetRepartition",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr); 1522*18c3aa7eSMark f2 = PETSC_TRUE; 1523*18c3aa7eSMark ierr = PetscOptionsBool("-pc_gamg_use_sa_esteig","Use eigen estimate from Smoothed aggregation for smoother","PCGAMGSetUseSAEstEig",f2,&f2,&flag);CHKERRQ(ierr); 1524*18c3aa7eSMark if (flag) pc_gamg->use_sa_esteig = f2 ? 1 : 0; 15251cc46a46SBarry Smith ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr); 1526a303c832SJed Brown ierr = PetscOptionsBool("-pc_gamg_asm_use_agg","Use aggregation aggregates for ASM smoother","PCGAMGASMSetUseAggs",pc_gamg->use_aggs_in_asm,&pc_gamg->use_aggs_in_asm,NULL);CHKERRQ(ierr); 1527cf8ae1d3SMark Adams ierr = PetscOptionsBool("-pc_gamg_use_parallel_coarse_grid_solver","Use parallel coarse grid solver (otherwise put last grid on one process)","PCGAMGSetUseParallelCoarseGridSolve",pc_gamg->use_parallel_coarse_grid_solver,&pc_gamg->use_parallel_coarse_grid_solver,NULL);CHKERRQ(ierr); 1528ce7c7f2fSMark Adams ierr = PetscOptionsBool("-pc_gamg_cpu_pin_coarse_grids","Pin coarse grids to the CPU","PCGAMGSetCpuPinCoarseGrids",pc_gamg->cpu_pin_coarse_grids,&pc_gamg->cpu_pin_coarse_grids,NULL);CHKERRQ(ierr); 1529a0095786SMark ierr = PetscOptionsEnum("-pc_gamg_coarse_grid_layout_type","compact: place reduced grids on processes in natural order; spread: distribute to whole machine for more memory bandwidth","PCGAMGSetCoarseGridLayoutType",LayoutTypes,(PetscEnum)pc_gamg->layout_type,(PetscEnum*)&pc_gamg->layout_type,NULL);CHKERRQ(ierr); 153094ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_gamg_process_eq_limit","Limit (goal) on number of equations per process on coarse grids","PCGAMGSetProcEqLim",pc_gamg->min_eq_proc,&pc_gamg->min_eq_proc,NULL);CHKERRQ(ierr); 1531*18c3aa7eSMark ierr = PetscOptionsInt("-pc_gamg_esteig_ksp_max_it","Number of iterations of eigen estimator","PCGAMGSetEstEigKSPMaxIt",pc_gamg->esteig_max_it,&pc_gamg->esteig_max_it,NULL);CHKERRQ(ierr); 153294ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_gamg_coarse_eq_limit","Limit on number of equations for the coarse grid","PCGAMGSetCoarseEqLim",pc_gamg->coarse_eq_limit,&pc_gamg->coarse_eq_limit,NULL);CHKERRQ(ierr); 1533a303c832SJed Brown ierr = PetscOptionsReal("-pc_gamg_threshold_scale","Scaling of threshold for each level not specified","PCGAMGSetThresholdScale",pc_gamg->threshold_scale,&pc_gamg->threshold_scale,NULL);CHKERRQ(ierr); 1534*18c3aa7eSMark n = PETSC_MG_MAXLEVELS; 1535c1eae691SMark Adams ierr = PetscOptionsRealArray("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&n,&flag);CHKERRQ(ierr); 1536*18c3aa7eSMark if (!flag || n < PETSC_MG_MAXLEVELS) { 1537efd3c5ceSMark Adams if (!flag) n = 1; 1538c1eae691SMark Adams i = n; 1539*18c3aa7eSMark do {pc_gamg->threshold[i] = pc_gamg->threshold[i-1]*pc_gamg->threshold_scale;} while (++i<PETSC_MG_MAXLEVELS); 1540c1eae691SMark Adams } 154194ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr); 1542*18c3aa7eSMark { 1543*18c3aa7eSMark PetscReal eminmax[2] = {0., 0.}; 1544*18c3aa7eSMark n = 2; 1545*18c3aa7eSMark ierr = PetscOptionsRealArray("-pc_gamg_eigenvalues","extreme eigenvalues for smoothed aggregation","PCGAMGSetEigenvalues",eminmax,&n,&flag);CHKERRQ(ierr); 1546*18c3aa7eSMark if (flag) { 1547*18c3aa7eSMark if (n != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"-pc_gamg_eigenvalues: must specify 2 parameters, min and max eigenvalues"); 1548*18c3aa7eSMark ierr = PCGAMGSetEigenvalues(pc, eminmax[1], eminmax[0]);CHKERRQ(ierr); 1549*18c3aa7eSMark } 1550*18c3aa7eSMark } 1551b7cbab4eSMark Adams /* set options for subtype */ 1552e55864a3SBarry Smith if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);} 1553*18c3aa7eSMark 155414a9496bSBarry Smith ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr); 155514a9496bSBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr); 1556676e1743SMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 15575b89ad90SMark F. Adams PetscFunctionReturn(0); 15585b89ad90SMark F. Adams } 15595b89ad90SMark F. Adams 15605b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 15615b89ad90SMark F. Adams /*MC 15621cc46a46SBarry Smith PCGAMG - Geometric algebraic multigrid (AMG) preconditioner 15635b89ad90SMark F. Adams 1564280d9858SJed Brown Options Database Keys: 1565cab9ed1eSBarry Smith + -pc_gamg_type <type> - one of agg, geo, or classical 1566cab9ed1eSBarry Smith . -pc_gamg_repartition <true,default=false> - repartition the degrees of freedom accross the coarse grids as they are determined 1567cab9ed1eSBarry Smith . -pc_gamg_reuse_interpolation <true,default=false> - when rebuilding the algebraic multigrid preconditioner reuse the previously computed interpolations 1568cab9ed1eSBarry Smith . -pc_gamg_asm_use_agg <true,default=false> - use the aggregates from the coasening process to defined the subdomains on each level for the PCASM smoother 1569cab9ed1eSBarry Smith . -pc_gamg_process_eq_limit <limit, default=50> - GAMG will reduce the number of MPI processes used directly on the coarse grids so that there are around <limit> 1570cab9ed1eSBarry Smith equations on each process that has degrees of freedom 1571cab9ed1eSBarry Smith . -pc_gamg_coarse_eq_limit <limit, default=50> - Set maximum number of equations on coarsest grid to aim for. 15726008e27bSRichard Tran Mills . -pc_gamg_threshold[] <thresh,default=0> - Before aggregating the graph GAMG will remove small values from the graph on each level 1573c1eae691SMark Adams - -pc_gamg_threshold_scale <scale,default=1> - Scaling of threshold on each coarser grid if not specified 1574cab9ed1eSBarry Smith 1575cab9ed1eSBarry Smith Options Database Keys for default Aggregation: 1576cab9ed1eSBarry Smith + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation 1577cab9ed1eSBarry Smith . -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation 1578cab9ed1eSBarry Smith - -pc_gamg_square_graph <n,default=1> - number of levels to square the graph before aggregating it 1579cab9ed1eSBarry Smith 1580db9745e2SBarry Smith Multigrid options: 1581db9745e2SBarry Smith + -pc_mg_cycles <v> - v or w, see PCMGSetCycleType() 1582db9745e2SBarry Smith . -pc_mg_distinct_smoothup - configure the up and down (pre and post) smoothers separately, see PCMGSetDistinctSmoothUp() 1583db9745e2SBarry Smith . -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade 1584cab9ed1eSBarry Smith - -pc_mg_levels <levels> - Number of levels of multigrid to use. 15855b89ad90SMark F. Adams 15861cc46a46SBarry Smith 158795452b02SPatrick Sanan Notes: 158895452b02SPatrick Sanan In order to obtain good performance for PCGAMG for vector valued problems you must 1589db9745e2SBarry Smith Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point 1590db9745e2SBarry Smith Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator 1591db9745e2SBarry Smith See the Users Manual Chapter 4 for more details 15921cc46a46SBarry Smith 15935b89ad90SMark F. Adams Level: intermediate 1594280d9858SJed Brown 15951cc46a46SBarry Smith .seealso: PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(), 1596*18c3aa7eSMark PCGAMGSetCoarseEqLim(), PCGAMGSetRepartition(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGASMSetUseAggs(), PCGAMGSetUseParallelCoarseGridSolve(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseSAEstEig(), PCGAMGSetEstEigKSPMaxIt(), PCGAMGSetEstEigKSPType() 15975b89ad90SMark F. Adams M*/ 1598b2573a8aSBarry Smith 15998cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc) 16005b89ad90SMark F. Adams { 1601c1eae691SMark Adams PetscErrorCode ierr,i; 16025b89ad90SMark F. Adams PC_GAMG *pc_gamg; 16035b89ad90SMark F. Adams PC_MG *mg; 16045b89ad90SMark F. Adams 16055b89ad90SMark F. Adams PetscFunctionBegin; 16061c1aac46SBarry Smith /* register AMG type */ 16071c1aac46SBarry Smith ierr = PCGAMGInitializePackage();CHKERRQ(ierr); 16081c1aac46SBarry Smith 16095b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 16101c1aac46SBarry Smith ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); 16115b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr); 16125b89ad90SMark F. Adams 16135b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 1614b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr); 161569aca0b8SBarry Smith ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr); 16165b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 16175b89ad90SMark F. Adams mg->innerctx = pc_gamg; 16185b89ad90SMark F. Adams 1619b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr); 16201ab5ffc9SJed Brown 16219d5b6da9SMark F. Adams pc_gamg->setup_count = 0; 16229d5b6da9SMark F. Adams /* these should be in subctx but repartitioning needs simple arrays */ 16239d5b6da9SMark F. Adams pc_gamg->data_sz = 0; 16249d5b6da9SMark F. Adams pc_gamg->data = 0; 16255b89ad90SMark F. Adams 16269d5b6da9SMark F. Adams /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 16275b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 16285b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 16295b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 16305b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 16315adeb434SBarry Smith mg->view = PCView_GAMG; 16325b89ad90SMark F. Adams 163397d33e41SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr); 163497d33e41SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr); 1635bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr); 1636bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr); 1637cab9ed1eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartition_C",PCGAMGSetRepartition_GAMG);CHKERRQ(ierr); 1638*18c3aa7eSMark ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEstEigKSPType_C",PCGAMGSetEstEigKSPType_GAMG);CHKERRQ(ierr); 1639*18c3aa7eSMark ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEstEigKSPMaxIt_C",PCGAMGSetEstEigKSPMaxIt_GAMG);CHKERRQ(ierr); 1640*18c3aa7eSMark ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetEigenvalues_C",PCGAMGSetEigenvalues_GAMG);CHKERRQ(ierr); 1641*18c3aa7eSMark ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseSAEstEig_C",PCGAMGSetUseSAEstEig_GAMG);CHKERRQ(ierr); 16421cc46a46SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr); 1643cab9ed1eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGASMSetUseAggs_C",PCGAMGASMSetUseAggs_GAMG);CHKERRQ(ierr); 1644171cca9aSMark Adams ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseParallelCoarseGridSolve_C",PCGAMGSetUseParallelCoarseGridSolve_GAMG);CHKERRQ(ierr); 1645ce7c7f2fSMark Adams ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCpuPinCoarseGrids_C",PCGAMGSetCpuPinCoarseGrids_GAMG);CHKERRQ(ierr); 1646ce7c7f2fSMark Adams ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseGridLayoutType_C",PCGAMGSetCoarseGridLayoutType_GAMG);CHKERRQ(ierr); 1647bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr); 1648c1eae691SMark Adams ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThresholdScale_C",PCGAMGSetThresholdScale_GAMG);CHKERRQ(ierr); 1649bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr); 1650c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr); 1651bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr); 16529d5b6da9SMark F. Adams pc_gamg->repart = PETSC_FALSE; 1653d3042614SMark Adams pc_gamg->reuse_prol = PETSC_FALSE; 16540c3bc534SBarry Smith pc_gamg->use_aggs_in_asm = PETSC_FALSE; 1655171cca9aSMark Adams pc_gamg->use_parallel_coarse_grid_solver = PETSC_FALSE; 1656a0095786SMark pc_gamg->cpu_pin_coarse_grids = PETSC_FALSE; 1657a0095786SMark pc_gamg->layout_type = PCGAMG_LAYOUT_SPREAD; 1658038f3aa4SMark F. Adams pc_gamg->min_eq_proc = 50; 165925a145a7SMark Adams pc_gamg->coarse_eq_limit = 50; 1660*18c3aa7eSMark for (i=0;i<PETSC_MG_MAXLEVELS;i++) pc_gamg->threshold[i] = 0.; 1661c1eae691SMark Adams pc_gamg->threshold_scale = 1.; 1662*18c3aa7eSMark pc_gamg->Nlevels = PETSC_MG_MAXLEVELS; 16639ab59c8bSMark Adams pc_gamg->current_level = 0; /* don't need to init really */ 1664*18c3aa7eSMark ierr = PetscStrcpy(pc_gamg->esteig_type,KSPGMRES);CHKERRQ(ierr); 1665*18c3aa7eSMark pc_gamg->esteig_max_it = 10; 1666*18c3aa7eSMark pc_gamg->use_sa_esteig = -1; 1667*18c3aa7eSMark pc_gamg->emin = 0; 1668*18c3aa7eSMark pc_gamg->emax = 0; 1669*18c3aa7eSMark 1670c238b0ebSToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 16719d5b6da9SMark F. Adams 1672bd94a7aaSJed Brown /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */ 1673bd94a7aaSJed Brown ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr); 16745b89ad90SMark F. Adams PetscFunctionReturn(0); 16755b89ad90SMark F. Adams } 16763e3471ccSMark Adams 16773e3471ccSMark Adams /*@C 16783e3471ccSMark Adams PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called 16798a690491SBarry Smith from PCInitializePackage(). 16803e3471ccSMark Adams 16813e3471ccSMark Adams Level: developer 16823e3471ccSMark Adams 16833e3471ccSMark Adams .seealso: PetscInitialize() 16843e3471ccSMark Adams @*/ 16853e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void) 16863e3471ccSMark Adams { 16873e3471ccSMark Adams PetscErrorCode ierr; 16883e3471ccSMark Adams 16893e3471ccSMark Adams PetscFunctionBegin; 16903e3471ccSMark Adams if (PCGAMGPackageInitialized) PetscFunctionReturn(0); 16913e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_TRUE; 16923e3471ccSMark Adams ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr); 16933e3471ccSMark Adams ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr); 16948e6d0c30SPeter Brune ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr); 16953e3471ccSMark Adams ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr); 1696c1c463dbSMark Adams 1697c1c463dbSMark Adams /* general events */ 1698fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr); 1699fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr); 1700fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr); 1701fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr); 1702c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr); 1703c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr); 1704fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr); 1705c1c463dbSMark Adams 17065b89ad90SMark F. Adams #if defined PETSC_GAMG_USE_LOG 17075b89ad90SMark F. Adams ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr); 17085b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr); 17095b89ad90SMark F. Adams /* PetscLogEventRegister(" G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */ 17105b89ad90SMark F. Adams /* PetscLogEventRegister(" G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */ 17115b89ad90SMark F. Adams /* PetscLogEventRegister(" G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */ 17125b89ad90SMark F. Adams ierr = PetscLogEventRegister(" MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr); 17135b89ad90SMark F. Adams ierr = PetscLogEventRegister(" geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr); 17145b89ad90SMark F. Adams ierr = PetscLogEventRegister(" geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr); 1715bb235841SBarry Smith ierr = PetscLogEventRegister(" search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr); 17165b89ad90SMark F. Adams ierr = PetscLogEventRegister(" SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr); 17175b89ad90SMark F. Adams ierr = PetscLogEventRegister(" SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr); 17185b89ad90SMark F. Adams ierr = PetscLogEventRegister(" SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr); 17195b89ad90SMark F. Adams ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr); 17205b89ad90SMark F. Adams ierr = PetscLogEventRegister(" repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr); 17215b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr); 17225b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr); 17235b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr); 17245b89ad90SMark F. Adams 17255b89ad90SMark F. Adams /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */ 17265b89ad90SMark F. Adams /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */ 17275b89ad90SMark F. Adams /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */ 17285b89ad90SMark F. Adams /* create timer stages */ 17295b89ad90SMark F. Adams #if defined GAMG_STAGES 17305b89ad90SMark F. Adams { 17315b89ad90SMark F. Adams char str[32]; 17325b89ad90SMark F. Adams PetscInt lidx; 17335b89ad90SMark F. Adams sprintf(str,"MG Level %d (finest)",0); 17345b89ad90SMark F. Adams ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr); 17355b89ad90SMark F. Adams for (lidx=1; lidx<9; lidx++) { 17365b89ad90SMark F. Adams sprintf(str,"MG Level %d",lidx); 17375b89ad90SMark F. Adams ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr); 17385b89ad90SMark F. Adams } 17395b89ad90SMark F. Adams } 17405b89ad90SMark F. Adams #endif 17415b89ad90SMark F. Adams #endif 17423e3471ccSMark Adams PetscFunctionReturn(0); 17433e3471ccSMark Adams } 17443e3471ccSMark Adams 17453e3471ccSMark Adams /*@C 17461c1aac46SBarry Smith PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is 17471c1aac46SBarry Smith called from PetscFinalize() automatically. 17483e3471ccSMark Adams 17493e3471ccSMark Adams Level: developer 17503e3471ccSMark Adams 17513e3471ccSMark Adams .seealso: PetscFinalize() 17523e3471ccSMark Adams @*/ 17533e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void) 17543e3471ccSMark Adams { 17553e3471ccSMark Adams PetscErrorCode ierr; 17563e3471ccSMark Adams 17573e3471ccSMark Adams PetscFunctionBegin; 17583e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_FALSE; 17593e3471ccSMark Adams ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr); 17603e3471ccSMark Adams PetscFunctionReturn(0); 17613e3471ccSMark Adams } 1762a36cf38bSToby Isaac 1763a36cf38bSToby Isaac /*@C 1764a36cf38bSToby Isaac PCGAMGRegister - Register a PCGAMG implementation. 1765a36cf38bSToby Isaac 1766a36cf38bSToby Isaac Input Parameters: 1767a36cf38bSToby Isaac + type - string that will be used as the name of the GAMG type. 1768a36cf38bSToby Isaac - create - function for creating the gamg context. 1769a36cf38bSToby Isaac 1770a36cf38bSToby Isaac Level: advanced 1771a36cf38bSToby Isaac 17721c1aac46SBarry Smith .seealso: PCGAMGType, PCGAMG, PCGAMGSetType() 1773a36cf38bSToby Isaac @*/ 1774a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC)) 1775a36cf38bSToby Isaac { 1776a36cf38bSToby Isaac PetscErrorCode ierr; 1777a36cf38bSToby Isaac 1778a36cf38bSToby Isaac PetscFunctionBegin; 1779a36cf38bSToby Isaac ierr = PCGAMGInitializePackage();CHKERRQ(ierr); 1780a36cf38bSToby Isaac ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr); 1781a36cf38bSToby Isaac PetscFunctionReturn(0); 1782a36cf38bSToby Isaac } 1783a36cf38bSToby Isaac 1784