xref: /libCEED/examples/fluids/src/setupdm.c (revision 7c5bba50effac2f172c1bac4d5ef42318aada6db)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Setup DM for Navier-Stokes example using PETSc
10 
11 #include <ceed.h>
12 #include <petscdmplex.h>
13 #include <petscds.h>
14 
15 #include "../navierstokes.h"
16 #include "../problems/stg_shur14.h"
17 
18 // Create mesh
19 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType mat_type, VecType vec_type, DM *dm) {
20   PetscFunctionBeginUser;
21   // Create DMPLEX
22   PetscCall(DMCreate(comm, dm));
23   PetscCall(DMSetType(*dm, DMPLEX));
24   {
25     PetscBool skip = PETSC_TRUE;
26     PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL));
27     PetscCall(DMSetMatrixPreallocateSkip(*dm, skip));
28   }
29   PetscCall(DMSetMatType(*dm, mat_type));
30   PetscCall(DMSetVecType(*dm, vec_type));
31 
32   // Set Tensor elements
33   PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"));
34   PetscCall(PetscOptionsSetValue(NULL, "-dm_sparse_localize", "0"));
35   // Set CL options
36   PetscCall(DMSetFromOptions(*dm));
37   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
38   PetscFunctionReturn(PETSC_SUCCESS);
39 }
40 
41 // Setup DM
42 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys) {
43   PetscInt num_comp_q = 5;
44   PetscFunctionBeginUser;
45 
46   PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, PETSC_DECIDE, q_extra, 1, &num_comp_q, dm));
47 
48   {  // Add strong boundary conditions to DM
49     DMLabel label;
50     PetscCall(DMGetLabel(dm, "Face Sets", &label));
51     PetscCall(DMPlexLabelComplete(dm, label));
52     // Set wall BCs
53     if (bc->num_wall > 0) {
54       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, bc->num_wall, bc->walls, 0, bc->num_comps, bc->wall_comps, NULL, NULL, NULL, NULL));
55     }
56     // Set symmetry BCs in the x direction
57     if (bc->num_symmetry[0] > 0) {
58       PetscInt comps[1] = {1};
59       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_x", label, bc->num_symmetry[0], bc->symmetries[0], 0, 1, comps, NULL, NULL, NULL, NULL));
60     }
61     // Set symmetry BCs in the y direction
62     if (bc->num_symmetry[1] > 0) {
63       PetscInt comps[1] = {2};
64       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_y", label, bc->num_symmetry[1], bc->symmetries[1], 0, 1, comps, NULL, NULL, NULL, NULL));
65     }
66     // Set symmetry BCs in the z direction
67     if (bc->num_symmetry[2] > 0) {
68       PetscInt comps[1] = {3};
69       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_z", label, bc->num_symmetry[2], bc->symmetries[2], 0, 1, comps, NULL, NULL, NULL, NULL));
70     }
71     {
72       PetscBool use_strongstg = PETSC_FALSE;
73       PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
74       if (use_strongstg) PetscCall(SetupStrongStg(dm, bc, problem, phys));
75     }
76   }
77 
78   PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm));
79 
80   // Empty name for conserved field (because there is only one field)
81   PetscSection section;
82   PetscCall(DMGetLocalSection(dm, &section));
83   PetscCall(PetscSectionSetFieldName(section, 0, ""));
84   switch (phys->state_var) {
85     case STATEVAR_CONSERVATIVE:
86       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density"));
87       PetscCall(PetscSectionSetComponentName(section, 0, 1, "MomentumX"));
88       PetscCall(PetscSectionSetComponentName(section, 0, 2, "MomentumY"));
89       PetscCall(PetscSectionSetComponentName(section, 0, 3, "MomentumZ"));
90       PetscCall(PetscSectionSetComponentName(section, 0, 4, "TotalEnergy"));
91       break;
92 
93     case STATEVAR_PRIMITIVE:
94       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure"));
95       PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX"));
96       PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY"));
97       PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ"));
98       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature"));
99       break;
100   }
101   PetscFunctionReturn(PETSC_SUCCESS);
102 }
103 
104 // Refine DM for high-order viz
105 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys) {
106   DM      dm_hierarchy[user->app_ctx->viz_refine + 1];
107   VecType vec_type;
108 
109   PetscFunctionBeginUser;
110   PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
111 
112   dm_hierarchy[0] = dm;
113   for (PetscInt i = 0, d = user->app_ctx->degree; i < user->app_ctx->viz_refine; i++) {
114     Mat interp_next;
115     PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1]));
116     PetscCall(DMClearDS(dm_hierarchy[i + 1]));
117     PetscCall(DMClearFields(dm_hierarchy[i + 1]));
118     PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i]));
119     d                = (d + 1) / 2;
120     PetscInt q_order = d + user->app_ctx->q_extra;
121     if (i + 1 == user->app_ctx->viz_refine) d = 1;
122     PetscCall(DMGetVecType(dm, &vec_type));
123     PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type));
124     PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, bc, phys));
125     PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL));
126     if (!i) user->interp_viz = interp_next;
127     else {
128       Mat C;
129       PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C));
130       PetscCall(MatDestroy(&interp_next));
131       PetscCall(MatDestroy(&user->interp_viz));
132       user->interp_viz = C;
133     }
134   }
135   for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) {
136     PetscCall(DMDestroy(&dm_hierarchy[i]));
137   }
138   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
139   PetscFunctionReturn(PETSC_SUCCESS);
140 }
141