xref: /libCEED/examples/fluids/src/setupdm.c (revision 0814d5a7c7108099a9a6106ea0d2f7f6b264fad8)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
377841947SLeila Ghaffari //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
577841947SLeila Ghaffari //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
777841947SLeila Ghaffari 
877841947SLeila Ghaffari /// @file
977841947SLeila Ghaffari /// Setup DM for Navier-Stokes example using PETSc
1077841947SLeila Ghaffari 
1149aac155SJeremy L Thompson #include <ceed.h>
1249aac155SJeremy L Thompson #include <petscdmplex.h>
13*0814d5a7SKenneth E. Jansen #include <petscds.h>
1449aac155SJeremy L Thompson 
1577841947SLeila Ghaffari #include "../navierstokes.h"
16363b60e1SJames Wright #include "../problems/stg_shur14.h"
1777841947SLeila Ghaffari 
181864f1c2SLeila Ghaffari // Create mesh
192b730f8bSJeremy L Thompson PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType mat_type, VecType vec_type, DM *dm) {
2077841947SLeila Ghaffari   PetscFunctionBeginUser;
211864f1c2SLeila Ghaffari   // Create DMPLEX
222b730f8bSJeremy L Thompson   PetscCall(DMCreate(comm, dm));
232b730f8bSJeremy L Thompson   PetscCall(DMSetType(*dm, DMPLEX));
24c32b0260SJed Brown   {
25c32b0260SJed Brown     PetscBool skip = PETSC_TRUE;
262b730f8bSJeremy L Thompson     PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL));
27c32b0260SJed Brown     PetscCall(DMSetMatrixPreallocateSkip(*dm, skip));
28c32b0260SJed Brown   }
292b730f8bSJeremy L Thompson   PetscCall(DMSetMatType(*dm, mat_type));
302b730f8bSJeremy L Thompson   PetscCall(DMSetVecType(*dm, vec_type));
314ea65e7bSJed Brown 
321864f1c2SLeila Ghaffari   // Set Tensor elements
332b730f8bSJeremy L Thompson   PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"));
342b730f8bSJeremy L Thompson   PetscCall(PetscOptionsSetValue(NULL, "-dm_sparse_localize", "0"));
351864f1c2SLeila Ghaffari   // Set CL options
362b730f8bSJeremy L Thompson   PetscCall(DMSetFromOptions(*dm));
372b730f8bSJeremy L Thompson   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
38ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3977841947SLeila Ghaffari }
4077841947SLeila Ghaffari 
4177841947SLeila Ghaffari // Setup DM
42*0814d5a7SKenneth E. Jansen PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys) {
43*0814d5a7SKenneth E. Jansen   PetscInt q_order = degree + q_extra;
4477841947SLeila Ghaffari   PetscFunctionBeginUser;
4577841947SLeila Ghaffari   {
46*0814d5a7SKenneth E. Jansen     PetscBool is_simplex = PETSC_TRUE;
47*0814d5a7SKenneth E. Jansen 
48*0814d5a7SKenneth E. Jansen     // Check if simplex or tensor-product mesh
49*0814d5a7SKenneth E. Jansen     PetscCall(DMPlexIsSimplex(dm, &is_simplex));
5077841947SLeila Ghaffari     // Configure the finite element space and boundary conditions
5177841947SLeila Ghaffari     PetscFE  fe;
5277841947SLeila Ghaffari     PetscInt num_comp_q = 5;
53496f2382SJed Brown     DMLabel  label;
54*0814d5a7SKenneth E. Jansen     PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q, is_simplex, degree, q_order, &fe));
552b730f8bSJeremy L Thompson     PetscCall(PetscObjectSetName((PetscObject)fe, "Q"));
562b730f8bSJeremy L Thompson     PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
572b730f8bSJeremy L Thompson     PetscCall(DMCreateDS(dm));
58*0814d5a7SKenneth E. Jansen     {  // Project coordinates to enrich quadrature space
59*0814d5a7SKenneth E. Jansen       DM             dm_coord;
60*0814d5a7SKenneth E. Jansen       PetscDS        ds_coord;
61*0814d5a7SKenneth E. Jansen       PetscFE        fe_coord_current, fe_coord_new, fe_coord_face_new;
62*0814d5a7SKenneth E. Jansen       PetscDualSpace fe_coord_dual_space;
63*0814d5a7SKenneth E. Jansen       PetscInt       fe_coord_order, num_comp_coord;
64*0814d5a7SKenneth E. Jansen 
65*0814d5a7SKenneth E. Jansen       PetscCall(DMGetCoordinateDM(dm, &dm_coord));
66*0814d5a7SKenneth E. Jansen       PetscCall(DMGetCoordinateDim(dm, &num_comp_coord));
67*0814d5a7SKenneth E. Jansen       PetscCall(DMGetRegionDS(dm_coord, NULL, NULL, &ds_coord, NULL));
68*0814d5a7SKenneth E. Jansen       PetscCall(PetscDSGetDiscretization(ds_coord, 0, (PetscObject *)&fe_coord_current));
69*0814d5a7SKenneth E. Jansen       PetscCall(PetscFEGetDualSpace(fe_coord_current, &fe_coord_dual_space));
70*0814d5a7SKenneth E. Jansen       PetscCall(PetscDualSpaceGetOrder(fe_coord_dual_space, &fe_coord_order));
71*0814d5a7SKenneth E. Jansen 
72*0814d5a7SKenneth E. Jansen       // Create FE for coordinates
73*0814d5a7SKenneth E. Jansen       PetscCheck(fe_coord_order == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER_INPUT,
74*0814d5a7SKenneth E. Jansen                  "Only linear mesh geometry supported. Recieved %d order", fe_coord_order);
75*0814d5a7SKenneth E. Jansen       PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_coord, is_simplex, fe_coord_order, q_order, &fe_coord_new));
76*0814d5a7SKenneth E. Jansen       PetscCall(PetscFEGetHeightSubspace(fe_coord_new, 1, &fe_coord_face_new));
77*0814d5a7SKenneth E. Jansen       PetscCall(DMProjectCoordinates(dm, fe_coord_new));
78*0814d5a7SKenneth E. Jansen       PetscCall(PetscFEDestroy(&fe_coord_new));
79*0814d5a7SKenneth E. Jansen     }
80*0814d5a7SKenneth E. Jansen 
812b730f8bSJeremy L Thompson     PetscCall(DMGetLabel(dm, "Face Sets", &label));
82ca69d878SAdeleke O. Bankole     PetscCall(DMPlexLabelComplete(dm, label));
832fe7aee7SLeila Ghaffari     // Set wall BCs
842fe7aee7SLeila Ghaffari     if (bc->num_wall > 0) {
858213d101SJames Wright       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, bc->num_wall, bc->walls, 0, bc->num_comps, bc->wall_comps, NULL, NULL, NULL, NULL));
862fe7aee7SLeila Ghaffari     }
872fe7aee7SLeila Ghaffari     // Set slip BCs in the x direction
882fe7aee7SLeila Ghaffari     if (bc->num_slip[0] > 0) {
892fe7aee7SLeila Ghaffari       PetscInt comps[1] = {1};
908213d101SJames Wright       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label, bc->num_slip[0], bc->slips[0], 0, 1, comps, NULL, NULL, NULL, NULL));
912fe7aee7SLeila Ghaffari     }
922fe7aee7SLeila Ghaffari     // Set slip BCs in the y direction
932fe7aee7SLeila Ghaffari     if (bc->num_slip[1] > 0) {
942fe7aee7SLeila Ghaffari       PetscInt comps[1] = {2};
958213d101SJames Wright       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label, bc->num_slip[1], bc->slips[1], 0, 1, comps, NULL, NULL, NULL, NULL));
962fe7aee7SLeila Ghaffari     }
972fe7aee7SLeila Ghaffari     // Set slip BCs in the z direction
982fe7aee7SLeila Ghaffari     if (bc->num_slip[2] > 0) {
992fe7aee7SLeila Ghaffari       PetscInt comps[1] = {3};
1008213d101SJames Wright       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label, bc->num_slip[2], bc->slips[2], 0, 1, comps, NULL, NULL, NULL, NULL));
1012fe7aee7SLeila Ghaffari     }
102363b60e1SJames Wright     {
103363b60e1SJames Wright       PetscBool use_strongstg = PETSC_FALSE;
1042b730f8bSJeremy L Thompson       PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
105f5452247SJames Wright       if (use_strongstg) PetscCall(SetupStrongSTG(dm, bc, problem, phys));
106363b60e1SJames Wright     }
107363b60e1SJames Wright 
108*0814d5a7SKenneth E. Jansen     if (!is_simplex) {
109*0814d5a7SKenneth E. Jansen       DM dm_coord;
110*0814d5a7SKenneth E. Jansen       PetscCall(DMGetCoordinateDM(dm, &dm_coord));
1112b730f8bSJeremy L Thompson       PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
112*0814d5a7SKenneth E. Jansen       PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL));
113*0814d5a7SKenneth E. Jansen     }
114*0814d5a7SKenneth E. Jansen 
1152b730f8bSJeremy L Thompson     PetscCall(PetscFEDestroy(&fe));
11677841947SLeila Ghaffari   }
117dc805cc4SLeila Ghaffari 
11877841947SLeila Ghaffari   // Empty name for conserved field (because there is only one field)
11977841947SLeila Ghaffari   PetscSection section;
1202b730f8bSJeremy L Thompson   PetscCall(DMGetLocalSection(dm, &section));
1212b730f8bSJeremy L Thompson   PetscCall(PetscSectionSetFieldName(section, 0, ""));
12297baf651SJames Wright   switch (phys->state_var) {
12397baf651SJames Wright     case STATEVAR_CONSERVATIVE:
1242b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density"));
1252b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "Momentum X"));
1262b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "Momentum Y"));
1272b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "Momentum Z"));
1282b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Energy Density"));
12997baf651SJames Wright       break;
13097baf651SJames Wright 
13197baf651SJames Wright     case STATEVAR_PRIMITIVE:
1322b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure"));
1332b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "Velocity X"));
1342b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "Velocity Y"));
1352b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "Velocity Z"));
1362b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature"));
13797baf651SJames Wright       break;
13877841947SLeila Ghaffari   }
139ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
14077841947SLeila Ghaffari }
14177841947SLeila Ghaffari 
14277841947SLeila Ghaffari // Refine DM for high-order viz
1432b730f8bSJeremy L Thompson PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys) {
14477841947SLeila Ghaffari   DM      dm_hierarchy[user->app_ctx->viz_refine + 1];
14577841947SLeila Ghaffari   VecType vec_type;
14677841947SLeila Ghaffari   PetscFunctionBeginUser;
14777841947SLeila Ghaffari 
1482b730f8bSJeremy L Thompson   PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
14977841947SLeila Ghaffari 
15077841947SLeila Ghaffari   dm_hierarchy[0] = dm;
1512b730f8bSJeremy L Thompson   for (PetscInt i = 0, d = user->app_ctx->degree; i < user->app_ctx->viz_refine; i++) {
15277841947SLeila Ghaffari     Mat interp_next;
1532b730f8bSJeremy L Thompson     PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1]));
1542b730f8bSJeremy L Thompson     PetscCall(DMClearDS(dm_hierarchy[i + 1]));
1552b730f8bSJeremy L Thompson     PetscCall(DMClearFields(dm_hierarchy[i + 1]));
1562b730f8bSJeremy L Thompson     PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i]));
15777841947SLeila Ghaffari     d                = (d + 1) / 2;
158*0814d5a7SKenneth E. Jansen     PetscInt q_order = d + user->app_ctx->q_extra;
15977841947SLeila Ghaffari     if (i + 1 == user->app_ctx->viz_refine) d = 1;
1602b730f8bSJeremy L Thompson     PetscCall(DMGetVecType(dm, &vec_type));
1612b730f8bSJeremy L Thompson     PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type));
162*0814d5a7SKenneth E. Jansen     PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, bc, phys));
1632b730f8bSJeremy L Thompson     PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL));
16477841947SLeila Ghaffari     if (!i) user->interp_viz = interp_next;
16577841947SLeila Ghaffari     else {
16677841947SLeila Ghaffari       Mat C;
1672b730f8bSJeremy L Thompson       PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C));
1682b730f8bSJeremy L Thompson       PetscCall(MatDestroy(&interp_next));
1692b730f8bSJeremy L Thompson       PetscCall(MatDestroy(&user->interp_viz));
17077841947SLeila Ghaffari       user->interp_viz = C;
17177841947SLeila Ghaffari     }
17277841947SLeila Ghaffari   }
17377841947SLeila Ghaffari   for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) {
1742b730f8bSJeremy L Thompson     PetscCall(DMDestroy(&dm_hierarchy[i]));
17577841947SLeila Ghaffari   }
17677841947SLeila Ghaffari   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
17777841947SLeila Ghaffari 
178ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
17977841947SLeila Ghaffari }
180