1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, 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 11a515125bSLeila Ghaffari #include "../navierstokes.h" 12866f23abSJames Wright #include "../problems/stg_shur14.h" 13a515125bSLeila Ghaffari 1405a512bdSLeila Ghaffari // Create mesh 15a9804fe9SJed Brown PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, 16a9804fe9SJed Brown MatType mat_type, VecType vec_type, 17a9804fe9SJed Brown DM *dm) { 18a515125bSLeila Ghaffari PetscErrorCode ierr; 19a515125bSLeila Ghaffari PetscFunctionBeginUser; 2005a512bdSLeila Ghaffari // Create DMPLEX 2105a512bdSLeila Ghaffari ierr = DMCreate(comm, dm); CHKERRQ(ierr); 2205a512bdSLeila Ghaffari ierr = DMSetType(*dm, DMPLEX); CHKERRQ(ierr); 23a9804fe9SJed Brown ierr = DMSetMatType(*dm, mat_type); CHKERRQ(ierr); 24a9804fe9SJed Brown ierr = DMSetVecType(*dm, vec_type); CHKERRQ(ierr); 25a9804fe9SJed Brown 2605a512bdSLeila Ghaffari // Set Tensor elements 2705a512bdSLeila Ghaffari ierr = PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"); CHKERRQ(ierr); 2805a512bdSLeila Ghaffari // Set CL options 2905a512bdSLeila Ghaffari ierr = DMSetFromOptions(*dm); CHKERRQ(ierr); 30a515125bSLeila Ghaffari ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr); 31a515125bSLeila Ghaffari PetscFunctionReturn(0); 32a515125bSLeila Ghaffari } 33a515125bSLeila Ghaffari 34a515125bSLeila Ghaffari // Setup DM 35a515125bSLeila Ghaffari PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, 36c848b8b7SJed Brown SimpleBC bc, Physics phys) { 37a515125bSLeila Ghaffari PetscErrorCode ierr; 38a515125bSLeila Ghaffari PetscFunctionBeginUser; 39a515125bSLeila Ghaffari { 40a515125bSLeila Ghaffari // Configure the finite element space and boundary conditions 41a515125bSLeila Ghaffari PetscFE fe; 42a515125bSLeila Ghaffari PetscInt num_comp_q = 5; 43047c2946SJed Brown DMLabel label; 44a515125bSLeila Ghaffari ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q, 45a515125bSLeila Ghaffari PETSC_FALSE, degree, PETSC_DECIDE, 46a515125bSLeila Ghaffari &fe); CHKERRQ(ierr); 47a515125bSLeila Ghaffari ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr); 48a515125bSLeila Ghaffari ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr); 49a515125bSLeila Ghaffari ierr = DMCreateDS(dm); CHKERRQ(ierr); 50047c2946SJed Brown ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr); 51002797a3SLeila Ghaffari // Set wall BCs 52002797a3SLeila Ghaffari if (bc->num_wall > 0) { 53002797a3SLeila Ghaffari ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 54002797a3SLeila Ghaffari bc->num_wall, bc->walls, 0, bc->num_comps, 55002797a3SLeila Ghaffari bc->wall_comps, (void(*)(void))problem->bc, 56c848b8b7SJed Brown NULL, problem->bc_ctx, NULL); CHKERRQ(ierr); 57002797a3SLeila Ghaffari } 58002797a3SLeila Ghaffari // Set slip BCs in the x direction 59002797a3SLeila Ghaffari if (bc->num_slip[0] > 0) { 60002797a3SLeila Ghaffari PetscInt comps[1] = {1}; 61002797a3SLeila Ghaffari ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label, 62002797a3SLeila Ghaffari bc->num_slip[0], bc->slips[0], 0, 1, comps, 63c848b8b7SJed Brown (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr); 64002797a3SLeila Ghaffari } 65002797a3SLeila Ghaffari // Set slip BCs in the y direction 66002797a3SLeila Ghaffari if (bc->num_slip[1] > 0) { 67002797a3SLeila Ghaffari PetscInt comps[1] = {2}; 68002797a3SLeila Ghaffari ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label, 69002797a3SLeila Ghaffari bc->num_slip[1], bc->slips[1], 0, 1, comps, 70c848b8b7SJed Brown (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr); 71002797a3SLeila Ghaffari } 72002797a3SLeila Ghaffari // Set slip BCs in the z direction 73002797a3SLeila Ghaffari if (bc->num_slip[2] > 0) { 74002797a3SLeila Ghaffari PetscInt comps[1] = {3}; 75002797a3SLeila Ghaffari ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label, 76002797a3SLeila Ghaffari bc->num_slip[2], bc->slips[2], 0, 1, comps, 77c848b8b7SJed Brown (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr); 78002797a3SLeila Ghaffari } 79866f23abSJames Wright { 80866f23abSJames Wright PetscBool use_strongstg = PETSC_FALSE; 81866f23abSJames Wright ierr = PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL); 82866f23abSJames Wright CHKERRQ(ierr); 83866f23abSJames Wright 84866f23abSJames Wright if (use_strongstg) { 85*8085925cSJames Wright ierr = SetupStrongSTG(dm, bc, problem); CHKERRQ(ierr); 86866f23abSJames Wright } 87866f23abSJames Wright } 88866f23abSJames Wright 89a515125bSLeila Ghaffari ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 90a515125bSLeila Ghaffari CHKERRQ(ierr); 91a515125bSLeila Ghaffari ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 92a515125bSLeila Ghaffari } 93a515125bSLeila Ghaffari { 94a515125bSLeila Ghaffari // Empty name for conserved field (because there is only one field) 95a515125bSLeila Ghaffari PetscSection section; 96a515125bSLeila Ghaffari ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 97a515125bSLeila Ghaffari ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr); 98a515125bSLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 0, "Density"); 99a515125bSLeila Ghaffari CHKERRQ(ierr); 100a515125bSLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 1, "Momentum X"); 101a515125bSLeila Ghaffari CHKERRQ(ierr); 102a515125bSLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 2, "Momentum Y"); 103a515125bSLeila Ghaffari CHKERRQ(ierr); 104a515125bSLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 3, "Momentum Z"); 105a515125bSLeila Ghaffari CHKERRQ(ierr); 106a515125bSLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 4, "Energy Density"); 107a515125bSLeila Ghaffari CHKERRQ(ierr); 108a515125bSLeila Ghaffari } 109a515125bSLeila Ghaffari PetscFunctionReturn(0); 110a515125bSLeila Ghaffari } 111a515125bSLeila Ghaffari 112a515125bSLeila Ghaffari // Refine DM for high-order viz 113a515125bSLeila Ghaffari PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, 114c848b8b7SJed Brown SimpleBC bc, Physics phys) { 115a515125bSLeila Ghaffari PetscErrorCode ierr; 116a515125bSLeila Ghaffari DM dm_hierarchy[user->app_ctx->viz_refine + 1]; 117a515125bSLeila Ghaffari VecType vec_type; 118a515125bSLeila Ghaffari PetscFunctionBeginUser; 119a515125bSLeila Ghaffari 120a515125bSLeila Ghaffari ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr); 121a515125bSLeila Ghaffari 122a515125bSLeila Ghaffari dm_hierarchy[0] = dm; 123a515125bSLeila Ghaffari for (PetscInt i = 0, d = user->app_ctx->degree; 124a515125bSLeila Ghaffari i < user->app_ctx->viz_refine; i++) { 125a515125bSLeila Ghaffari Mat interp_next; 126a515125bSLeila Ghaffari ierr = DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i+1]); 127a515125bSLeila Ghaffari CHKERRQ(ierr); 128a515125bSLeila Ghaffari ierr = DMClearDS(dm_hierarchy[i+1]); CHKERRQ(ierr); 129a515125bSLeila Ghaffari ierr = DMClearFields(dm_hierarchy[i+1]); CHKERRQ(ierr); 130a515125bSLeila Ghaffari ierr = DMSetCoarseDM(dm_hierarchy[i+1], dm_hierarchy[i]); CHKERRQ(ierr); 131a515125bSLeila Ghaffari d = (d + 1) / 2; 132a515125bSLeila Ghaffari if (i + 1 == user->app_ctx->viz_refine) d = 1; 133a515125bSLeila Ghaffari ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr); 134a515125bSLeila Ghaffari ierr = DMSetVecType(dm_hierarchy[i+1], vec_type); CHKERRQ(ierr); 135c848b8b7SJed Brown ierr = SetUpDM(dm_hierarchy[i+1], problem, d, bc, phys); 136a515125bSLeila Ghaffari CHKERRQ(ierr); 137a515125bSLeila Ghaffari ierr = DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i+1], &interp_next, 138a515125bSLeila Ghaffari NULL); CHKERRQ(ierr); 139a515125bSLeila Ghaffari if (!i) user->interp_viz = interp_next; 140a515125bSLeila Ghaffari else { 141a515125bSLeila Ghaffari Mat C; 142a515125bSLeila Ghaffari ierr = MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, 143a515125bSLeila Ghaffari PETSC_DECIDE, &C); CHKERRQ(ierr); 144a515125bSLeila Ghaffari ierr = MatDestroy(&interp_next); CHKERRQ(ierr); 145a515125bSLeila Ghaffari ierr = MatDestroy(&user->interp_viz); CHKERRQ(ierr); 146a515125bSLeila Ghaffari user->interp_viz = C; 147a515125bSLeila Ghaffari } 148a515125bSLeila Ghaffari } 149a515125bSLeila Ghaffari for (PetscInt i=1; i<user->app_ctx->viz_refine; i++) { 150a515125bSLeila Ghaffari ierr = DMDestroy(&dm_hierarchy[i]); CHKERRQ(ierr); 151a515125bSLeila Ghaffari } 152a515125bSLeila Ghaffari user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine]; 153a515125bSLeila Ghaffari 154a515125bSLeila Ghaffari PetscFunctionReturn(0); 155a515125bSLeila Ghaffari } 156