1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3a515125bSLeila Ghaffari 4a515125bSLeila Ghaffari /// @file 5a515125bSLeila Ghaffari /// Setup DM for Navier-Stokes example using PETSc 6a515125bSLeila Ghaffari 7e419654dSJeremy L Thompson #include <ceed.h> 8e419654dSJeremy L Thompson #include <petscdmplex.h> 967263decSKenneth E. Jansen #include <petscds.h> 10e419654dSJeremy L Thompson 11149fb536SJames Wright #include <navierstokes.h> 12866f23abSJames Wright #include "../problems/stg_shur14.h" 13a515125bSLeila Ghaffari 1405a512bdSLeila Ghaffari // Create mesh 15991aef52SJames Wright PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType mat_type, VecType vec_type, DM *dm) { 16a515125bSLeila Ghaffari PetscFunctionBeginUser; 1705a512bdSLeila Ghaffari // Create DMPLEX 182b916ea7SJeremy L Thompson PetscCall(DMCreate(comm, dm)); 192b916ea7SJeremy L Thompson PetscCall(DMSetType(*dm, DMPLEX)); 20b150a244SJed Brown { 21b150a244SJed Brown PetscBool skip = PETSC_TRUE; 222b916ea7SJeremy L Thompson PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL)); 23b150a244SJed Brown PetscCall(DMSetMatrixPreallocateSkip(*dm, skip)); 24b150a244SJed Brown } 252b916ea7SJeremy L Thompson PetscCall(DMSetMatType(*dm, mat_type)); 262b916ea7SJeremy L Thompson PetscCall(DMSetVecType(*dm, vec_type)); 27a9804fe9SJed Brown 2805a512bdSLeila Ghaffari // Set Tensor elements 2959572072SJames Wright PetscCall(HoneeOptionsSetValueDefault(NULL, "-dm_plex_simplex", "0")); 30*e747eef9SJames Wright PetscCall(HoneeOptionsSetValueDefault(NULL, "-dm_localize", "0")); // Delay localization until DMSetupByOrderEnd_FEM 31*e747eef9SJames Wright PetscCall(DMSetSparseLocalize(*dm, PETSC_FALSE)); 32*e747eef9SJames Wright PetscCall(DMSetBlockingType(*dm, DM_BLOCKING_FIELD_NODE)); 33ebfabadfSJames Wright 3405a512bdSLeila Ghaffari // Set CL options 352b916ea7SJeremy L Thompson PetscCall(DMSetFromOptions(*dm)); 362b916ea7SJeremy L Thompson PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 37d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 38a515125bSLeila Ghaffari } 39a515125bSLeila Ghaffari 40a515125bSLeila Ghaffari // Setup DM 41991aef52SJames Wright PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys) { 42a515125bSLeila Ghaffari PetscInt num_comp_q = 5; 43da4ca0cfSJames Wright PetscFunctionBeginUser; 44da4ca0cfSJames Wright 45b2e5b5b3SJames Wright PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, PETSC_DECIDE, q_extra, 1, &num_comp_q, dm)); 46da4ca0cfSJames Wright 47da4ca0cfSJames Wright { // Add strong boundary conditions to DM 48047c2946SJed Brown DMLabel label; 492b916ea7SJeremy L Thompson PetscCall(DMGetLabel(dm, "Face Sets", &label)); 505892ddc8SJames Wright PetscCheck(label, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DM does not have 'Face Sets' label"); 51c5e9980aSAdeleke O. Bankole PetscCall(DMPlexLabelComplete(dm, label)); 52487a3b6eSJames Wright 53487a3b6eSJames Wright for (PetscInt i = 0; i < problem->num_bc_defs; i++) { 54487a3b6eSJames Wright BCDefinition bc_def = problem->bc_defs[i]; 55487a3b6eSJames Wright PetscInt num_essential_comps, num_label_values; 56487a3b6eSJames Wright const PetscInt *essential_comps, *label_values; 57487a3b6eSJames Wright const char *name; 58487a3b6eSJames Wright 59487a3b6eSJames Wright PetscCall(BCDefinitionGetEssential(bc_def, &num_essential_comps, &essential_comps)); 60487a3b6eSJames Wright if (essential_comps > 0) { 61487a3b6eSJames Wright PetscCall(BCDefinitionGetInfo(bc_def, &name, &num_label_values, &label_values)); 62487a3b6eSJames Wright PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, name, label, num_label_values, label_values, 0, num_essential_comps, essential_comps, NULL, NULL, 63487a3b6eSJames Wright NULL, NULL)); 64002797a3SLeila Ghaffari } 65002797a3SLeila Ghaffari } 66866f23abSJames Wright { 67866f23abSJames Wright PetscBool use_strongstg = PETSC_FALSE; 682b916ea7SJeremy L Thompson PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL)); 6942454adaSJames Wright if (use_strongstg) PetscCall(SetupStrongStg(dm, bc, problem, phys)); 70866f23abSJames Wright } 7167263decSKenneth E. Jansen } 7267263decSKenneth E. Jansen 73da4ca0cfSJames Wright PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm)); 74cbe60e31SLeila Ghaffari 75a515125bSLeila Ghaffari // Empty name for conserved field (because there is only one field) 76a515125bSLeila Ghaffari PetscSection section; 772b916ea7SJeremy L Thompson PetscCall(DMGetLocalSection(dm, §ion)); 782b916ea7SJeremy L Thompson PetscCall(PetscSectionSetFieldName(section, 0, "")); 793636f6a4SJames Wright switch (phys->state_var) { 803636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 812b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density")); 822b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 1, "MomentumX")); 832b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 2, "MomentumY")); 842b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 3, "MomentumZ")); 85525869caSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "TotalEnergy")); 863636f6a4SJames Wright break; 873636f6a4SJames Wright 883636f6a4SJames Wright case STATEVAR_PRIMITIVE: 892b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure")); 902b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX")); 912b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY")); 922b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ")); 932b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature")); 943636f6a4SJames Wright break; 959b103f75SJames Wright 969b103f75SJames Wright case STATEVAR_ENTROPY: 979b103f75SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "EntropyDensity")); 989b103f75SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "EntropyMomentumX")); 999b103f75SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "EntropyMomentumY")); 1009b103f75SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "EntropyMomentumZ")); 1019b103f75SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "EntropyTotalEnergy")); 1029b103f75SJames Wright break; 103a515125bSLeila Ghaffari } 104d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 105a515125bSLeila Ghaffari } 106a515125bSLeila Ghaffari 107a515125bSLeila Ghaffari // Refine DM for high-order viz 108991aef52SJames 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