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")); 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 420814d5a7SKenneth E. Jansen PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys) { 4377841947SLeila Ghaffari PetscInt num_comp_q = 5; 44*d68a66c4SJames Wright PetscFunctionBeginUser; 45*d68a66c4SJames Wright 46*d68a66c4SJames Wright PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &num_comp_q, dm)); 47*d68a66c4SJames Wright 48*d68a66c4SJames Wright { // Add strong boundary conditions to DM 49496f2382SJed Brown DMLabel label; 502b730f8bSJeremy L Thompson PetscCall(DMGetLabel(dm, "Face Sets", &label)); 51ca69d878SAdeleke O. Bankole PetscCall(DMPlexLabelComplete(dm, label)); 522fe7aee7SLeila Ghaffari // Set wall BCs 532fe7aee7SLeila Ghaffari if (bc->num_wall > 0) { 548213d101SJames 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)); 552fe7aee7SLeila Ghaffari } 562fe7aee7SLeila Ghaffari // Set slip BCs in the x direction 572fe7aee7SLeila Ghaffari if (bc->num_slip[0] > 0) { 582fe7aee7SLeila Ghaffari PetscInt comps[1] = {1}; 598213d101SJames Wright PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label, bc->num_slip[0], bc->slips[0], 0, 1, comps, NULL, NULL, NULL, NULL)); 602fe7aee7SLeila Ghaffari } 612fe7aee7SLeila Ghaffari // Set slip BCs in the y direction 622fe7aee7SLeila Ghaffari if (bc->num_slip[1] > 0) { 632fe7aee7SLeila Ghaffari PetscInt comps[1] = {2}; 648213d101SJames Wright PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label, bc->num_slip[1], bc->slips[1], 0, 1, comps, NULL, NULL, NULL, NULL)); 652fe7aee7SLeila Ghaffari } 662fe7aee7SLeila Ghaffari // Set slip BCs in the z direction 672fe7aee7SLeila Ghaffari if (bc->num_slip[2] > 0) { 682fe7aee7SLeila Ghaffari PetscInt comps[1] = {3}; 698213d101SJames Wright PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label, bc->num_slip[2], bc->slips[2], 0, 1, comps, NULL, NULL, NULL, NULL)); 702fe7aee7SLeila Ghaffari } 71363b60e1SJames Wright { 72363b60e1SJames Wright PetscBool use_strongstg = PETSC_FALSE; 732b730f8bSJeremy L Thompson PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL)); 74f5452247SJames Wright if (use_strongstg) PetscCall(SetupStrongSTG(dm, bc, problem, phys)); 75363b60e1SJames Wright } 760814d5a7SKenneth E. Jansen } 770814d5a7SKenneth E. Jansen 78*d68a66c4SJames Wright PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm)); 79dc805cc4SLeila Ghaffari 8077841947SLeila Ghaffari // Empty name for conserved field (because there is only one field) 8177841947SLeila Ghaffari PetscSection section; 822b730f8bSJeremy L Thompson PetscCall(DMGetLocalSection(dm, §ion)); 832b730f8bSJeremy L Thompson PetscCall(PetscSectionSetFieldName(section, 0, "")); 8497baf651SJames Wright switch (phys->state_var) { 8597baf651SJames Wright case STATEVAR_CONSERVATIVE: 862b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density")); 872b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 1, "MomentumX")); 882b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 2, "MomentumY")); 892b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 3, "MomentumZ")); 90116622c9SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "TotalEnergy")); 9197baf651SJames Wright break; 9297baf651SJames Wright 9397baf651SJames Wright case STATEVAR_PRIMITIVE: 942b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure")); 952b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX")); 962b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY")); 972b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ")); 982b730f8bSJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature")); 9997baf651SJames Wright break; 10077841947SLeila Ghaffari } 101ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 10277841947SLeila Ghaffari } 10377841947SLeila Ghaffari 10477841947SLeila Ghaffari // Refine DM for high-order viz 1052b730f8bSJeremy L Thompson PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys) { 10677841947SLeila Ghaffari DM dm_hierarchy[user->app_ctx->viz_refine + 1]; 10777841947SLeila Ghaffari VecType vec_type; 10877841947SLeila Ghaffari PetscFunctionBeginUser; 10977841947SLeila Ghaffari 1102b730f8bSJeremy L Thompson PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 11177841947SLeila Ghaffari 11277841947SLeila Ghaffari dm_hierarchy[0] = dm; 1132b730f8bSJeremy L Thompson for (PetscInt i = 0, d = user->app_ctx->degree; i < user->app_ctx->viz_refine; i++) { 11477841947SLeila Ghaffari Mat interp_next; 1152b730f8bSJeremy L Thompson PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1])); 1162b730f8bSJeremy L Thompson PetscCall(DMClearDS(dm_hierarchy[i + 1])); 1172b730f8bSJeremy L Thompson PetscCall(DMClearFields(dm_hierarchy[i + 1])); 1182b730f8bSJeremy L Thompson PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i])); 11977841947SLeila Ghaffari d = (d + 1) / 2; 1200814d5a7SKenneth E. Jansen PetscInt q_order = d + user->app_ctx->q_extra; 12177841947SLeila Ghaffari if (i + 1 == user->app_ctx->viz_refine) d = 1; 1222b730f8bSJeremy L Thompson PetscCall(DMGetVecType(dm, &vec_type)); 1232b730f8bSJeremy L Thompson PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type)); 1240814d5a7SKenneth E. Jansen PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, bc, phys)); 1252b730f8bSJeremy L Thompson PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL)); 12677841947SLeila Ghaffari if (!i) user->interp_viz = interp_next; 12777841947SLeila Ghaffari else { 12877841947SLeila Ghaffari Mat C; 1292b730f8bSJeremy L Thompson PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C)); 1302b730f8bSJeremy L Thompson PetscCall(MatDestroy(&interp_next)); 1312b730f8bSJeremy L Thompson PetscCall(MatDestroy(&user->interp_viz)); 13277841947SLeila Ghaffari user->interp_viz = C; 13377841947SLeila Ghaffari } 13477841947SLeila Ghaffari } 13577841947SLeila Ghaffari for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) { 1362b730f8bSJeremy L Thompson PetscCall(DMDestroy(&dm_hierarchy[i])); 13777841947SLeila Ghaffari } 13877841947SLeila Ghaffari user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine]; 13977841947SLeila Ghaffari 140ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 14177841947SLeila Ghaffari } 142