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*/ 6af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 75b42dca8SJed Brown #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /* Hack to access same_local_solves */ 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 23b4fbaa2aSMark F. Adams #define GAMG_MAXLEVELS 30 24b4fbaa2aSMark F. Adams 25b8fd24d8SMark F. Adams /* #define GAMG_STAGES */ 260cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 27b4fbaa2aSMark F. Adams static PetscLogStage gamg_stages[GAMG_MAXLEVELS]; 28b4fbaa2aSMark F. Adams #endif 29f96513f1SMatthew G Knepley 30140e18c1SBarry Smith static PetscFunctionList GAMGList = 0; 313e3471ccSMark Adams static PetscBool PCGAMGPackageInitialized; 329d5b6da9SMark F. Adams 33d3d6bff4SMark F. Adams /* ----------------------------------------------------------------------------- */ 34d3d6bff4SMark F. Adams #undef __FUNCT__ 35d3d6bff4SMark F. Adams #define __FUNCT__ "PCReset_GAMG" 36d3d6bff4SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 37d3d6bff4SMark F. Adams { 38d3d6bff4SMark F. Adams PetscErrorCode ierr; 39d3d6bff4SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 40d3d6bff4SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 41d3d6bff4SMark F. Adams 42d3d6bff4SMark F. Adams PetscFunctionBegin; 4362294041SBarry Smith if (pc_gamg->data) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"This should not happen, cleaned up in SetUp\n"); 441c1aac46SBarry Smith pc_gamg->data_sz = 0; 45878e152fSMark F. Adams ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr); 46a2f3521dSMark F. Adams PetscFunctionReturn(0); 47a2f3521dSMark F. Adams } 48a2f3521dSMark F. Adams 495b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 505b89ad90SMark F. Adams /* 51c238b0ebSToby Isaac PCGAMGCreateLevel_GAMG: create coarse op with RAP. repartition and/or reduce number 52a147abb0SMark F. Adams of active processors. 535b89ad90SMark F. Adams 545b89ad90SMark F. Adams Input Parameter: 55a2f3521dSMark F. Adams . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 56a2f3521dSMark F. Adams 'pc_gamg->data_sz' are changed via repartitioning/reduction. 579d5b6da9SMark F. Adams . Amat_fine - matrix on this fine (k) level 58c5bfad50SMark F. Adams . cr_bs - coarse block size 593530afc2SMark F. Adams In/Output Parameter: 60a2f3521dSMark F. Adams . a_P_inout - prolongation operator to the next level (k-->k-1) 61afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 6211e60469SMark F. Adams Output Parameter: 633530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 645b89ad90SMark F. Adams */ 655cb416c2SMark F. Adams 665b89ad90SMark F. Adams #undef __FUNCT__ 67c238b0ebSToby Isaac #define __FUNCT__ "PCGAMGCreateLevel_GAMG" 6862294041SBarry Smith 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) 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 853ae0bb68SMark Adams /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/ 860298fd71SBarry Smith ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr); 873ae0bb68SMark Adams if (pc_gamg->data_cell_rows>0) { 883ae0bb68SMark Adams ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows; 8973911c69SBarry Smith } else { 903ae0bb68SMark Adams PetscInt bs; 913ae0bb68SMark Adams ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr); 923ae0bb68SMark Adams ncrs = ncrs_eq/bs; 933ae0bb68SMark Adams } 94a2f3521dSMark F. Adams 95c5df96a5SBarry Smith /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */ 96a2f3521dSMark F. Adams { 97472110cdSMark F. Adams PetscInt ncrs_eq_glob; 980298fd71SBarry Smith ierr = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr); 99a90e85d9SMark Adams new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 1007f66b68fSMark Adams if (!new_size) new_size = 1; /* not likely, posible? */ 101c5df96a5SBarry Smith else if (new_size >= nactive) new_size = nactive; /* no change, rare */ 102a2f3521dSMark F. Adams } 103f852f58cSMark F. Adams 1043cb8563fSToby Isaac if (Pcolumnperm) *Pcolumnperm = NULL; 1053cb8563fSToby Isaac 106a90e85d9SMark Adams if (!pc_gamg->repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 1072fa5cd67SKarl Rupp else { 1083ae0bb68SMark Adams PetscInt *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old; 109885364a3SMark Adams IS is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices; 110e33ef3b1SMark F. Adams 11171959b99SBarry Smith nloc_old = ncrs_eq/cr_bs; 11271959b99SBarry 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); 1130cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 1140cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 115b4fbaa2aSMark F. Adams #endif 116a2f3521dSMark F. Adams /* make 'is_eq_newproc' */ 117785e854fSJed Brown ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr); 118a90e85d9SMark Adams if (pc_gamg->repart) { 119a2f3521dSMark F. Adams /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */ 1205a9b9e01SMark F. Adams Mat adj; 1215a9b9e01SMark F. Adams 122302440fdSBarry Smith ierr = PetscInfo3(pc,"Repartition: size (active): %D --> %D, neq = %D\n",*a_nactive_proc,new_size,ncrs_eq);CHKERRQ(ierr); 1235a9b9e01SMark F. Adams 124a2f3521dSMark F. Adams /* get 'adj' */ 125c5bfad50SMark F. Adams if (cr_bs == 1) { 126038e3b61SMark F. Adams ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 127806fa848SBarry Smith } else { 128a2f3521dSMark F. Adams /* make a scalar matrix to partition (no Stokes here) */ 129eb07cef2SMark F. Adams Mat tMat; 130a2f3521dSMark F. Adams PetscInt Istart_crs,Iend_crs,ncols,jj,Ii; 131b4fbaa2aSMark F. Adams const PetscScalar *vals; 132b4fbaa2aSMark F. Adams const PetscInt *idx; 133a2f3521dSMark F. Adams PetscInt *d_nnz, *o_nnz, M, N; 1349057884aSMark F. Adams static PetscInt llev = 0; 135d9558ea9SBarry Smith MatType mtype; 136b4fbaa2aSMark F. Adams 137e632b94dSBarry Smith ierr = PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);CHKERRQ(ierr); 138a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr); 139a2f3521dSMark F. Adams ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr); 140c5bfad50SMark F. Adams for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) { 14158471d46SMark F. Adams ierr = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr); 142c5bfad50SMark F. Adams d_nnz[jj] = ncols/cr_bs; 143c5bfad50SMark F. Adams o_nnz[jj] = ncols/cr_bs; 14458471d46SMark F. Adams ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr); 1453ae0bb68SMark Adams if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs; 1463ae0bb68SMark Adams if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs; 14758471d46SMark F. Adams } 1486876a03eSMark F. Adams 149d9558ea9SBarry Smith ierr = MatGetType(Amat_fine,&mtype);CHKERRQ(ierr); 1503b4367a7SBarry Smith ierr = MatCreate(comm, &tMat);CHKERRQ(ierr); 1513ae0bb68SMark Adams ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 152d9558ea9SBarry Smith ierr = MatSetType(tMat,mtype);CHKERRQ(ierr); 153a2f3521dSMark F. Adams ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr); 154a2f3521dSMark F. Adams ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 155e632b94dSBarry Smith ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr); 156eb07cef2SMark F. Adams 157a2f3521dSMark F. Adams for (ii = Istart_crs; ii < Iend_crs; ii++) { 158c5bfad50SMark F. Adams PetscInt dest_row = ii/cr_bs; 15922063be5SMark F. Adams ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 160eb07cef2SMark F. Adams for (jj = 0; jj < ncols; jj++) { 161c5bfad50SMark F. Adams PetscInt dest_col = idx[jj]/cr_bs; 162eb07cef2SMark F. Adams PetscScalar v = 1.0; 163eb07cef2SMark F. Adams ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr); 164eb07cef2SMark F. Adams } 16522063be5SMark F. Adams ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 166eb07cef2SMark F. Adams } 167eb07cef2SMark F. Adams ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168eb07cef2SMark F. Adams ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 169eb07cef2SMark F. Adams 170b4fbaa2aSMark F. Adams if (llev++ == -1) { 171b4fbaa2aSMark F. Adams PetscViewer viewer; char fname[32]; 1728caf3d72SBarry Smith ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr); 1733b4367a7SBarry Smith PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer); 174b4fbaa2aSMark F. Adams ierr = MatView(tMat, viewer);CHKERRQ(ierr); 1753bf036e2SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 176b4fbaa2aSMark F. Adams } 177eb07cef2SMark F. Adams ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 178eb07cef2SMark F. Adams ierr = MatDestroy(&tMat);CHKERRQ(ierr); 179a2f3521dSMark F. Adams } /* create 'adj' */ 180f150b916SMark F. Adams 181a2f3521dSMark F. Adams { /* partition: get newproc_idx */ 1825a9b9e01SMark F. Adams char prefix[256]; 1835a9b9e01SMark F. Adams const char *pcpre; 184b4fbaa2aSMark F. Adams const PetscInt *is_idx; 185b4fbaa2aSMark F. Adams MatPartitioning mpart; 186a4b7d37bSMark F. Adams IS proc_is; 187a2f3521dSMark F. Adams PetscInt targetPE; 1882f03bc48SMark F. Adams 1893b4367a7SBarry Smith ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr); 1905ef31b24SMark F. Adams ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr); 1919d5b6da9SMark F. Adams ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr); 1928caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr); 19359a0be82SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 19411e60469SMark F. Adams ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 195c5df96a5SBarry Smith ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr); 196a4b7d37bSMark F. Adams ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr); 19711e60469SMark F. Adams ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 1985a9b9e01SMark F. Adams 1995ef31b24SMark F. Adams /* collect IS info */ 200785e854fSJed Brown ierr = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr); 201a4b7d37bSMark F. Adams ierr = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr); 202a2f3521dSMark F. Adams targetPE = 1; /* bring to "front" of machine */ 203c5df96a5SBarry Smith /*targetPE = size/new_size;*/ /* spread partitioning across machine */ 204a2f3521dSMark F. Adams for (kk = jj = 0 ; kk < nloc_old ; kk++) { 205c5bfad50SMark F. Adams for (ii = 0 ; ii < cr_bs ; ii++, jj++) { 206a2f3521dSMark F. Adams newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */ 207eb07cef2SMark F. Adams } 2085ef31b24SMark F. Adams } 209a4b7d37bSMark F. Adams ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr); 210a4b7d37bSMark F. Adams ierr = ISDestroy(&proc_is);CHKERRQ(ierr); 2115ef31b24SMark F. Adams } 2125ef31b24SMark F. Adams ierr = MatDestroy(&adj);CHKERRQ(ierr); 2135a9b9e01SMark F. Adams 2143b4367a7SBarry Smith ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr); 2158263b398SMark F. Adams ierr = PetscFree(newproc_idx);CHKERRQ(ierr); 216806fa848SBarry Smith } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */ 217a2f3521dSMark F. Adams PetscInt rfactor,targetPE; 21862294041SBarry Smith 2195a9b9e01SMark F. Adams /* find factor */ 220c5df96a5SBarry Smith if (new_size == 1) rfactor = size; /* easy */ 2215a9b9e01SMark F. Adams else { 2225a9b9e01SMark F. Adams PetscReal best_fact = 0.; 2235a9b9e01SMark F. Adams jj = -1; 224c5df96a5SBarry Smith for (kk = 1 ; kk <= size ; kk++) { 2252a82295dSMark Adams if (!(size%kk)) { /* a candidate */ 226c5df96a5SBarry Smith PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size; 2275a9b9e01SMark F. Adams if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */ 2285a9b9e01SMark F. Adams if (fact > best_fact) { 2295a9b9e01SMark F. Adams best_fact = fact; jj = kk; 2305a9b9e01SMark F. Adams } 2315a9b9e01SMark F. Adams } 2325a9b9e01SMark F. Adams } 2335a9b9e01SMark F. Adams if (jj != -1) rfactor = jj; 234a2f3521dSMark F. Adams else rfactor = 1; /* does this happen .. a prime */ 2355a9b9e01SMark F. Adams } 236c5df96a5SBarry Smith new_size = size/rfactor; 2375a9b9e01SMark F. Adams 238c5df96a5SBarry Smith if (new_size==nactive) { 239a2f3521dSMark F. Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */ 2405a9b9e01SMark F. Adams ierr = PetscFree(counts);CHKERRQ(ierr); 241302440fdSBarry Smith ierr = PetscInfo2(pc,"Aggregate processors noop: new_size=%D, neq(loc)=%D\n",new_size,ncrs_eq);CHKERRQ(ierr); 242fe682e87SBarry Smith #if defined PETSC_GAMG_USE_LOG 243fe682e87SBarry Smith ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 244fe682e87SBarry Smith #endif 2455a9b9e01SMark F. Adams PetscFunctionReturn(0); 2465a9b9e01SMark F. Adams } 2475a9b9e01SMark F. Adams 248302440fdSBarry Smith ierr = PetscInfo1(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);CHKERRQ(ierr); 249c5df96a5SBarry Smith targetPE = rank/rfactor; 2503b4367a7SBarry Smith ierr = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr); 251a2f3521dSMark F. Adams } /* end simple 'is_eq_newproc' */ 252e33ef3b1SMark F. Adams 25311e60469SMark F. Adams /* 254a2f3521dSMark F. Adams Create an index set from the is_eq_newproc index set to indicate the mapping TO 25511e60469SMark F. Adams */ 256a2f3521dSMark F. Adams ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr); 2577700e67bSMark Adams is_eq_num_prim = is_eq_num; 25811e60469SMark F. Adams /* 259a2f3521dSMark F. Adams Determine how many equations/vertices are assigned to each processor 26011e60469SMark F. Adams */ 261c5df96a5SBarry Smith ierr = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr); 262c5df96a5SBarry Smith ncrs_eq_new = counts[rank]; 263a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_newproc);CHKERRQ(ierr); 2643ae0bb68SMark Adams ncrs_new = ncrs_eq_new/cr_bs; /* eqs */ 265a2f3521dSMark F. Adams 266a2f3521dSMark F. Adams ierr = PetscFree(counts);CHKERRQ(ierr); 2670cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 2680cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 269b4fbaa2aSMark F. Adams #endif 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 278a2f3521dSMark F. Adams /* move data (for primal equations only) */ 27922063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 2803b4367a7SBarry Smith ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr); 2813ae0bb68SMark Adams ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr); 282c0dedaeaSBarry Smith ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */ 28311e60469SMark F. Adams /* 2849d5b6da9SMark F. Adams There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 285c5bfad50SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'cr_bs'. 28611e60469SMark F. Adams */ 287854ce69bSBarry Smith ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr); 288a2f3521dSMark F. Adams ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 2893ae0bb68SMark Adams for (ii=0,jj=0; ii<ncrs; ii++) { 290c5bfad50SMark F. Adams PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */ 291a2f3521dSMark F. Adams for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk; 29211e60469SMark F. Adams } 293a2f3521dSMark F. Adams ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 2943ae0bb68SMark Adams ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr); 29592a756f0SMark F. Adams ierr = PetscFree(tidx);CHKERRQ(ierr); 29611e60469SMark F. Adams /* 29711e60469SMark F. Adams Create a vector to contain the original vertex information for each element 29811e60469SMark F. Adams */ 2993ae0bb68SMark Adams ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr); 3009d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols; jj++) { 3013ae0bb68SMark Adams const PetscInt stride0=ncrs*pc_gamg->data_cell_rows; 3023ae0bb68SMark Adams for (ii=0; ii<ncrs; ii++) { 3039d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows; kk++) { 304a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj; 305c8b0795cSMark F. Adams PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 306676e1743SMark F. Adams ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr); 307d3d6bff4SMark F. Adams } 308038e3b61SMark F. Adams } 309eb07cef2SMark F. Adams } 310eb07cef2SMark F. Adams ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr); 311eb07cef2SMark F. Adams ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr); 31211e60469SMark F. Adams /* 31311e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 31411e60469SMark F. Adams to the correct processor 31511e60469SMark F. Adams */ 3160298fd71SBarry Smith ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr); 31711e60469SMark F. Adams ierr = ISDestroy(&isscat);CHKERRQ(ierr); 31811e60469SMark F. Adams ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 31911e60469SMark F. Adams ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 32011e60469SMark F. Adams ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 32111e60469SMark F. Adams ierr = VecDestroy(&src_crd);CHKERRQ(ierr); 32211e60469SMark F. Adams /* 32311e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 32411e60469SMark F. Adams */ 325c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 326578f55a3SPeter Brune ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr); 3272fa5cd67SKarl Rupp 3283ae0bb68SMark Adams pc_gamg->data_sz = node_data_sz*ncrs_new; 3293ae0bb68SMark Adams strideNew = ncrs_new*ndata_rows; 3302fa5cd67SKarl Rupp 33111e60469SMark F. Adams ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr); 3329d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols; jj++) { 3333ae0bb68SMark Adams for (ii=0; ii<ncrs_new; ii++) { 3349d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows; kk++) { 335a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj; 336c8b0795cSMark F. Adams pc_gamg->data[ix] = PetscRealPart(array[jx]); 337d3d6bff4SMark F. Adams } 338038e3b61SMark F. Adams } 339038e3b61SMark F. Adams } 34011e60469SMark F. Adams ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr); 34111e60469SMark F. Adams ierr = VecDestroy(&dest_crd);CHKERRQ(ierr); 342885364a3SMark Adams } 343a2f3521dSMark F. Adams /* move A and P (columns) with new layout */ 3440cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3450cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 346ed3f9983SMark F. Adams #endif 347a2f3521dSMark F. Adams 34811e60469SMark F. Adams /* 34911e60469SMark F. Adams Invert for MatGetSubMatrix 35011e60469SMark F. Adams */ 351a2f3521dSMark F. Adams ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr); 352a2f3521dSMark F. Adams ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */ 353c5bfad50SMark F. Adams ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr); 354a2f3521dSMark F. Adams if (is_eq_num != is_eq_num_prim) { 355a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */ 356a2f3521dSMark F. Adams } 3573cb8563fSToby Isaac if (Pcolumnperm) { 3583cb8563fSToby Isaac ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr); 3593cb8563fSToby Isaac *Pcolumnperm = new_eq_indices; 3603cb8563fSToby Isaac } 361a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr); 3620cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3630cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 3640cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 365ed3f9983SMark F. Adams #endif 366a2f3521dSMark F. Adams /* 'a_Amat_crs' output */ 367a2f3521dSMark F. Adams { 368a2f3521dSMark F. Adams Mat mat; 369806fa848SBarry Smith ierr = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr); 370a2f3521dSMark F. Adams *a_Amat_crs = mat; 371a2f3521dSMark F. Adams } 372038e3b61SMark F. Adams ierr = MatDestroy(&Cmat);CHKERRQ(ierr); 373a2f3521dSMark F. Adams 3740cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3750cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 376ed3f9983SMark F. Adams #endif 37711e60469SMark F. Adams /* prolongator */ 37811e60469SMark F. Adams { 37911e60469SMark F. Adams IS findices; 380a2f3521dSMark F. Adams PetscInt Istart,Iend; 381a2f3521dSMark F. Adams Mat Pnew; 38262294041SBarry Smith 383a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr); 3840cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3850cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 386ed3f9983SMark F. Adams #endif 3873b4367a7SBarry Smith ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr); 388c5bfad50SMark F. Adams ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr); 389806fa848SBarry Smith ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr); 39011e60469SMark F. Adams ierr = ISDestroy(&findices);CHKERRQ(ierr); 391c5bfad50SMark F. Adams 3920cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3930cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 394ed3f9983SMark F. Adams #endif 3953530afc2SMark F. Adams ierr = MatDestroy(a_P_inout);CHKERRQ(ierr); 396a2f3521dSMark F. Adams 397a2f3521dSMark F. Adams /* output - repartitioned */ 398a2f3521dSMark F. Adams *a_P_inout = Pnew; 399e33ef3b1SMark F. Adams } 400a2f3521dSMark F. Adams ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr); 4015b89ad90SMark F. Adams 402c5df96a5SBarry Smith *a_nactive_proc = new_size; /* output */ 403a2f3521dSMark F. Adams } 4045b89ad90SMark F. Adams PetscFunctionReturn(0); 4055b89ad90SMark F. Adams } 4065b89ad90SMark F. Adams 4075b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 4085b89ad90SMark F. Adams /* 4095b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 4105b89ad90SMark F. Adams by setting data structures and options. 4115b89ad90SMark F. Adams 4125b89ad90SMark F. Adams Input Parameter: 4135b89ad90SMark F. Adams . pc - the preconditioner context 4145b89ad90SMark F. Adams 4155b89ad90SMark F. Adams */ 4165b89ad90SMark F. Adams #undef __FUNCT__ 4175b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG" 4189d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc) 4195b89ad90SMark F. Adams { 4205b89ad90SMark F. Adams PetscErrorCode ierr; 4219d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 4225b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 4232adcac29SMark F. Adams Mat Pmat = pc->pmat; 424a2f3521dSMark F. Adams PetscInt fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS]; 4253b4367a7SBarry Smith MPI_Comm comm; 426c5df96a5SBarry Smith PetscMPIInt rank,size,nactivepe; 427587fa25dSMark F. Adams Mat Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS]; 428e696ed0bSMark F. Adams IS *ASMLocalIDsArr[GAMG_MAXLEVELS]; 429a2f3521dSMark F. Adams PetscLogDouble nnz0=0.,nnztot=0.; 430569f4572SMark Adams MatInfo info; 4315ef31b24SMark F. Adams 4325b89ad90SMark F. Adams PetscFunctionBegin; 4333b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 4343b4367a7SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 4353b4367a7SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 436dfd5c07aSMark F. Adams 43784d3f75bSMark F. Adams if (pc_gamg->setup_count++ > 0) { 4381c1aac46SBarry Smith if ((PetscBool)(!pc_gamg->reuse_prol)) { 439878e152fSMark F. Adams /* reset everything */ 440878e152fSMark F. Adams ierr = PCReset_MG(pc);CHKERRQ(ierr); 441878e152fSMark F. Adams pc->setupcalled = 0; 442806fa848SBarry Smith } else { 44384d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 44403a628feSMark F. Adams /* just do Galerkin grids */ 44558471d46SMark F. Adams Mat B,dA,dB; 44658471d46SMark F. Adams 44771959b99SBarry Smith if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet"); 4489d5b6da9SMark F. Adams if (pc_gamg->Nlevels > 1) { 44958471d46SMark F. Adams /* currently only handle case where mat and pmat are the same on coarser levels */ 45023ee1639SBarry Smith ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 45158471d46SMark F. Adams /* (re)set to get dirty flag */ 45223ee1639SBarry Smith ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr); 45358471d46SMark F. Adams 4542fb0b348SMark F. Adams for (level=pc_gamg->Nlevels-2; level>=0; level--) { 45503a628feSMark F. Adams /* the first time through the matrix structure has changed from repartitioning */ 4560a97e771SToby Isaac if (pc_gamg->setup_count==2) { 45703a628feSMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 458084a8fe3SJed Brown ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr); 4592fa5cd67SKarl Rupp 46003a628feSMark F. Adams mglevels[level]->A = B; 461806fa848SBarry Smith } else { 46223ee1639SBarry Smith ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr); 46358471d46SMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 46403a628feSMark F. Adams } 46523ee1639SBarry Smith ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr); 46658471d46SMark F. Adams dB = B; 46758471d46SMark F. Adams } 4685f8cf99dSMark F. Adams } 469d5280255SMark F. Adams 470d5280255SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 47158471d46SMark F. Adams PetscFunctionReturn(0); 472eb07cef2SMark F. Adams } 473878e152fSMark F. Adams } 474f6536408SMark F. Adams 475878e152fSMark F. Adams if (!pc_gamg->data) { 476878e152fSMark F. Adams if (pc_gamg->orig_data) { 477878e152fSMark F. Adams ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 4780298fd71SBarry Smith ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr); 4792fa5cd67SKarl Rupp 480878e152fSMark F. Adams pc_gamg->data_sz = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols; 481878e152fSMark F. Adams pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 482878e152fSMark F. Adams pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 4832fa5cd67SKarl Rupp 484785e854fSJed Brown ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr); 485878e152fSMark F. Adams for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 486806fa848SBarry Smith } else { 4871ab5ffc9SJed Brown if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data"); 4887700e67bSMark Adams ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr); 4899d5b6da9SMark F. Adams } 490878e152fSMark F. Adams } 491878e152fSMark F. Adams 492878e152fSMark F. Adams /* cache original data for reuse */ 4931c1aac46SBarry Smith if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) { 494785e854fSJed Brown ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr); 495878e152fSMark F. Adams for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 496878e152fSMark F. Adams pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 497878e152fSMark F. Adams pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 498878e152fSMark F. Adams } 499038e3b61SMark F. Adams 500302f38e8SMark F. Adams /* get basic dims */ 501302f38e8SMark F. Adams ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 502a2f3521dSMark F. Adams ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr); 50384d3f75bSMark F. Adams 504569f4572SMark Adams ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */ 505569f4572SMark Adams nnz0 = info.nz_used; 506569f4572SMark Adams nnztot = info.nz_used; 50762294041SBarry 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); 508569f4572SMark Adams 509a2f3521dSMark F. Adams /* Get A_i and R_i */ 51062294041SBarry Smith for (level=0, Aarr[0]=Pmat, nactivepe = size; level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); level++) { 5119ab59c8bSMark Adams pc_gamg->current_level = level; 5125b89ad90SMark F. Adams level1 = level + 1; 5130cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 5140cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 515a2f3521dSMark F. Adams #if (defined GAMG_STAGES) 516a2f3521dSMark F. Adams ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr); 517b4fbaa2aSMark F. Adams #endif 518a2f3521dSMark F. Adams #endif 519c8b0795cSMark F. Adams { /* construct prolongator */ 520725b86d8SJed Brown Mat Gmat; 5210cbbd2e1SMark F. Adams PetscCoarsenData *agg_lists; 5227700e67bSMark Adams Mat Prol11; 523c8b0795cSMark F. Adams 5247700e67bSMark Adams ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr); 5251ab5ffc9SJed Brown ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr); 5267700e67bSMark Adams ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr); 527c8b0795cSMark F. Adams 528a2f3521dSMark F. Adams /* could have failed to create new level */ 529a2f3521dSMark F. Adams if (Prol11) { 5309d5b6da9SMark F. Adams /* get new block size of coarse matrices */ 5310298fd71SBarry Smith ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr); 532a2f3521dSMark F. Adams 533fd1112cbSBarry Smith if (pc_gamg->ops->optprolongator) { 534c8b0795cSMark F. Adams /* smooth */ 535fd1112cbSBarry Smith ierr = pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);CHKERRQ(ierr); 536c8b0795cSMark F. Adams } 537c8b0795cSMark F. Adams 5387700e67bSMark Adams Parr[level1] = Prol11; 5390298fd71SBarry Smith } else Parr[level1] = NULL; 540ffc955d6SMark F. Adams 541ffc955d6SMark F. Adams if (pc_gamg->use_aggs_in_gasm) { 5421b18a24aSMark Adams PetscInt bs; 5431b18a24aSMark Adams ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr); 544806fa848SBarry Smith ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr); 545ffc955d6SMark F. Adams } 546ffc955d6SMark F. Adams 547a2f3521dSMark F. Adams ierr = MatDestroy(&Gmat);CHKERRQ(ierr); 54841b27cdeSMark F. Adams ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr); 549a2f3521dSMark F. Adams } /* construct prolongator scope */ 5500cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 5510cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 552c8b0795cSMark F. Adams #endif 5537f66b68fSMark Adams if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 554c8b0795cSMark F. Adams if (!Parr[level1]) { 555569f4572SMark Adams ierr = PetscInfo1(pc,"Stop gridding, level %D\n",level);CHKERRQ(ierr); 55662294041SBarry Smith #if defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES 557a90e85d9SMark Adams ierr = PetscLogStagePop();CHKERRQ(ierr); 558a90e85d9SMark Adams #endif 559c8b0795cSMark F. Adams break; 560c8b0795cSMark F. Adams } 5610cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 5620cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 563b4fbaa2aSMark F. Adams #endif 564a2f3521dSMark F. Adams 5651c1aac46SBarry Smith ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs,&Parr[level1], &Aarr[level1], &nactivepe, NULL);CHKERRQ(ierr); 566a2f3521dSMark F. Adams 5670cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 5680cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 569b4fbaa2aSMark F. Adams #endif 570a2f3521dSMark F. Adams ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr); 571569f4572SMark Adams ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr); 572569f4572SMark Adams nnztot += info.nz_used; 5731d5b2942SMark 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); 574569f4572SMark Adams 5750cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 576b4fbaa2aSMark F. Adams ierr = PetscLogStagePop();CHKERRQ(ierr); 577b4fbaa2aSMark F. Adams #endif 578a90e85d9SMark Adams /* stop if one node or one proc -- could pull back for singular problems */ 5799ab59c8bSMark Adams if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2) ) { 5809ab59c8bSMark Adams ierr = PetscInfo2(pc,"HARD stop of coarsening on level %D. Grid too small: %D block nodes\n",level,M/bs);CHKERRQ(ierr); 581a90e85d9SMark Adams level++; 582a90e85d9SMark Adams break; 583a90e85d9SMark Adams } 584c8b0795cSMark F. Adams } /* levels */ 585c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 586c8b0795cSMark F. Adams 587569f4572SMark Adams ierr = PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);CHKERRQ(ierr); 5889d5b6da9SMark F. Adams pc_gamg->Nlevels = level + 1; 5895b89ad90SMark F. Adams fine_level = level; 5900298fd71SBarry Smith ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr); 5915b89ad90SMark F. Adams 59262294041SBarry Smith if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */ 593d5280255SMark F. Adams /* set default smoothers & set operators */ 59462294041SBarry Smith for (lidx = 1, level = pc_gamg->Nlevels-2; lidx <= fine_level; lidx++, level--) { 595ffc955d6SMark F. Adams KSP smoother; 596ffc955d6SMark F. Adams PC subpc; 597a2f3521dSMark F. Adams 5989d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 599f6536408SMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 600ffc955d6SMark F. Adams 601a2f3521dSMark F. Adams ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 602a2f3521dSMark F. Adams /* set ops */ 60323ee1639SBarry Smith ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr); 604a2f3521dSMark F. Adams ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr); 605a2f3521dSMark F. Adams 606a2f3521dSMark F. Adams /* set defaults */ 6076c9de887SHong Zhang ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr); 608a2f3521dSMark F. Adams 6091b18a24aSMark Adams /* set blocks for GASM smoother that uses the 'aggregates' */ 610ffc955d6SMark F. Adams if (pc_gamg->use_aggs_in_gasm) { 6112d3561bbSSatish Balay PetscInt sz; 6122d3561bbSSatish Balay IS *is; 613a2f3521dSMark F. Adams 6142d3561bbSSatish Balay sz = nASMBlocksArr[level]; 6152d3561bbSSatish Balay is = ASMLocalIDsArr[level]; 616ffc955d6SMark F. Adams ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr); 6171b18a24aSMark Adams ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr); 6187f66b68fSMark Adams if (!sz) { 619ffc955d6SMark F. Adams IS is; 620ffc955d6SMark F. Adams PetscInt my0,kk; 621ffc955d6SMark F. Adams ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr); 622ffc955d6SMark F. Adams ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr); 6230298fd71SBarry Smith ierr = PCGASMSetSubdomains(subpc, 1, &is, NULL);CHKERRQ(ierr); 624a94c3b12SMark F. Adams ierr = ISDestroy(&is);CHKERRQ(ierr); 625806fa848SBarry Smith } else { 626a94c3b12SMark F. Adams PetscInt kk; 6270298fd71SBarry Smith ierr = PCGASMSetSubdomains(subpc, sz, is, NULL);CHKERRQ(ierr); 628a94c3b12SMark F. Adams for (kk=0; kk<sz; kk++) { 629a94c3b12SMark F. Adams ierr = ISDestroy(&is[kk]);CHKERRQ(ierr); 630a94c3b12SMark F. Adams } 631ffc955d6SMark F. Adams ierr = PetscFree(is);CHKERRQ(ierr); 632ffc955d6SMark F. Adams } 6330298fd71SBarry Smith ASMLocalIDsArr[level] = NULL; 634ffc955d6SMark F. Adams nASMBlocksArr[level] = 0; 635ffc955d6SMark F. Adams ierr = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr); 636806fa848SBarry Smith } else { 637890ffe84SMark Adams ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr); 638ffc955d6SMark F. Adams } 639d5280255SMark F. Adams } 640d5280255SMark F. Adams { 641d5280255SMark F. Adams /* coarse grid */ 642d5280255SMark F. Adams KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 643d5280255SMark F. Adams Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0; 644d5280255SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 64523ee1639SBarry Smith ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr); 646d5280255SMark F. Adams ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 647d5280255SMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 648d5280255SMark F. Adams ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr); 649d5280255SMark F. Adams ierr = PCSetUp(subpc);CHKERRQ(ierr); 65071959b99SBarry Smith ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); 65171959b99SBarry Smith if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii); 652d5280255SMark F. Adams ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 653d5280255SMark F. Adams ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr); 6549dbfc187SHong Zhang ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 6552fb0b348SMark F. Adams ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr); 65608e36f19SMark Adams ierr = KSPSetType(k2[0], KSPPREONLY);CHKERRQ(ierr); 6575b42dca8SJed Brown /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in 6585b42dca8SJed Brown * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to 6595b42dca8SJed Brown * view every subdomain as though they were different. */ 6605b42dca8SJed Brown ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE; 661d5280255SMark F. Adams } 662d5280255SMark F. Adams 663d5280255SMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 664d5280255SMark F. Adams ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 665e55864a3SBarry Smith ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr); 666d5280255SMark F. Adams ierr = PetscOptionsEnd();CHKERRQ(ierr); 6671c1aac46SBarry Smith if (!mg->galerkin) SETERRQ(comm,PETSC_ERR_USER,"PCGAMG must use Galerkin for coarse operators."); 6681c1aac46SBarry Smith if (mg->galerkin == 1) mg->galerkin = 2; 669d5280255SMark F. Adams 670d5280255SMark F. Adams /* clean up */ 671d5280255SMark F. Adams for (level=1; level<pc_gamg->Nlevels; level++) { 672587fa25dSMark F. Adams ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr); 673587fa25dSMark F. Adams ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr); 6745b89ad90SMark F. Adams } 6750cbbd2e1SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 676806fa848SBarry Smith } else { 6775f8cf99dSMark F. Adams KSP smoother; 678302440fdSBarry Smith ierr = PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");CHKERRQ(ierr); 6799d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr); 68023ee1639SBarry Smith ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr); 6815f8cf99dSMark F. Adams ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr); 6829d5b6da9SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 6835f8cf99dSMark F. Adams } 6845b89ad90SMark F. Adams PetscFunctionReturn(0); 6855b89ad90SMark F. Adams } 6865b89ad90SMark F. Adams 687eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */ 6885b89ad90SMark F. Adams /* 6895b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 6905b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 6915b89ad90SMark F. Adams 6925b89ad90SMark F. Adams Input Parameter: 6935b89ad90SMark F. Adams . pc - the preconditioner context 6945b89ad90SMark F. Adams 6955b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 6965b89ad90SMark F. Adams */ 6975b89ad90SMark F. Adams #undef __FUNCT__ 6985b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG" 6995b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc) 7005b89ad90SMark F. Adams { 7015b89ad90SMark F. Adams PetscErrorCode ierr; 7025b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 7035b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 7045b89ad90SMark F. Adams 7055b89ad90SMark F. Adams PetscFunctionBegin; 7065b89ad90SMark F. Adams ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 7079b8ffb57SJed Brown if (pc_gamg->ops->destroy) { 7089b8ffb57SJed Brown ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr); 7099b8ffb57SJed Brown } 71014a9496bSBarry Smith ierr = PetscRandomDestroy(&pc_gamg->random);CHKERRQ(ierr); 7111ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr); 7121ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr); 7135b89ad90SMark F. Adams ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 7145b89ad90SMark F. Adams ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 7155b89ad90SMark F. Adams PetscFunctionReturn(0); 7165b89ad90SMark F. Adams } 7175b89ad90SMark F. Adams 718676e1743SMark F. Adams #undef __FUNCT__ 719676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim" 720676e1743SMark F. Adams /*@ 7211cc46a46SBarry Smith PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via processor reduction. 722676e1743SMark F. Adams 7231cc46a46SBarry Smith Logically Collective on PC 724676e1743SMark F. Adams 725676e1743SMark F. Adams Input Parameters: 7261cc46a46SBarry Smith + pc - the preconditioner context 7271cc46a46SBarry Smith - n - the number of equations 728676e1743SMark F. Adams 729676e1743SMark F. Adams 730676e1743SMark F. Adams Options Database Key: 7311cc46a46SBarry Smith . -pc_gamg_process_eq_limit <limit> 732676e1743SMark F. Adams 733676e1743SMark F. Adams Level: intermediate 734676e1743SMark F. Adams 7351c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 736676e1743SMark F. Adams 7371c1aac46SBarry Smith .seealso: PCGAMGSetCoarseEqLim() 738676e1743SMark F. Adams @*/ 739676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 740676e1743SMark F. Adams { 741676e1743SMark F. Adams PetscErrorCode ierr; 742676e1743SMark F. Adams 743676e1743SMark F. Adams PetscFunctionBegin; 744676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 745676e1743SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 746676e1743SMark F. Adams PetscFunctionReturn(0); 747676e1743SMark F. Adams } 748676e1743SMark F. Adams 749676e1743SMark F. Adams #undef __FUNCT__ 750676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG" 7511e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 752676e1743SMark F. Adams { 753c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 754c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 755676e1743SMark F. Adams 756676e1743SMark F. Adams PetscFunctionBegin; 7579d5b6da9SMark F. Adams if (n>0) pc_gamg->min_eq_proc = n; 758676e1743SMark F. Adams PetscFunctionReturn(0); 759676e1743SMark F. Adams } 760676e1743SMark F. Adams 761676e1743SMark F. Adams #undef __FUNCT__ 762389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim" 763389730f3SMark F. Adams /*@ 764389730f3SMark F. Adams PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids. 765389730f3SMark F. Adams 766389730f3SMark F. Adams Collective on PC 767389730f3SMark F. Adams 768389730f3SMark F. Adams Input Parameters: 7691cc46a46SBarry Smith + pc - the preconditioner context 7701cc46a46SBarry Smith - n - maximum number of equations to aim for 771389730f3SMark F. Adams 772389730f3SMark F. Adams Options Database Key: 7731cc46a46SBarry Smith . -pc_gamg_coarse_eq_limit <limit> 774389730f3SMark F. Adams 775389730f3SMark F. Adams Level: intermediate 776389730f3SMark F. Adams 7771c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 778389730f3SMark F. Adams 7791c1aac46SBarry Smith .seealso: PCGAMGSetProcEqLim() 780389730f3SMark F. Adams @*/ 781389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 782389730f3SMark F. Adams { 783389730f3SMark F. Adams PetscErrorCode ierr; 784389730f3SMark F. Adams 785389730f3SMark F. Adams PetscFunctionBegin; 786389730f3SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 787389730f3SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 788389730f3SMark F. Adams PetscFunctionReturn(0); 789389730f3SMark F. Adams } 790389730f3SMark F. Adams 791389730f3SMark F. Adams #undef __FUNCT__ 792389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG" 7931e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 794389730f3SMark F. Adams { 795389730f3SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 796389730f3SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 797389730f3SMark F. Adams 798389730f3SMark F. Adams PetscFunctionBegin; 7999d5b6da9SMark F. Adams if (n>0) pc_gamg->coarse_eq_limit = n; 800389730f3SMark F. Adams PetscFunctionReturn(0); 801389730f3SMark F. Adams } 802389730f3SMark F. Adams 803389730f3SMark F. Adams #undef __FUNCT__ 8048263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning" 805676e1743SMark F. Adams /*@ 8068263b398SMark F. Adams PCGAMGSetRepartitioning - Repartition the coarse grids 807676e1743SMark F. Adams 808676e1743SMark F. Adams Collective on PC 809676e1743SMark F. Adams 810676e1743SMark F. Adams Input Parameters: 8111cc46a46SBarry Smith + pc - the preconditioner context 8121cc46a46SBarry Smith - n - PETSC_TRUE or PETSC_FALSE 813676e1743SMark F. Adams 814676e1743SMark F. Adams Options Database Key: 8151cc46a46SBarry Smith . -pc_gamg_repartition <true,false> 816676e1743SMark F. Adams 817676e1743SMark F. Adams Level: intermediate 818676e1743SMark F. Adams 8191c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 820676e1743SMark F. Adams 821676e1743SMark F. Adams .seealso: () 822676e1743SMark F. Adams @*/ 8238263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n) 824676e1743SMark F. Adams { 825676e1743SMark F. Adams PetscErrorCode ierr; 826676e1743SMark F. Adams 827676e1743SMark F. Adams PetscFunctionBegin; 828676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 8298263b398SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 830676e1743SMark F. Adams PetscFunctionReturn(0); 831676e1743SMark F. Adams } 832676e1743SMark F. Adams 833676e1743SMark F. Adams #undef __FUNCT__ 8348263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG" 8351e6b0712SBarry Smith static PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n) 836676e1743SMark F. Adams { 837c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 838c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 839676e1743SMark F. Adams 840676e1743SMark F. Adams PetscFunctionBegin; 8419d5b6da9SMark F. Adams pc_gamg->repart = n; 842676e1743SMark F. Adams PetscFunctionReturn(0); 843676e1743SMark F. Adams } 844676e1743SMark F. Adams 845676e1743SMark F. Adams #undef __FUNCT__ 8461cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation" 847dfd5c07aSMark F. Adams /*@ 8481cc46a46SBarry Smith PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding preconditioner 849dfd5c07aSMark F. Adams 850dfd5c07aSMark F. Adams Collective on PC 851dfd5c07aSMark F. Adams 852dfd5c07aSMark F. Adams Input Parameters: 8531cc46a46SBarry Smith + pc - the preconditioner context 8541cc46a46SBarry Smith - n - PETSC_TRUE or PETSC_FALSE 855dfd5c07aSMark F. Adams 856dfd5c07aSMark F. Adams Options Database Key: 8571cc46a46SBarry Smith . -pc_gamg_reuse_interpolation <true,false> 858dfd5c07aSMark F. Adams 859dfd5c07aSMark F. Adams Level: intermediate 860dfd5c07aSMark F. Adams 8611c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 862dfd5c07aSMark F. Adams 863dfd5c07aSMark F. Adams .seealso: () 864dfd5c07aSMark F. Adams @*/ 8651cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n) 866dfd5c07aSMark F. Adams { 867dfd5c07aSMark F. Adams PetscErrorCode ierr; 868dfd5c07aSMark F. Adams 869dfd5c07aSMark F. Adams PetscFunctionBegin; 870dfd5c07aSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 8711cc46a46SBarry Smith ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 872dfd5c07aSMark F. Adams PetscFunctionReturn(0); 873dfd5c07aSMark F. Adams } 874dfd5c07aSMark F. Adams 875dfd5c07aSMark F. Adams #undef __FUNCT__ 8761cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation_GAMG" 8771cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n) 878dfd5c07aSMark F. Adams { 879dfd5c07aSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 880dfd5c07aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 881dfd5c07aSMark F. Adams 882dfd5c07aSMark F. Adams PetscFunctionBegin; 883dfd5c07aSMark F. Adams pc_gamg->reuse_prol = n; 884dfd5c07aSMark F. Adams PetscFunctionReturn(0); 885dfd5c07aSMark F. Adams } 886dfd5c07aSMark F. Adams 887dfd5c07aSMark F. Adams #undef __FUNCT__ 888ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs" 889ffc955d6SMark F. Adams /*@ 890ffc955d6SMark F. Adams PCGAMGSetUseASMAggs - 891ffc955d6SMark F. Adams 892ffc955d6SMark F. Adams Collective on PC 893ffc955d6SMark F. Adams 894ffc955d6SMark F. Adams Input Parameters: 895ffc955d6SMark F. Adams . pc - the preconditioner context 896ffc955d6SMark F. Adams 897ffc955d6SMark F. Adams 898ffc955d6SMark F. Adams Options Database Key: 899ffc955d6SMark F. Adams . -pc_gamg_use_agg_gasm 900ffc955d6SMark F. Adams 901ffc955d6SMark F. Adams Level: intermediate 902ffc955d6SMark F. Adams 9031c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 904ffc955d6SMark F. Adams 905ffc955d6SMark F. Adams .seealso: () 906ffc955d6SMark F. Adams @*/ 907ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n) 908ffc955d6SMark F. Adams { 909ffc955d6SMark F. Adams PetscErrorCode ierr; 910ffc955d6SMark F. Adams 911ffc955d6SMark F. Adams PetscFunctionBegin; 912ffc955d6SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 913ffc955d6SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 914ffc955d6SMark F. Adams PetscFunctionReturn(0); 915ffc955d6SMark F. Adams } 916ffc955d6SMark F. Adams 917ffc955d6SMark F. Adams #undef __FUNCT__ 918ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG" 9191e6b0712SBarry Smith static PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n) 920ffc955d6SMark F. Adams { 921ffc955d6SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 922ffc955d6SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 923ffc955d6SMark F. Adams 924ffc955d6SMark F. Adams PetscFunctionBegin; 925ffc955d6SMark F. Adams pc_gamg->use_aggs_in_gasm = n; 926ffc955d6SMark F. Adams PetscFunctionReturn(0); 927ffc955d6SMark F. Adams } 928ffc955d6SMark F. Adams 929ffc955d6SMark F. Adams #undef __FUNCT__ 9304ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels" 9314ef23d27SMark F. Adams /*@ 9321cc46a46SBarry Smith PCGAMGSetNlevels - Sets the maximum number of levels PCGAMG will use 9334ef23d27SMark F. Adams 9344ef23d27SMark F. Adams Not collective on PC 9354ef23d27SMark F. Adams 9364ef23d27SMark F. Adams Input Parameters: 9371cc46a46SBarry Smith + pc - the preconditioner 9381cc46a46SBarry Smith - n - the maximum number of levels to use 9394ef23d27SMark F. Adams 9404ef23d27SMark F. Adams Options Database Key: 9414ef23d27SMark F. Adams . -pc_mg_levels 9424ef23d27SMark F. Adams 9434ef23d27SMark F. Adams Level: intermediate 9444ef23d27SMark F. Adams 9451c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 9464ef23d27SMark F. Adams 9474ef23d27SMark F. Adams .seealso: () 9484ef23d27SMark F. Adams @*/ 9494ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 9504ef23d27SMark F. Adams { 9514ef23d27SMark F. Adams PetscErrorCode ierr; 9524ef23d27SMark F. Adams 9534ef23d27SMark F. Adams PetscFunctionBegin; 9544ef23d27SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9554ef23d27SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 9564ef23d27SMark F. Adams PetscFunctionReturn(0); 9574ef23d27SMark F. Adams } 9584ef23d27SMark F. Adams 9594ef23d27SMark F. Adams #undef __FUNCT__ 9604ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG" 9611e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 9624ef23d27SMark F. Adams { 9634ef23d27SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 9644ef23d27SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9654ef23d27SMark F. Adams 9664ef23d27SMark F. Adams PetscFunctionBegin; 9679d5b6da9SMark F. Adams pc_gamg->Nlevels = n; 9684ef23d27SMark F. Adams PetscFunctionReturn(0); 9694ef23d27SMark F. Adams } 9704ef23d27SMark F. Adams 9714ef23d27SMark F. Adams #undef __FUNCT__ 9723542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold" 9733542efc5SMark F. Adams /*@ 9743542efc5SMark F. Adams PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 9753542efc5SMark F. Adams 9763542efc5SMark F. Adams Not collective on PC 9773542efc5SMark F. Adams 9783542efc5SMark F. Adams Input Parameters: 9791cc46a46SBarry Smith + pc - the preconditioner context 980b001cb0fSBarry Smith - threshold - the threshold value, 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph 9813542efc5SMark F. Adams 9823542efc5SMark F. Adams Options Database Key: 9831cc46a46SBarry Smith . -pc_gamg_threshold <threshold> 9843542efc5SMark F. Adams 9853542efc5SMark F. Adams Level: intermediate 9863542efc5SMark F. Adams 9871c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 9883542efc5SMark F. Adams 9893542efc5SMark F. Adams .seealso: () 9903542efc5SMark F. Adams @*/ 9913542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n) 9923542efc5SMark F. Adams { 9933542efc5SMark F. Adams PetscErrorCode ierr; 9943542efc5SMark F. Adams 9953542efc5SMark F. Adams PetscFunctionBegin; 9963542efc5SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9973542efc5SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr); 9983542efc5SMark F. Adams PetscFunctionReturn(0); 9993542efc5SMark F. Adams } 10003542efc5SMark F. Adams 10013542efc5SMark F. Adams #undef __FUNCT__ 10023542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG" 10031e6b0712SBarry Smith static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n) 10043542efc5SMark F. Adams { 1005c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1006c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 10073542efc5SMark F. Adams 10083542efc5SMark F. Adams PetscFunctionBegin; 10099d5b6da9SMark F. Adams pc_gamg->threshold = n; 10103542efc5SMark F. Adams PetscFunctionReturn(0); 10113542efc5SMark F. Adams } 10123542efc5SMark F. Adams 10133542efc5SMark F. Adams #undef __FUNCT__ 10149d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType" 1015676e1743SMark F. Adams /*@ 1016c60c7ad4SBarry Smith PCGAMGSetType - Set solution method 1017676e1743SMark F. Adams 1018676e1743SMark F. Adams Collective on PC 1019676e1743SMark F. Adams 1020676e1743SMark F. Adams Input Parameters: 1021c60c7ad4SBarry Smith + pc - the preconditioner context 1022c60c7ad4SBarry Smith - type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL 1023676e1743SMark F. Adams 1024676e1743SMark F. Adams Options Database Key: 1025c60c7ad4SBarry Smith . -pc_gamg_type <agg,geo,classical> 1026676e1743SMark F. Adams 1027676e1743SMark F. Adams Level: intermediate 1028676e1743SMark F. Adams 10291c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 1030676e1743SMark F. Adams 1031c60c7ad4SBarry Smith .seealso: PCGAMGGetType(), PCGAMG 1032676e1743SMark F. Adams @*/ 103319fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) 1034676e1743SMark F. Adams { 1035676e1743SMark F. Adams PetscErrorCode ierr; 1036676e1743SMark F. Adams 1037676e1743SMark F. Adams PetscFunctionBegin; 1038676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1039806fa848SBarry Smith ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr); 1040676e1743SMark F. Adams PetscFunctionReturn(0); 1041676e1743SMark F. Adams } 1042676e1743SMark F. Adams 1043676e1743SMark F. Adams #undef __FUNCT__ 1044c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType" 1045c60c7ad4SBarry Smith /*@ 1046c60c7ad4SBarry Smith PCGAMGGetType - Get solution method 1047c60c7ad4SBarry Smith 1048c60c7ad4SBarry Smith Collective on PC 1049c60c7ad4SBarry Smith 1050c60c7ad4SBarry Smith Input Parameter: 1051c60c7ad4SBarry Smith . pc - the preconditioner context 1052c60c7ad4SBarry Smith 1053c60c7ad4SBarry Smith Output Parameter: 1054c60c7ad4SBarry Smith . type - the type of algorithm used 1055c60c7ad4SBarry Smith 1056c60c7ad4SBarry Smith Level: intermediate 1057c60c7ad4SBarry Smith 10581c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 1059c60c7ad4SBarry Smith 10601c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType 1061c60c7ad4SBarry Smith @*/ 1062c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type) 1063c60c7ad4SBarry Smith { 1064c60c7ad4SBarry Smith PetscErrorCode ierr; 1065c60c7ad4SBarry Smith 1066c60c7ad4SBarry Smith PetscFunctionBegin; 1067c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1068c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr); 1069c60c7ad4SBarry Smith PetscFunctionReturn(0); 1070c60c7ad4SBarry Smith } 1071c60c7ad4SBarry Smith 1072c60c7ad4SBarry Smith #undef __FUNCT__ 1073c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType_GAMG" 1074c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type) 1075c60c7ad4SBarry Smith { 1076c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1077c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1078c60c7ad4SBarry Smith 1079c60c7ad4SBarry Smith PetscFunctionBegin; 1080c60c7ad4SBarry Smith *type = pc_gamg->type; 1081c60c7ad4SBarry Smith PetscFunctionReturn(0); 1082c60c7ad4SBarry Smith } 1083c60c7ad4SBarry Smith 1084c60c7ad4SBarry Smith #undef __FUNCT__ 10859d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG" 10861e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) 1087676e1743SMark F. Adams { 10889d5b6da9SMark F. Adams PetscErrorCode ierr,(*r)(PC); 10891ab5ffc9SJed Brown PC_MG *mg = (PC_MG*)pc->data; 10901ab5ffc9SJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1091676e1743SMark F. Adams 1092676e1743SMark F. Adams PetscFunctionBegin; 1093c60c7ad4SBarry Smith pc_gamg->type = type; 10941c9cd337SJed Brown ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr); 10959d5b6da9SMark F. Adams if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 10961ab5ffc9SJed Brown if (pc_gamg->ops->destroy) { 10971ab5ffc9SJed Brown ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr); 10981ab5ffc9SJed Brown ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr); 1099e616c208SToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 11003ae0bb68SMark Adams /* cleaning up common data in pc_gamg - this should disapear someday */ 11013ae0bb68SMark Adams pc_gamg->data_cell_cols = 0; 11023ae0bb68SMark Adams pc_gamg->data_cell_rows = 0; 11033ae0bb68SMark Adams pc_gamg->orig_data_cell_cols = 0; 11043ae0bb68SMark Adams pc_gamg->orig_data_cell_rows = 0; 11053ae0bb68SMark Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 11063ae0bb68SMark Adams pc_gamg->data_sz = 0; 11071ab5ffc9SJed Brown } 11081ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr); 11091ab5ffc9SJed Brown ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr); 11109d5b6da9SMark F. Adams ierr = (*r)(pc);CHKERRQ(ierr); 1111676e1743SMark F. Adams PetscFunctionReturn(0); 1112676e1743SMark F. Adams } 1113676e1743SMark F. Adams 11145b89ad90SMark F. Adams #undef __FUNCT__ 11155adeb434SBarry Smith #define __FUNCT__ "PCView_GAMG" 11165adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer) 11175adeb434SBarry Smith { 11185adeb434SBarry Smith PetscErrorCode ierr; 11195adeb434SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 11205adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 11215adeb434SBarry Smith 11225adeb434SBarry Smith PetscFunctionBegin; 11235adeb434SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," GAMG specific options\n");CHKERRQ(ierr); 1124b001cb0fSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Threshold for dropping small values from graph %g\n",(double)pc_gamg->threshold);CHKERRQ(ierr); 11255adeb434SBarry Smith if (pc_gamg->ops->view) { 11265adeb434SBarry Smith ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr); 11275adeb434SBarry Smith } 11285adeb434SBarry Smith PetscFunctionReturn(0); 11295adeb434SBarry Smith } 11305adeb434SBarry Smith 11315adeb434SBarry Smith #undef __FUNCT__ 11325b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG" 11334416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc) 11345b89ad90SMark F. Adams { 1135676e1743SMark F. Adams PetscErrorCode ierr; 1136676e1743SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1137676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1138676e1743SMark F. Adams PetscBool flag; 11393b4367a7SBarry Smith MPI_Comm comm; 114014a9496bSBarry Smith char prefix[256]; 114114a9496bSBarry Smith const char *pcpre; 11425b89ad90SMark F. Adams 11435b89ad90SMark F. Adams PetscFunctionBegin; 11443b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1145e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr); 1146676e1743SMark F. Adams { 1147bd94a7aaSJed Brown char tname[256]; 11481a1c1e04SBarry Smith ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr); 1149bd94a7aaSJed Brown if (flag) { 1150bd94a7aaSJed Brown ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr); 11511ab5ffc9SJed Brown } 115294ae4db5SBarry Smith ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGRepartitioning",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr); 11531cc46a46SBarry Smith ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr); 115494ae4db5SBarry Smith ierr = PetscOptionsBool("-pc_gamg_use_agg_gasm","Use aggregation agragates for GASM smoother","PCGAMGUseASMAggs",pc_gamg->use_aggs_in_gasm,&pc_gamg->use_aggs_in_gasm,NULL);CHKERRQ(ierr); 115594ae4db5SBarry 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); 115694ae4db5SBarry 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); 115794ae4db5SBarry Smith ierr = PetscOptionsReal("-pc_gamg_threshold","Relative threshold to use for dropping edges in aggregation graph","PCGAMGSetThreshold",pc_gamg->threshold,&pc_gamg->threshold,&flag);CHKERRQ(ierr); 115894ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr); 1159b7cbab4eSMark Adams 1160b7cbab4eSMark Adams /* set options for subtype */ 1161e55864a3SBarry Smith if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);} 1162676e1743SMark F. Adams } 116314a9496bSBarry Smith ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr); 116414a9496bSBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr); 116514a9496bSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)pc_gamg->random,prefix);CHKERRQ(ierr); 116614a9496bSBarry Smith ierr = PetscRandomSetFromOptions(pc_gamg->random);CHKERRQ(ierr); 1167676e1743SMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 11685b89ad90SMark F. Adams PetscFunctionReturn(0); 11695b89ad90SMark F. Adams } 11705b89ad90SMark F. Adams 11715b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 11725b89ad90SMark F. Adams /*MC 11731cc46a46SBarry Smith PCGAMG - Geometric algebraic multigrid (AMG) preconditioner 11745b89ad90SMark F. Adams 1175280d9858SJed Brown Options Database Keys: 11765b89ad90SMark F. Adams Multigrid options(inherited) 11771cc46a46SBarry Smith + -pc_mg_cycles <v>: v or w (PCMGSetCycleType()) 1178280d9858SJed Brown . -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp) 1179280d9858SJed Brown . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown) 11808c1c2452SJed Brown - -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade 11815b89ad90SMark F. Adams 11821cc46a46SBarry Smith 11831cc46a46SBarry Smith Notes: In order to obtain good performance for PCGAMG for vector valued problems you must 11841cc46a46SBarry Smith $ Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point 11851cc46a46SBarry Smith $ Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator 11861cc46a46SBarry Smith $ See the Users Manual Chapter 4 for more details 11871cc46a46SBarry Smith 11885b89ad90SMark F. Adams Level: intermediate 1189280d9858SJed Brown 11901cc46a46SBarry Smith Concepts: algebraic multigrid 11915b89ad90SMark F. Adams 11921cc46a46SBarry Smith .seealso: PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(), 11931cc46a46SBarry Smith PCGAMGSetCoarseEqLim(), PCGAMGSetRepartitioning(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseASMAggs(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType() 11945b89ad90SMark F. Adams M*/ 1195b2573a8aSBarry Smith 11965b89ad90SMark F. Adams #undef __FUNCT__ 11975b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG" 11988cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc) 11995b89ad90SMark F. Adams { 12005b89ad90SMark F. Adams PetscErrorCode ierr; 12015b89ad90SMark F. Adams PC_GAMG *pc_gamg; 12025b89ad90SMark F. Adams PC_MG *mg; 12035b89ad90SMark F. Adams 12045b89ad90SMark F. Adams PetscFunctionBegin; 12051c1aac46SBarry Smith /* register AMG type */ 12061c1aac46SBarry Smith ierr = PCGAMGInitializePackage();CHKERRQ(ierr); 12071c1aac46SBarry Smith 12085b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 12091c1aac46SBarry Smith ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); 12105b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr); 12115b89ad90SMark F. Adams 12125b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 1213b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr); 12145b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 12151c1aac46SBarry Smith mg->galerkin = 2; /* Use Galerkin, but it is computed externally from PCMG by GAMG code */ 12165b89ad90SMark F. Adams mg->innerctx = pc_gamg; 12175b89ad90SMark F. Adams 1218b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr); 12191ab5ffc9SJed Brown 12209d5b6da9SMark F. Adams pc_gamg->setup_count = 0; 12219d5b6da9SMark F. Adams /* these should be in subctx but repartitioning needs simple arrays */ 12229d5b6da9SMark F. Adams pc_gamg->data_sz = 0; 12239d5b6da9SMark F. Adams pc_gamg->data = 0; 12245b89ad90SMark F. Adams 12259d5b6da9SMark F. Adams /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 12265b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 12275b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 12285b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 12295b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 12305adeb434SBarry Smith mg->view = PCView_GAMG; 12315b89ad90SMark F. Adams 1232bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr); 1233bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr); 1234bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartitioning_C",PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr); 12351cc46a46SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr); 1236bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseASMAggs_C",PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr); 1237bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr); 1238bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr); 1239c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr); 1240bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr); 12419d5b6da9SMark F. Adams pc_gamg->repart = PETSC_FALSE; 1242d3042614SMark Adams pc_gamg->reuse_prol = PETSC_FALSE; 1243ffc955d6SMark F. Adams pc_gamg->use_aggs_in_gasm = PETSC_FALSE; 1244038f3aa4SMark F. Adams pc_gamg->min_eq_proc = 50; 124525a145a7SMark Adams pc_gamg->coarse_eq_limit = 50; 1246d3042614SMark Adams pc_gamg->threshold = 0.; 12479d5b6da9SMark F. Adams pc_gamg->Nlevels = GAMG_MAXLEVELS; 12489ab59c8bSMark Adams pc_gamg->current_level = 0; /* don't need to init really */ 1249c238b0ebSToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 12509d5b6da9SMark F. Adams 125114a9496bSBarry Smith ierr = PetscRandomCreate(PetscObjectComm((PetscObject)pc),&pc_gamg->random);CHKERRQ(ierr); 125214a9496bSBarry Smith 1253bd94a7aaSJed Brown /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */ 1254bd94a7aaSJed Brown ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr); 12555b89ad90SMark F. Adams PetscFunctionReturn(0); 12565b89ad90SMark F. Adams } 12573e3471ccSMark Adams 12583e3471ccSMark Adams #undef __FUNCT__ 12593e3471ccSMark Adams #define __FUNCT__ "PCGAMGInitializePackage" 12603e3471ccSMark Adams /*@C 12613e3471ccSMark Adams PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called 12623e3471ccSMark Adams from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG() 12633e3471ccSMark Adams when using static libraries. 12643e3471ccSMark Adams 12653e3471ccSMark Adams Level: developer 12663e3471ccSMark Adams 12673e3471ccSMark Adams .keywords: PC, PCGAMG, initialize, package 12683e3471ccSMark Adams .seealso: PetscInitialize() 12693e3471ccSMark Adams @*/ 12703e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void) 12713e3471ccSMark Adams { 12723e3471ccSMark Adams PetscErrorCode ierr; 12733e3471ccSMark Adams 12743e3471ccSMark Adams PetscFunctionBegin; 12753e3471ccSMark Adams if (PCGAMGPackageInitialized) PetscFunctionReturn(0); 12763e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_TRUE; 12773e3471ccSMark Adams ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr); 12783e3471ccSMark Adams ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr); 12798e6d0c30SPeter Brune ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr); 12803e3471ccSMark Adams ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr); 1281c1c463dbSMark Adams 1282c1c463dbSMark Adams /* general events */ 1283fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr); 1284fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr); 1285fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr); 1286fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr); 1287c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr); 1288c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr); 1289fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr); 1290c1c463dbSMark Adams 12915b89ad90SMark F. Adams #if defined PETSC_GAMG_USE_LOG 12925b89ad90SMark F. Adams ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr); 12935b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr); 12945b89ad90SMark F. Adams /* PetscLogEventRegister(" G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */ 12955b89ad90SMark F. Adams /* PetscLogEventRegister(" G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */ 12965b89ad90SMark F. Adams /* PetscLogEventRegister(" G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */ 12975b89ad90SMark F. Adams ierr = PetscLogEventRegister(" MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr); 12985b89ad90SMark F. Adams ierr = PetscLogEventRegister(" geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr); 12995b89ad90SMark F. Adams ierr = PetscLogEventRegister(" geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr); 1300*bb235841SBarry Smith ierr = PetscLogEventRegister(" search-set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr); 13015b89ad90SMark F. Adams ierr = PetscLogEventRegister(" SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr); 13025b89ad90SMark F. Adams ierr = PetscLogEventRegister(" SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr); 13035b89ad90SMark F. Adams ierr = PetscLogEventRegister(" SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr); 13045b89ad90SMark F. Adams ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr); 13055b89ad90SMark F. Adams ierr = PetscLogEventRegister(" repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr); 13065b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr); 13075b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr); 13085b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr); 13095b89ad90SMark F. Adams 13105b89ad90SMark F. Adams /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */ 13115b89ad90SMark F. Adams /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */ 13125b89ad90SMark F. Adams /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */ 13135b89ad90SMark F. Adams /* create timer stages */ 13145b89ad90SMark F. Adams #if defined GAMG_STAGES 13155b89ad90SMark F. Adams { 13165b89ad90SMark F. Adams char str[32]; 13175b89ad90SMark F. Adams PetscInt lidx; 13185b89ad90SMark F. Adams sprintf(str,"MG Level %d (finest)",0); 13195b89ad90SMark F. Adams ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr); 13205b89ad90SMark F. Adams for (lidx=1; lidx<9; lidx++) { 13215b89ad90SMark F. Adams sprintf(str,"MG Level %d",lidx); 13225b89ad90SMark F. Adams ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr); 13235b89ad90SMark F. Adams } 13245b89ad90SMark F. Adams } 13255b89ad90SMark F. Adams #endif 13265b89ad90SMark F. Adams #endif 13273e3471ccSMark Adams PetscFunctionReturn(0); 13283e3471ccSMark Adams } 13293e3471ccSMark Adams 13303e3471ccSMark Adams #undef __FUNCT__ 13313e3471ccSMark Adams #define __FUNCT__ "PCGAMGFinalizePackage" 13323e3471ccSMark Adams /*@C 13331c1aac46SBarry Smith PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is 13341c1aac46SBarry Smith called from PetscFinalize() automatically. 13353e3471ccSMark Adams 13363e3471ccSMark Adams Level: developer 13373e3471ccSMark Adams 13383e3471ccSMark Adams .keywords: Petsc, destroy, package 13393e3471ccSMark Adams .seealso: PetscFinalize() 13403e3471ccSMark Adams @*/ 13413e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void) 13423e3471ccSMark Adams { 13433e3471ccSMark Adams PetscErrorCode ierr; 13443e3471ccSMark Adams 13453e3471ccSMark Adams PetscFunctionBegin; 13463e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_FALSE; 13473e3471ccSMark Adams ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr); 13483e3471ccSMark Adams PetscFunctionReturn(0); 13493e3471ccSMark Adams } 1350a36cf38bSToby Isaac 1351a36cf38bSToby Isaac #undef __FUNCT__ 1352a36cf38bSToby Isaac #define __FUNCT__ "PCGAMGRegister" 1353a36cf38bSToby Isaac /*@C 1354a36cf38bSToby Isaac PCGAMGRegister - Register a PCGAMG implementation. 1355a36cf38bSToby Isaac 1356a36cf38bSToby Isaac Input Parameters: 1357a36cf38bSToby Isaac + type - string that will be used as the name of the GAMG type. 1358a36cf38bSToby Isaac - create - function for creating the gamg context. 1359a36cf38bSToby Isaac 1360a36cf38bSToby Isaac Level: advanced 1361a36cf38bSToby Isaac 13621c1aac46SBarry Smith .seealso: PCGAMGType, PCGAMG, PCGAMGSetType() 1363a36cf38bSToby Isaac @*/ 1364a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC)) 1365a36cf38bSToby Isaac { 1366a36cf38bSToby Isaac PetscErrorCode ierr; 1367a36cf38bSToby Isaac 1368a36cf38bSToby Isaac PetscFunctionBegin; 1369a36cf38bSToby Isaac ierr = PCGAMGInitializePackage();CHKERRQ(ierr); 1370a36cf38bSToby Isaac ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr); 1371a36cf38bSToby Isaac PetscFunctionReturn(0); 1372a36cf38bSToby Isaac } 1373a36cf38bSToby Isaac 1374