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 5ea615d4cSJames Wright /// Setup DM for HONEE 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) { 16eb9b4fe1SJames Wright PetscBool isCGNS; 17eb9b4fe1SJames Wright char filename[PETSC_MAX_PATH_LEN] = ""; 18eb9b4fe1SJames Wright 19a515125bSLeila Ghaffari PetscFunctionBeginUser; 20eb9b4fe1SJames Wright PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_filename", filename, sizeof(filename), NULL)); 21eb9b4fe1SJames Wright PetscCall(HoneeCheckFilenameExtension(comm, filename, ".cgns", &isCGNS)); 22eb9b4fe1SJames Wright 23eb9b4fe1SJames Wright if (isCGNS) { 24eb9b4fe1SJames Wright // Must create from file to keep the sfNatural field in tact 25eb9b4fe1SJames Wright PetscCall(DMPlexCreateFromFile(comm, filename, "HONEE", PETSC_TRUE, dm)); 26eb9b4fe1SJames Wright PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_filename", "")); 27eb9b4fe1SJames Wright } else { 282b916ea7SJeremy L Thompson PetscCall(DMCreate(comm, dm)); 292b916ea7SJeremy L Thompson PetscCall(DMSetType(*dm, DMPLEX)); 30eb9b4fe1SJames Wright } 31eb9b4fe1SJames Wright 32b150a244SJed Brown { 33b150a244SJed Brown PetscBool skip = PETSC_TRUE; 342b916ea7SJeremy L Thompson PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL)); 35b150a244SJed Brown PetscCall(DMSetMatrixPreallocateSkip(*dm, skip)); 36b150a244SJed Brown } 372b916ea7SJeremy L Thompson PetscCall(DMSetMatType(*dm, mat_type)); 382b916ea7SJeremy L Thompson PetscCall(DMSetVecType(*dm, vec_type)); 3959572072SJames Wright PetscCall(HoneeOptionsSetValueDefault(NULL, "-dm_plex_simplex", "0")); 40e747eef9SJames Wright PetscCall(HoneeOptionsSetValueDefault(NULL, "-dm_localize", "0")); // Delay localization until DMSetupByOrderEnd_FEM 41*517ae6cdSJames Wright { // Use Mat graph generation for partitioning. See https://gitlab.com/petsc/petsc/-/merge_requests/8336 42*517ae6cdSJames Wright PetscBool distribute = PETSC_FALSE; 43*517ae6cdSJames Wright PetscMPIInt num_ranks; 44*517ae6cdSJames Wright PetscCall(DMPlexDistributeGetDefault(*dm, &distribute)); 45*517ae6cdSJames Wright PetscCallMPI(MPI_Comm_size(comm, &num_ranks)); 46*517ae6cdSJames Wright if (distribute && num_ranks > 1) PetscCall(HoneeOptionsSetValueDefault(NULL, "-dm_plex_csr_alg", "mat")); 47*517ae6cdSJames Wright } 48e747eef9SJames Wright PetscCall(DMSetSparseLocalize(*dm, PETSC_FALSE)); 49e747eef9SJames Wright PetscCall(DMSetBlockingType(*dm, DM_BLOCKING_FIELD_NODE)); 50e7934809SJames Wright PetscCall(DMPlexSetPartitionBalance(*dm, PETSC_TRUE)); 51ebfabadfSJames Wright 5205a512bdSLeila Ghaffari // Set CL options 532b916ea7SJeremy L Thompson PetscCall(DMSetFromOptions(*dm)); 542b916ea7SJeremy L Thompson PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 55eb9b4fe1SJames Wright { // Force isoperiodic point SF to be created to update sfNatural. 56eb9b4fe1SJames Wright // Needs to be done before removing the field corresponding to sfNatural, but after IsoperiodicFaces have been set 57c4a0f6c7SJames Wright PetscInt num_fields; 58c4a0f6c7SJames Wright PetscCall(DMGetNumFields(*dm, &num_fields)); 59c4a0f6c7SJames Wright if (num_fields) { 60eb9b4fe1SJames Wright PetscSection dummy; 61eb9b4fe1SJames Wright PetscCall(DMGetGlobalSection(*dm, &dummy)); 62eb9b4fe1SJames Wright } 63c4a0f6c7SJames Wright } 64d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 65a515125bSLeila Ghaffari } 66a515125bSLeila Ghaffari 67a515125bSLeila Ghaffari // Setup DM 68d3c60affSJames Wright PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, Physics phys) { 69a515125bSLeila Ghaffari PetscInt num_comp_q = 5; 70da4ca0cfSJames Wright 71eb9b4fe1SJames Wright PetscFunctionBeginUser; 72eb9b4fe1SJames Wright PetscCall(DMClearFields(dm)); 73eb9b4fe1SJames Wright PetscCall(DMSetLocalSection(dm, NULL)); 74b2e5b5b3SJames Wright PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, PETSC_DECIDE, q_extra, 1, &num_comp_q, dm)); 75da4ca0cfSJames Wright 76d6cac220SJames Wright { 77d6cac220SJames Wright PetscBool use_strongstg = PETSC_FALSE; 78d6cac220SJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL)); 79d6cac220SJames Wright if (use_strongstg) PetscCall(SetupStrongStg(dm, problem, phys)); 80d6cac220SJames Wright } 81d6cac220SJames Wright 82da4ca0cfSJames Wright { // Add strong boundary conditions to DM 83047c2946SJed Brown DMLabel label; 842b916ea7SJeremy L Thompson PetscCall(DMGetLabel(dm, "Face Sets", &label)); 855892ddc8SJames Wright PetscCheck(label, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DM does not have 'Face Sets' label"); 86c5e9980aSAdeleke O. Bankole PetscCall(DMPlexLabelComplete(dm, label)); 87487a3b6eSJames Wright 88487a3b6eSJames Wright for (PetscInt i = 0; i < problem->num_bc_defs; i++) { 89487a3b6eSJames Wright BCDefinition bc_def = problem->bc_defs[i]; 90487a3b6eSJames Wright PetscInt num_essential_comps, num_label_values; 91487a3b6eSJames Wright const PetscInt *essential_comps, *label_values; 92487a3b6eSJames Wright const char *name; 93487a3b6eSJames Wright 9409240e0aSJames Wright PetscCall(BCDefinitionSetDM(problem->bc_defs[i], dm)); 95487a3b6eSJames Wright PetscCall(BCDefinitionGetEssential(bc_def, &num_essential_comps, &essential_comps)); 96487a3b6eSJames Wright if (essential_comps > 0) { 97487a3b6eSJames Wright PetscCall(BCDefinitionGetInfo(bc_def, &name, &num_label_values, &label_values)); 98487a3b6eSJames Wright PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, name, label, num_label_values, label_values, 0, num_essential_comps, essential_comps, NULL, NULL, 99487a3b6eSJames Wright NULL, NULL)); 100002797a3SLeila Ghaffari } 101002797a3SLeila Ghaffari } 10267263decSKenneth E. Jansen } 10367263decSKenneth E. Jansen 104da4ca0cfSJames Wright PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm)); 105cbe60e31SLeila Ghaffari 106a515125bSLeila Ghaffari // Empty name for conserved field (because there is only one field) 107a515125bSLeila Ghaffari PetscSection section; 1082b916ea7SJeremy L Thompson PetscCall(DMGetLocalSection(dm, §ion)); 1092b916ea7SJeremy L Thompson PetscCall(PetscSectionSetFieldName(section, 0, "")); 1103636f6a4SJames Wright switch (phys->state_var) { 1113636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 1122b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density")); 1132b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 1, "MomentumX")); 1142b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 2, "MomentumY")); 1152b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 3, "MomentumZ")); 116525869caSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "TotalEnergy")); 1173636f6a4SJames Wright break; 1183636f6a4SJames Wright 1193636f6a4SJames Wright case STATEVAR_PRIMITIVE: 1202b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure")); 1212b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX")); 1222b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY")); 1232b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ")); 1242b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature")); 1253636f6a4SJames Wright break; 1269b103f75SJames Wright 1279b103f75SJames Wright case STATEVAR_ENTROPY: 1289b103f75SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "EntropyDensity")); 1299b103f75SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "EntropyMomentumX")); 1309b103f75SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "EntropyMomentumY")); 1319b103f75SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "EntropyMomentumZ")); 1329b103f75SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "EntropyTotalEnergy")); 1339b103f75SJames Wright break; 134a515125bSLeila Ghaffari } 135d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 136a515125bSLeila Ghaffari } 137a515125bSLeila Ghaffari 138a515125bSLeila Ghaffari // Refine DM for high-order viz 139d3c60affSJames Wright PetscErrorCode VizRefineDM(DM dm, Honee honee, ProblemData problem, Physics phys) { 1400c373b74SJames Wright DM dm_hierarchy[honee->app_ctx->viz_refine + 1]; 141a515125bSLeila Ghaffari VecType vec_type; 142a515125bSLeila Ghaffari 14306f41313SJames Wright PetscFunctionBeginUser; 1442b916ea7SJeremy L Thompson PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 145a515125bSLeila Ghaffari 146a515125bSLeila Ghaffari dm_hierarchy[0] = dm; 1470c373b74SJames Wright for (PetscInt i = 0, d = honee->app_ctx->degree; i < honee->app_ctx->viz_refine; i++) { 148a515125bSLeila Ghaffari Mat interp_next; 1492b916ea7SJeremy L Thompson PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1])); 1502b916ea7SJeremy L Thompson PetscCall(DMClearDS(dm_hierarchy[i + 1])); 1512b916ea7SJeremy L Thompson PetscCall(DMClearFields(dm_hierarchy[i + 1])); 1522b916ea7SJeremy L Thompson PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i])); 153a515125bSLeila Ghaffari d = (d + 1) / 2; 1540c373b74SJames Wright PetscInt q_order = d + honee->app_ctx->q_extra; 1550c373b74SJames Wright if (i + 1 == honee->app_ctx->viz_refine) d = 1; 1562b916ea7SJeremy L Thompson PetscCall(DMGetVecType(dm, &vec_type)); 1572b916ea7SJeremy L Thompson PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type)); 158d3c60affSJames Wright PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, phys)); 1592b916ea7SJeremy L Thompson PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL)); 1600c373b74SJames Wright if (!i) honee->interp_viz = interp_next; 161a515125bSLeila Ghaffari else { 162a515125bSLeila Ghaffari Mat C; 1630c373b74SJames Wright PetscCall(MatMatMult(interp_next, honee->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C)); 1642b916ea7SJeremy L Thompson PetscCall(MatDestroy(&interp_next)); 1650c373b74SJames Wright PetscCall(MatDestroy(&honee->interp_viz)); 1660c373b74SJames Wright honee->interp_viz = C; 167a515125bSLeila Ghaffari } 168a515125bSLeila Ghaffari } 1690c373b74SJames Wright for (PetscInt i = 1; i < honee->app_ctx->viz_refine; i++) { 1702b916ea7SJeremy L Thompson PetscCall(DMDestroy(&dm_hierarchy[i])); 171a515125bSLeila Ghaffari } 1720c373b74SJames Wright honee->dm_viz = dm_hierarchy[honee->app_ctx->viz_refine]; 173d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 174a515125bSLeila Ghaffari } 175