1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 4 /// @file 5 /// Setup DM for Navier-Stokes example using PETSc 6 7 #include <ceed.h> 8 #include <petscdmplex.h> 9 #include <petscds.h> 10 11 #include <navierstokes.h> 12 #include "../problems/stg_shur14.h" 13 14 // Create mesh 15 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType mat_type, VecType vec_type, DM *dm) { 16 PetscBool isCGNS; 17 char filename[PETSC_MAX_PATH_LEN] = ""; 18 19 PetscFunctionBeginUser; 20 PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_filename", filename, sizeof(filename), NULL)); 21 PetscCall(HoneeCheckFilenameExtension(comm, filename, ".cgns", &isCGNS)); 22 23 if (isCGNS) { 24 // Must create from file to keep the sfNatural field in tact 25 PetscCall(DMPlexCreateFromFile(comm, filename, "HONEE", PETSC_TRUE, dm)); 26 PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_filename", "")); 27 } else { 28 PetscCall(DMCreate(comm, dm)); 29 PetscCall(DMSetType(*dm, DMPLEX)); 30 } 31 32 { 33 PetscBool skip = PETSC_TRUE; 34 PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL)); 35 PetscCall(DMSetMatrixPreallocateSkip(*dm, skip)); 36 } 37 PetscCall(DMSetMatType(*dm, mat_type)); 38 PetscCall(DMSetVecType(*dm, vec_type)); 39 PetscCall(HoneeOptionsSetValueDefault(NULL, "-dm_plex_simplex", "0")); 40 PetscCall(HoneeOptionsSetValueDefault(NULL, "-dm_localize", "0")); // Delay localization until DMSetupByOrderEnd_FEM 41 PetscCall(DMSetSparseLocalize(*dm, PETSC_FALSE)); 42 PetscCall(DMSetBlockingType(*dm, DM_BLOCKING_FIELD_NODE)); 43 44 // Set CL options 45 PetscCall(DMSetFromOptions(*dm)); 46 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 47 { // Force isoperiodic point SF to be created to update sfNatural. 48 // Needs to be done before removing the field corresponding to sfNatural, but after IsoperiodicFaces have been set 49 PetscSection dummy; 50 PetscCall(DMGetGlobalSection(*dm, &dummy)); 51 } 52 PetscFunctionReturn(PETSC_SUCCESS); 53 } 54 55 // Setup DM 56 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys) { 57 PetscInt num_comp_q = 5; 58 59 PetscFunctionBeginUser; 60 PetscCall(DMClearFields(dm)); 61 PetscCall(DMSetLocalSection(dm, NULL)); 62 PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, PETSC_DECIDE, q_extra, 1, &num_comp_q, dm)); 63 64 { // Add strong boundary conditions to DM 65 DMLabel label; 66 PetscCall(DMGetLabel(dm, "Face Sets", &label)); 67 PetscCheck(label, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DM does not have 'Face Sets' label"); 68 PetscCall(DMPlexLabelComplete(dm, label)); 69 70 for (PetscInt i = 0; i < problem->num_bc_defs; i++) { 71 BCDefinition bc_def = problem->bc_defs[i]; 72 PetscInt num_essential_comps, num_label_values; 73 const PetscInt *essential_comps, *label_values; 74 const char *name; 75 76 PetscCall(BCDefinitionGetEssential(bc_def, &num_essential_comps, &essential_comps)); 77 if (essential_comps > 0) { 78 PetscCall(BCDefinitionGetInfo(bc_def, &name, &num_label_values, &label_values)); 79 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, name, label, num_label_values, label_values, 0, num_essential_comps, essential_comps, NULL, NULL, 80 NULL, NULL)); 81 } 82 } 83 { 84 PetscBool use_strongstg = PETSC_FALSE; 85 PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL)); 86 if (use_strongstg) PetscCall(SetupStrongStg(dm, bc, problem, phys)); 87 } 88 } 89 90 PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm)); 91 92 // Empty name for conserved field (because there is only one field) 93 PetscSection section; 94 PetscCall(DMGetLocalSection(dm, §ion)); 95 PetscCall(PetscSectionSetFieldName(section, 0, "")); 96 switch (phys->state_var) { 97 case STATEVAR_CONSERVATIVE: 98 PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density")); 99 PetscCall(PetscSectionSetComponentName(section, 0, 1, "MomentumX")); 100 PetscCall(PetscSectionSetComponentName(section, 0, 2, "MomentumY")); 101 PetscCall(PetscSectionSetComponentName(section, 0, 3, "MomentumZ")); 102 PetscCall(PetscSectionSetComponentName(section, 0, 4, "TotalEnergy")); 103 break; 104 105 case STATEVAR_PRIMITIVE: 106 PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure")); 107 PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX")); 108 PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY")); 109 PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ")); 110 PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature")); 111 break; 112 113 case STATEVAR_ENTROPY: 114 PetscCall(PetscSectionSetComponentName(section, 0, 0, "EntropyDensity")); 115 PetscCall(PetscSectionSetComponentName(section, 0, 1, "EntropyMomentumX")); 116 PetscCall(PetscSectionSetComponentName(section, 0, 2, "EntropyMomentumY")); 117 PetscCall(PetscSectionSetComponentName(section, 0, 3, "EntropyMomentumZ")); 118 PetscCall(PetscSectionSetComponentName(section, 0, 4, "EntropyTotalEnergy")); 119 break; 120 } 121 PetscFunctionReturn(PETSC_SUCCESS); 122 } 123 124 // Refine DM for high-order viz 125 PetscErrorCode VizRefineDM(DM dm, Honee honee, ProblemData problem, SimpleBC bc, Physics phys) { 126 DM dm_hierarchy[honee->app_ctx->viz_refine + 1]; 127 VecType vec_type; 128 129 PetscFunctionBeginUser; 130 PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 131 132 dm_hierarchy[0] = dm; 133 for (PetscInt i = 0, d = honee->app_ctx->degree; i < honee->app_ctx->viz_refine; i++) { 134 Mat interp_next; 135 PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1])); 136 PetscCall(DMClearDS(dm_hierarchy[i + 1])); 137 PetscCall(DMClearFields(dm_hierarchy[i + 1])); 138 PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i])); 139 d = (d + 1) / 2; 140 PetscInt q_order = d + honee->app_ctx->q_extra; 141 if (i + 1 == honee->app_ctx->viz_refine) d = 1; 142 PetscCall(DMGetVecType(dm, &vec_type)); 143 PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type)); 144 PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, bc, phys)); 145 PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL)); 146 if (!i) honee->interp_viz = interp_next; 147 else { 148 Mat C; 149 PetscCall(MatMatMult(interp_next, honee->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C)); 150 PetscCall(MatDestroy(&interp_next)); 151 PetscCall(MatDestroy(&honee->interp_viz)); 152 honee->interp_viz = C; 153 } 154 } 155 for (PetscInt i = 1; i < honee->app_ctx->viz_refine; i++) { 156 PetscCall(DMDestroy(&dm_hierarchy[i])); 157 } 158 honee->dm_viz = dm_hierarchy[honee->app_ctx->viz_refine]; 159 PetscFunctionReturn(PETSC_SUCCESS); 160 } 161