xref: /honee/src/setupdm.c (revision cac8aa24e3e327ca01b405ba38750eb1df7ccf42)
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"
12866f23abSJames Wright #include "../problems/stg_shur14.h"
13a515125bSLeila Ghaffari 
1405a512bdSLeila Ghaffari // Create mesh
15a9804fe9SJed Brown PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem,
16a9804fe9SJed Brown                         MatType mat_type, VecType vec_type,
17a9804fe9SJed Brown                         DM *dm) {
18a515125bSLeila Ghaffari   PetscErrorCode   ierr;
19a515125bSLeila Ghaffari   PetscFunctionBeginUser;
2005a512bdSLeila Ghaffari   // Create DMPLEX
2105a512bdSLeila Ghaffari   ierr = DMCreate(comm, dm); CHKERRQ(ierr);
2205a512bdSLeila Ghaffari   ierr = DMSetType(*dm, DMPLEX); CHKERRQ(ierr);
23a9804fe9SJed Brown   ierr = DMSetMatType(*dm, mat_type); CHKERRQ(ierr);
24a9804fe9SJed Brown   ierr = DMSetVecType(*dm, vec_type); CHKERRQ(ierr);
25a9804fe9SJed Brown 
2605a512bdSLeila Ghaffari   // Set Tensor elements
2705a512bdSLeila Ghaffari   ierr = PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"); CHKERRQ(ierr);
28*cac8aa24SJed Brown   ierr = PetscOptionsSetValue(NULL, "-dm_sparse_localize", "0"); CHKERRQ(ierr);
2905a512bdSLeila Ghaffari   // Set CL options
3005a512bdSLeila Ghaffari   ierr = DMSetFromOptions(*dm); CHKERRQ(ierr);
31a515125bSLeila Ghaffari   ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr);
32a515125bSLeila Ghaffari   PetscFunctionReturn(0);
33a515125bSLeila Ghaffari }
34a515125bSLeila Ghaffari 
35a515125bSLeila Ghaffari // Setup DM
36a515125bSLeila Ghaffari PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree,
37c848b8b7SJed Brown                        SimpleBC bc, Physics phys) {
38a515125bSLeila Ghaffari   PetscErrorCode ierr;
39a515125bSLeila Ghaffari   PetscFunctionBeginUser;
40a515125bSLeila Ghaffari   {
41a515125bSLeila Ghaffari     // Configure the finite element space and boundary conditions
42a515125bSLeila Ghaffari     PetscFE  fe;
43a515125bSLeila Ghaffari     PetscInt num_comp_q = 5;
44047c2946SJed Brown     DMLabel label;
45a515125bSLeila Ghaffari     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q,
46a515125bSLeila Ghaffari                                  PETSC_FALSE, degree, PETSC_DECIDE,
47a515125bSLeila Ghaffari                                  &fe); CHKERRQ(ierr);
48a515125bSLeila Ghaffari     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
49a515125bSLeila Ghaffari     ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr);
50a515125bSLeila Ghaffari     ierr = DMCreateDS(dm); CHKERRQ(ierr);
51047c2946SJed Brown     ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr);
52002797a3SLeila Ghaffari     // Set wall BCs
53002797a3SLeila Ghaffari     if (bc->num_wall > 0) {
54002797a3SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label,
55002797a3SLeila Ghaffari                            bc->num_wall, bc->walls, 0, bc->num_comps,
56002797a3SLeila Ghaffari                            bc->wall_comps, (void(*)(void))problem->bc,
57c848b8b7SJed Brown                            NULL, problem->bc_ctx, NULL);  CHKERRQ(ierr);
58002797a3SLeila Ghaffari     }
59002797a3SLeila Ghaffari     // Set slip BCs in the x direction
60002797a3SLeila Ghaffari     if (bc->num_slip[0] > 0) {
61002797a3SLeila Ghaffari       PetscInt comps[1] = {1};
62002797a3SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label,
63002797a3SLeila Ghaffari                            bc->num_slip[0], bc->slips[0], 0, 1, comps,
64c848b8b7SJed Brown                            (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr);
65002797a3SLeila Ghaffari     }
66002797a3SLeila Ghaffari     // Set slip BCs in the y direction
67002797a3SLeila Ghaffari     if (bc->num_slip[1] > 0) {
68002797a3SLeila Ghaffari       PetscInt comps[1] = {2};
69002797a3SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label,
70002797a3SLeila Ghaffari                            bc->num_slip[1], bc->slips[1], 0, 1, comps,
71c848b8b7SJed Brown                            (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr);
72002797a3SLeila Ghaffari     }
73002797a3SLeila Ghaffari     // Set slip BCs in the z direction
74002797a3SLeila Ghaffari     if (bc->num_slip[2] > 0) {
75002797a3SLeila Ghaffari       PetscInt comps[1] = {3};
76002797a3SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label,
77002797a3SLeila Ghaffari                            bc->num_slip[2], bc->slips[2], 0, 1, comps,
78c848b8b7SJed Brown                            (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr);
79002797a3SLeila Ghaffari     }
80866f23abSJames Wright     {
81866f23abSJames Wright       PetscBool use_strongstg = PETSC_FALSE;
82866f23abSJames Wright       ierr = PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL);
83866f23abSJames Wright       CHKERRQ(ierr);
84866f23abSJames Wright 
85866f23abSJames Wright       if (use_strongstg) {
868085925cSJames Wright         ierr = SetupStrongSTG(dm, bc, problem); CHKERRQ(ierr);
87866f23abSJames Wright       }
88866f23abSJames Wright     }
89866f23abSJames Wright 
90a515125bSLeila Ghaffari     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
91a515125bSLeila Ghaffari     CHKERRQ(ierr);
92a515125bSLeila Ghaffari     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
93a515125bSLeila Ghaffari   }
94a515125bSLeila Ghaffari   {
95a515125bSLeila Ghaffari     // Empty name for conserved field (because there is only one field)
96a515125bSLeila Ghaffari     PetscSection section;
97a515125bSLeila Ghaffari     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
98a515125bSLeila Ghaffari     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
99a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
100a515125bSLeila Ghaffari     CHKERRQ(ierr);
101a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 1, "Momentum X");
102a515125bSLeila Ghaffari     CHKERRQ(ierr);
103a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 2, "Momentum Y");
104a515125bSLeila Ghaffari     CHKERRQ(ierr);
105a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 3, "Momentum Z");
106a515125bSLeila Ghaffari     CHKERRQ(ierr);
107a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 4, "Energy Density");
108a515125bSLeila Ghaffari     CHKERRQ(ierr);
109a515125bSLeila Ghaffari   }
110a515125bSLeila Ghaffari   PetscFunctionReturn(0);
111a515125bSLeila Ghaffari }
112a515125bSLeila Ghaffari 
113a515125bSLeila Ghaffari // Refine DM for high-order viz
114a515125bSLeila Ghaffari PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem,
115c848b8b7SJed Brown                            SimpleBC bc, Physics phys) {
116a515125bSLeila Ghaffari   PetscErrorCode ierr;
117a515125bSLeila Ghaffari   DM             dm_hierarchy[user->app_ctx->viz_refine + 1];
118a515125bSLeila Ghaffari   VecType        vec_type;
119a515125bSLeila Ghaffari   PetscFunctionBeginUser;
120a515125bSLeila Ghaffari 
121a515125bSLeila Ghaffari   ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
122a515125bSLeila Ghaffari 
123a515125bSLeila Ghaffari   dm_hierarchy[0] = dm;
124a515125bSLeila Ghaffari   for (PetscInt i = 0, d = user->app_ctx->degree;
125a515125bSLeila Ghaffari        i < user->app_ctx->viz_refine; i++) {
126a515125bSLeila Ghaffari     Mat interp_next;
127a515125bSLeila Ghaffari     ierr = DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i+1]);
128a515125bSLeila Ghaffari     CHKERRQ(ierr);
129a515125bSLeila Ghaffari     ierr = DMClearDS(dm_hierarchy[i+1]); CHKERRQ(ierr);
130a515125bSLeila Ghaffari     ierr = DMClearFields(dm_hierarchy[i+1]); CHKERRQ(ierr);
131a515125bSLeila Ghaffari     ierr = DMSetCoarseDM(dm_hierarchy[i+1], dm_hierarchy[i]); CHKERRQ(ierr);
132a515125bSLeila Ghaffari     d = (d + 1) / 2;
133a515125bSLeila Ghaffari     if (i + 1 == user->app_ctx->viz_refine) d = 1;
134a515125bSLeila Ghaffari     ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr);
135a515125bSLeila Ghaffari     ierr = DMSetVecType(dm_hierarchy[i+1], vec_type); CHKERRQ(ierr);
136c848b8b7SJed Brown     ierr = SetUpDM(dm_hierarchy[i+1], problem, d, bc, phys);
137a515125bSLeila Ghaffari     CHKERRQ(ierr);
138a515125bSLeila Ghaffari     ierr = DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i+1], &interp_next,
139a515125bSLeila Ghaffari                                  NULL); CHKERRQ(ierr);
140a515125bSLeila Ghaffari     if (!i) user->interp_viz = interp_next;
141a515125bSLeila Ghaffari     else {
142a515125bSLeila Ghaffari       Mat C;
143a515125bSLeila Ghaffari       ierr = MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX,
144a515125bSLeila Ghaffari                         PETSC_DECIDE, &C); CHKERRQ(ierr);
145a515125bSLeila Ghaffari       ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
146a515125bSLeila Ghaffari       ierr = MatDestroy(&user->interp_viz); CHKERRQ(ierr);
147a515125bSLeila Ghaffari       user->interp_viz = C;
148a515125bSLeila Ghaffari     }
149a515125bSLeila Ghaffari   }
150a515125bSLeila Ghaffari   for (PetscInt i=1; i<user->app_ctx->viz_refine; i++) {
151a515125bSLeila Ghaffari     ierr = DMDestroy(&dm_hierarchy[i]); CHKERRQ(ierr);
152a515125bSLeila Ghaffari   }
153a515125bSLeila Ghaffari   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
154a515125bSLeila Ghaffari 
155a515125bSLeila Ghaffari   PetscFunctionReturn(0);
156a515125bSLeila Ghaffari }
157