xref: /petsc/src/ksp/pc/tests/ex10.c (revision 834855d6effb0d027771461c8e947ee1ce5a1e17)
142e5ec60SJeff-Hadley #include <petscksp.h>
242e5ec60SJeff-Hadley #include <petscpc.h>
342e5ec60SJeff-Hadley #include <petscviewer.h>
442e5ec60SJeff-Hadley 
542e5ec60SJeff-Hadley typedef struct {
642e5ec60SJeff-Hadley   PetscInt  num_levels;
742e5ec60SJeff-Hadley   PetscInt *n_per_level;
842e5ec60SJeff-Hadley   Mat       stiff;
942e5ec60SJeff-Hadley   Mat      *ProlongationOps;
1042e5ec60SJeff-Hadley   PetscBT  *CFMarkers;
1142e5ec60SJeff-Hadley   KSP       kspHypre;
1242e5ec60SJeff-Hadley } *DataCompression;
1342e5ec60SJeff-Hadley 
1442e5ec60SJeff-Hadley PetscErrorCode Create1dLaplacian(PetscInt, Mat *);
1542e5ec60SJeff-Hadley PetscErrorCode DataCompExportMats(DataCompression);
1642e5ec60SJeff-Hadley PetscErrorCode DataCompDestroy(DataCompression);
1742e5ec60SJeff-Hadley 
main(int Argc,char ** Args)1842e5ec60SJeff-Hadley int main(int Argc, char **Args)
1942e5ec60SJeff-Hadley {
2042e5ec60SJeff-Hadley   PetscInt        n_nodes = 33;
2142e5ec60SJeff-Hadley   Vec             x, b;
2242e5ec60SJeff-Hadley   PC              pcHypre;
2342e5ec60SJeff-Hadley   DataCompression data_comp;
2442e5ec60SJeff-Hadley 
2542e5ec60SJeff-Hadley   PetscFunctionBeginUser;
2642e5ec60SJeff-Hadley   PetscCall(PetscInitialize(&Argc, &Args, NULL, NULL));
2742e5ec60SJeff-Hadley 
2842e5ec60SJeff-Hadley   PetscCall(PetscNew(&data_comp));
2942e5ec60SJeff-Hadley 
3042e5ec60SJeff-Hadley   // Creating stiffness matrix
3142e5ec60SJeff-Hadley   PetscCall(Create1dLaplacian(n_nodes, &data_comp->stiff));
3242e5ec60SJeff-Hadley   PetscCall(PetscObjectSetName((PetscObject)data_comp->stiff, "Stiffness"));
3342e5ec60SJeff-Hadley 
3442e5ec60SJeff-Hadley   // Set-up BoomerAMG PC to get Prolongation Operators and Coarse/Fine splittings
3542e5ec60SJeff-Hadley   PetscCall(KSPCreate(PETSC_COMM_WORLD, &data_comp->kspHypre));
3642e5ec60SJeff-Hadley   PetscCall(KSPSetType(data_comp->kspHypre, KSPRICHARDSON));
3742e5ec60SJeff-Hadley   PetscCall(KSPGetPC(data_comp->kspHypre, &pcHypre));
3842e5ec60SJeff-Hadley   PetscCall(PCSetType(pcHypre, PCHYPRE));
3942e5ec60SJeff-Hadley   PetscCall(PCHYPRESetType(pcHypre, "boomeramg"));
4042e5ec60SJeff-Hadley   PetscCall(PCSetFromOptions(pcHypre));
4142e5ec60SJeff-Hadley   PetscCall(PCSetOperators(pcHypre, data_comp->stiff, data_comp->stiff));
4242e5ec60SJeff-Hadley   PetscCall(PCSetUp(pcHypre));
4342e5ec60SJeff-Hadley 
4442e5ec60SJeff-Hadley   PetscCall(MatCreateVecs(data_comp->stiff, &x, &b));
4542e5ec60SJeff-Hadley   PetscCall(PCApply(pcHypre, x, b));
4642e5ec60SJeff-Hadley   PetscCall(VecDestroy(&x));
4742e5ec60SJeff-Hadley   PetscCall(VecDestroy(&b));
4842e5ec60SJeff-Hadley 
4942e5ec60SJeff-Hadley   //Viewing the PC and Extracting the Prolongation Operators and CFMarkers
5042e5ec60SJeff-Hadley   PetscCall(PCView(pcHypre, NULL));
5142e5ec60SJeff-Hadley   PetscCall(PCGetInterpolations(pcHypre, &data_comp->num_levels, &data_comp->ProlongationOps));
5242e5ec60SJeff-Hadley   PetscCall(PCHYPREGetCFMarkers(pcHypre, &data_comp->n_per_level, &data_comp->CFMarkers));
5342e5ec60SJeff-Hadley 
5442e5ec60SJeff-Hadley   PetscCall(DataCompExportMats(data_comp));
5542e5ec60SJeff-Hadley   PetscCall(DataCompDestroy(data_comp));
5642e5ec60SJeff-Hadley   PetscCall(PetscFinalize());
5742e5ec60SJeff-Hadley   return 0;
5842e5ec60SJeff-Hadley }
5942e5ec60SJeff-Hadley 
Create1dLaplacian(PetscInt n,Mat * mat)6042e5ec60SJeff-Hadley PetscErrorCode Create1dLaplacian(PetscInt n, Mat *mat)
6142e5ec60SJeff-Hadley {
6242e5ec60SJeff-Hadley   PetscFunctionBeginUser;
6342e5ec60SJeff-Hadley   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, n, n, 3, NULL, mat));
6442e5ec60SJeff-Hadley   PetscCall(MatSetValue(*mat, n - 1, n - 1, 2.0, INSERT_VALUES));
6542e5ec60SJeff-Hadley   for (PetscInt i = 0; i < n - 1; i++) {
6642e5ec60SJeff-Hadley     PetscCall(MatSetValue(*mat, i, i, 2.0, INSERT_VALUES));
6742e5ec60SJeff-Hadley     PetscCall(MatSetValue(*mat, i + 1, i, -1.0, INSERT_VALUES));
6842e5ec60SJeff-Hadley     PetscCall(MatSetValue(*mat, i, i + 1, -1.0, INSERT_VALUES));
6942e5ec60SJeff-Hadley   }
7042e5ec60SJeff-Hadley   PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY));
7142e5ec60SJeff-Hadley   PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY));
7242e5ec60SJeff-Hadley   PetscFunctionReturn(PETSC_SUCCESS);
7342e5ec60SJeff-Hadley }
7442e5ec60SJeff-Hadley 
DataCompExportMats(DataCompression data_comp)7542e5ec60SJeff-Hadley PetscErrorCode DataCompExportMats(DataCompression data_comp)
7642e5ec60SJeff-Hadley {
7742e5ec60SJeff-Hadley   PetscFunctionBeginUser;
7842e5ec60SJeff-Hadley   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Num levels: %" PetscInt_FMT "\n", data_comp->num_levels));
7942e5ec60SJeff-Hadley   PetscCall(PetscPrintf(PETSC_COMM_WORLD, " -- Nodes per level --\n"));
80*3a7d0413SPierre Jolivet   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]));
8142e5ec60SJeff-Hadley 
8242e5ec60SJeff-Hadley   for (PetscInt i = 0; i < data_comp->num_levels - 1; i++) {
8342e5ec60SJeff-Hadley     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Prolongation Operator - Level %" PetscInt_FMT "\n", i));
8442e5ec60SJeff-Hadley     PetscCall(PetscObjectSetName((PetscObject)data_comp->ProlongationOps[i], "P"));
8542e5ec60SJeff-Hadley     PetscCall(MatView(data_comp->ProlongationOps[i], PETSC_VIEWER_STDOUT_WORLD));
8642e5ec60SJeff-Hadley     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
8742e5ec60SJeff-Hadley   }
8842e5ec60SJeff-Hadley 
8942e5ec60SJeff-Hadley   for (PetscInt i = 0; i < data_comp->num_levels - 1; i++) {
9042e5ec60SJeff-Hadley     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Coarse/Fine splitting - Level %" PetscInt_FMT "\n", i + 1));
9142e5ec60SJeff-Hadley     PetscCall(PetscBTView(data_comp->n_per_level[i + 1], data_comp->CFMarkers[i], PETSC_VIEWER_STDOUT_WORLD));
9242e5ec60SJeff-Hadley   }
9342e5ec60SJeff-Hadley   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Stiffness matrix, sparse format:\n"));
9442e5ec60SJeff-Hadley   PetscCall(MatViewFromOptions(data_comp->stiff, NULL, "-mat_view_stiff"));
9542e5ec60SJeff-Hadley   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Finished calling the Viewer functions\n"));
9642e5ec60SJeff-Hadley   PetscFunctionReturn(PETSC_SUCCESS);
9742e5ec60SJeff-Hadley }
9842e5ec60SJeff-Hadley 
DataCompDestroy(DataCompression data_comp)9942e5ec60SJeff-Hadley PetscErrorCode DataCompDestroy(DataCompression data_comp)
10042e5ec60SJeff-Hadley {
10142e5ec60SJeff-Hadley   PetscFunctionBeginUser;
10242e5ec60SJeff-Hadley   if (data_comp == NULL) PetscFunctionReturn(PETSC_SUCCESS);
10342e5ec60SJeff-Hadley   PetscCall(MatDestroy(&data_comp->stiff));
10442e5ec60SJeff-Hadley   PetscCall(KSPDestroy(&data_comp->kspHypre));
10542e5ec60SJeff-Hadley   for (PetscInt i = 0; i < data_comp->num_levels - 1; i++) {
10642e5ec60SJeff-Hadley     PetscCall(MatDestroy(&data_comp->ProlongationOps[i]));
10742e5ec60SJeff-Hadley     PetscCall(PetscBTDestroy(&data_comp->CFMarkers[i]));
10842e5ec60SJeff-Hadley   }
10942e5ec60SJeff-Hadley   PetscCall(PetscFree(data_comp->ProlongationOps));
11042e5ec60SJeff-Hadley   PetscCall(PetscFree(data_comp->n_per_level));
11142e5ec60SJeff-Hadley   PetscCall(PetscFree(data_comp->CFMarkers));
11242e5ec60SJeff-Hadley   PetscCall(PetscFree(data_comp));
11342e5ec60SJeff-Hadley   PetscFunctionReturn(PETSC_SUCCESS);
11442e5ec60SJeff-Hadley }
11542e5ec60SJeff-Hadley 
11642e5ec60SJeff-Hadley /*TEST
11742e5ec60SJeff-Hadley 
11842e5ec60SJeff-Hadley    test:
11942e5ec60SJeff-Hadley       requires: hypre !complex !defined(PETSC_HAVE_HYPRE_DEVICE)
12042e5ec60SJeff-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
12142e5ec60SJeff-Hadley 
12242e5ec60SJeff-Hadley TEST*/
123