xref: /libCEED/examples/fluids/src/setupdm.c (revision 731c13d707f61db93aa32d87ad4bd1d287319f81)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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>
130814d5a7SKenneth 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
19*731c13d7SJames Wright 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"));
352127b05cSJames Wright   PetscCall(PetscOptionsSetValue(NULL, "-dm_localize", "0"));  // Localization done in DMSetupByOrderEnd_FEM
3691c97f41SJames Wright   PetscCall(PetscOptionsSetValue(NULL, "-dm_blocking_type", "field_node"));
3791c97f41SJames Wright 
381864f1c2SLeila Ghaffari   // Set CL options
392b730f8bSJeremy L Thompson   PetscCall(DMSetFromOptions(*dm));
402b730f8bSJeremy L Thompson   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
41ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4277841947SLeila Ghaffari }
4377841947SLeila Ghaffari 
4477841947SLeila Ghaffari // Setup DM
45*731c13d7SJames Wright PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys) {
4677841947SLeila Ghaffari   PetscInt num_comp_q = 5;
47d68a66c4SJames Wright   PetscFunctionBeginUser;
48d68a66c4SJames Wright 
49a0b9cdb5SJames Wright   PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, PETSC_DECIDE, q_extra, 1, &num_comp_q, dm));
50d68a66c4SJames Wright 
51d68a66c4SJames Wright   {  // Add strong boundary conditions to DM
52496f2382SJed Brown     DMLabel label;
532b730f8bSJeremy L Thompson     PetscCall(DMGetLabel(dm, "Face Sets", &label));
54ca69d878SAdeleke O. Bankole     PetscCall(DMPlexLabelComplete(dm, label));
552fe7aee7SLeila Ghaffari     // Set wall BCs
562fe7aee7SLeila Ghaffari     if (bc->num_wall > 0) {
578213d101SJames 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));
582fe7aee7SLeila Ghaffari     }
597c5bba50SJames Wright     // Set symmetry BCs in the x direction
607c5bba50SJames Wright     if (bc->num_symmetry[0] > 0) {
612fe7aee7SLeila Ghaffari       PetscInt comps[1] = {1};
627c5bba50SJames Wright       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_x", label, bc->num_symmetry[0], bc->symmetries[0], 0, 1, comps, NULL, NULL, NULL, NULL));
632fe7aee7SLeila Ghaffari     }
647c5bba50SJames Wright     // Set symmetry BCs in the y direction
657c5bba50SJames Wright     if (bc->num_symmetry[1] > 0) {
662fe7aee7SLeila Ghaffari       PetscInt comps[1] = {2};
677c5bba50SJames Wright       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_y", label, bc->num_symmetry[1], bc->symmetries[1], 0, 1, comps, NULL, NULL, NULL, NULL));
682fe7aee7SLeila Ghaffari     }
697c5bba50SJames Wright     // Set symmetry BCs in the z direction
707c5bba50SJames Wright     if (bc->num_symmetry[2] > 0) {
712fe7aee7SLeila Ghaffari       PetscInt comps[1] = {3};
727c5bba50SJames Wright       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_z", label, bc->num_symmetry[2], bc->symmetries[2], 0, 1, comps, NULL, NULL, NULL, NULL));
732fe7aee7SLeila Ghaffari     }
74363b60e1SJames Wright     {
75363b60e1SJames Wright       PetscBool use_strongstg = PETSC_FALSE;
762b730f8bSJeremy L Thompson       PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
77cbef7084SJames Wright       if (use_strongstg) PetscCall(SetupStrongStg(dm, bc, problem, phys));
78363b60e1SJames Wright     }
790814d5a7SKenneth E. Jansen   }
800814d5a7SKenneth E. Jansen 
81d68a66c4SJames Wright   PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm));
82dc805cc4SLeila Ghaffari 
8377841947SLeila Ghaffari   // Empty name for conserved field (because there is only one field)
8477841947SLeila Ghaffari   PetscSection section;
852b730f8bSJeremy L Thompson   PetscCall(DMGetLocalSection(dm, &section));
862b730f8bSJeremy L Thompson   PetscCall(PetscSectionSetFieldName(section, 0, ""));
8797baf651SJames Wright   switch (phys->state_var) {
8897baf651SJames Wright     case STATEVAR_CONSERVATIVE:
892b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density"));
902b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "MomentumX"));
912b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "MomentumY"));
922b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "MomentumZ"));
93116622c9SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 4, "TotalEnergy"));
9497baf651SJames Wright       break;
9597baf651SJames Wright 
9697baf651SJames Wright     case STATEVAR_PRIMITIVE:
972b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure"));
982b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX"));
992b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY"));
1002b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ"));
1012b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature"));
10297baf651SJames Wright       break;
10377841947SLeila Ghaffari   }
104ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
10577841947SLeila Ghaffari }
10677841947SLeila Ghaffari 
10777841947SLeila Ghaffari // Refine DM for high-order viz
108*731c13d7SJames Wright PetscErrorCode VizRefineDM(DM dm, User user, ProblemData problem, SimpleBC bc, Physics phys) {
10977841947SLeila Ghaffari   DM      dm_hierarchy[user->app_ctx->viz_refine + 1];
11077841947SLeila Ghaffari   VecType vec_type;
11177841947SLeila Ghaffari 
112f17d818dSJames Wright   PetscFunctionBeginUser;
1132b730f8bSJeremy L Thompson   PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
11477841947SLeila Ghaffari 
11577841947SLeila Ghaffari   dm_hierarchy[0] = dm;
1162b730f8bSJeremy L Thompson   for (PetscInt i = 0, d = user->app_ctx->degree; i < user->app_ctx->viz_refine; i++) {
11777841947SLeila Ghaffari     Mat interp_next;
1182b730f8bSJeremy L Thompson     PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1]));
1192b730f8bSJeremy L Thompson     PetscCall(DMClearDS(dm_hierarchy[i + 1]));
1202b730f8bSJeremy L Thompson     PetscCall(DMClearFields(dm_hierarchy[i + 1]));
1212b730f8bSJeremy L Thompson     PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i]));
12277841947SLeila Ghaffari     d                = (d + 1) / 2;
1230814d5a7SKenneth E. Jansen     PetscInt q_order = d + user->app_ctx->q_extra;
12477841947SLeila Ghaffari     if (i + 1 == user->app_ctx->viz_refine) d = 1;
1252b730f8bSJeremy L Thompson     PetscCall(DMGetVecType(dm, &vec_type));
1262b730f8bSJeremy L Thompson     PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type));
1270814d5a7SKenneth E. Jansen     PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, bc, phys));
1282b730f8bSJeremy L Thompson     PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL));
12977841947SLeila Ghaffari     if (!i) user->interp_viz = interp_next;
13077841947SLeila Ghaffari     else {
13177841947SLeila Ghaffari       Mat C;
1322b730f8bSJeremy L Thompson       PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C));
1332b730f8bSJeremy L Thompson       PetscCall(MatDestroy(&interp_next));
1342b730f8bSJeremy L Thompson       PetscCall(MatDestroy(&user->interp_viz));
13577841947SLeila Ghaffari       user->interp_viz = C;
13677841947SLeila Ghaffari     }
13777841947SLeila Ghaffari   }
13877841947SLeila Ghaffari   for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) {
1392b730f8bSJeremy L Thompson     PetscCall(DMDestroy(&dm_hierarchy[i]));
14077841947SLeila Ghaffari   }
14177841947SLeila Ghaffari   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
142ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
14377841947SLeila Ghaffari }
144