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; 43a2f3521dSMark F. Adams if (pc_gamg->data) { /* this should not happen, cleaned up in SetUp */ 44ce94432eSBarry Smith PetscPrintf(PetscObjectComm((PetscObject)pc),"***[%d]%s this should not happen, cleaned up in SetUp\n",0,__FUNCT__); 459d5b6da9SMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 4658471d46SMark F. Adams } 471c1aac46SBarry Smith pc_gamg->data_sz = 0; 48878e152fSMark F. Adams ierr = PetscFree(pc_gamg->orig_data);CHKERRQ(ierr); 49a2f3521dSMark F. Adams PetscFunctionReturn(0); 50a2f3521dSMark F. Adams } 51a2f3521dSMark F. Adams 525b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 535b89ad90SMark F. Adams /* 54c238b0ebSToby Isaac PCGAMGCreateLevel_GAMG: create coarse op with RAP. repartition and/or reduce number 55a147abb0SMark F. Adams of active processors. 565b89ad90SMark F. Adams 575b89ad90SMark F. Adams Input Parameter: 58a2f3521dSMark F. Adams . pc - parameters + side effect: coarse data in 'pc_gamg->data' and 59a2f3521dSMark F. Adams 'pc_gamg->data_sz' are changed via repartitioning/reduction. 609d5b6da9SMark F. Adams . Amat_fine - matrix on this fine (k) level 61c5bfad50SMark F. Adams . cr_bs - coarse block size 623530afc2SMark F. Adams In/Output Parameter: 63a2f3521dSMark F. Adams . a_P_inout - prolongation operator to the next level (k-->k-1) 64afc97cdcSMark F. Adams . a_nactive_proc - number of active procs 6511e60469SMark F. Adams Output Parameter: 663530afc2SMark F. Adams . a_Amat_crs - coarse matrix that is created (k-1) 675b89ad90SMark F. Adams */ 685cb416c2SMark F. Adams 695b89ad90SMark F. Adams #undef __FUNCT__ 70c238b0ebSToby Isaac #define __FUNCT__ "PCGAMGCreateLevel_GAMG" 71b34066adSToby Isaac static PetscErrorCode PCGAMGCreateLevel_GAMG(PC pc,Mat Amat_fine,PetscInt cr_bs, 723cb8563fSToby Isaac Mat *a_P_inout,Mat *a_Amat_crs,PetscMPIInt *a_nactive_proc, 733cb8563fSToby Isaac IS * Pcolumnperm) 745b89ad90SMark F. Adams { 75a2f3521dSMark F. Adams PetscErrorCode ierr; 769d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 77486a8d0bSJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 78a2f3521dSMark F. Adams Mat Cmat,Pold=*a_P_inout; 793b4367a7SBarry Smith MPI_Comm comm; 80c5df96a5SBarry Smith PetscMPIInt rank,size,new_size,nactive=*a_nactive_proc; 813ae0bb68SMark Adams PetscInt ncrs_eq,ncrs,f_bs; 825b89ad90SMark F. Adams 835b89ad90SMark F. Adams PetscFunctionBegin; 843b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat_fine,&comm);CHKERRQ(ierr); 853b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 863b4367a7SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 87c5bfad50SMark F. Adams ierr = MatGetBlockSize(Amat_fine, &f_bs);CHKERRQ(ierr); 889d5b6da9SMark F. Adams ierr = MatPtAP(Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Cmat);CHKERRQ(ierr); 89038e3b61SMark F. Adams 903ae0bb68SMark Adams /* set 'ncrs' (nodes), 'ncrs_eq' (equations)*/ 910298fd71SBarry Smith ierr = MatGetLocalSize(Cmat, &ncrs_eq, NULL);CHKERRQ(ierr); 923ae0bb68SMark Adams if (pc_gamg->data_cell_rows>0) { 933ae0bb68SMark Adams ncrs = pc_gamg->data_sz/pc_gamg->data_cell_cols/pc_gamg->data_cell_rows; 9473911c69SBarry Smith } else { 953ae0bb68SMark Adams PetscInt bs; 963ae0bb68SMark Adams ierr = MatGetBlockSize(Cmat, &bs);CHKERRQ(ierr); 973ae0bb68SMark Adams ncrs = ncrs_eq/bs; 983ae0bb68SMark Adams } 99a2f3521dSMark F. Adams 100c5df96a5SBarry Smith /* get number of PEs to make active 'new_size', reduce, can be any integer 1-P */ 101a2f3521dSMark F. Adams { 102472110cdSMark F. Adams PetscInt ncrs_eq_glob; 1030298fd71SBarry Smith ierr = MatGetSize(Cmat, &ncrs_eq_glob, NULL);CHKERRQ(ierr); 104a90e85d9SMark Adams new_size = (PetscMPIInt)((float)ncrs_eq_glob/(float)pc_gamg->min_eq_proc + 0.5); /* hardwire min. number of eq/proc */ 1057f66b68fSMark Adams if (!new_size) new_size = 1; /* not likely, posible? */ 106c5df96a5SBarry Smith else if (new_size >= nactive) new_size = nactive; /* no change, rare */ 107a2f3521dSMark F. Adams } 108f852f58cSMark F. Adams 1093cb8563fSToby Isaac if (Pcolumnperm) *Pcolumnperm = NULL; 1103cb8563fSToby Isaac 111a90e85d9SMark Adams if (!pc_gamg->repart && new_size==nactive) *a_Amat_crs = Cmat; /* output - no repartitioning or reduction - could bail here */ 1122fa5cd67SKarl Rupp else { 1133ae0bb68SMark Adams PetscInt *counts,*newproc_idx,ii,jj,kk,strideNew,*tidx,ncrs_new,ncrs_eq_new,nloc_old; 114885364a3SMark Adams IS is_eq_newproc,is_eq_num,is_eq_num_prim,new_eq_indices; 115e33ef3b1SMark F. Adams 11671959b99SBarry Smith nloc_old = ncrs_eq/cr_bs; 11771959b99SBarry 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); 1180cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 1190cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 120b4fbaa2aSMark F. Adams #endif 121a2f3521dSMark F. Adams /* make 'is_eq_newproc' */ 122785e854fSJed Brown ierr = PetscMalloc1(size, &counts);CHKERRQ(ierr); 123a90e85d9SMark Adams if (pc_gamg->repart) { 124a2f3521dSMark F. Adams /* Repartition Cmat_{k} and move colums of P^{k}_{k-1} and coordinates of primal part accordingly */ 1255a9b9e01SMark F. Adams Mat adj; 1265a9b9e01SMark F. Adams 127302440fdSBarry Smith ierr = PetscInfo3(pc,"Repartition: size (active): %D --> %D, neq = %D\n",*a_nactive_proc,new_size,ncrs_eq);CHKERRQ(ierr); 1285a9b9e01SMark F. Adams 129a2f3521dSMark F. Adams /* get 'adj' */ 130c5bfad50SMark F. Adams if (cr_bs == 1) { 131038e3b61SMark F. Adams ierr = MatConvert(Cmat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 132806fa848SBarry Smith } else { 133a2f3521dSMark F. Adams /* make a scalar matrix to partition (no Stokes here) */ 134eb07cef2SMark F. Adams Mat tMat; 135a2f3521dSMark F. Adams PetscInt Istart_crs,Iend_crs,ncols,jj,Ii; 136b4fbaa2aSMark F. Adams const PetscScalar *vals; 137b4fbaa2aSMark F. Adams const PetscInt *idx; 138a2f3521dSMark F. Adams PetscInt *d_nnz, *o_nnz, M, N; 1399057884aSMark F. Adams static PetscInt llev = 0; 140d9558ea9SBarry Smith MatType mtype; 141b4fbaa2aSMark F. Adams 142e632b94dSBarry Smith ierr = PetscMalloc2(ncrs, &d_nnz,ncrs, &o_nnz);CHKERRQ(ierr); 143a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(Cmat, &Istart_crs, &Iend_crs);CHKERRQ(ierr); 144a2f3521dSMark F. Adams ierr = MatGetSize(Cmat, &M, &N);CHKERRQ(ierr); 145c5bfad50SMark F. Adams for (Ii = Istart_crs, jj = 0; Ii < Iend_crs; Ii += cr_bs, jj++) { 14658471d46SMark F. Adams ierr = MatGetRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr); 147c5bfad50SMark F. Adams d_nnz[jj] = ncols/cr_bs; 148c5bfad50SMark F. Adams o_nnz[jj] = ncols/cr_bs; 14958471d46SMark F. Adams ierr = MatRestoreRow(Cmat,Ii,&ncols,0,0);CHKERRQ(ierr); 1503ae0bb68SMark Adams if (d_nnz[jj] > ncrs) d_nnz[jj] = ncrs; 1513ae0bb68SMark Adams if (o_nnz[jj] > (M/cr_bs-ncrs)) o_nnz[jj] = M/cr_bs-ncrs; 15258471d46SMark F. Adams } 1536876a03eSMark F. Adams 154d9558ea9SBarry Smith ierr = MatGetType(Amat_fine,&mtype);CHKERRQ(ierr); 1553b4367a7SBarry Smith ierr = MatCreate(comm, &tMat);CHKERRQ(ierr); 1563ae0bb68SMark Adams ierr = MatSetSizes(tMat, ncrs, ncrs,PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 157d9558ea9SBarry Smith ierr = MatSetType(tMat,mtype);CHKERRQ(ierr); 158a2f3521dSMark F. Adams ierr = MatSeqAIJSetPreallocation(tMat,0,d_nnz);CHKERRQ(ierr); 159a2f3521dSMark F. Adams ierr = MatMPIAIJSetPreallocation(tMat,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 160e632b94dSBarry Smith ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr); 161eb07cef2SMark F. Adams 162a2f3521dSMark F. Adams for (ii = Istart_crs; ii < Iend_crs; ii++) { 163c5bfad50SMark F. Adams PetscInt dest_row = ii/cr_bs; 16422063be5SMark F. Adams ierr = MatGetRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 165eb07cef2SMark F. Adams for (jj = 0; jj < ncols; jj++) { 166c5bfad50SMark F. Adams PetscInt dest_col = idx[jj]/cr_bs; 167eb07cef2SMark F. Adams PetscScalar v = 1.0; 168eb07cef2SMark F. Adams ierr = MatSetValues(tMat,1,&dest_row,1,&dest_col,&v,ADD_VALUES);CHKERRQ(ierr); 169eb07cef2SMark F. Adams } 17022063be5SMark F. Adams ierr = MatRestoreRow(Cmat,ii,&ncols,&idx,&vals);CHKERRQ(ierr); 171eb07cef2SMark F. Adams } 172eb07cef2SMark F. Adams ierr = MatAssemblyBegin(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 173eb07cef2SMark F. Adams ierr = MatAssemblyEnd(tMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 174eb07cef2SMark F. Adams 175b4fbaa2aSMark F. Adams if (llev++ == -1) { 176b4fbaa2aSMark F. Adams PetscViewer viewer; char fname[32]; 1778caf3d72SBarry Smith ierr = PetscSNPrintf(fname,sizeof(fname),"part_mat_%D.mat",llev);CHKERRQ(ierr); 1783b4367a7SBarry Smith PetscViewerBinaryOpen(comm,fname,FILE_MODE_WRITE,&viewer); 179b4fbaa2aSMark F. Adams ierr = MatView(tMat, viewer);CHKERRQ(ierr); 1803bf036e2SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 181b4fbaa2aSMark F. Adams } 182b4fbaa2aSMark F. Adams 183eb07cef2SMark F. Adams ierr = MatConvert(tMat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj);CHKERRQ(ierr); 184eb07cef2SMark F. Adams 185eb07cef2SMark F. Adams ierr = MatDestroy(&tMat);CHKERRQ(ierr); 186a2f3521dSMark F. Adams } /* create 'adj' */ 187f150b916SMark F. Adams 188a2f3521dSMark F. Adams { /* partition: get newproc_idx */ 1895a9b9e01SMark F. Adams char prefix[256]; 1905a9b9e01SMark F. Adams const char *pcpre; 191b4fbaa2aSMark F. Adams const PetscInt *is_idx; 192b4fbaa2aSMark F. Adams MatPartitioning mpart; 193a4b7d37bSMark F. Adams IS proc_is; 194a2f3521dSMark F. Adams PetscInt targetPE; 1952f03bc48SMark F. Adams 1963b4367a7SBarry Smith ierr = MatPartitioningCreate(comm, &mpart);CHKERRQ(ierr); 1975ef31b24SMark F. Adams ierr = MatPartitioningSetAdjacency(mpart, adj);CHKERRQ(ierr); 1989d5b6da9SMark F. Adams ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr); 1998caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr); 20059a0be82SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 20111e60469SMark F. Adams ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 202c5df96a5SBarry Smith ierr = MatPartitioningSetNParts(mpart, new_size);CHKERRQ(ierr); 203a4b7d37bSMark F. Adams ierr = MatPartitioningApply(mpart, &proc_is);CHKERRQ(ierr); 20411e60469SMark F. Adams ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 2055a9b9e01SMark F. Adams 2065ef31b24SMark F. Adams /* collect IS info */ 207785e854fSJed Brown ierr = PetscMalloc1(ncrs_eq, &newproc_idx);CHKERRQ(ierr); 208a4b7d37bSMark F. Adams ierr = ISGetIndices(proc_is, &is_idx);CHKERRQ(ierr); 209a2f3521dSMark F. Adams targetPE = 1; /* bring to "front" of machine */ 210c5df96a5SBarry Smith /*targetPE = size/new_size;*/ /* spread partitioning across machine */ 211a2f3521dSMark F. Adams for (kk = jj = 0 ; kk < nloc_old ; kk++) { 212c5bfad50SMark F. Adams for (ii = 0 ; ii < cr_bs ; ii++, jj++) { 213a2f3521dSMark F. Adams newproc_idx[jj] = is_idx[kk] * targetPE; /* distribution */ 214eb07cef2SMark F. Adams } 2155ef31b24SMark F. Adams } 216a4b7d37bSMark F. Adams ierr = ISRestoreIndices(proc_is, &is_idx);CHKERRQ(ierr); 217a4b7d37bSMark F. Adams ierr = ISDestroy(&proc_is);CHKERRQ(ierr); 2185ef31b24SMark F. Adams } 2195ef31b24SMark F. Adams ierr = MatDestroy(&adj);CHKERRQ(ierr); 2205a9b9e01SMark F. Adams 2213b4367a7SBarry Smith ierr = ISCreateGeneral(comm, ncrs_eq, newproc_idx, PETSC_COPY_VALUES, &is_eq_newproc);CHKERRQ(ierr); 2228263b398SMark F. Adams ierr = PetscFree(newproc_idx);CHKERRQ(ierr); 223806fa848SBarry Smith } else { /* simple aggreagtion of parts -- 'is_eq_newproc' */ 224a2f3521dSMark F. Adams 225a2f3521dSMark F. Adams PetscInt rfactor,targetPE; 2265a9b9e01SMark F. Adams /* find factor */ 227c5df96a5SBarry Smith if (new_size == 1) rfactor = size; /* easy */ 2285a9b9e01SMark F. Adams else { 2295a9b9e01SMark F. Adams PetscReal best_fact = 0.; 2305a9b9e01SMark F. Adams jj = -1; 231c5df96a5SBarry Smith for (kk = 1 ; kk <= size ; kk++) { 2322a82295dSMark Adams if (!(size%kk)) { /* a candidate */ 233c5df96a5SBarry Smith PetscReal nactpe = (PetscReal)size/(PetscReal)kk, fact = nactpe/(PetscReal)new_size; 2345a9b9e01SMark F. Adams if (fact > 1.0) fact = 1./fact; /* keep fact < 1 */ 2355a9b9e01SMark F. Adams if (fact > best_fact) { 2365a9b9e01SMark F. Adams best_fact = fact; jj = kk; 2375a9b9e01SMark F. Adams } 2385a9b9e01SMark F. Adams } 2395a9b9e01SMark F. Adams } 2405a9b9e01SMark F. Adams if (jj != -1) rfactor = jj; 241a2f3521dSMark F. Adams else rfactor = 1; /* does this happen .. a prime */ 2425a9b9e01SMark F. Adams } 243c5df96a5SBarry Smith new_size = size/rfactor; 2445a9b9e01SMark F. Adams 245c5df96a5SBarry Smith if (new_size==nactive) { 246a2f3521dSMark F. Adams *a_Amat_crs = Cmat; /* output - no repartitioning or reduction, bail out because nested here */ 2475a9b9e01SMark F. Adams ierr = PetscFree(counts);CHKERRQ(ierr); 248302440fdSBarry Smith ierr = PetscInfo2(pc,"Aggregate processors noop: new_size=%D, neq(loc)=%D\n",new_size,ncrs_eq);CHKERRQ(ierr); 2495a9b9e01SMark F. Adams PetscFunctionReturn(0); 2505a9b9e01SMark F. Adams } 2515a9b9e01SMark F. Adams 252302440fdSBarry Smith ierr = PetscInfo1(pc,"Number of equations (loc) %D with simple aggregation\n",ncrs_eq);CHKERRQ(ierr); 253c5df96a5SBarry Smith targetPE = rank/rfactor; 2543b4367a7SBarry Smith ierr = ISCreateStride(comm, ncrs_eq, targetPE, 0, &is_eq_newproc);CHKERRQ(ierr); 255a2f3521dSMark F. Adams } /* end simple 'is_eq_newproc' */ 256e33ef3b1SMark F. Adams 25711e60469SMark F. Adams /* 258a2f3521dSMark F. Adams Create an index set from the is_eq_newproc index set to indicate the mapping TO 25911e60469SMark F. Adams */ 260a2f3521dSMark F. Adams ierr = ISPartitioningToNumbering(is_eq_newproc, &is_eq_num);CHKERRQ(ierr); 2617700e67bSMark Adams is_eq_num_prim = is_eq_num; 26211e60469SMark F. Adams /* 263a2f3521dSMark F. Adams Determine how many equations/vertices are assigned to each processor 26411e60469SMark F. Adams */ 265c5df96a5SBarry Smith ierr = ISPartitioningCount(is_eq_newproc, size, counts);CHKERRQ(ierr); 266c5df96a5SBarry Smith ncrs_eq_new = counts[rank]; 267a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_newproc);CHKERRQ(ierr); 2683ae0bb68SMark Adams ncrs_new = ncrs_eq_new/cr_bs; /* eqs */ 269a2f3521dSMark F. Adams 270a2f3521dSMark F. Adams ierr = PetscFree(counts);CHKERRQ(ierr); 2710cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 2720cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET12],0,0,0,0);CHKERRQ(ierr); 273b4fbaa2aSMark F. Adams #endif 274885364a3SMark 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 */ 275885364a3SMark Adams { 276885364a3SMark Adams Vec src_crd, dest_crd; 277885364a3SMark 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; 278885364a3SMark Adams VecScatter vecscat; 279885364a3SMark Adams PetscScalar *array; 280885364a3SMark Adams IS isscat; 281a2f3521dSMark F. Adams 282a2f3521dSMark F. Adams /* move data (for primal equations only) */ 28322063be5SMark F. Adams /* Create a vector to contain the newly ordered element information */ 2843b4367a7SBarry Smith ierr = VecCreate(comm, &dest_crd);CHKERRQ(ierr); 2853ae0bb68SMark Adams ierr = VecSetSizes(dest_crd, node_data_sz*ncrs_new, PETSC_DECIDE);CHKERRQ(ierr); 286c0dedaeaSBarry Smith ierr = VecSetType(dest_crd,VECSTANDARD);CHKERRQ(ierr); /* this is needed! */ 28711e60469SMark F. Adams /* 2889d5b6da9SMark F. Adams There are 'ndata_rows*ndata_cols' data items per node, (one can think of the vectors of having 289c5bfad50SMark F. Adams a block size of ...). Note, ISs are expanded into equation space by 'cr_bs'. 29011e60469SMark F. Adams */ 291854ce69bSBarry Smith ierr = PetscMalloc1(ncrs*node_data_sz, &tidx);CHKERRQ(ierr); 292a2f3521dSMark F. Adams ierr = ISGetIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 2933ae0bb68SMark Adams for (ii=0,jj=0; ii<ncrs; ii++) { 294c5bfad50SMark F. Adams PetscInt id = idx[ii*cr_bs]/cr_bs; /* get node back */ 295a2f3521dSMark F. Adams for (kk=0; kk<node_data_sz; kk++, jj++) tidx[jj] = id*node_data_sz + kk; 29611e60469SMark F. Adams } 297a2f3521dSMark F. Adams ierr = ISRestoreIndices(is_eq_num_prim, &idx);CHKERRQ(ierr); 2983ae0bb68SMark Adams ierr = ISCreateGeneral(comm, node_data_sz*ncrs, tidx, PETSC_COPY_VALUES, &isscat);CHKERRQ(ierr); 29992a756f0SMark F. Adams ierr = PetscFree(tidx);CHKERRQ(ierr); 30011e60469SMark F. Adams /* 30111e60469SMark F. Adams Create a vector to contain the original vertex information for each element 30211e60469SMark F. Adams */ 3033ae0bb68SMark Adams ierr = VecCreateSeq(PETSC_COMM_SELF, node_data_sz*ncrs, &src_crd);CHKERRQ(ierr); 3049d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols; jj++) { 3053ae0bb68SMark Adams const PetscInt stride0=ncrs*pc_gamg->data_cell_rows; 3063ae0bb68SMark Adams for (ii=0; ii<ncrs; ii++) { 3079d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows; kk++) { 308a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*stride0, jx = ii*node_data_sz + kk*ndata_cols + jj; 309c8b0795cSMark F. Adams PetscScalar tt = (PetscScalar)pc_gamg->data[ix]; 310676e1743SMark F. Adams ierr = VecSetValues(src_crd, 1, &jx, &tt, INSERT_VALUES);CHKERRQ(ierr); 311d3d6bff4SMark F. Adams } 312038e3b61SMark F. Adams } 313eb07cef2SMark F. Adams } 314eb07cef2SMark F. Adams ierr = VecAssemblyBegin(src_crd);CHKERRQ(ierr); 315eb07cef2SMark F. Adams ierr = VecAssemblyEnd(src_crd);CHKERRQ(ierr); 31611e60469SMark F. Adams /* 31711e60469SMark F. Adams Scatter the element vertex information (still in the original vertex ordering) 31811e60469SMark F. Adams to the correct processor 31911e60469SMark F. Adams */ 3200298fd71SBarry Smith ierr = VecScatterCreate(src_crd, NULL, dest_crd, isscat, &vecscat);CHKERRQ(ierr); 32111e60469SMark F. Adams ierr = ISDestroy(&isscat);CHKERRQ(ierr); 32211e60469SMark F. Adams ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 32311e60469SMark F. Adams ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 32411e60469SMark F. Adams ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 32511e60469SMark F. Adams ierr = VecDestroy(&src_crd);CHKERRQ(ierr); 32611e60469SMark F. Adams /* 32711e60469SMark F. Adams Put the element vertex data into a new allocation of the gdata->ele 32811e60469SMark F. Adams */ 329c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 330578f55a3SPeter Brune ierr = PetscMalloc1(node_data_sz*ncrs_new, &pc_gamg->data);CHKERRQ(ierr); 3312fa5cd67SKarl Rupp 3323ae0bb68SMark Adams pc_gamg->data_sz = node_data_sz*ncrs_new; 3333ae0bb68SMark Adams strideNew = ncrs_new*ndata_rows; 3342fa5cd67SKarl Rupp 33511e60469SMark F. Adams ierr = VecGetArray(dest_crd, &array);CHKERRQ(ierr); 3369d5b6da9SMark F. Adams for (jj=0; jj<ndata_cols; jj++) { 3373ae0bb68SMark Adams for (ii=0; ii<ncrs_new; ii++) { 3389d5b6da9SMark F. Adams for (kk=0; kk<ndata_rows; kk++) { 339a2f3521dSMark F. Adams PetscInt ix = ii*ndata_rows + kk + jj*strideNew, jx = ii*node_data_sz + kk*ndata_cols + jj; 340c8b0795cSMark F. Adams pc_gamg->data[ix] = PetscRealPart(array[jx]); 341d3d6bff4SMark F. Adams } 342038e3b61SMark F. Adams } 343038e3b61SMark F. Adams } 34411e60469SMark F. Adams ierr = VecRestoreArray(dest_crd, &array);CHKERRQ(ierr); 34511e60469SMark F. Adams ierr = VecDestroy(&dest_crd);CHKERRQ(ierr); 346885364a3SMark Adams } 347a2f3521dSMark F. Adams /* move A and P (columns) with new layout */ 3480cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3490cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 350ed3f9983SMark F. Adams #endif 351a2f3521dSMark F. Adams 35211e60469SMark F. Adams /* 35311e60469SMark F. Adams Invert for MatGetSubMatrix 35411e60469SMark F. Adams */ 355a2f3521dSMark F. Adams ierr = ISInvertPermutation(is_eq_num, ncrs_eq_new, &new_eq_indices);CHKERRQ(ierr); 356a2f3521dSMark F. Adams ierr = ISSort(new_eq_indices);CHKERRQ(ierr); /* is this needed? */ 357c5bfad50SMark F. Adams ierr = ISSetBlockSize(new_eq_indices, cr_bs);CHKERRQ(ierr); 358a2f3521dSMark F. Adams if (is_eq_num != is_eq_num_prim) { 359a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_num_prim);CHKERRQ(ierr); /* could be same as 'is_eq_num' */ 360a2f3521dSMark F. Adams } 3613cb8563fSToby Isaac if (Pcolumnperm) { 3623cb8563fSToby Isaac ierr = PetscObjectReference((PetscObject)new_eq_indices);CHKERRQ(ierr); 3633cb8563fSToby Isaac *Pcolumnperm = new_eq_indices; 3643cb8563fSToby Isaac } 365a2f3521dSMark F. Adams ierr = ISDestroy(&is_eq_num);CHKERRQ(ierr); 3660cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3670cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET13],0,0,0,0);CHKERRQ(ierr); 3680cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 369ed3f9983SMark F. Adams #endif 370a2f3521dSMark F. Adams /* 'a_Amat_crs' output */ 371a2f3521dSMark F. Adams { 372a2f3521dSMark F. Adams Mat mat; 373806fa848SBarry Smith ierr = MatGetSubMatrix(Cmat, new_eq_indices, new_eq_indices, MAT_INITIAL_MATRIX, &mat);CHKERRQ(ierr); 374a2f3521dSMark F. Adams *a_Amat_crs = mat; 375c5bfad50SMark F. Adams 376c5bfad50SMark F. Adams if (!PETSC_TRUE) { 377c5bfad50SMark F. Adams PetscInt cbs, rbs; 378c5bfad50SMark F. Adams ierr = MatGetBlockSizes(Cmat, &rbs, &cbs);CHKERRQ(ierr); 379c5df96a5SBarry Smith ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Old Mat rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr); 380c5bfad50SMark F. Adams ierr = MatGetBlockSizes(mat, &rbs, &cbs);CHKERRQ(ierr); 381c5df96a5SBarry Smith ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s New Mat rbs=%d cbs=%d cr_bs=%d\n",rank,__FUNCT__,rbs,cbs,cr_bs);CHKERRQ(ierr); 382c5bfad50SMark F. Adams } 383a2f3521dSMark F. Adams } 384038e3b61SMark F. Adams ierr = MatDestroy(&Cmat);CHKERRQ(ierr); 385a2f3521dSMark F. Adams 3860cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3870cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET14],0,0,0,0);CHKERRQ(ierr); 388ed3f9983SMark F. Adams #endif 38911e60469SMark F. Adams /* prolongator */ 39011e60469SMark F. Adams { 39111e60469SMark F. Adams IS findices; 392a2f3521dSMark F. Adams PetscInt Istart,Iend; 393a2f3521dSMark F. Adams Mat Pnew; 394a2f3521dSMark F. Adams ierr = MatGetOwnershipRange(Pold, &Istart, &Iend);CHKERRQ(ierr); 3950cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 3960cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 397ed3f9983SMark F. Adams #endif 3983b4367a7SBarry Smith ierr = ISCreateStride(comm,Iend-Istart,Istart,1,&findices);CHKERRQ(ierr); 399c5bfad50SMark F. Adams ierr = ISSetBlockSize(findices,f_bs);CHKERRQ(ierr); 400806fa848SBarry Smith ierr = MatGetSubMatrix(Pold, findices, new_eq_indices, MAT_INITIAL_MATRIX, &Pnew);CHKERRQ(ierr); 40111e60469SMark F. Adams ierr = ISDestroy(&findices);CHKERRQ(ierr); 402c5bfad50SMark F. Adams 403c5bfad50SMark F. Adams if (!PETSC_TRUE) { 404c5bfad50SMark F. Adams PetscInt cbs, rbs; 405c5bfad50SMark F. Adams ierr = MatGetBlockSizes(Pold, &rbs, &cbs);CHKERRQ(ierr); 406c5df96a5SBarry Smith ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pold rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr); 407c5bfad50SMark F. Adams ierr = MatGetBlockSizes(Pnew, &rbs, &cbs);CHKERRQ(ierr); 408c5df96a5SBarry Smith ierr = PetscPrintf(MPI_COMM_SELF,"[%d]%s Pnew rbs=%d cbs=%d\n",rank,__FUNCT__,rbs,cbs);CHKERRQ(ierr); 409c5bfad50SMark F. Adams } 4100cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 4110cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET15],0,0,0,0);CHKERRQ(ierr); 412ed3f9983SMark F. Adams #endif 4133530afc2SMark F. Adams ierr = MatDestroy(a_P_inout);CHKERRQ(ierr); 414a2f3521dSMark F. Adams 415a2f3521dSMark F. Adams /* output - repartitioned */ 416a2f3521dSMark F. Adams *a_P_inout = Pnew; 417e33ef3b1SMark F. Adams } 418a2f3521dSMark F. Adams ierr = ISDestroy(&new_eq_indices);CHKERRQ(ierr); 4195b89ad90SMark F. Adams 420c5df96a5SBarry Smith *a_nactive_proc = new_size; /* output */ 421a2f3521dSMark F. Adams } 4225a9b9e01SMark F. Adams 423a2f3521dSMark F. Adams /* outout matrix data */ 424c8b0795cSMark F. Adams if (!PETSC_TRUE) { 425c8b0795cSMark F. Adams PetscViewer viewer; char fname[32]; static int llev=0; Cmat = *a_Amat_crs; 4267f66b68fSMark Adams if (!llev) { 427c8b0795cSMark F. Adams sprintf(fname,"Cmat_%d.m",llev++); 4283b4367a7SBarry Smith PetscViewerASCIIOpen(comm,fname,&viewer); 429*6a9046bcSBarry Smith ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 430c8b0795cSMark F. Adams ierr = MatView(Amat_fine, viewer);CHKERRQ(ierr); 431*6a9046bcSBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 432c8b0795cSMark F. Adams ierr = PetscViewerDestroy(&viewer); 433c8b0795cSMark F. Adams } 434c8b0795cSMark F. Adams sprintf(fname,"Cmat_%d.m",llev++); 4353b4367a7SBarry Smith PetscViewerASCIIOpen(comm,fname,&viewer); 436*6a9046bcSBarry Smith ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 437c8b0795cSMark F. Adams ierr = MatView(Cmat, viewer);CHKERRQ(ierr); 438*6a9046bcSBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 439c8b0795cSMark F. Adams ierr = PetscViewerDestroy(&viewer); 440c8b0795cSMark F. Adams } 4415b89ad90SMark F. Adams PetscFunctionReturn(0); 4425b89ad90SMark F. Adams } 4435b89ad90SMark F. Adams 4445b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 4455b89ad90SMark F. Adams /* 4465b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 4475b89ad90SMark F. Adams by setting data structures and options. 4485b89ad90SMark F. Adams 4495b89ad90SMark F. Adams Input Parameter: 4505b89ad90SMark F. Adams . pc - the preconditioner context 4515b89ad90SMark F. Adams 4525b89ad90SMark F. Adams */ 4535b89ad90SMark F. Adams #undef __FUNCT__ 4545b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG" 4559d5b6da9SMark F. Adams PetscErrorCode PCSetUp_GAMG(PC pc) 4565b89ad90SMark F. Adams { 4575b89ad90SMark F. Adams PetscErrorCode ierr; 4589d5b6da9SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 4595b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 4602adcac29SMark F. Adams Mat Pmat = pc->pmat; 461a2f3521dSMark F. Adams PetscInt fine_level,level,level1,bs,M,qq,lidx,nASMBlocksArr[GAMG_MAXLEVELS]; 4623b4367a7SBarry Smith MPI_Comm comm; 463c5df96a5SBarry Smith PetscMPIInt rank,size,nactivepe; 464587fa25dSMark F. Adams Mat Aarr[GAMG_MAXLEVELS],Parr[GAMG_MAXLEVELS]; 465e696ed0bSMark F. Adams IS *ASMLocalIDsArr[GAMG_MAXLEVELS]; 466a2f3521dSMark F. Adams PetscLogDouble nnz0=0.,nnztot=0.; 467569f4572SMark Adams MatInfo info; 4685ef31b24SMark F. Adams 4695b89ad90SMark F. Adams PetscFunctionBegin; 4703b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 4713b4367a7SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 4723b4367a7SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 473dfd5c07aSMark F. Adams 47484d3f75bSMark F. Adams if (pc_gamg->setup_count++ > 0) { 4751c1aac46SBarry Smith if ((PetscBool)(!pc_gamg->reuse_prol)) { 476878e152fSMark F. Adams /* reset everything */ 477878e152fSMark F. Adams ierr = PCReset_MG(pc);CHKERRQ(ierr); 478878e152fSMark F. Adams pc->setupcalled = 0; 479806fa848SBarry Smith } else { 48084d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 48103a628feSMark F. Adams /* just do Galerkin grids */ 48258471d46SMark F. Adams Mat B,dA,dB; 48358471d46SMark F. Adams 48471959b99SBarry Smith if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet"); 4859d5b6da9SMark F. Adams if (pc_gamg->Nlevels > 1) { 48658471d46SMark F. Adams /* currently only handle case where mat and pmat are the same on coarser levels */ 48723ee1639SBarry Smith ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 48858471d46SMark F. Adams /* (re)set to get dirty flag */ 48923ee1639SBarry Smith ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB);CHKERRQ(ierr); 49058471d46SMark F. Adams 4912fb0b348SMark F. Adams for (level=pc_gamg->Nlevels-2; level>=0; level--) { 49203a628feSMark F. Adams /* the first time through the matrix structure has changed from repartitioning */ 4930a97e771SToby Isaac if (pc_gamg->setup_count==2) { 49403a628feSMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 495084a8fe3SJed Brown ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr); 4962fa5cd67SKarl Rupp 49703a628feSMark F. Adams mglevels[level]->A = B; 498806fa848SBarry Smith } else { 49923ee1639SBarry Smith ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B);CHKERRQ(ierr); 50058471d46SMark F. Adams ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 50103a628feSMark F. Adams } 50223ee1639SBarry Smith ierr = KSPSetOperators(mglevels[level]->smoothd,B,B);CHKERRQ(ierr); 50358471d46SMark F. Adams dB = B; 50458471d46SMark F. Adams } 5055f8cf99dSMark F. Adams } 506d5280255SMark F. Adams 507d5280255SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 508d5280255SMark F. Adams 50958471d46SMark F. Adams PetscFunctionReturn(0); 510eb07cef2SMark F. Adams } 511878e152fSMark F. Adams } 512f6536408SMark F. Adams 513878e152fSMark F. Adams if (!pc_gamg->data) { 514878e152fSMark F. Adams if (pc_gamg->orig_data) { 515878e152fSMark F. Adams ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 5160298fd71SBarry Smith ierr = MatGetLocalSize(Pmat, &qq, NULL);CHKERRQ(ierr); 5172fa5cd67SKarl Rupp 518878e152fSMark F. Adams pc_gamg->data_sz = (qq/bs)*pc_gamg->orig_data_cell_rows*pc_gamg->orig_data_cell_cols; 519878e152fSMark F. Adams pc_gamg->data_cell_rows = pc_gamg->orig_data_cell_rows; 520878e152fSMark F. Adams pc_gamg->data_cell_cols = pc_gamg->orig_data_cell_cols; 5212fa5cd67SKarl Rupp 522785e854fSJed Brown ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->data);CHKERRQ(ierr); 523878e152fSMark F. Adams for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->data[qq] = pc_gamg->orig_data[qq]; 524806fa848SBarry Smith } else { 5251ab5ffc9SJed Brown if (!pc_gamg->ops->createdefaultdata) SETERRQ(comm,PETSC_ERR_PLIB,"'createdefaultdata' not set(?) need to support NULL data"); 5267700e67bSMark Adams ierr = pc_gamg->ops->createdefaultdata(pc,Pmat);CHKERRQ(ierr); 5279d5b6da9SMark F. Adams } 528878e152fSMark F. Adams } 529878e152fSMark F. Adams 530878e152fSMark F. Adams /* cache original data for reuse */ 5311c1aac46SBarry Smith if (!pc_gamg->orig_data && (PetscBool)(!pc_gamg->reuse_prol)) { 532785e854fSJed Brown ierr = PetscMalloc1(pc_gamg->data_sz, &pc_gamg->orig_data);CHKERRQ(ierr); 533878e152fSMark F. Adams for (qq=0; qq<pc_gamg->data_sz; qq++) pc_gamg->orig_data[qq] = pc_gamg->data[qq]; 534878e152fSMark F. Adams pc_gamg->orig_data_cell_rows = pc_gamg->data_cell_rows; 535878e152fSMark F. Adams pc_gamg->orig_data_cell_cols = pc_gamg->data_cell_cols; 536878e152fSMark F. Adams } 537038e3b61SMark F. Adams 538302f38e8SMark F. Adams /* get basic dims */ 539302f38e8SMark F. Adams ierr = MatGetBlockSize(Pmat, &bs);CHKERRQ(ierr); 540a2f3521dSMark F. Adams ierr = MatGetSize(Pmat, &M, &qq);CHKERRQ(ierr); 54184d3f75bSMark F. Adams 542569f4572SMark Adams ierr = MatGetInfo(Pmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); /* global reduction */ 543569f4572SMark Adams nnz0 = info.nz_used; 544569f4572SMark Adams nnztot = info.nz_used; 5451d5b2942SMark Adams ierr = PetscInfo6(pc,"level %d) N=%D, n data rows=%d, n data cols=%d, nnz/row (ave)=%d, np=%d\n", 546569f4572SMark Adams 0,M,pc_gamg->data_cell_rows,pc_gamg->data_cell_cols, 547569f4572SMark Adams (int)(nnz0/(PetscReal)M+0.5),size); 548569f4572SMark Adams CHKERRQ(ierr); 549569f4572SMark Adams 550a2f3521dSMark F. Adams /* Get A_i and R_i */ 551c5df96a5SBarry Smith for (level=0, Aarr[0]=Pmat, nactivepe = size; /* hard wired stopping logic */ 5527f66b68fSMark Adams level < (pc_gamg->Nlevels-1) && (!level || M>pc_gamg->coarse_eq_limit); 5530205a208SMark F. Adams level++) { 5549ab59c8bSMark Adams pc_gamg->current_level = level; 5555b89ad90SMark F. Adams level1 = level + 1; 5560cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 5570cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 558a2f3521dSMark F. Adams #if (defined GAMG_STAGES) 559a2f3521dSMark F. Adams ierr = PetscLogStagePush(gamg_stages[level]);CHKERRQ(ierr); 560b4fbaa2aSMark F. Adams #endif 561a2f3521dSMark F. Adams #endif 562c8b0795cSMark F. Adams { /* construct prolongator */ 563725b86d8SJed Brown Mat Gmat; 5640cbbd2e1SMark F. Adams PetscCoarsenData *agg_lists; 5657700e67bSMark Adams Mat Prol11; 566c8b0795cSMark F. Adams 5677700e67bSMark Adams ierr = pc_gamg->ops->graph(pc,Aarr[level], &Gmat);CHKERRQ(ierr); 5681ab5ffc9SJed Brown ierr = pc_gamg->ops->coarsen(pc, &Gmat, &agg_lists);CHKERRQ(ierr); 5697700e67bSMark Adams ierr = pc_gamg->ops->prolongator(pc,Aarr[level],Gmat,agg_lists,&Prol11);CHKERRQ(ierr); 570c8b0795cSMark F. Adams 571a2f3521dSMark F. Adams /* could have failed to create new level */ 572a2f3521dSMark F. Adams if (Prol11) { 5739d5b6da9SMark F. Adams /* get new block size of coarse matrices */ 5740298fd71SBarry Smith ierr = MatGetBlockSizes(Prol11, NULL, &bs);CHKERRQ(ierr); 575a2f3521dSMark F. Adams 576fd1112cbSBarry Smith if (pc_gamg->ops->optprolongator) { 577c8b0795cSMark F. Adams /* smooth */ 578fd1112cbSBarry Smith ierr = pc_gamg->ops->optprolongator(pc, Aarr[level], &Prol11);CHKERRQ(ierr); 579c8b0795cSMark F. Adams } 580c8b0795cSMark F. Adams 5817700e67bSMark Adams Parr[level1] = Prol11; 5820298fd71SBarry Smith } else Parr[level1] = NULL; 583ffc955d6SMark F. Adams 584ffc955d6SMark F. Adams if (pc_gamg->use_aggs_in_gasm) { 5851b18a24aSMark Adams PetscInt bs; 5861b18a24aSMark Adams ierr = MatGetBlockSizes(Prol11, &bs, NULL);CHKERRQ(ierr); 587806fa848SBarry Smith ierr = PetscCDGetASMBlocks(agg_lists, bs, &nASMBlocksArr[level], &ASMLocalIDsArr[level]);CHKERRQ(ierr); 588ffc955d6SMark F. Adams } 589ffc955d6SMark F. Adams 590a2f3521dSMark F. Adams ierr = MatDestroy(&Gmat);CHKERRQ(ierr); 59141b27cdeSMark F. Adams ierr = PetscCDDestroy(agg_lists);CHKERRQ(ierr); 592a2f3521dSMark F. Adams } /* construct prolongator scope */ 5930cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 5940cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET1],0,0,0,0);CHKERRQ(ierr); 595c8b0795cSMark F. Adams #endif 5967f66b68fSMark Adams if (!level) Aarr[0] = Pmat; /* use Pmat for finest level setup */ 597c8b0795cSMark F. Adams if (!Parr[level1]) { 598569f4572SMark Adams ierr = PetscInfo1(pc,"Stop gridding, level %D\n",level);CHKERRQ(ierr); 599a90e85d9SMark Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 600a90e85d9SMark Adams ierr = PetscLogStagePop();CHKERRQ(ierr); 601a90e85d9SMark Adams #endif 602c8b0795cSMark F. Adams break; 603c8b0795cSMark F. Adams } 6040cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 6050cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 606b4fbaa2aSMark F. Adams #endif 607a2f3521dSMark F. Adams 6081c1aac46SBarry Smith ierr = pc_gamg->ops->createlevel(pc, Aarr[level], bs,&Parr[level1], &Aarr[level1], &nactivepe, NULL);CHKERRQ(ierr); 609a2f3521dSMark F. Adams 6100cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 6110cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET2],0,0,0,0);CHKERRQ(ierr); 612b4fbaa2aSMark F. Adams #endif 613a2f3521dSMark F. Adams ierr = MatGetSize(Aarr[level1], &M, &qq);CHKERRQ(ierr); 614569f4572SMark Adams 615569f4572SMark Adams ierr = MatGetInfo(Aarr[level1], MAT_GLOBAL_SUM, &info);CHKERRQ(ierr); 616569f4572SMark Adams nnztot += info.nz_used; 6171d5b2942SMark 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); 618569f4572SMark Adams 6190cbbd2e1SMark F. Adams #if (defined PETSC_GAMG_USE_LOG && defined GAMG_STAGES) 620b4fbaa2aSMark F. Adams ierr = PetscLogStagePop();CHKERRQ(ierr); 621b4fbaa2aSMark F. Adams #endif 622a90e85d9SMark Adams /* stop if one node or one proc -- could pull back for singular problems */ 6239ab59c8bSMark Adams if ( (pc_gamg->data_cell_cols && M/pc_gamg->data_cell_cols < 2) || (!pc_gamg->data_cell_cols && M/bs < 2) ) { 6249ab59c8bSMark Adams ierr = PetscInfo2(pc,"HARD stop of coarsening on level %D. Grid too small: %D block nodes\n",level,M/bs);CHKERRQ(ierr); 625a90e85d9SMark Adams level++; 626a90e85d9SMark Adams break; 627a90e85d9SMark Adams } 628c8b0795cSMark F. Adams } /* levels */ 629c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 630c8b0795cSMark F. Adams 631569f4572SMark Adams ierr = PetscInfo2(pc,"%D levels, grid complexity = %g\n",level+1,nnztot/nnz0);CHKERRQ(ierr); 6329d5b6da9SMark F. Adams pc_gamg->Nlevels = level + 1; 6335b89ad90SMark F. Adams fine_level = level; 6340298fd71SBarry Smith ierr = PCMGSetLevels(pc,pc_gamg->Nlevels,NULL);CHKERRQ(ierr); 6355b89ad90SMark F. Adams 63684d3f75bSMark F. Adams /* simple setup */ 63784d3f75bSMark F. Adams if (!PETSC_TRUE) { 63884d3f75bSMark F. Adams PC_MG_Levels **mglevels = mg->levels; 63984d3f75bSMark F. Adams for (lidx=0,level=pc_gamg->Nlevels-1; 64084d3f75bSMark F. Adams lidx<fine_level; 64184d3f75bSMark F. Adams lidx++, level--) { 64284d3f75bSMark F. Adams ierr = PCMGSetInterpolation(pc, lidx+1, Parr[level]);CHKERRQ(ierr); 64323ee1639SBarry Smith ierr = KSPSetOperators(mglevels[lidx]->smoothd, Aarr[level], Aarr[level]);CHKERRQ(ierr); 64484d3f75bSMark F. Adams ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr); 64584d3f75bSMark F. Adams ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr); 64684d3f75bSMark F. Adams } 64723ee1639SBarry Smith ierr = KSPSetOperators(mglevels[fine_level]->smoothd, Aarr[0], Aarr[0]);CHKERRQ(ierr); 64884d3f75bSMark F. Adams 64984d3f75bSMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 650806fa848SBarry Smith } else if (pc_gamg->Nlevels > 1) { /* don't setup MG if one level */ 651d5280255SMark F. Adams /* set default smoothers & set operators */ 6529d5b6da9SMark F. Adams for (lidx = 1, level = pc_gamg->Nlevels-2; 653587fa25dSMark F. Adams lidx <= fine_level; 654587fa25dSMark F. Adams lidx++, level--) { 655ffc955d6SMark F. Adams KSP smoother; 656ffc955d6SMark F. Adams PC subpc; 657a2f3521dSMark F. Adams 6589d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 659f6536408SMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 660ffc955d6SMark F. Adams 661a2f3521dSMark F. Adams ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 662a2f3521dSMark F. Adams /* set ops */ 66323ee1639SBarry Smith ierr = KSPSetOperators(smoother, Aarr[level], Aarr[level]);CHKERRQ(ierr); 664a2f3521dSMark F. Adams ierr = PCMGSetInterpolation(pc, lidx, Parr[level+1]);CHKERRQ(ierr); 665a2f3521dSMark F. Adams 666a2f3521dSMark F. Adams /* set defaults */ 6676c9de887SHong Zhang ierr = KSPSetType(smoother, KSPCHEBYSHEV);CHKERRQ(ierr); 668a2f3521dSMark F. Adams 6691b18a24aSMark Adams /* set blocks for GASM smoother that uses the 'aggregates' */ 670ffc955d6SMark F. Adams if (pc_gamg->use_aggs_in_gasm) { 6712d3561bbSSatish Balay PetscInt sz; 6722d3561bbSSatish Balay IS *is; 673a2f3521dSMark F. Adams 6742d3561bbSSatish Balay sz = nASMBlocksArr[level]; 6752d3561bbSSatish Balay is = ASMLocalIDsArr[level]; 676ffc955d6SMark F. Adams ierr = PCSetType(subpc, PCGASM);CHKERRQ(ierr); 6771b18a24aSMark Adams ierr = PCGASMSetOverlap(subpc, 0);CHKERRQ(ierr); 6787f66b68fSMark Adams if (!sz) { 679ffc955d6SMark F. Adams IS is; 680ffc955d6SMark F. Adams PetscInt my0,kk; 681ffc955d6SMark F. Adams ierr = MatGetOwnershipRange(Aarr[level], &my0, &kk);CHKERRQ(ierr); 682ffc955d6SMark F. Adams ierr = ISCreateGeneral(PETSC_COMM_SELF, 1, &my0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr); 6830298fd71SBarry Smith ierr = PCGASMSetSubdomains(subpc, 1, &is, NULL);CHKERRQ(ierr); 684a94c3b12SMark F. Adams ierr = ISDestroy(&is);CHKERRQ(ierr); 685806fa848SBarry Smith } else { 686a94c3b12SMark F. Adams PetscInt kk; 6870298fd71SBarry Smith ierr = PCGASMSetSubdomains(subpc, sz, is, NULL);CHKERRQ(ierr); 688a94c3b12SMark F. Adams for (kk=0; kk<sz; kk++) { 689a94c3b12SMark F. Adams ierr = ISDestroy(&is[kk]);CHKERRQ(ierr); 690a94c3b12SMark F. Adams } 691ffc955d6SMark F. Adams ierr = PetscFree(is);CHKERRQ(ierr); 692ffc955d6SMark F. Adams } 6930298fd71SBarry Smith ASMLocalIDsArr[level] = NULL; 694ffc955d6SMark F. Adams nASMBlocksArr[level] = 0; 695ffc955d6SMark F. Adams ierr = PCGASMSetType(subpc, PC_GASM_BASIC);CHKERRQ(ierr); 696806fa848SBarry Smith } else { 697890ffe84SMark Adams ierr = PCSetType(subpc, PCSOR);CHKERRQ(ierr); 698ffc955d6SMark F. Adams } 699d5280255SMark F. Adams } 700d5280255SMark F. Adams { 701d5280255SMark F. Adams /* coarse grid */ 702d5280255SMark F. Adams KSP smoother,*k2; PC subpc,pc2; PetscInt ii,first; 703d5280255SMark F. Adams Mat Lmat = Aarr[(level=pc_gamg->Nlevels-1)]; lidx = 0; 704d5280255SMark F. Adams ierr = PCMGGetSmoother(pc, lidx, &smoother);CHKERRQ(ierr); 70523ee1639SBarry Smith ierr = KSPSetOperators(smoother, Lmat, Lmat);CHKERRQ(ierr); 706d5280255SMark F. Adams ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr); 707d5280255SMark F. Adams ierr = KSPGetPC(smoother, &subpc);CHKERRQ(ierr); 708d5280255SMark F. Adams ierr = PCSetType(subpc, PCBJACOBI);CHKERRQ(ierr); 709d5280255SMark F. Adams ierr = PCSetUp(subpc);CHKERRQ(ierr); 71071959b99SBarry Smith ierr = PCBJacobiGetSubKSP(subpc,&ii,&first,&k2);CHKERRQ(ierr); 71171959b99SBarry Smith if (ii != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D is not one",ii); 712d5280255SMark F. Adams ierr = KSPGetPC(k2[0],&pc2);CHKERRQ(ierr); 713d5280255SMark F. Adams ierr = PCSetType(pc2, PCLU);CHKERRQ(ierr); 7149dbfc187SHong Zhang ierr = PCFactorSetShiftType(pc2,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 7152fb0b348SMark F. Adams ierr = KSPSetTolerances(k2[0],PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr); 71608e36f19SMark Adams ierr = KSPSetType(k2[0], KSPPREONLY);CHKERRQ(ierr); 7175b42dca8SJed Brown /* This flag gets reset by PCBJacobiGetSubKSP(), but our BJacobi really does the same algorithm everywhere (and in 7185b42dca8SJed Brown * fact, all but one process will have zero dofs), so we reset the flag to avoid having PCView_BJacobi attempt to 7195b42dca8SJed Brown * view every subdomain as though they were different. */ 7205b42dca8SJed Brown ((PC_BJacobi*)subpc->data)->same_local_solves = PETSC_TRUE; 721d5280255SMark F. Adams } 722d5280255SMark F. Adams 723d5280255SMark F. Adams /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 724d5280255SMark F. Adams ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 725e55864a3SBarry Smith ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr); 726d5280255SMark F. Adams ierr = PetscOptionsEnd();CHKERRQ(ierr); 7271c1aac46SBarry Smith if (!mg->galerkin) SETERRQ(comm,PETSC_ERR_USER,"PCGAMG must use Galerkin for coarse operators."); 7281c1aac46SBarry Smith if (mg->galerkin == 1) mg->galerkin = 2; 729d5280255SMark F. Adams 730d5280255SMark F. Adams /* clean up */ 731d5280255SMark F. Adams for (level=1; level<pc_gamg->Nlevels; level++) { 732587fa25dSMark F. Adams ierr = MatDestroy(&Parr[level]);CHKERRQ(ierr); 733587fa25dSMark F. Adams ierr = MatDestroy(&Aarr[level]);CHKERRQ(ierr); 7345b89ad90SMark F. Adams } 7350cbbd2e1SMark F. Adams 7360cbbd2e1SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 737806fa848SBarry Smith } else { 7385f8cf99dSMark F. Adams KSP smoother; 739302440fdSBarry Smith ierr = PetscInfo(pc,"One level solver used (system is seen as DD). Using default solver.\n");CHKERRQ(ierr); 7409d5b6da9SMark F. Adams ierr = PCMGGetSmoother(pc, 0, &smoother);CHKERRQ(ierr); 74123ee1639SBarry Smith ierr = KSPSetOperators(smoother, Aarr[0], Aarr[0]);CHKERRQ(ierr); 7425f8cf99dSMark F. Adams ierr = KSPSetType(smoother, KSPPREONLY);CHKERRQ(ierr); 7439d5b6da9SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 7445f8cf99dSMark F. Adams } 7455b89ad90SMark F. Adams PetscFunctionReturn(0); 7465b89ad90SMark F. Adams } 7475b89ad90SMark F. Adams 748eb07cef2SMark F. Adams /* ------------------------------------------------------------------------- */ 7495b89ad90SMark F. Adams /* 7505b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 7515b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 7525b89ad90SMark F. Adams 7535b89ad90SMark F. Adams Input Parameter: 7545b89ad90SMark F. Adams . pc - the preconditioner context 7555b89ad90SMark F. Adams 7565b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 7575b89ad90SMark F. Adams */ 7585b89ad90SMark F. Adams #undef __FUNCT__ 7595b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG" 7605b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc) 7615b89ad90SMark F. Adams { 7625b89ad90SMark F. Adams PetscErrorCode ierr; 7635b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 7645b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 7655b89ad90SMark F. Adams 7665b89ad90SMark F. Adams PetscFunctionBegin; 7675b89ad90SMark F. Adams ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 7689b8ffb57SJed Brown if (pc_gamg->ops->destroy) { 7699b8ffb57SJed Brown ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr); 7709b8ffb57SJed Brown } 77114a9496bSBarry Smith ierr = PetscRandomDestroy(&pc_gamg->random);CHKERRQ(ierr); 7721ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->ops);CHKERRQ(ierr); 7731ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr); 7745b89ad90SMark F. Adams ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 7755b89ad90SMark F. Adams ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 7765b89ad90SMark F. Adams PetscFunctionReturn(0); 7775b89ad90SMark F. Adams } 7785b89ad90SMark F. Adams 779676e1743SMark F. Adams 780676e1743SMark F. Adams #undef __FUNCT__ 781676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim" 782676e1743SMark F. Adams /*@ 7831cc46a46SBarry Smith PCGAMGSetProcEqLim - Set number of equations to aim for on coarse grids via processor reduction. 784676e1743SMark F. Adams 7851cc46a46SBarry Smith Logically Collective on PC 786676e1743SMark F. Adams 787676e1743SMark F. Adams Input Parameters: 7881cc46a46SBarry Smith + pc - the preconditioner context 7891cc46a46SBarry Smith - n - the number of equations 790676e1743SMark F. Adams 791676e1743SMark F. Adams 792676e1743SMark F. Adams Options Database Key: 7931cc46a46SBarry Smith . -pc_gamg_process_eq_limit <limit> 794676e1743SMark F. Adams 795676e1743SMark F. Adams Level: intermediate 796676e1743SMark F. Adams 7971c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 798676e1743SMark F. Adams 7991c1aac46SBarry Smith .seealso: PCGAMGSetCoarseEqLim() 800676e1743SMark F. Adams @*/ 801676e1743SMark F. Adams PetscErrorCode PCGAMGSetProcEqLim(PC pc, PetscInt n) 802676e1743SMark F. Adams { 803676e1743SMark F. Adams PetscErrorCode ierr; 804676e1743SMark F. Adams 805676e1743SMark F. Adams PetscFunctionBegin; 806676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 807676e1743SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetProcEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 808676e1743SMark F. Adams PetscFunctionReturn(0); 809676e1743SMark F. Adams } 810676e1743SMark F. Adams 811676e1743SMark F. Adams #undef __FUNCT__ 812676e1743SMark F. Adams #define __FUNCT__ "PCGAMGSetProcEqLim_GAMG" 8131e6b0712SBarry Smith static PetscErrorCode PCGAMGSetProcEqLim_GAMG(PC pc, PetscInt n) 814676e1743SMark F. Adams { 815c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 816c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 817676e1743SMark F. Adams 818676e1743SMark F. Adams PetscFunctionBegin; 8199d5b6da9SMark F. Adams if (n>0) pc_gamg->min_eq_proc = n; 820676e1743SMark F. Adams PetscFunctionReturn(0); 821676e1743SMark F. Adams } 822676e1743SMark F. Adams 823676e1743SMark F. Adams #undef __FUNCT__ 824389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim" 825389730f3SMark F. Adams /*@ 826389730f3SMark F. Adams PCGAMGSetCoarseEqLim - Set max number of equations on coarse grids. 827389730f3SMark F. Adams 828389730f3SMark F. Adams Collective on PC 829389730f3SMark F. Adams 830389730f3SMark F. Adams Input Parameters: 8311cc46a46SBarry Smith + pc - the preconditioner context 8321cc46a46SBarry Smith - n - maximum number of equations to aim for 833389730f3SMark F. Adams 834389730f3SMark F. Adams Options Database Key: 8351cc46a46SBarry Smith . -pc_gamg_coarse_eq_limit <limit> 836389730f3SMark F. Adams 837389730f3SMark F. Adams Level: intermediate 838389730f3SMark F. Adams 8391c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 840389730f3SMark F. Adams 8411c1aac46SBarry Smith .seealso: PCGAMGSetProcEqLim() 842389730f3SMark F. Adams @*/ 843389730f3SMark F. Adams PetscErrorCode PCGAMGSetCoarseEqLim(PC pc, PetscInt n) 844389730f3SMark F. Adams { 845389730f3SMark F. Adams PetscErrorCode ierr; 846389730f3SMark F. Adams 847389730f3SMark F. Adams PetscFunctionBegin; 848389730f3SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 849389730f3SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetCoarseEqLim_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 850389730f3SMark F. Adams PetscFunctionReturn(0); 851389730f3SMark F. Adams } 852389730f3SMark F. Adams 853389730f3SMark F. Adams #undef __FUNCT__ 854389730f3SMark F. Adams #define __FUNCT__ "PCGAMGSetCoarseEqLim_GAMG" 8551e6b0712SBarry Smith static PetscErrorCode PCGAMGSetCoarseEqLim_GAMG(PC pc, PetscInt n) 856389730f3SMark F. Adams { 857389730f3SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 858389730f3SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 859389730f3SMark F. Adams 860389730f3SMark F. Adams PetscFunctionBegin; 8619d5b6da9SMark F. Adams if (n>0) pc_gamg->coarse_eq_limit = n; 862389730f3SMark F. Adams PetscFunctionReturn(0); 863389730f3SMark F. Adams } 864389730f3SMark F. Adams 865389730f3SMark F. Adams #undef __FUNCT__ 8668263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning" 867676e1743SMark F. Adams /*@ 8688263b398SMark F. Adams PCGAMGSetRepartitioning - Repartition the coarse grids 869676e1743SMark F. Adams 870676e1743SMark F. Adams Collective on PC 871676e1743SMark F. Adams 872676e1743SMark F. Adams Input Parameters: 8731cc46a46SBarry Smith + pc - the preconditioner context 8741cc46a46SBarry Smith - n - PETSC_TRUE or PETSC_FALSE 875676e1743SMark F. Adams 876676e1743SMark F. Adams Options Database Key: 8771cc46a46SBarry Smith . -pc_gamg_repartition <true,false> 878676e1743SMark F. Adams 879676e1743SMark F. Adams Level: intermediate 880676e1743SMark F. Adams 8811c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 882676e1743SMark F. Adams 883676e1743SMark F. Adams .seealso: () 884676e1743SMark F. Adams @*/ 8858263b398SMark F. Adams PetscErrorCode PCGAMGSetRepartitioning(PC pc, PetscBool n) 886676e1743SMark F. Adams { 887676e1743SMark F. Adams PetscErrorCode ierr; 888676e1743SMark F. Adams 889676e1743SMark F. Adams PetscFunctionBegin; 890676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 8918263b398SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetRepartitioning_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 892676e1743SMark F. Adams PetscFunctionReturn(0); 893676e1743SMark F. Adams } 894676e1743SMark F. Adams 895676e1743SMark F. Adams #undef __FUNCT__ 8968263b398SMark F. Adams #define __FUNCT__ "PCGAMGSetRepartitioning_GAMG" 8971e6b0712SBarry Smith static PetscErrorCode PCGAMGSetRepartitioning_GAMG(PC pc, PetscBool n) 898676e1743SMark F. Adams { 899c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 900c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 901676e1743SMark F. Adams 902676e1743SMark F. Adams PetscFunctionBegin; 9039d5b6da9SMark F. Adams pc_gamg->repart = n; 904676e1743SMark F. Adams PetscFunctionReturn(0); 905676e1743SMark F. Adams } 906676e1743SMark F. Adams 907676e1743SMark F. Adams #undef __FUNCT__ 9081cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation" 909dfd5c07aSMark F. Adams /*@ 9101cc46a46SBarry Smith PCGAMGSetReuseInterpolation - Reuse prolongation when rebuilding preconditioner 911dfd5c07aSMark F. Adams 912dfd5c07aSMark F. Adams Collective on PC 913dfd5c07aSMark F. Adams 914dfd5c07aSMark F. Adams Input Parameters: 9151cc46a46SBarry Smith + pc - the preconditioner context 9161cc46a46SBarry Smith - n - PETSC_TRUE or PETSC_FALSE 917dfd5c07aSMark F. Adams 918dfd5c07aSMark F. Adams Options Database Key: 9191cc46a46SBarry Smith . -pc_gamg_reuse_interpolation <true,false> 920dfd5c07aSMark F. Adams 921dfd5c07aSMark F. Adams Level: intermediate 922dfd5c07aSMark F. Adams 9231c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 924dfd5c07aSMark F. Adams 925dfd5c07aSMark F. Adams .seealso: () 926dfd5c07aSMark F. Adams @*/ 9271cc46a46SBarry Smith PetscErrorCode PCGAMGSetReuseInterpolation(PC pc, PetscBool n) 928dfd5c07aSMark F. Adams { 929dfd5c07aSMark F. Adams PetscErrorCode ierr; 930dfd5c07aSMark F. Adams 931dfd5c07aSMark F. Adams PetscFunctionBegin; 932dfd5c07aSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9331cc46a46SBarry Smith ierr = PetscTryMethod(pc,"PCGAMGSetReuseInterpolation_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 934dfd5c07aSMark F. Adams PetscFunctionReturn(0); 935dfd5c07aSMark F. Adams } 936dfd5c07aSMark F. Adams 937dfd5c07aSMark F. Adams #undef __FUNCT__ 9381cc46a46SBarry Smith #define __FUNCT__ "PCGAMGSetReuseInterpolation_GAMG" 9391cc46a46SBarry Smith static PetscErrorCode PCGAMGSetReuseInterpolation_GAMG(PC pc, PetscBool n) 940dfd5c07aSMark F. Adams { 941dfd5c07aSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 942dfd5c07aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 943dfd5c07aSMark F. Adams 944dfd5c07aSMark F. Adams PetscFunctionBegin; 945dfd5c07aSMark F. Adams pc_gamg->reuse_prol = n; 946dfd5c07aSMark F. Adams PetscFunctionReturn(0); 947dfd5c07aSMark F. Adams } 948dfd5c07aSMark F. Adams 949dfd5c07aSMark F. Adams #undef __FUNCT__ 950ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs" 951ffc955d6SMark F. Adams /*@ 952ffc955d6SMark F. Adams PCGAMGSetUseASMAggs - 953ffc955d6SMark F. Adams 954ffc955d6SMark F. Adams Collective on PC 955ffc955d6SMark F. Adams 956ffc955d6SMark F. Adams Input Parameters: 957ffc955d6SMark F. Adams . pc - the preconditioner context 958ffc955d6SMark F. Adams 959ffc955d6SMark F. Adams 960ffc955d6SMark F. Adams Options Database Key: 961ffc955d6SMark F. Adams . -pc_gamg_use_agg_gasm 962ffc955d6SMark F. Adams 963ffc955d6SMark F. Adams Level: intermediate 964ffc955d6SMark F. Adams 9651c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 966ffc955d6SMark F. Adams 967ffc955d6SMark F. Adams .seealso: () 968ffc955d6SMark F. Adams @*/ 969ffc955d6SMark F. Adams PetscErrorCode PCGAMGSetUseASMAggs(PC pc, PetscBool n) 970ffc955d6SMark F. Adams { 971ffc955d6SMark F. Adams PetscErrorCode ierr; 972ffc955d6SMark F. Adams 973ffc955d6SMark F. Adams PetscFunctionBegin; 974ffc955d6SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 975ffc955d6SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetUseASMAggs_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 976ffc955d6SMark F. Adams PetscFunctionReturn(0); 977ffc955d6SMark F. Adams } 978ffc955d6SMark F. Adams 979ffc955d6SMark F. Adams #undef __FUNCT__ 980ffc955d6SMark F. Adams #define __FUNCT__ "PCGAMGSetUseASMAggs_GAMG" 9811e6b0712SBarry Smith static PetscErrorCode PCGAMGSetUseASMAggs_GAMG(PC pc, PetscBool n) 982ffc955d6SMark F. Adams { 983ffc955d6SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 984ffc955d6SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 985ffc955d6SMark F. Adams 986ffc955d6SMark F. Adams PetscFunctionBegin; 987ffc955d6SMark F. Adams pc_gamg->use_aggs_in_gasm = n; 988ffc955d6SMark F. Adams PetscFunctionReturn(0); 989ffc955d6SMark F. Adams } 990ffc955d6SMark F. Adams 991ffc955d6SMark F. Adams #undef __FUNCT__ 9924ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels" 9934ef23d27SMark F. Adams /*@ 9941cc46a46SBarry Smith PCGAMGSetNlevels - Sets the maximum number of levels PCGAMG will use 9954ef23d27SMark F. Adams 9964ef23d27SMark F. Adams Not collective on PC 9974ef23d27SMark F. Adams 9984ef23d27SMark F. Adams Input Parameters: 9991cc46a46SBarry Smith + pc - the preconditioner 10001cc46a46SBarry Smith - n - the maximum number of levels to use 10014ef23d27SMark F. Adams 10024ef23d27SMark F. Adams Options Database Key: 10034ef23d27SMark F. Adams . -pc_mg_levels 10044ef23d27SMark F. Adams 10054ef23d27SMark F. Adams Level: intermediate 10064ef23d27SMark F. Adams 10071c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 10084ef23d27SMark F. Adams 10094ef23d27SMark F. Adams .seealso: () 10104ef23d27SMark F. Adams @*/ 10114ef23d27SMark F. Adams PetscErrorCode PCGAMGSetNlevels(PC pc, PetscInt n) 10124ef23d27SMark F. Adams { 10134ef23d27SMark F. Adams PetscErrorCode ierr; 10144ef23d27SMark F. Adams 10154ef23d27SMark F. Adams PetscFunctionBegin; 10164ef23d27SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10174ef23d27SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetNlevels_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 10184ef23d27SMark F. Adams PetscFunctionReturn(0); 10194ef23d27SMark F. Adams } 10204ef23d27SMark F. Adams 10214ef23d27SMark F. Adams #undef __FUNCT__ 10224ef23d27SMark F. Adams #define __FUNCT__ "PCGAMGSetNlevels_GAMG" 10231e6b0712SBarry Smith static PetscErrorCode PCGAMGSetNlevels_GAMG(PC pc, PetscInt n) 10244ef23d27SMark F. Adams { 10254ef23d27SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 10264ef23d27SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 10274ef23d27SMark F. Adams 10284ef23d27SMark F. Adams PetscFunctionBegin; 10299d5b6da9SMark F. Adams pc_gamg->Nlevels = n; 10304ef23d27SMark F. Adams PetscFunctionReturn(0); 10314ef23d27SMark F. Adams } 10324ef23d27SMark F. Adams 10334ef23d27SMark F. Adams #undef __FUNCT__ 10343542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold" 10353542efc5SMark F. Adams /*@ 10363542efc5SMark F. Adams PCGAMGSetThreshold - Relative threshold to use for dropping edges in aggregation graph 10373542efc5SMark F. Adams 10383542efc5SMark F. Adams Not collective on PC 10393542efc5SMark F. Adams 10403542efc5SMark F. Adams Input Parameters: 10411cc46a46SBarry Smith + pc - the preconditioner context 1042b001cb0fSBarry Smith - threshold - the threshold value, 0.0 means keep all nonzero entries in the graph; negative means keep even zero entries in the graph 10433542efc5SMark F. Adams 10443542efc5SMark F. Adams Options Database Key: 10451cc46a46SBarry Smith . -pc_gamg_threshold <threshold> 10463542efc5SMark F. Adams 10473542efc5SMark F. Adams Level: intermediate 10483542efc5SMark F. Adams 10491c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 10503542efc5SMark F. Adams 10513542efc5SMark F. Adams .seealso: () 10523542efc5SMark F. Adams @*/ 10533542efc5SMark F. Adams PetscErrorCode PCGAMGSetThreshold(PC pc, PetscReal n) 10543542efc5SMark F. Adams { 10553542efc5SMark F. Adams PetscErrorCode ierr; 10563542efc5SMark F. Adams 10573542efc5SMark F. Adams PetscFunctionBegin; 10583542efc5SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10593542efc5SMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetThreshold_C",(PC,PetscReal),(pc,n));CHKERRQ(ierr); 10603542efc5SMark F. Adams PetscFunctionReturn(0); 10613542efc5SMark F. Adams } 10623542efc5SMark F. Adams 10633542efc5SMark F. Adams #undef __FUNCT__ 10643542efc5SMark F. Adams #define __FUNCT__ "PCGAMGSetThreshold_GAMG" 10651e6b0712SBarry Smith static PetscErrorCode PCGAMGSetThreshold_GAMG(PC pc, PetscReal n) 10663542efc5SMark F. Adams { 1067c20e4228SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1068c20e4228SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 10693542efc5SMark F. Adams 10703542efc5SMark F. Adams PetscFunctionBegin; 10719d5b6da9SMark F. Adams pc_gamg->threshold = n; 10723542efc5SMark F. Adams PetscFunctionReturn(0); 10733542efc5SMark F. Adams } 10743542efc5SMark F. Adams 10753542efc5SMark F. Adams #undef __FUNCT__ 10769d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType" 1077676e1743SMark F. Adams /*@ 1078c60c7ad4SBarry Smith PCGAMGSetType - Set solution method 1079676e1743SMark F. Adams 1080676e1743SMark F. Adams Collective on PC 1081676e1743SMark F. Adams 1082676e1743SMark F. Adams Input Parameters: 1083c60c7ad4SBarry Smith + pc - the preconditioner context 1084c60c7ad4SBarry Smith - type - PCGAMGAGG, PCGAMGGEO, or PCGAMGCLASSICAL 1085676e1743SMark F. Adams 1086676e1743SMark F. Adams Options Database Key: 1087c60c7ad4SBarry Smith . -pc_gamg_type <agg,geo,classical> 1088676e1743SMark F. Adams 1089676e1743SMark F. Adams Level: intermediate 1090676e1743SMark F. Adams 10911c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 1092676e1743SMark F. Adams 1093c60c7ad4SBarry Smith .seealso: PCGAMGGetType(), PCGAMG 1094676e1743SMark F. Adams @*/ 109519fd82e9SBarry Smith PetscErrorCode PCGAMGSetType(PC pc, PCGAMGType type) 1096676e1743SMark F. Adams { 1097676e1743SMark F. Adams PetscErrorCode ierr; 1098676e1743SMark F. Adams 1099676e1743SMark F. Adams PetscFunctionBegin; 1100676e1743SMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1101806fa848SBarry Smith ierr = PetscTryMethod(pc,"PCGAMGSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr); 1102676e1743SMark F. Adams PetscFunctionReturn(0); 1103676e1743SMark F. Adams } 1104676e1743SMark F. Adams 1105676e1743SMark F. Adams #undef __FUNCT__ 1106c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType" 1107c60c7ad4SBarry Smith /*@ 1108c60c7ad4SBarry Smith PCGAMGGetType - Get solution method 1109c60c7ad4SBarry Smith 1110c60c7ad4SBarry Smith Collective on PC 1111c60c7ad4SBarry Smith 1112c60c7ad4SBarry Smith Input Parameter: 1113c60c7ad4SBarry Smith . pc - the preconditioner context 1114c60c7ad4SBarry Smith 1115c60c7ad4SBarry Smith Output Parameter: 1116c60c7ad4SBarry Smith . type - the type of algorithm used 1117c60c7ad4SBarry Smith 1118c60c7ad4SBarry Smith Level: intermediate 1119c60c7ad4SBarry Smith 11201c1aac46SBarry Smith Concepts: Unstructured multigrid preconditioner 1121c60c7ad4SBarry Smith 11221c1aac46SBarry Smith .seealso: PCGAMGSetType(), PCGAMGType 1123c60c7ad4SBarry Smith @*/ 1124c60c7ad4SBarry Smith PetscErrorCode PCGAMGGetType(PC pc, PCGAMGType *type) 1125c60c7ad4SBarry Smith { 1126c60c7ad4SBarry Smith PetscErrorCode ierr; 1127c60c7ad4SBarry Smith 1128c60c7ad4SBarry Smith PetscFunctionBegin; 1129c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1130c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCGAMGGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr); 1131c60c7ad4SBarry Smith PetscFunctionReturn(0); 1132c60c7ad4SBarry Smith } 1133c60c7ad4SBarry Smith 1134c60c7ad4SBarry Smith #undef __FUNCT__ 1135c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGGetType_GAMG" 1136c60c7ad4SBarry Smith static PetscErrorCode PCGAMGGetType_GAMG(PC pc, PCGAMGType *type) 1137c60c7ad4SBarry Smith { 1138c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1139c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1140c60c7ad4SBarry Smith 1141c60c7ad4SBarry Smith PetscFunctionBegin; 1142c60c7ad4SBarry Smith *type = pc_gamg->type; 1143c60c7ad4SBarry Smith PetscFunctionReturn(0); 1144c60c7ad4SBarry Smith } 1145c60c7ad4SBarry Smith 1146c60c7ad4SBarry Smith #undef __FUNCT__ 11479d5b6da9SMark F. Adams #define __FUNCT__ "PCGAMGSetType_GAMG" 11481e6b0712SBarry Smith static PetscErrorCode PCGAMGSetType_GAMG(PC pc, PCGAMGType type) 1149676e1743SMark F. Adams { 11509d5b6da9SMark F. Adams PetscErrorCode ierr,(*r)(PC); 11511ab5ffc9SJed Brown PC_MG *mg = (PC_MG*)pc->data; 11521ab5ffc9SJed Brown PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1153676e1743SMark F. Adams 1154676e1743SMark F. Adams PetscFunctionBegin; 1155c60c7ad4SBarry Smith pc_gamg->type = type; 11561c9cd337SJed Brown ierr = PetscFunctionListFind(GAMGList,type,&r);CHKERRQ(ierr); 11579d5b6da9SMark F. Adams if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown GAMG type %s given",type); 11581ab5ffc9SJed Brown if (pc_gamg->ops->destroy) { 11591ab5ffc9SJed Brown ierr = (*pc_gamg->ops->destroy)(pc);CHKERRQ(ierr); 11601ab5ffc9SJed Brown ierr = PetscMemzero(pc_gamg->ops,sizeof(struct _PCGAMGOps));CHKERRQ(ierr); 1161e616c208SToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 11623ae0bb68SMark Adams /* cleaning up common data in pc_gamg - this should disapear someday */ 11633ae0bb68SMark Adams pc_gamg->data_cell_cols = 0; 11643ae0bb68SMark Adams pc_gamg->data_cell_rows = 0; 11653ae0bb68SMark Adams pc_gamg->orig_data_cell_cols = 0; 11663ae0bb68SMark Adams pc_gamg->orig_data_cell_rows = 0; 11673ae0bb68SMark Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 11683ae0bb68SMark Adams pc_gamg->data_sz = 0; 11691ab5ffc9SJed Brown } 11701ab5ffc9SJed Brown ierr = PetscFree(pc_gamg->gamg_type_name);CHKERRQ(ierr); 11711ab5ffc9SJed Brown ierr = PetscStrallocpy(type,&pc_gamg->gamg_type_name);CHKERRQ(ierr); 11729d5b6da9SMark F. Adams ierr = (*r)(pc);CHKERRQ(ierr); 1173676e1743SMark F. Adams PetscFunctionReturn(0); 1174676e1743SMark F. Adams } 1175676e1743SMark F. Adams 11765b89ad90SMark F. Adams #undef __FUNCT__ 11775adeb434SBarry Smith #define __FUNCT__ "PCView_GAMG" 11785adeb434SBarry Smith static PetscErrorCode PCView_GAMG(PC pc,PetscViewer viewer) 11795adeb434SBarry Smith { 11805adeb434SBarry Smith PetscErrorCode ierr; 11815adeb434SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 11825adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 11835adeb434SBarry Smith 11845adeb434SBarry Smith PetscFunctionBegin; 11855adeb434SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," GAMG specific options\n");CHKERRQ(ierr); 1186b001cb0fSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Threshold for dropping small values from graph %g\n",(double)pc_gamg->threshold);CHKERRQ(ierr); 11875adeb434SBarry Smith if (pc_gamg->ops->view) { 11885adeb434SBarry Smith ierr = (*pc_gamg->ops->view)(pc,viewer);CHKERRQ(ierr); 11895adeb434SBarry Smith } 11905adeb434SBarry Smith PetscFunctionReturn(0); 11915adeb434SBarry Smith } 11925adeb434SBarry Smith 11935adeb434SBarry Smith #undef __FUNCT__ 11945b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG" 11954416b707SBarry Smith PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems *PetscOptionsObject,PC pc) 11965b89ad90SMark F. Adams { 1197676e1743SMark F. Adams PetscErrorCode ierr; 1198676e1743SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1199676e1743SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1200676e1743SMark F. Adams PetscBool flag; 12013b4367a7SBarry Smith MPI_Comm comm; 120214a9496bSBarry Smith char prefix[256]; 120314a9496bSBarry Smith const char *pcpre; 12045b89ad90SMark F. Adams 12055b89ad90SMark F. Adams PetscFunctionBegin; 12063b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1207e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"GAMG options");CHKERRQ(ierr); 1208676e1743SMark F. Adams { 1209bd94a7aaSJed Brown char tname[256]; 12101a1c1e04SBarry Smith ierr = PetscOptionsFList("-pc_gamg_type","Type of AMG method","PCGAMGSetType",GAMGList, pc_gamg->gamg_type_name, tname, sizeof(tname), &flag);CHKERRQ(ierr); 1211bd94a7aaSJed Brown if (flag) { 1212bd94a7aaSJed Brown ierr = PCGAMGSetType(pc,tname);CHKERRQ(ierr); 12131ab5ffc9SJed Brown } 121494ae4db5SBarry Smith ierr = PetscOptionsBool("-pc_gamg_repartition","Repartion coarse grids","PCGAMGRepartitioning",pc_gamg->repart,&pc_gamg->repart,NULL);CHKERRQ(ierr); 12151cc46a46SBarry Smith ierr = PetscOptionsBool("-pc_gamg_reuse_interpolation","Reuse prolongation operator","PCGAMGReuseInterpolation",pc_gamg->reuse_prol,&pc_gamg->reuse_prol,NULL);CHKERRQ(ierr); 121694ae4db5SBarry 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); 121794ae4db5SBarry 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); 121894ae4db5SBarry 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); 121994ae4db5SBarry 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); 122094ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Set number of MG levels","PCGAMGSetNlevels",pc_gamg->Nlevels,&pc_gamg->Nlevels,NULL);CHKERRQ(ierr); 1221b7cbab4eSMark Adams 1222b7cbab4eSMark Adams /* set options for subtype */ 1223e55864a3SBarry Smith if (pc_gamg->ops->setfromoptions) {ierr = (*pc_gamg->ops->setfromoptions)(PetscOptionsObject,pc);CHKERRQ(ierr);} 1224676e1743SMark F. Adams } 122514a9496bSBarry Smith ierr = PCGetOptionsPrefix(pc, &pcpre);CHKERRQ(ierr); 122614a9496bSBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%spc_gamg_",pcpre ? pcpre : "");CHKERRQ(ierr); 122714a9496bSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)pc_gamg->random,prefix);CHKERRQ(ierr); 122814a9496bSBarry Smith ierr = PetscRandomSetFromOptions(pc_gamg->random);CHKERRQ(ierr); 1229676e1743SMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 12305b89ad90SMark F. Adams PetscFunctionReturn(0); 12315b89ad90SMark F. Adams } 12325b89ad90SMark F. Adams 12335b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 12345b89ad90SMark F. Adams /*MC 12351cc46a46SBarry Smith PCGAMG - Geometric algebraic multigrid (AMG) preconditioner 12365b89ad90SMark F. Adams 1237280d9858SJed Brown Options Database Keys: 12385b89ad90SMark F. Adams Multigrid options(inherited) 12391cc46a46SBarry Smith + -pc_mg_cycles <v>: v or w (PCMGSetCycleType()) 1240280d9858SJed Brown . -pc_mg_smoothup <1>: Number of post-smoothing steps (PCMGSetNumberSmoothUp) 1241280d9858SJed Brown . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (PCMGSetNumberSmoothDown) 12428c1c2452SJed Brown - -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade 12435b89ad90SMark F. Adams 12441cc46a46SBarry Smith 12451cc46a46SBarry Smith Notes: In order to obtain good performance for PCGAMG for vector valued problems you must 12461cc46a46SBarry Smith $ Call MatSetBlockSize() to indicate the number of degrees of freedom per grid point 12471cc46a46SBarry Smith $ Call MatSetNearNullSpace() (or PCSetCoordinates() if solving the equations of elasticity) to indicate the near null space of the operator 12481cc46a46SBarry Smith $ See the Users Manual Chapter 4 for more details 12491cc46a46SBarry Smith 12505b89ad90SMark F. Adams Level: intermediate 1251280d9858SJed Brown 12521cc46a46SBarry Smith Concepts: algebraic multigrid 12535b89ad90SMark F. Adams 12541cc46a46SBarry Smith .seealso: PCCreate(), PCSetType(), MatSetBlockSize(), PCMGType, PCSetCoordinates(), MatSetNearNullSpace(), PCGAMGSetType(), PCGAMGAGG, PCGAMGGEO, PCGAMGCLASSICAL, PCGAMGSetProcEqLim(), 12551cc46a46SBarry Smith PCGAMGSetCoarseEqLim(), PCGAMGSetRepartitioning(), PCGAMGRegister(), PCGAMGSetReuseInterpolation(), PCGAMGSetUseASMAggs(), PCGAMGSetNlevels(), PCGAMGSetThreshold(), PCGAMGGetType() 12565b89ad90SMark F. Adams M*/ 1257b2573a8aSBarry Smith 12585b89ad90SMark F. Adams #undef __FUNCT__ 12595b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG" 12608cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GAMG(PC pc) 12615b89ad90SMark F. Adams { 12625b89ad90SMark F. Adams PetscErrorCode ierr; 12635b89ad90SMark F. Adams PC_GAMG *pc_gamg; 12645b89ad90SMark F. Adams PC_MG *mg; 12655b89ad90SMark F. Adams 12665b89ad90SMark F. Adams PetscFunctionBegin; 12671c1aac46SBarry Smith /* register AMG type */ 12681c1aac46SBarry Smith ierr = PCGAMGInitializePackage();CHKERRQ(ierr); 12691c1aac46SBarry Smith 12705b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 12711c1aac46SBarry Smith ierr = PCSetType(pc, PCMG);CHKERRQ(ierr); 12725b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName((PetscObject)pc, PCGAMG);CHKERRQ(ierr); 12735b89ad90SMark F. Adams 12745b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 1275b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg);CHKERRQ(ierr); 12765b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 12771c1aac46SBarry Smith mg->galerkin = 2; /* Use Galerkin, but it is computed externally from PCMG by GAMG code */ 12785b89ad90SMark F. Adams mg->innerctx = pc_gamg; 12795b89ad90SMark F. Adams 1280b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg->ops);CHKERRQ(ierr); 12811ab5ffc9SJed Brown 12829d5b6da9SMark F. Adams pc_gamg->setup_count = 0; 12839d5b6da9SMark F. Adams /* these should be in subctx but repartitioning needs simple arrays */ 12849d5b6da9SMark F. Adams pc_gamg->data_sz = 0; 12859d5b6da9SMark F. Adams pc_gamg->data = 0; 12865b89ad90SMark F. Adams 12879d5b6da9SMark F. Adams /* overwrite the pointers of PCMG by the functions of base class PCGAMG */ 12885b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 12895b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 12905b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 12915b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 12925adeb434SBarry Smith mg->view = PCView_GAMG; 12935b89ad90SMark F. Adams 1294bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetProcEqLim_C",PCGAMGSetProcEqLim_GAMG);CHKERRQ(ierr); 1295bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetCoarseEqLim_C",PCGAMGSetCoarseEqLim_GAMG);CHKERRQ(ierr); 1296bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetRepartitioning_C",PCGAMGSetRepartitioning_GAMG);CHKERRQ(ierr); 12971cc46a46SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetReuseInterpolation_C",PCGAMGSetReuseInterpolation_GAMG);CHKERRQ(ierr); 1298bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetUseASMAggs_C",PCGAMGSetUseASMAggs_GAMG);CHKERRQ(ierr); 1299bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetThreshold_C",PCGAMGSetThreshold_GAMG);CHKERRQ(ierr); 1300bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetType_C",PCGAMGSetType_GAMG);CHKERRQ(ierr); 1301c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGGetType_C",PCGAMGGetType_GAMG);CHKERRQ(ierr); 1302bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNlevels_C",PCGAMGSetNlevels_GAMG);CHKERRQ(ierr); 13039d5b6da9SMark F. Adams pc_gamg->repart = PETSC_FALSE; 1304d3042614SMark Adams pc_gamg->reuse_prol = PETSC_FALSE; 1305ffc955d6SMark F. Adams pc_gamg->use_aggs_in_gasm = PETSC_FALSE; 1306038f3aa4SMark F. Adams pc_gamg->min_eq_proc = 50; 130725a145a7SMark Adams pc_gamg->coarse_eq_limit = 50; 1308d3042614SMark Adams pc_gamg->threshold = 0.; 13099d5b6da9SMark F. Adams pc_gamg->Nlevels = GAMG_MAXLEVELS; 13109ab59c8bSMark Adams pc_gamg->current_level = 0; /* don't need to init really */ 1311c238b0ebSToby Isaac pc_gamg->ops->createlevel = PCGAMGCreateLevel_GAMG; 13129d5b6da9SMark F. Adams 131314a9496bSBarry Smith ierr = PetscRandomCreate(PetscObjectComm((PetscObject)pc),&pc_gamg->random);CHKERRQ(ierr); 131414a9496bSBarry Smith 1315bd94a7aaSJed Brown /* PCSetUp_GAMG assumes that the type has been set, so set it to the default now */ 1316bd94a7aaSJed Brown ierr = PCGAMGSetType(pc,PCGAMGAGG);CHKERRQ(ierr); 13175b89ad90SMark F. Adams PetscFunctionReturn(0); 13185b89ad90SMark F. Adams } 13193e3471ccSMark Adams 13203e3471ccSMark Adams #undef __FUNCT__ 13213e3471ccSMark Adams #define __FUNCT__ "PCGAMGInitializePackage" 13223e3471ccSMark Adams /*@C 13233e3471ccSMark Adams PCGAMGInitializePackage - This function initializes everything in the PCGAMG package. It is called 13243e3471ccSMark Adams from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to PCCreate_GAMG() 13253e3471ccSMark Adams when using static libraries. 13263e3471ccSMark Adams 13273e3471ccSMark Adams Level: developer 13283e3471ccSMark Adams 13293e3471ccSMark Adams .keywords: PC, PCGAMG, initialize, package 13303e3471ccSMark Adams .seealso: PetscInitialize() 13313e3471ccSMark Adams @*/ 13323e3471ccSMark Adams PetscErrorCode PCGAMGInitializePackage(void) 13333e3471ccSMark Adams { 13343e3471ccSMark Adams PetscErrorCode ierr; 13353e3471ccSMark Adams 13363e3471ccSMark Adams PetscFunctionBegin; 13373e3471ccSMark Adams if (PCGAMGPackageInitialized) PetscFunctionReturn(0); 13383e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_TRUE; 13393e3471ccSMark Adams ierr = PetscFunctionListAdd(&GAMGList,PCGAMGGEO,PCCreateGAMG_GEO);CHKERRQ(ierr); 13403e3471ccSMark Adams ierr = PetscFunctionListAdd(&GAMGList,PCGAMGAGG,PCCreateGAMG_AGG);CHKERRQ(ierr); 13418e6d0c30SPeter Brune ierr = PetscFunctionListAdd(&GAMGList,PCGAMGCLASSICAL,PCCreateGAMG_Classical);CHKERRQ(ierr); 13423e3471ccSMark Adams ierr = PetscRegisterFinalize(PCGAMGFinalizePackage);CHKERRQ(ierr); 1343c1c463dbSMark Adams 1344c1c463dbSMark Adams /* general events */ 1345fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGGraph_AGG", 0, &PC_GAMGGraph_AGG);CHKERRQ(ierr); 1346fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGGraph_GEO", PC_CLASSID, &PC_GAMGGraph_GEO);CHKERRQ(ierr); 1347fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGCoarse_AGG", PC_CLASSID, &PC_GAMGCoarsen_AGG);CHKERRQ(ierr); 1348fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGCoarse_GEO", PC_CLASSID, &PC_GAMGCoarsen_GEO);CHKERRQ(ierr); 1349c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGProl_AGG", PC_CLASSID, &PC_GAMGProlongator_AGG);CHKERRQ(ierr); 1350c1c463dbSMark Adams ierr = PetscLogEventRegister("PCGAMGProl_GEO", PC_CLASSID, &PC_GAMGProlongator_GEO);CHKERRQ(ierr); 1351fd1112cbSBarry Smith ierr = PetscLogEventRegister("PCGAMGPOpt_AGG", PC_CLASSID, &PC_GAMGOptProlongator_AGG);CHKERRQ(ierr); 1352c1c463dbSMark Adams 13535b89ad90SMark F. Adams #if defined PETSC_GAMG_USE_LOG 13545b89ad90SMark F. Adams ierr = PetscLogEventRegister("GAMG: createProl", PC_CLASSID, &petsc_gamg_setup_events[SET1]);CHKERRQ(ierr); 13555b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Graph", PC_CLASSID, &petsc_gamg_setup_events[GRAPH]);CHKERRQ(ierr); 13565b89ad90SMark F. Adams /* PetscLogEventRegister(" G.Mat", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_MAT]); */ 13575b89ad90SMark F. Adams /* PetscLogEventRegister(" G.Filter", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_FILTER]); */ 13585b89ad90SMark F. Adams /* PetscLogEventRegister(" G.Square", PC_CLASSID, &petsc_gamg_setup_events[GRAPH_SQR]); */ 13595b89ad90SMark F. Adams ierr = PetscLogEventRegister(" MIS/Agg", PC_CLASSID, &petsc_gamg_setup_events[SET4]);CHKERRQ(ierr); 13605b89ad90SMark F. Adams ierr = PetscLogEventRegister(" geo: growSupp", PC_CLASSID, &petsc_gamg_setup_events[SET5]);CHKERRQ(ierr); 13615b89ad90SMark F. Adams ierr = PetscLogEventRegister(" geo: triangle", PC_CLASSID, &petsc_gamg_setup_events[SET6]);CHKERRQ(ierr); 13625b89ad90SMark F. Adams ierr = PetscLogEventRegister(" search&set", PC_CLASSID, &petsc_gamg_setup_events[FIND_V]);CHKERRQ(ierr); 13635b89ad90SMark F. Adams ierr = PetscLogEventRegister(" SA: col data", PC_CLASSID, &petsc_gamg_setup_events[SET7]);CHKERRQ(ierr); 13645b89ad90SMark F. Adams ierr = PetscLogEventRegister(" SA: frmProl0", PC_CLASSID, &petsc_gamg_setup_events[SET8]);CHKERRQ(ierr); 13655b89ad90SMark F. Adams ierr = PetscLogEventRegister(" SA: smooth", PC_CLASSID, &petsc_gamg_setup_events[SET9]);CHKERRQ(ierr); 13665b89ad90SMark F. Adams ierr = PetscLogEventRegister("GAMG: partLevel", PC_CLASSID, &petsc_gamg_setup_events[SET2]);CHKERRQ(ierr); 13675b89ad90SMark F. Adams ierr = PetscLogEventRegister(" repartition", PC_CLASSID, &petsc_gamg_setup_events[SET12]);CHKERRQ(ierr); 13685b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Invert-Sort", PC_CLASSID, &petsc_gamg_setup_events[SET13]);CHKERRQ(ierr); 13695b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Move A", PC_CLASSID, &petsc_gamg_setup_events[SET14]);CHKERRQ(ierr); 13705b89ad90SMark F. Adams ierr = PetscLogEventRegister(" Move P", PC_CLASSID, &petsc_gamg_setup_events[SET15]);CHKERRQ(ierr); 13715b89ad90SMark F. Adams 13725b89ad90SMark F. Adams /* PetscLogEventRegister(" PL move data", PC_CLASSID, &petsc_gamg_setup_events[SET13]); */ 13735b89ad90SMark F. Adams /* PetscLogEventRegister("GAMG: fix", PC_CLASSID, &petsc_gamg_setup_events[SET10]); */ 13745b89ad90SMark F. Adams /* PetscLogEventRegister("GAMG: set levels", PC_CLASSID, &petsc_gamg_setup_events[SET11]); */ 13755b89ad90SMark F. Adams /* create timer stages */ 13765b89ad90SMark F. Adams #if defined GAMG_STAGES 13775b89ad90SMark F. Adams { 13785b89ad90SMark F. Adams char str[32]; 13795b89ad90SMark F. Adams PetscInt lidx; 13805b89ad90SMark F. Adams sprintf(str,"MG Level %d (finest)",0); 13815b89ad90SMark F. Adams ierr = PetscLogStageRegister(str, &gamg_stages[0]);CHKERRQ(ierr); 13825b89ad90SMark F. Adams for (lidx=1; lidx<9; lidx++) { 13835b89ad90SMark F. Adams sprintf(str,"MG Level %d",lidx); 13845b89ad90SMark F. Adams ierr = PetscLogStageRegister(str, &gamg_stages[lidx]);CHKERRQ(ierr); 13855b89ad90SMark F. Adams } 13865b89ad90SMark F. Adams } 13875b89ad90SMark F. Adams #endif 13885b89ad90SMark F. Adams #endif 13893e3471ccSMark Adams PetscFunctionReturn(0); 13903e3471ccSMark Adams } 13913e3471ccSMark Adams 13923e3471ccSMark Adams #undef __FUNCT__ 13933e3471ccSMark Adams #define __FUNCT__ "PCGAMGFinalizePackage" 13943e3471ccSMark Adams /*@C 13951c1aac46SBarry Smith PCGAMGFinalizePackage - This function frees everything from the PCGAMG package. It is 13961c1aac46SBarry Smith called from PetscFinalize() automatically. 13973e3471ccSMark Adams 13983e3471ccSMark Adams Level: developer 13993e3471ccSMark Adams 14003e3471ccSMark Adams .keywords: Petsc, destroy, package 14013e3471ccSMark Adams .seealso: PetscFinalize() 14023e3471ccSMark Adams @*/ 14033e3471ccSMark Adams PetscErrorCode PCGAMGFinalizePackage(void) 14043e3471ccSMark Adams { 14053e3471ccSMark Adams PetscErrorCode ierr; 14063e3471ccSMark Adams 14073e3471ccSMark Adams PetscFunctionBegin; 14083e3471ccSMark Adams PCGAMGPackageInitialized = PETSC_FALSE; 14093e3471ccSMark Adams ierr = PetscFunctionListDestroy(&GAMGList);CHKERRQ(ierr); 14103e3471ccSMark Adams PetscFunctionReturn(0); 14113e3471ccSMark Adams } 1412a36cf38bSToby Isaac 1413a36cf38bSToby Isaac #undef __FUNCT__ 1414a36cf38bSToby Isaac #define __FUNCT__ "PCGAMGRegister" 1415a36cf38bSToby Isaac /*@C 1416a36cf38bSToby Isaac PCGAMGRegister - Register a PCGAMG implementation. 1417a36cf38bSToby Isaac 1418a36cf38bSToby Isaac Input Parameters: 1419a36cf38bSToby Isaac + type - string that will be used as the name of the GAMG type. 1420a36cf38bSToby Isaac - create - function for creating the gamg context. 1421a36cf38bSToby Isaac 1422a36cf38bSToby Isaac Level: advanced 1423a36cf38bSToby Isaac 14241c1aac46SBarry Smith .seealso: PCGAMGType, PCGAMG, PCGAMGSetType() 1425a36cf38bSToby Isaac @*/ 1426a36cf38bSToby Isaac PetscErrorCode PCGAMGRegister(PCGAMGType type, PetscErrorCode (*create)(PC)) 1427a36cf38bSToby Isaac { 1428a36cf38bSToby Isaac PetscErrorCode ierr; 1429a36cf38bSToby Isaac 1430a36cf38bSToby Isaac PetscFunctionBegin; 1431a36cf38bSToby Isaac ierr = PCGAMGInitializePackage();CHKERRQ(ierr); 1432a36cf38bSToby Isaac ierr = PetscFunctionListAdd(&GAMGList,type,create);CHKERRQ(ierr); 1433a36cf38bSToby Isaac PetscFunctionReturn(0); 1434a36cf38bSToby Isaac } 1435a36cf38bSToby Isaac 1436