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> 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 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")); 35*91c97f41SJames Wright PetscCall(PetscOptionsSetValue(NULL, "-dm_blocking_type", "field_node")); 36*91c97f41SJames Wright 371864f1c2SLeila Ghaffari // Set CL options 382b730f8bSJeremy L Thompson PetscCall(DMSetFromOptions(*dm)); 392b730f8bSJeremy L Thompson PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 40ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4177841947SLeila Ghaffari } 4277841947SLeila Ghaffari 4377841947SLeila Ghaffari // Setup DM 440814d5a7SKenneth E. Jansen PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys) { 4577841947SLeila Ghaffari PetscInt num_comp_q = 5; 46d68a66c4SJames Wright PetscFunctionBeginUser; 47d68a66c4SJames Wright 48a0b9cdb5SJames Wright PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, PETSC_DECIDE, q_extra, 1, &num_comp_q, dm)); 49d68a66c4SJames Wright 50d68a66c4SJames Wright { // Add strong boundary conditions to DM 51496f2382SJed Brown DMLabel label; 522b730f8bSJeremy L Thompson PetscCall(DMGetLabel(dm, "Face Sets", &label)); 53ca69d878SAdeleke O. Bankole PetscCall(DMPlexLabelComplete(dm, label)); 542fe7aee7SLeila Ghaffari // Set wall BCs 552fe7aee7SLeila Ghaffari if (bc->num_wall > 0) { 568213d101SJames 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)); 572fe7aee7SLeila Ghaffari } 587c5bba50SJames Wright // Set symmetry BCs in the x direction 597c5bba50SJames Wright if (bc->num_symmetry[0] > 0) { 602fe7aee7SLeila Ghaffari PetscInt comps[1] = {1}; 617c5bba50SJames Wright PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_x", label, bc->num_symmetry[0], bc->symmetries[0], 0, 1, comps, NULL, NULL, NULL, NULL)); 622fe7aee7SLeila Ghaffari } 637c5bba50SJames Wright // Set symmetry BCs in the y direction 647c5bba50SJames Wright if (bc->num_symmetry[1] > 0) { 652fe7aee7SLeila Ghaffari PetscInt comps[1] = {2}; 667c5bba50SJames Wright PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_y", label, bc->num_symmetry[1], bc->symmetries[1], 0, 1, comps, NULL, NULL, NULL, NULL)); 672fe7aee7SLeila Ghaffari } 687c5bba50SJames Wright // Set symmetry BCs in the z direction 697c5bba50SJames Wright if (bc->num_symmetry[2] > 0) { 702fe7aee7SLeila Ghaffari PetscInt comps[1] = {3}; 717c5bba50SJames Wright PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_z", label, bc->num_symmetry[2], bc->symmetries[2], 0, 1, comps, NULL, NULL, NULL, NULL)); 722fe7aee7SLeila Ghaffari } 73363b60e1SJames Wright { 74363b60e1SJames Wright PetscBool use_strongstg = PETSC_FALSE; 752b730f8bSJeremy L Thompson PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL)); 76cbef7084SJames Wright if (use_strongstg) PetscCall(SetupStrongStg(dm, bc, problem, phys)); 77363b60e1SJames Wright } 780814d5a7SKenneth E. Jansen } 790814d5a7SKenneth E. Jansen 80d68a66c4SJames Wright PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm)); 81dc805cc4SLeila Ghaffari 8277841947SLeila Ghaffari // Empty name for conserved field (because there is only one field) 8377841947SLeila Ghaffari PetscSection section; 842b730f8bSJeremy L Thompson PetscCall(DMGetLocalSection(dm, §ion)); 852b730f8bSJeremy L Thompson PetscCall(PetscSectionSetFieldName(section, 0, "")); 8697baf651SJames Wright switch (phys->state_var) { 8797baf651SJames Wright case STATEVAR_CONSERVATIVE: 882b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density")); 892b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 1, "MomentumX")); 902b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 2, "MomentumY")); 912b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 3, "MomentumZ")); 92116622c9SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "TotalEnergy")); 9397baf651SJames Wright break; 9497baf651SJames Wright 9597baf651SJames Wright case STATEVAR_PRIMITIVE: 962b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure")); 972b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX")); 982b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY")); 992b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ")); 1002b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature")); 10197baf651SJames Wright break; 10277841947SLeila Ghaffari } 103ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 10477841947SLeila Ghaffari } 10577841947SLeila Ghaffari 10677841947SLeila Ghaffari // Refine DM for high-order viz 1072b730f8bSJeremy L Thompson PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys) { 10877841947SLeila Ghaffari DM dm_hierarchy[user->app_ctx->viz_refine + 1]; 10977841947SLeila Ghaffari VecType vec_type; 11077841947SLeila Ghaffari 111f17d818dSJames Wright PetscFunctionBeginUser; 1122b730f8bSJeremy L Thompson PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 11377841947SLeila Ghaffari 11477841947SLeila Ghaffari dm_hierarchy[0] = dm; 1152b730f8bSJeremy L Thompson for (PetscInt i = 0, d = user->app_ctx->degree; i < user->app_ctx->viz_refine; i++) { 11677841947SLeila Ghaffari Mat interp_next; 1172b730f8bSJeremy L Thompson PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1])); 1182b730f8bSJeremy L Thompson PetscCall(DMClearDS(dm_hierarchy[i + 1])); 1192b730f8bSJeremy L Thompson PetscCall(DMClearFields(dm_hierarchy[i + 1])); 1202b730f8bSJeremy L Thompson PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i])); 12177841947SLeila Ghaffari d = (d + 1) / 2; 1220814d5a7SKenneth E. Jansen PetscInt q_order = d + user->app_ctx->q_extra; 12377841947SLeila Ghaffari if (i + 1 == user->app_ctx->viz_refine) d = 1; 1242b730f8bSJeremy L Thompson PetscCall(DMGetVecType(dm, &vec_type)); 1252b730f8bSJeremy L Thompson PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type)); 1260814d5a7SKenneth E. Jansen PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, bc, phys)); 1272b730f8bSJeremy L Thompson PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL)); 12877841947SLeila Ghaffari if (!i) user->interp_viz = interp_next; 12977841947SLeila Ghaffari else { 13077841947SLeila Ghaffari Mat C; 1312b730f8bSJeremy L Thompson PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C)); 1322b730f8bSJeremy L Thompson PetscCall(MatDestroy(&interp_next)); 1332b730f8bSJeremy L Thompson PetscCall(MatDestroy(&user->interp_viz)); 13477841947SLeila Ghaffari user->interp_viz = C; 13577841947SLeila Ghaffari } 13677841947SLeila Ghaffari } 13777841947SLeila Ghaffari for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) { 1382b730f8bSJeremy L Thompson PetscCall(DMDestroy(&dm_hierarchy[i])); 13977841947SLeila Ghaffari } 14077841947SLeila Ghaffari user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine]; 141ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 14277841947SLeila Ghaffari } 143