xref: /libCEED/examples/fluids/src/setupdm.c (revision 13964f0727a62e5421e6d3b433e838b96a9ce891)
1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3 // reserved. See files LICENSE and NOTICE for details.
4 //
5 // This file is part of CEED, a collection of benchmarks, miniapps, software
6 // libraries and APIs for efficient high-order finite element and spectral
7 // element discretizations for exascale applications. For more information and
8 // source code availability see http://github.com/ceed.
9 //
10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11 // a collaborative effort of two U.S. Department of Energy organizations (Office
12 // of Science and the National Nuclear Security Administration) responsible for
13 // the planning and preparation of a capable exascale ecosystem, including
14 // software, applications, hardware, advanced system engineering and early
15 // testbed platforms, in support of the nation's exascale computing imperative.
16 
17 /// @file
18 /// Setup DM for Navier-Stokes example using PETSc
19 
20 #include "../navierstokes.h"
21 
22 // Create mesh
23 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, DM *dm) {
24   PetscErrorCode   ierr;
25   PetscFunctionBeginUser;
26   // Create DMPLEX
27   ierr = DMCreate(comm, dm); CHKERRQ(ierr);
28   ierr = DMSetType(*dm, DMPLEX); CHKERRQ(ierr);
29   // Set Tensor elements
30   ierr = PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"); CHKERRQ(ierr);
31   // Set CL options
32   ierr = DMSetFromOptions(*dm); CHKERRQ(ierr);
33   ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr);
34   PetscFunctionReturn(0);
35 }
36 
37 // Setup DM
38 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree,
39                        SimpleBC bc, Physics phys, void *setup_ctx) {
40   PetscErrorCode ierr;
41   PetscFunctionBeginUser;
42   {
43     // Configure the finite element space and boundary conditions
44     PetscFE  fe;
45     PetscInt num_comp_q = 5;
46     DMLabel label;
47     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q,
48                                  PETSC_FALSE, degree, PETSC_DECIDE,
49                                  &fe); CHKERRQ(ierr);
50     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
51     ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr);
52     ierr = DMCreateDS(dm); CHKERRQ(ierr);
53     {
54       /* create FE field for coordinates */
55       PetscFE fe_coords;
56       PetscInt num_comp_coord;
57       ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr);
58       ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_coord,
59                                    PETSC_FALSE, 1, 1, &fe_coords); CHKERRQ(ierr);
60       ierr = DMProjectCoordinates(dm, fe_coords); CHKERRQ(ierr);
61       ierr = PetscFEDestroy(&fe_coords); CHKERRQ(ierr);
62     }
63     ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr);
64     // Set wall BCs
65     if (bc->num_wall > 0) {
66       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label,
67                            bc->num_wall, bc->walls, 0, bc->num_comps,
68                            bc->wall_comps, (void(*)(void))problem->bc,
69                            NULL, setup_ctx, NULL);  CHKERRQ(ierr);
70     }
71     // Set slip BCs in the x direction
72     if (bc->num_slip[0] > 0) {
73       PetscInt comps[1] = {1};
74       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label,
75                            bc->num_slip[0], bc->slips[0], 0, 1, comps,
76                            (void(*)(void))NULL, NULL, setup_ctx, NULL); CHKERRQ(ierr);
77     }
78     // Set slip BCs in the y direction
79     if (bc->num_slip[1] > 0) {
80       PetscInt comps[1] = {2};
81       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label,
82                            bc->num_slip[1], bc->slips[1], 0, 1, comps,
83                            (void(*)(void))NULL, NULL, setup_ctx, NULL); CHKERRQ(ierr);
84     }
85     // Set slip BCs in the z direction
86     if (bc->num_slip[2] > 0) {
87       PetscInt comps[1] = {3};
88       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label,
89                            bc->num_slip[2], bc->slips[2], 0, 1, comps,
90                            (void(*)(void))NULL, NULL, setup_ctx, NULL); CHKERRQ(ierr);
91     }
92     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
93     CHKERRQ(ierr);
94     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
95   }
96   {
97     // Empty name for conserved field (because there is only one field)
98     PetscSection section;
99     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
100     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
101     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
102     CHKERRQ(ierr);
103     ierr = PetscSectionSetComponentName(section, 0, 1, "Momentum X");
104     CHKERRQ(ierr);
105     ierr = PetscSectionSetComponentName(section, 0, 2, "Momentum Y");
106     CHKERRQ(ierr);
107     ierr = PetscSectionSetComponentName(section, 0, 3, "Momentum Z");
108     CHKERRQ(ierr);
109     ierr = PetscSectionSetComponentName(section, 0, 4, "Energy Density");
110     CHKERRQ(ierr);
111   }
112   PetscFunctionReturn(0);
113 }
114 
115 // Refine DM for high-order viz
116 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem,
117                            SimpleBC bc, Physics phys, void *setup_ctx) {
118   PetscErrorCode ierr;
119   DM             dm_hierarchy[user->app_ctx->viz_refine + 1];
120   VecType        vec_type;
121   PetscFunctionBeginUser;
122 
123   ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
124 
125   dm_hierarchy[0] = dm;
126   for (PetscInt i = 0, d = user->app_ctx->degree;
127        i < user->app_ctx->viz_refine; i++) {
128     Mat interp_next;
129     ierr = DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i+1]);
130     CHKERRQ(ierr);
131     ierr = DMClearDS(dm_hierarchy[i+1]); CHKERRQ(ierr);
132     ierr = DMClearFields(dm_hierarchy[i+1]); CHKERRQ(ierr);
133     ierr = DMSetCoarseDM(dm_hierarchy[i+1], dm_hierarchy[i]); CHKERRQ(ierr);
134     d = (d + 1) / 2;
135     if (i + 1 == user->app_ctx->viz_refine) d = 1;
136     ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr);
137     ierr = DMSetVecType(dm_hierarchy[i+1], vec_type); CHKERRQ(ierr);
138     ierr = SetUpDM(dm_hierarchy[i+1], problem, d, bc, phys, setup_ctx);
139     CHKERRQ(ierr);
140     ierr = DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i+1], &interp_next,
141                                  NULL); CHKERRQ(ierr);
142     if (!i) user->interp_viz = interp_next;
143     else {
144       Mat C;
145       ierr = MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX,
146                         PETSC_DECIDE, &C); CHKERRQ(ierr);
147       ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
148       ierr = MatDestroy(&user->interp_viz); CHKERRQ(ierr);
149       user->interp_viz = C;
150     }
151   }
152   for (PetscInt i=1; i<user->app_ctx->viz_refine; i++) {
153     ierr = DMDestroy(&dm_hierarchy[i]); CHKERRQ(ierr);
154   }
155   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
156 
157   PetscFunctionReturn(0);
158 }
159