xref: /honee/src/setupdm.c (revision c5e9980ad559552e5e02639665d11e7187691949)
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
152b916ea7SJeremy L Thompson PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType mat_type, VecType vec_type, DM *dm) {
16a515125bSLeila Ghaffari   PetscFunctionBeginUser;
1705a512bdSLeila Ghaffari   // Create DMPLEX
182b916ea7SJeremy L Thompson   PetscCall(DMCreate(comm, dm));
192b916ea7SJeremy L Thompson   PetscCall(DMSetType(*dm, DMPLEX));
20b150a244SJed Brown   {
21b150a244SJed Brown     PetscBool skip = PETSC_TRUE;
222b916ea7SJeremy L Thompson     PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL));
23b150a244SJed Brown     PetscCall(DMSetMatrixPreallocateSkip(*dm, skip));
24b150a244SJed Brown   }
252b916ea7SJeremy L Thompson   PetscCall(DMSetMatType(*dm, mat_type));
262b916ea7SJeremy L Thompson   PetscCall(DMSetVecType(*dm, vec_type));
27a9804fe9SJed Brown 
2805a512bdSLeila Ghaffari   // Set Tensor elements
292b916ea7SJeremy L Thompson   PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"));
302b916ea7SJeremy L Thompson   PetscCall(PetscOptionsSetValue(NULL, "-dm_sparse_localize", "0"));
3105a512bdSLeila Ghaffari   // Set CL options
322b916ea7SJeremy L Thompson   PetscCall(DMSetFromOptions(*dm));
332b916ea7SJeremy L Thompson   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
34a515125bSLeila Ghaffari   PetscFunctionReturn(0);
35a515125bSLeila Ghaffari }
36a515125bSLeila Ghaffari 
37a515125bSLeila Ghaffari // Setup DM
382b916ea7SJeremy L Thompson PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, SimpleBC bc, Physics phys) {
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;
452b916ea7SJeremy L Thompson     PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q, PETSC_FALSE, degree, PETSC_DECIDE, &fe));
462b916ea7SJeremy L Thompson     PetscCall(PetscObjectSetName((PetscObject)fe, "Q"));
472b916ea7SJeremy L Thompson     PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
482b916ea7SJeremy L Thompson     PetscCall(DMCreateDS(dm));
492b916ea7SJeremy L Thompson     PetscCall(DMGetLabel(dm, "Face Sets", &label));
50*c5e9980aSAdeleke O. Bankole     PetscCall(DMPlexLabelComplete(dm, label));
51002797a3SLeila Ghaffari     // Set wall BCs
52002797a3SLeila Ghaffari     if (bc->num_wall > 0) {
532b916ea7SJeremy L Thompson       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, bc->num_wall, bc->walls, 0, bc->num_comps, bc->wall_comps,
542b916ea7SJeremy L Thompson                               (void (*)(void))problem->bc, NULL, problem->bc_ctx, NULL));
55002797a3SLeila Ghaffari     }
56002797a3SLeila Ghaffari     // Set slip BCs in the x direction
57002797a3SLeila Ghaffari     if (bc->num_slip[0] > 0) {
58002797a3SLeila Ghaffari       PetscInt comps[1] = {1};
592b916ea7SJeremy L Thompson       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label, bc->num_slip[0], bc->slips[0], 0, 1, comps, (void (*)(void))NULL, NULL,
602b916ea7SJeremy L Thompson                               problem->bc_ctx, NULL));
61002797a3SLeila Ghaffari     }
62002797a3SLeila Ghaffari     // Set slip BCs in the y direction
63002797a3SLeila Ghaffari     if (bc->num_slip[1] > 0) {
64002797a3SLeila Ghaffari       PetscInt comps[1] = {2};
652b916ea7SJeremy L Thompson       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label, bc->num_slip[1], bc->slips[1], 0, 1, comps, (void (*)(void))NULL, NULL,
662b916ea7SJeremy L Thompson                               problem->bc_ctx, NULL));
67002797a3SLeila Ghaffari     }
68002797a3SLeila Ghaffari     // Set slip BCs in the z direction
69002797a3SLeila Ghaffari     if (bc->num_slip[2] > 0) {
70002797a3SLeila Ghaffari       PetscInt comps[1] = {3};
712b916ea7SJeremy L Thompson       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label, bc->num_slip[2], bc->slips[2], 0, 1, comps, (void (*)(void))NULL, NULL,
722b916ea7SJeremy L Thompson                               problem->bc_ctx, NULL));
73002797a3SLeila Ghaffari     }
74866f23abSJames Wright     {
75866f23abSJames Wright       PetscBool use_strongstg = PETSC_FALSE;
762b916ea7SJeremy L Thompson       PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
77866f23abSJames Wright 
78866f23abSJames Wright       if (use_strongstg) {
792b916ea7SJeremy L Thompson         PetscCall(SetupStrongSTG(dm, bc, problem, phys));
80866f23abSJames Wright       }
81866f23abSJames Wright     }
82866f23abSJames Wright 
832b916ea7SJeremy L Thompson     PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
842b916ea7SJeremy L Thompson     PetscCall(PetscFEDestroy(&fe));
85a515125bSLeila Ghaffari   }
86cbe60e31SLeila Ghaffari 
87a515125bSLeila Ghaffari   // Empty name for conserved field (because there is only one field)
88a515125bSLeila Ghaffari   PetscSection section;
892b916ea7SJeremy L Thompson   PetscCall(DMGetLocalSection(dm, &section));
902b916ea7SJeremy L Thompson   PetscCall(PetscSectionSetFieldName(section, 0, ""));
913636f6a4SJames Wright   switch (phys->state_var) {
923636f6a4SJames Wright     case STATEVAR_CONSERVATIVE:
932b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density"));
942b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "Momentum X"));
952b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "Momentum Y"));
962b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "Momentum Z"));
972b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Energy Density"));
983636f6a4SJames Wright       break;
993636f6a4SJames Wright 
1003636f6a4SJames Wright     case STATEVAR_PRIMITIVE:
1012b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure"));
1022b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "Velocity X"));
1032b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "Velocity Y"));
1042b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "Velocity Z"));
1052b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature"));
1063636f6a4SJames Wright       break;
107a515125bSLeila Ghaffari   }
108a515125bSLeila Ghaffari   PetscFunctionReturn(0);
109a515125bSLeila Ghaffari }
110a515125bSLeila Ghaffari 
111a515125bSLeila Ghaffari // Refine DM for high-order viz
1122b916ea7SJeremy L Thompson PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys) {
113a515125bSLeila Ghaffari   DM      dm_hierarchy[user->app_ctx->viz_refine + 1];
114a515125bSLeila Ghaffari   VecType vec_type;
115a515125bSLeila Ghaffari   PetscFunctionBeginUser;
116a515125bSLeila Ghaffari 
1172b916ea7SJeremy L Thompson   PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
118a515125bSLeila Ghaffari 
119a515125bSLeila Ghaffari   dm_hierarchy[0] = dm;
1202b916ea7SJeremy L Thompson   for (PetscInt i = 0, d = user->app_ctx->degree; i < user->app_ctx->viz_refine; i++) {
121a515125bSLeila Ghaffari     Mat interp_next;
1222b916ea7SJeremy L Thompson     PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1]));
1232b916ea7SJeremy L Thompson     PetscCall(DMClearDS(dm_hierarchy[i + 1]));
1242b916ea7SJeremy L Thompson     PetscCall(DMClearFields(dm_hierarchy[i + 1]));
1252b916ea7SJeremy L Thompson     PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i]));
126a515125bSLeila Ghaffari     d = (d + 1) / 2;
127a515125bSLeila Ghaffari     if (i + 1 == user->app_ctx->viz_refine) d = 1;
1282b916ea7SJeremy L Thompson     PetscCall(DMGetVecType(dm, &vec_type));
1292b916ea7SJeremy L Thompson     PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type));
1302b916ea7SJeremy L Thompson     PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, bc, phys));
1312b916ea7SJeremy L Thompson     PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL));
132a515125bSLeila Ghaffari     if (!i) user->interp_viz = interp_next;
133a515125bSLeila Ghaffari     else {
134a515125bSLeila Ghaffari       Mat C;
1352b916ea7SJeremy L Thompson       PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C));
1362b916ea7SJeremy L Thompson       PetscCall(MatDestroy(&interp_next));
1372b916ea7SJeremy L Thompson       PetscCall(MatDestroy(&user->interp_viz));
138a515125bSLeila Ghaffari       user->interp_viz = C;
139a515125bSLeila Ghaffari     }
140a515125bSLeila Ghaffari   }
141a515125bSLeila Ghaffari   for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) {
1422b916ea7SJeremy L Thompson     PetscCall(DMDestroy(&dm_hierarchy[i]));
143a515125bSLeila Ghaffari   }
144a515125bSLeila Ghaffari   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
145a515125bSLeila Ghaffari 
146a515125bSLeila Ghaffari   PetscFunctionReturn(0);
147a515125bSLeila Ghaffari }
148