xref: /honee/src/setupdm.c (revision b150a24444ad01c0b386e5a3a34c2afd101f684f)
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);
23*b150a244SJed Brown   {
24*b150a244SJed Brown     PetscBool skip = PETSC_TRUE;
25*b150a244SJed Brown     PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip,
26*b150a244SJed Brown                                   NULL));
27*b150a244SJed Brown     PetscCall(DMSetMatrixPreallocateSkip(*dm, skip));
28*b150a244SJed Brown   }
29a9804fe9SJed Brown   ierr = DMSetMatType(*dm, mat_type); CHKERRQ(ierr);
30a9804fe9SJed Brown   ierr = DMSetVecType(*dm, vec_type); CHKERRQ(ierr);
31a9804fe9SJed Brown 
3205a512bdSLeila Ghaffari   // Set Tensor elements
3305a512bdSLeila Ghaffari   ierr = PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"); CHKERRQ(ierr);
34cac8aa24SJed Brown   ierr = PetscOptionsSetValue(NULL, "-dm_sparse_localize", "0"); CHKERRQ(ierr);
3505a512bdSLeila Ghaffari   // Set CL options
3605a512bdSLeila Ghaffari   ierr = DMSetFromOptions(*dm); CHKERRQ(ierr);
37a515125bSLeila Ghaffari   ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr);
38a515125bSLeila Ghaffari   PetscFunctionReturn(0);
39a515125bSLeila Ghaffari }
40a515125bSLeila Ghaffari 
41a515125bSLeila Ghaffari // Setup DM
42a515125bSLeila Ghaffari PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree,
43c848b8b7SJed Brown                        SimpleBC bc, Physics phys) {
44a515125bSLeila Ghaffari   PetscErrorCode ierr;
45a515125bSLeila Ghaffari   PetscFunctionBeginUser;
46a515125bSLeila Ghaffari   {
47a515125bSLeila Ghaffari     // Configure the finite element space and boundary conditions
48a515125bSLeila Ghaffari     PetscFE  fe;
49a515125bSLeila Ghaffari     PetscInt num_comp_q = 5;
50047c2946SJed Brown     DMLabel label;
51a515125bSLeila Ghaffari     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q,
52a515125bSLeila Ghaffari                                  PETSC_FALSE, degree, PETSC_DECIDE,
53a515125bSLeila Ghaffari                                  &fe); CHKERRQ(ierr);
54a515125bSLeila Ghaffari     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
55a515125bSLeila Ghaffari     ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr);
56a515125bSLeila Ghaffari     ierr = DMCreateDS(dm); CHKERRQ(ierr);
57047c2946SJed Brown     ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr);
58002797a3SLeila Ghaffari     // Set wall BCs
59002797a3SLeila Ghaffari     if (bc->num_wall > 0) {
60002797a3SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label,
61002797a3SLeila Ghaffari                            bc->num_wall, bc->walls, 0, bc->num_comps,
62002797a3SLeila Ghaffari                            bc->wall_comps, (void(*)(void))problem->bc,
63c848b8b7SJed Brown                            NULL, problem->bc_ctx, NULL);  CHKERRQ(ierr);
64002797a3SLeila Ghaffari     }
65002797a3SLeila Ghaffari     // Set slip BCs in the x direction
66002797a3SLeila Ghaffari     if (bc->num_slip[0] > 0) {
67002797a3SLeila Ghaffari       PetscInt comps[1] = {1};
68002797a3SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label,
69002797a3SLeila Ghaffari                            bc->num_slip[0], bc->slips[0], 0, 1, comps,
70c848b8b7SJed Brown                            (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr);
71002797a3SLeila Ghaffari     }
72002797a3SLeila Ghaffari     // Set slip BCs in the y direction
73002797a3SLeila Ghaffari     if (bc->num_slip[1] > 0) {
74002797a3SLeila Ghaffari       PetscInt comps[1] = {2};
75002797a3SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label,
76002797a3SLeila Ghaffari                            bc->num_slip[1], bc->slips[1], 0, 1, comps,
77c848b8b7SJed Brown                            (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr);
78002797a3SLeila Ghaffari     }
79002797a3SLeila Ghaffari     // Set slip BCs in the z direction
80002797a3SLeila Ghaffari     if (bc->num_slip[2] > 0) {
81002797a3SLeila Ghaffari       PetscInt comps[1] = {3};
82002797a3SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label,
83002797a3SLeila Ghaffari                            bc->num_slip[2], bc->slips[2], 0, 1, comps,
84c848b8b7SJed Brown                            (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr);
85002797a3SLeila Ghaffari     }
86866f23abSJames Wright     {
87866f23abSJames Wright       PetscBool use_strongstg = PETSC_FALSE;
88866f23abSJames Wright       ierr = PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL);
89866f23abSJames Wright       CHKERRQ(ierr);
90866f23abSJames Wright 
91866f23abSJames Wright       if (use_strongstg) {
928085925cSJames Wright         ierr = SetupStrongSTG(dm, bc, problem); CHKERRQ(ierr);
93866f23abSJames Wright       }
94866f23abSJames Wright     }
95866f23abSJames Wright 
96a515125bSLeila Ghaffari     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
97a515125bSLeila Ghaffari     CHKERRQ(ierr);
98a515125bSLeila Ghaffari     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
99a515125bSLeila Ghaffari   }
100a515125bSLeila Ghaffari   {
101a515125bSLeila Ghaffari     // Empty name for conserved field (because there is only one field)
102a515125bSLeila Ghaffari     PetscSection section;
103a515125bSLeila Ghaffari     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
104a515125bSLeila Ghaffari     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
105a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
106a515125bSLeila Ghaffari     CHKERRQ(ierr);
107a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 1, "Momentum X");
108a515125bSLeila Ghaffari     CHKERRQ(ierr);
109a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 2, "Momentum Y");
110a515125bSLeila Ghaffari     CHKERRQ(ierr);
111a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 3, "Momentum Z");
112a515125bSLeila Ghaffari     CHKERRQ(ierr);
113a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 4, "Energy Density");
114a515125bSLeila Ghaffari     CHKERRQ(ierr);
115a515125bSLeila Ghaffari   }
116a515125bSLeila Ghaffari   PetscFunctionReturn(0);
117a515125bSLeila Ghaffari }
118a515125bSLeila Ghaffari 
119a515125bSLeila Ghaffari // Refine DM for high-order viz
120a515125bSLeila Ghaffari PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem,
121c848b8b7SJed Brown                            SimpleBC bc, Physics phys) {
122a515125bSLeila Ghaffari   PetscErrorCode ierr;
123a515125bSLeila Ghaffari   DM             dm_hierarchy[user->app_ctx->viz_refine + 1];
124a515125bSLeila Ghaffari   VecType        vec_type;
125a515125bSLeila Ghaffari   PetscFunctionBeginUser;
126a515125bSLeila Ghaffari 
127a515125bSLeila Ghaffari   ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
128a515125bSLeila Ghaffari 
129a515125bSLeila Ghaffari   dm_hierarchy[0] = dm;
130a515125bSLeila Ghaffari   for (PetscInt i = 0, d = user->app_ctx->degree;
131a515125bSLeila Ghaffari        i < user->app_ctx->viz_refine; i++) {
132a515125bSLeila Ghaffari     Mat interp_next;
133a515125bSLeila Ghaffari     ierr = DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i+1]);
134a515125bSLeila Ghaffari     CHKERRQ(ierr);
135a515125bSLeila Ghaffari     ierr = DMClearDS(dm_hierarchy[i+1]); CHKERRQ(ierr);
136a515125bSLeila Ghaffari     ierr = DMClearFields(dm_hierarchy[i+1]); CHKERRQ(ierr);
137a515125bSLeila Ghaffari     ierr = DMSetCoarseDM(dm_hierarchy[i+1], dm_hierarchy[i]); CHKERRQ(ierr);
138a515125bSLeila Ghaffari     d = (d + 1) / 2;
139a515125bSLeila Ghaffari     if (i + 1 == user->app_ctx->viz_refine) d = 1;
140a515125bSLeila Ghaffari     ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr);
141a515125bSLeila Ghaffari     ierr = DMSetVecType(dm_hierarchy[i+1], vec_type); CHKERRQ(ierr);
142c848b8b7SJed Brown     ierr = SetUpDM(dm_hierarchy[i+1], problem, d, bc, phys);
143a515125bSLeila Ghaffari     CHKERRQ(ierr);
144a515125bSLeila Ghaffari     ierr = DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i+1], &interp_next,
145a515125bSLeila Ghaffari                                  NULL); CHKERRQ(ierr);
146a515125bSLeila Ghaffari     if (!i) user->interp_viz = interp_next;
147a515125bSLeila Ghaffari     else {
148a515125bSLeila Ghaffari       Mat C;
149a515125bSLeila Ghaffari       ierr = MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX,
150a515125bSLeila Ghaffari                         PETSC_DECIDE, &C); CHKERRQ(ierr);
151a515125bSLeila Ghaffari       ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
152a515125bSLeila Ghaffari       ierr = MatDestroy(&user->interp_viz); CHKERRQ(ierr);
153a515125bSLeila Ghaffari       user->interp_viz = C;
154a515125bSLeila Ghaffari     }
155a515125bSLeila Ghaffari   }
156a515125bSLeila Ghaffari   for (PetscInt i=1; i<user->app_ctx->viz_refine; i++) {
157a515125bSLeila Ghaffari     ierr = DMDestroy(&dm_hierarchy[i]); CHKERRQ(ierr);
158a515125bSLeila Ghaffari   }
159a515125bSLeila Ghaffari   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
160a515125bSLeila Ghaffari 
161a515125bSLeila Ghaffari   PetscFunctionReturn(0);
162a515125bSLeila Ghaffari }
163