1*42e5ec60SJeff-Hadley #include <petscksp.h> 2*42e5ec60SJeff-Hadley #include <petscpc.h> 3*42e5ec60SJeff-Hadley #include <petscviewer.h> 4*42e5ec60SJeff-Hadley 5*42e5ec60SJeff-Hadley typedef struct { 6*42e5ec60SJeff-Hadley PetscInt num_levels; 7*42e5ec60SJeff-Hadley PetscInt *n_per_level; 8*42e5ec60SJeff-Hadley Mat stiff; 9*42e5ec60SJeff-Hadley Mat *ProlongationOps; 10*42e5ec60SJeff-Hadley PetscBT *CFMarkers; 11*42e5ec60SJeff-Hadley KSP kspHypre; 12*42e5ec60SJeff-Hadley } *DataCompression; 13*42e5ec60SJeff-Hadley 14*42e5ec60SJeff-Hadley PetscErrorCode Create1dLaplacian(PetscInt, Mat *); 15*42e5ec60SJeff-Hadley PetscErrorCode DataCompExportMats(DataCompression); 16*42e5ec60SJeff-Hadley PetscErrorCode DataCompDestroy(DataCompression); 17*42e5ec60SJeff-Hadley 18*42e5ec60SJeff-Hadley int main(int Argc, char **Args) 19*42e5ec60SJeff-Hadley { 20*42e5ec60SJeff-Hadley PetscInt n_nodes = 33; 21*42e5ec60SJeff-Hadley Vec x, b; 22*42e5ec60SJeff-Hadley PC pcHypre; 23*42e5ec60SJeff-Hadley DataCompression data_comp; 24*42e5ec60SJeff-Hadley 25*42e5ec60SJeff-Hadley PetscFunctionBeginUser; 26*42e5ec60SJeff-Hadley PetscCall(PetscInitialize(&Argc, &Args, NULL, NULL)); 27*42e5ec60SJeff-Hadley 28*42e5ec60SJeff-Hadley PetscCall(PetscNew(&data_comp)); 29*42e5ec60SJeff-Hadley 30*42e5ec60SJeff-Hadley // Creating stiffness matrix 31*42e5ec60SJeff-Hadley PetscCall(Create1dLaplacian(n_nodes, &data_comp->stiff)); 32*42e5ec60SJeff-Hadley PetscCall(PetscObjectSetName((PetscObject)data_comp->stiff, "Stiffness")); 33*42e5ec60SJeff-Hadley 34*42e5ec60SJeff-Hadley // Set-up BoomerAMG PC to get Prolongation Operators and Coarse/Fine splittings 35*42e5ec60SJeff-Hadley PetscCall(KSPCreate(PETSC_COMM_WORLD, &data_comp->kspHypre)); 36*42e5ec60SJeff-Hadley PetscCall(KSPSetType(data_comp->kspHypre, KSPRICHARDSON)); 37*42e5ec60SJeff-Hadley PetscCall(KSPGetPC(data_comp->kspHypre, &pcHypre)); 38*42e5ec60SJeff-Hadley PetscCall(PCSetType(pcHypre, PCHYPRE)); 39*42e5ec60SJeff-Hadley PetscCall(PCHYPRESetType(pcHypre, "boomeramg")); 40*42e5ec60SJeff-Hadley PetscCall(PCSetFromOptions(pcHypre)); 41*42e5ec60SJeff-Hadley PetscCall(PCSetOperators(pcHypre, data_comp->stiff, data_comp->stiff)); 42*42e5ec60SJeff-Hadley PetscCall(PCSetUp(pcHypre)); 43*42e5ec60SJeff-Hadley 44*42e5ec60SJeff-Hadley PetscCall(MatCreateVecs(data_comp->stiff, &x, &b)); 45*42e5ec60SJeff-Hadley PetscCall(PCApply(pcHypre, x, b)); 46*42e5ec60SJeff-Hadley PetscCall(VecDestroy(&x)); 47*42e5ec60SJeff-Hadley PetscCall(VecDestroy(&b)); 48*42e5ec60SJeff-Hadley 49*42e5ec60SJeff-Hadley //Viewing the PC and Extracting the Prolongation Operators and CFMarkers 50*42e5ec60SJeff-Hadley PetscCall(PCView(pcHypre, NULL)); 51*42e5ec60SJeff-Hadley PetscCall(PCGetInterpolations(pcHypre, &data_comp->num_levels, &data_comp->ProlongationOps)); 52*42e5ec60SJeff-Hadley PetscCall(PCHYPREGetCFMarkers(pcHypre, &data_comp->n_per_level, &data_comp->CFMarkers)); 53*42e5ec60SJeff-Hadley 54*42e5ec60SJeff-Hadley PetscCall(DataCompExportMats(data_comp)); 55*42e5ec60SJeff-Hadley PetscCall(DataCompDestroy(data_comp)); 56*42e5ec60SJeff-Hadley PetscCall(PetscFinalize()); 57*42e5ec60SJeff-Hadley return 0; 58*42e5ec60SJeff-Hadley } 59*42e5ec60SJeff-Hadley 60*42e5ec60SJeff-Hadley PetscErrorCode Create1dLaplacian(PetscInt n, Mat *mat) 61*42e5ec60SJeff-Hadley { 62*42e5ec60SJeff-Hadley PetscFunctionBeginUser; 63*42e5ec60SJeff-Hadley PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, n, n, 3, NULL, mat)); 64*42e5ec60SJeff-Hadley PetscCall(MatSetValue(*mat, n - 1, n - 1, 2.0, INSERT_VALUES)); 65*42e5ec60SJeff-Hadley for (PetscInt i = 0; i < n - 1; i++) { 66*42e5ec60SJeff-Hadley PetscCall(MatSetValue(*mat, i, i, 2.0, INSERT_VALUES)); 67*42e5ec60SJeff-Hadley PetscCall(MatSetValue(*mat, i + 1, i, -1.0, INSERT_VALUES)); 68*42e5ec60SJeff-Hadley PetscCall(MatSetValue(*mat, i, i + 1, -1.0, INSERT_VALUES)); 69*42e5ec60SJeff-Hadley } 70*42e5ec60SJeff-Hadley PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY)); 71*42e5ec60SJeff-Hadley PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY)); 72*42e5ec60SJeff-Hadley PetscFunctionReturn(PETSC_SUCCESS); 73*42e5ec60SJeff-Hadley } 74*42e5ec60SJeff-Hadley 75*42e5ec60SJeff-Hadley PetscErrorCode DataCompExportMats(DataCompression data_comp) 76*42e5ec60SJeff-Hadley { 77*42e5ec60SJeff-Hadley PetscFunctionBeginUser; 78*42e5ec60SJeff-Hadley PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Num levels: %" PetscInt_FMT "\n", data_comp->num_levels)); 79*42e5ec60SJeff-Hadley PetscCall(PetscPrintf(PETSC_COMM_WORLD, " -- Nodes per level --\n")); 80*42e5ec60SJeff-Hadley for (PetscInt i = 0; i < data_comp->num_levels; i++) { PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Level %" PetscInt_FMT ": %" PetscInt_FMT "\n", i, data_comp->n_per_level[i])); } 81*42e5ec60SJeff-Hadley 82*42e5ec60SJeff-Hadley for (PetscInt i = 0; i < data_comp->num_levels - 1; i++) { 83*42e5ec60SJeff-Hadley PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Prolongation Operator - Level %" PetscInt_FMT "\n", i)); 84*42e5ec60SJeff-Hadley PetscCall(PetscObjectSetName((PetscObject)data_comp->ProlongationOps[i], "P")); 85*42e5ec60SJeff-Hadley PetscCall(MatView(data_comp->ProlongationOps[i], PETSC_VIEWER_STDOUT_WORLD)); 86*42e5ec60SJeff-Hadley PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n")); 87*42e5ec60SJeff-Hadley } 88*42e5ec60SJeff-Hadley 89*42e5ec60SJeff-Hadley for (PetscInt i = 0; i < data_comp->num_levels - 1; i++) { 90*42e5ec60SJeff-Hadley PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Coarse/Fine splitting - Level %" PetscInt_FMT "\n", i + 1)); 91*42e5ec60SJeff-Hadley PetscCall(PetscBTView(data_comp->n_per_level[i + 1], data_comp->CFMarkers[i], PETSC_VIEWER_STDOUT_WORLD)); 92*42e5ec60SJeff-Hadley } 93*42e5ec60SJeff-Hadley PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Stiffness matrix, sparse format:\n")); 94*42e5ec60SJeff-Hadley PetscCall(MatViewFromOptions(data_comp->stiff, NULL, "-mat_view_stiff")); 95*42e5ec60SJeff-Hadley PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Finished calling the Viewer functions\n")); 96*42e5ec60SJeff-Hadley PetscFunctionReturn(PETSC_SUCCESS); 97*42e5ec60SJeff-Hadley } 98*42e5ec60SJeff-Hadley 99*42e5ec60SJeff-Hadley PetscErrorCode DataCompDestroy(DataCompression data_comp) 100*42e5ec60SJeff-Hadley { 101*42e5ec60SJeff-Hadley PetscFunctionBeginUser; 102*42e5ec60SJeff-Hadley if (data_comp == NULL) PetscFunctionReturn(PETSC_SUCCESS); 103*42e5ec60SJeff-Hadley PetscCall(MatDestroy(&data_comp->stiff)); 104*42e5ec60SJeff-Hadley PetscCall(KSPDestroy(&data_comp->kspHypre)); 105*42e5ec60SJeff-Hadley for (PetscInt i = 0; i < data_comp->num_levels - 1; i++) { 106*42e5ec60SJeff-Hadley PetscCall(MatDestroy(&data_comp->ProlongationOps[i])); 107*42e5ec60SJeff-Hadley PetscCall(PetscBTDestroy(&data_comp->CFMarkers[i])); 108*42e5ec60SJeff-Hadley } 109*42e5ec60SJeff-Hadley PetscCall(PetscFree(data_comp->ProlongationOps)); 110*42e5ec60SJeff-Hadley PetscCall(PetscFree(data_comp->n_per_level)); 111*42e5ec60SJeff-Hadley PetscCall(PetscFree(data_comp->CFMarkers)); 112*42e5ec60SJeff-Hadley PetscCall(PetscFree(data_comp)); 113*42e5ec60SJeff-Hadley PetscFunctionReturn(PETSC_SUCCESS); 114*42e5ec60SJeff-Hadley } 115*42e5ec60SJeff-Hadley 116*42e5ec60SJeff-Hadley /*TEST 117*42e5ec60SJeff-Hadley 118*42e5ec60SJeff-Hadley test: 119*42e5ec60SJeff-Hadley requires: hypre !complex !defined(PETSC_HAVE_HYPRE_DEVICE) 120*42e5ec60SJeff-Hadley args: -pc_hypre_boomeramg_coarsen_type modifiedRuge-Stueben -pc_hypre_boomeramg_interp_type classical -pc_hypre_boomeramg_strong_threshold 0.25 pc_hypre_boomeramg_numfunctions 1 -pc_hypre_boomeramg_max_row_sum 1.0 -mat_view_stiff 121*42e5ec60SJeff-Hadley 122*42e5ec60SJeff-Hadley TEST*/ 123