xref: /honee/src/setupdm.c (revision a9804fe94170f53922be4c5dcd1d89307c34dd7d)
1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Setup DM for Navier-Stokes example using PETSc
10a515125bSLeila Ghaffari 
11a515125bSLeila Ghaffari #include "../navierstokes.h"
12a515125bSLeila Ghaffari 
1305a512bdSLeila Ghaffari // Create mesh
14*a9804fe9SJed Brown PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem,
15*a9804fe9SJed Brown                         MatType mat_type, VecType vec_type,
16*a9804fe9SJed Brown                         DM *dm) {
17a515125bSLeila Ghaffari   PetscErrorCode   ierr;
18a515125bSLeila Ghaffari   PetscFunctionBeginUser;
1905a512bdSLeila Ghaffari   // Create DMPLEX
2005a512bdSLeila Ghaffari   ierr = DMCreate(comm, dm); CHKERRQ(ierr);
2105a512bdSLeila Ghaffari   ierr = DMSetType(*dm, DMPLEX); CHKERRQ(ierr);
22*a9804fe9SJed Brown   ierr = DMSetMatType(*dm, mat_type); CHKERRQ(ierr);
23*a9804fe9SJed Brown   ierr = DMSetVecType(*dm, vec_type); CHKERRQ(ierr);
24*a9804fe9SJed Brown 
2505a512bdSLeila Ghaffari   // Set Tensor elements
2605a512bdSLeila Ghaffari   ierr = PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"); CHKERRQ(ierr);
2705a512bdSLeila Ghaffari   // Set CL options
2805a512bdSLeila Ghaffari   ierr = DMSetFromOptions(*dm); CHKERRQ(ierr);
29a515125bSLeila Ghaffari   ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr);
30a515125bSLeila Ghaffari   PetscFunctionReturn(0);
31a515125bSLeila Ghaffari }
32a515125bSLeila Ghaffari 
33a515125bSLeila Ghaffari // Setup DM
34a515125bSLeila Ghaffari PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree,
35c848b8b7SJed Brown                        SimpleBC bc, Physics phys) {
36a515125bSLeila Ghaffari   PetscErrorCode ierr;
37a515125bSLeila Ghaffari   PetscFunctionBeginUser;
38a515125bSLeila Ghaffari   {
39a515125bSLeila Ghaffari     // Configure the finite element space and boundary conditions
40a515125bSLeila Ghaffari     PetscFE  fe;
41a515125bSLeila Ghaffari     PetscInt num_comp_q = 5;
42047c2946SJed Brown     DMLabel label;
43a515125bSLeila Ghaffari     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q,
44a515125bSLeila Ghaffari                                  PETSC_FALSE, degree, PETSC_DECIDE,
45a515125bSLeila Ghaffari                                  &fe); CHKERRQ(ierr);
46a515125bSLeila Ghaffari     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
47a515125bSLeila Ghaffari     ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr);
48a515125bSLeila Ghaffari     ierr = DMCreateDS(dm); CHKERRQ(ierr);
49047c2946SJed Brown     ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr);
50002797a3SLeila Ghaffari     // Set wall BCs
51002797a3SLeila Ghaffari     if (bc->num_wall > 0) {
52002797a3SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label,
53002797a3SLeila Ghaffari                            bc->num_wall, bc->walls, 0, bc->num_comps,
54002797a3SLeila Ghaffari                            bc->wall_comps, (void(*)(void))problem->bc,
55c848b8b7SJed Brown                            NULL, problem->bc_ctx, NULL);  CHKERRQ(ierr);
56002797a3SLeila Ghaffari     }
57002797a3SLeila Ghaffari     // Set slip BCs in the x direction
58002797a3SLeila Ghaffari     if (bc->num_slip[0] > 0) {
59002797a3SLeila Ghaffari       PetscInt comps[1] = {1};
60002797a3SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label,
61002797a3SLeila Ghaffari                            bc->num_slip[0], bc->slips[0], 0, 1, comps,
62c848b8b7SJed Brown                            (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr);
63002797a3SLeila Ghaffari     }
64002797a3SLeila Ghaffari     // Set slip BCs in the y direction
65002797a3SLeila Ghaffari     if (bc->num_slip[1] > 0) {
66002797a3SLeila Ghaffari       PetscInt comps[1] = {2};
67002797a3SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label,
68002797a3SLeila Ghaffari                            bc->num_slip[1], bc->slips[1], 0, 1, comps,
69c848b8b7SJed Brown                            (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr);
70002797a3SLeila Ghaffari     }
71002797a3SLeila Ghaffari     // Set slip BCs in the z direction
72002797a3SLeila Ghaffari     if (bc->num_slip[2] > 0) {
73002797a3SLeila Ghaffari       PetscInt comps[1] = {3};
74002797a3SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label,
75002797a3SLeila Ghaffari                            bc->num_slip[2], bc->slips[2], 0, 1, comps,
76c848b8b7SJed Brown                            (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr);
77002797a3SLeila Ghaffari     }
78a515125bSLeila Ghaffari     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
79a515125bSLeila Ghaffari     CHKERRQ(ierr);
80a515125bSLeila Ghaffari     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
81a515125bSLeila Ghaffari   }
82a515125bSLeila Ghaffari   {
83a515125bSLeila Ghaffari     // Empty name for conserved field (because there is only one field)
84a515125bSLeila Ghaffari     PetscSection section;
85a515125bSLeila Ghaffari     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
86a515125bSLeila Ghaffari     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
87a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
88a515125bSLeila Ghaffari     CHKERRQ(ierr);
89a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 1, "Momentum X");
90a515125bSLeila Ghaffari     CHKERRQ(ierr);
91a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 2, "Momentum Y");
92a515125bSLeila Ghaffari     CHKERRQ(ierr);
93a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 3, "Momentum Z");
94a515125bSLeila Ghaffari     CHKERRQ(ierr);
95a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 4, "Energy Density");
96a515125bSLeila Ghaffari     CHKERRQ(ierr);
97a515125bSLeila Ghaffari   }
98a515125bSLeila Ghaffari   PetscFunctionReturn(0);
99a515125bSLeila Ghaffari }
100a515125bSLeila Ghaffari 
101a515125bSLeila Ghaffari // Refine DM for high-order viz
102a515125bSLeila Ghaffari PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem,
103c848b8b7SJed Brown                            SimpleBC bc, Physics phys) {
104a515125bSLeila Ghaffari   PetscErrorCode ierr;
105a515125bSLeila Ghaffari   DM             dm_hierarchy[user->app_ctx->viz_refine + 1];
106a515125bSLeila Ghaffari   VecType        vec_type;
107a515125bSLeila Ghaffari   PetscFunctionBeginUser;
108a515125bSLeila Ghaffari 
109a515125bSLeila Ghaffari   ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
110a515125bSLeila Ghaffari 
111a515125bSLeila Ghaffari   dm_hierarchy[0] = dm;
112a515125bSLeila Ghaffari   for (PetscInt i = 0, d = user->app_ctx->degree;
113a515125bSLeila Ghaffari        i < user->app_ctx->viz_refine; i++) {
114a515125bSLeila Ghaffari     Mat interp_next;
115a515125bSLeila Ghaffari     ierr = DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i+1]);
116a515125bSLeila Ghaffari     CHKERRQ(ierr);
117a515125bSLeila Ghaffari     ierr = DMClearDS(dm_hierarchy[i+1]); CHKERRQ(ierr);
118a515125bSLeila Ghaffari     ierr = DMClearFields(dm_hierarchy[i+1]); CHKERRQ(ierr);
119a515125bSLeila Ghaffari     ierr = DMSetCoarseDM(dm_hierarchy[i+1], dm_hierarchy[i]); CHKERRQ(ierr);
120a515125bSLeila Ghaffari     d = (d + 1) / 2;
121a515125bSLeila Ghaffari     if (i + 1 == user->app_ctx->viz_refine) d = 1;
122a515125bSLeila Ghaffari     ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr);
123a515125bSLeila Ghaffari     ierr = DMSetVecType(dm_hierarchy[i+1], vec_type); CHKERRQ(ierr);
124c848b8b7SJed Brown     ierr = SetUpDM(dm_hierarchy[i+1], problem, d, bc, phys);
125a515125bSLeila Ghaffari     CHKERRQ(ierr);
126a515125bSLeila Ghaffari     ierr = DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i+1], &interp_next,
127a515125bSLeila Ghaffari                                  NULL); CHKERRQ(ierr);
128a515125bSLeila Ghaffari     if (!i) user->interp_viz = interp_next;
129a515125bSLeila Ghaffari     else {
130a515125bSLeila Ghaffari       Mat C;
131a515125bSLeila Ghaffari       ierr = MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX,
132a515125bSLeila Ghaffari                         PETSC_DECIDE, &C); CHKERRQ(ierr);
133a515125bSLeila Ghaffari       ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
134a515125bSLeila Ghaffari       ierr = MatDestroy(&user->interp_viz); CHKERRQ(ierr);
135a515125bSLeila Ghaffari       user->interp_viz = C;
136a515125bSLeila Ghaffari     }
137a515125bSLeila Ghaffari   }
138a515125bSLeila Ghaffari   for (PetscInt i=1; i<user->app_ctx->viz_refine; i++) {
139a515125bSLeila Ghaffari     ierr = DMDestroy(&dm_hierarchy[i]); CHKERRQ(ierr);
140a515125bSLeila Ghaffari   }
141a515125bSLeila Ghaffari   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
142a515125bSLeila Ghaffari 
143a515125bSLeila Ghaffari   PetscFunctionReturn(0);
144a515125bSLeila Ghaffari }
145