1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, 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 11e419654dSJeremy L Thompson #include <ceed.h> 12e419654dSJeremy L Thompson #include <petscdmplex.h> 1367263decSKenneth E. Jansen #include <petscds.h> 14e419654dSJeremy L Thompson 15a515125bSLeila Ghaffari #include "../navierstokes.h" 16866f23abSJames Wright #include "../problems/stg_shur14.h" 17a515125bSLeila Ghaffari 1805a512bdSLeila Ghaffari // Create mesh 19*991aef52SJames Wright PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType mat_type, VecType vec_type, DM *dm) { 20a515125bSLeila Ghaffari PetscFunctionBeginUser; 2105a512bdSLeila Ghaffari // Create DMPLEX 222b916ea7SJeremy L Thompson PetscCall(DMCreate(comm, dm)); 232b916ea7SJeremy L Thompson PetscCall(DMSetType(*dm, DMPLEX)); 24b150a244SJed Brown { 25b150a244SJed Brown PetscBool skip = PETSC_TRUE; 262b916ea7SJeremy L Thompson PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL)); 27b150a244SJed Brown PetscCall(DMSetMatrixPreallocateSkip(*dm, skip)); 28b150a244SJed Brown } 292b916ea7SJeremy L Thompson PetscCall(DMSetMatType(*dm, mat_type)); 302b916ea7SJeremy L Thompson PetscCall(DMSetVecType(*dm, vec_type)); 31a9804fe9SJed Brown 3205a512bdSLeila Ghaffari // Set Tensor elements 332b916ea7SJeremy L Thompson PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0")); 342b916ea7SJeremy L Thompson PetscCall(PetscOptionsSetValue(NULL, "-dm_sparse_localize", "0")); 359bfedf88SJames Wright PetscCall(PetscOptionsSetValue(NULL, "-dm_localize", "0")); // Localization done in DMSetupByOrderEnd_FEM 36ebfabadfSJames Wright PetscCall(PetscOptionsSetValue(NULL, "-dm_blocking_type", "field_node")); 37ebfabadfSJames Wright 3805a512bdSLeila Ghaffari // Set CL options 392b916ea7SJeremy L Thompson PetscCall(DMSetFromOptions(*dm)); 402b916ea7SJeremy L Thompson PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 41d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 42a515125bSLeila Ghaffari } 43a515125bSLeila Ghaffari 44a515125bSLeila Ghaffari // Setup DM 45*991aef52SJames Wright PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys) { 46a515125bSLeila Ghaffari PetscInt num_comp_q = 5; 47da4ca0cfSJames Wright PetscFunctionBeginUser; 48da4ca0cfSJames Wright 49b2e5b5b3SJames Wright PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, PETSC_DECIDE, q_extra, 1, &num_comp_q, dm)); 50da4ca0cfSJames Wright 51da4ca0cfSJames Wright { // Add strong boundary conditions to DM 52047c2946SJed Brown DMLabel label; 532b916ea7SJeremy L Thompson PetscCall(DMGetLabel(dm, "Face Sets", &label)); 54c5e9980aSAdeleke O. Bankole PetscCall(DMPlexLabelComplete(dm, label)); 55002797a3SLeila Ghaffari // Set wall BCs 56002797a3SLeila Ghaffari if (bc->num_wall > 0) { 57bedfd28dSJames 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)); 58002797a3SLeila Ghaffari } 59fce2147eSJames Wright // Set symmetry BCs in the x direction 60fce2147eSJames Wright if (bc->num_symmetry[0] > 0) { 61002797a3SLeila Ghaffari PetscInt comps[1] = {1}; 62fce2147eSJames Wright PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_x", label, bc->num_symmetry[0], bc->symmetries[0], 0, 1, comps, NULL, NULL, NULL, NULL)); 63002797a3SLeila Ghaffari } 64fce2147eSJames Wright // Set symmetry BCs in the y direction 65fce2147eSJames Wright if (bc->num_symmetry[1] > 0) { 66002797a3SLeila Ghaffari PetscInt comps[1] = {2}; 67fce2147eSJames Wright PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_y", label, bc->num_symmetry[1], bc->symmetries[1], 0, 1, comps, NULL, NULL, NULL, NULL)); 68002797a3SLeila Ghaffari } 69fce2147eSJames Wright // Set symmetry BCs in the z direction 70fce2147eSJames Wright if (bc->num_symmetry[2] > 0) { 71002797a3SLeila Ghaffari PetscInt comps[1] = {3}; 72fce2147eSJames Wright PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_z", label, bc->num_symmetry[2], bc->symmetries[2], 0, 1, comps, NULL, NULL, NULL, NULL)); 73002797a3SLeila Ghaffari } 74866f23abSJames Wright { 75866f23abSJames Wright PetscBool use_strongstg = PETSC_FALSE; 762b916ea7SJeremy L Thompson PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL)); 7742454adaSJames Wright if (use_strongstg) PetscCall(SetupStrongStg(dm, bc, problem, phys)); 78866f23abSJames Wright } 7967263decSKenneth E. Jansen } 8067263decSKenneth E. Jansen 81da4ca0cfSJames Wright PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm)); 82cbe60e31SLeila Ghaffari 83a515125bSLeila Ghaffari // Empty name for conserved field (because there is only one field) 84a515125bSLeila Ghaffari PetscSection section; 852b916ea7SJeremy L Thompson PetscCall(DMGetLocalSection(dm, §ion)); 862b916ea7SJeremy L Thompson PetscCall(PetscSectionSetFieldName(section, 0, "")); 873636f6a4SJames Wright switch (phys->state_var) { 883636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 892b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density")); 902b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 1, "MomentumX")); 912b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 2, "MomentumY")); 922b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 3, "MomentumZ")); 93525869caSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "TotalEnergy")); 943636f6a4SJames Wright break; 953636f6a4SJames Wright 963636f6a4SJames Wright case STATEVAR_PRIMITIVE: 972b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure")); 982b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX")); 992b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY")); 1002b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ")); 1012b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature")); 1023636f6a4SJames Wright break; 103a515125bSLeila Ghaffari } 104d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 105a515125bSLeila Ghaffari } 106a515125bSLeila Ghaffari 107a515125bSLeila Ghaffari // Refine DM for high-order viz 108*991aef52SJames Wright PetscErrorCode VizRefineDM(DM dm, User user, ProblemData problem, SimpleBC bc, Physics phys) { 109a515125bSLeila Ghaffari DM dm_hierarchy[user->app_ctx->viz_refine + 1]; 110a515125bSLeila Ghaffari VecType vec_type; 111a515125bSLeila Ghaffari 11206f41313SJames Wright PetscFunctionBeginUser; 1132b916ea7SJeremy L Thompson PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 114a515125bSLeila Ghaffari 115a515125bSLeila Ghaffari dm_hierarchy[0] = dm; 1162b916ea7SJeremy L Thompson for (PetscInt i = 0, d = user->app_ctx->degree; i < user->app_ctx->viz_refine; i++) { 117a515125bSLeila Ghaffari Mat interp_next; 1182b916ea7SJeremy L Thompson PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1])); 1192b916ea7SJeremy L Thompson PetscCall(DMClearDS(dm_hierarchy[i + 1])); 1202b916ea7SJeremy L Thompson PetscCall(DMClearFields(dm_hierarchy[i + 1])); 1212b916ea7SJeremy L Thompson PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i])); 122a515125bSLeila Ghaffari d = (d + 1) / 2; 12367263decSKenneth E. Jansen PetscInt q_order = d + user->app_ctx->q_extra; 124a515125bSLeila Ghaffari if (i + 1 == user->app_ctx->viz_refine) d = 1; 1252b916ea7SJeremy L Thompson PetscCall(DMGetVecType(dm, &vec_type)); 1262b916ea7SJeremy L Thompson PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type)); 12767263decSKenneth E. Jansen PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, bc, phys)); 1282b916ea7SJeremy L Thompson PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL)); 129a515125bSLeila Ghaffari if (!i) user->interp_viz = interp_next; 130a515125bSLeila Ghaffari else { 131a515125bSLeila Ghaffari Mat C; 1322b916ea7SJeremy L Thompson PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C)); 1332b916ea7SJeremy L Thompson PetscCall(MatDestroy(&interp_next)); 1342b916ea7SJeremy L Thompson PetscCall(MatDestroy(&user->interp_viz)); 135a515125bSLeila Ghaffari user->interp_viz = C; 136a515125bSLeila Ghaffari } 137a515125bSLeila Ghaffari } 138a515125bSLeila Ghaffari for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) { 1392b916ea7SJeremy L Thompson PetscCall(DMDestroy(&dm_hierarchy[i])); 140a515125bSLeila Ghaffari } 141a515125bSLeila Ghaffari user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine]; 142d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 143a515125bSLeila Ghaffari } 144