xref: /libCEED/examples/fluids/src/setupdm.c (revision 019b76820d7ff306c177822c4e76ffe5939c204b)
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 "../navierstokes.h"
12 
13 // Create mesh
14 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, DM *dm) {
15   PetscErrorCode   ierr;
16   PetscFunctionBeginUser;
17   // Create DMPLEX
18   ierr = DMCreate(comm, dm); CHKERRQ(ierr);
19   ierr = DMSetType(*dm, DMPLEX); CHKERRQ(ierr);
20   // Set Tensor elements
21   ierr = PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"); CHKERRQ(ierr);
22   // Set CL options
23   ierr = DMSetFromOptions(*dm); CHKERRQ(ierr);
24   ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr);
25   PetscFunctionReturn(0);
26 }
27 
28 // Setup DM
29 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree,
30                        SimpleBC bc, Physics phys, void *setup_ctx) {
31   PetscErrorCode ierr;
32   PetscFunctionBeginUser;
33   {
34     // Configure the finite element space and boundary conditions
35     PetscFE  fe;
36     PetscInt num_comp_q = 5;
37     DMLabel label;
38     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q,
39                                  PETSC_FALSE, degree, PETSC_DECIDE,
40                                  &fe); CHKERRQ(ierr);
41     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
42     ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr);
43     ierr = DMCreateDS(dm); CHKERRQ(ierr);
44     ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr);
45     // Set wall BCs
46     if (bc->num_wall > 0) {
47       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label,
48                            bc->num_wall, bc->walls, 0, bc->num_comps,
49                            bc->wall_comps, (void(*)(void))problem->bc,
50                            NULL, setup_ctx, NULL);  CHKERRQ(ierr);
51     }
52     // Set slip BCs in the x direction
53     if (bc->num_slip[0] > 0) {
54       PetscInt comps[1] = {1};
55       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label,
56                            bc->num_slip[0], bc->slips[0], 0, 1, comps,
57                            (void(*)(void))NULL, NULL, setup_ctx, NULL); CHKERRQ(ierr);
58     }
59     // Set slip BCs in the y direction
60     if (bc->num_slip[1] > 0) {
61       PetscInt comps[1] = {2};
62       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label,
63                            bc->num_slip[1], bc->slips[1], 0, 1, comps,
64                            (void(*)(void))NULL, NULL, setup_ctx, NULL); CHKERRQ(ierr);
65     }
66     // Set slip BCs in the z direction
67     if (bc->num_slip[2] > 0) {
68       PetscInt comps[1] = {3};
69       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label,
70                            bc->num_slip[2], bc->slips[2], 0, 1, comps,
71                            (void(*)(void))NULL, NULL, setup_ctx, NULL); CHKERRQ(ierr);
72     }
73     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
74     CHKERRQ(ierr);
75     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
76   }
77   {
78     // Empty name for conserved field (because there is only one field)
79     PetscSection section;
80     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
81     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
82     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
83     CHKERRQ(ierr);
84     ierr = PetscSectionSetComponentName(section, 0, 1, "Momentum X");
85     CHKERRQ(ierr);
86     ierr = PetscSectionSetComponentName(section, 0, 2, "Momentum Y");
87     CHKERRQ(ierr);
88     ierr = PetscSectionSetComponentName(section, 0, 3, "Momentum Z");
89     CHKERRQ(ierr);
90     ierr = PetscSectionSetComponentName(section, 0, 4, "Energy Density");
91     CHKERRQ(ierr);
92   }
93   PetscFunctionReturn(0);
94 }
95 
96 // Refine DM for high-order viz
97 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem,
98                            SimpleBC bc, Physics phys, void *setup_ctx) {
99   PetscErrorCode ierr;
100   DM             dm_hierarchy[user->app_ctx->viz_refine + 1];
101   VecType        vec_type;
102   PetscFunctionBeginUser;
103 
104   ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
105 
106   dm_hierarchy[0] = dm;
107   for (PetscInt i = 0, d = user->app_ctx->degree;
108        i < user->app_ctx->viz_refine; i++) {
109     Mat interp_next;
110     ierr = DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i+1]);
111     CHKERRQ(ierr);
112     ierr = DMClearDS(dm_hierarchy[i+1]); CHKERRQ(ierr);
113     ierr = DMClearFields(dm_hierarchy[i+1]); CHKERRQ(ierr);
114     ierr = DMSetCoarseDM(dm_hierarchy[i+1], dm_hierarchy[i]); CHKERRQ(ierr);
115     d = (d + 1) / 2;
116     if (i + 1 == user->app_ctx->viz_refine) d = 1;
117     ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr);
118     ierr = DMSetVecType(dm_hierarchy[i+1], vec_type); CHKERRQ(ierr);
119     ierr = SetUpDM(dm_hierarchy[i+1], problem, d, bc, phys, setup_ctx);
120     CHKERRQ(ierr);
121     ierr = DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i+1], &interp_next,
122                                  NULL); CHKERRQ(ierr);
123     if (!i) user->interp_viz = interp_next;
124     else {
125       Mat C;
126       ierr = MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX,
127                         PETSC_DECIDE, &C); CHKERRQ(ierr);
128       ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
129       ierr = MatDestroy(&user->interp_viz); CHKERRQ(ierr);
130       user->interp_viz = C;
131     }
132   }
133   for (PetscInt i=1; i<user->app_ctx->viz_refine; i++) {
134     ierr = DMDestroy(&dm_hierarchy[i]); CHKERRQ(ierr);
135   }
136   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
137 
138   PetscFunctionReturn(0);
139 }
140