13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 377841947SLeila Ghaffari // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 577841947SLeila Ghaffari // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 777841947SLeila Ghaffari 877841947SLeila Ghaffari /// @file 977841947SLeila Ghaffari /// Setup DM for Navier-Stokes example using PETSc 1077841947SLeila Ghaffari 1177841947SLeila Ghaffari #include "../navierstokes.h" 12363b60e1SJames Wright #include "../problems/stg_shur14.h" 1377841947SLeila Ghaffari 141864f1c2SLeila Ghaffari // Create mesh 154ea65e7bSJed Brown PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, 164ea65e7bSJed Brown MatType mat_type, VecType vec_type, 174ea65e7bSJed Brown DM *dm) { 1877841947SLeila Ghaffari PetscErrorCode ierr; 1977841947SLeila Ghaffari PetscFunctionBeginUser; 201864f1c2SLeila Ghaffari // Create DMPLEX 211864f1c2SLeila Ghaffari ierr = DMCreate(comm, dm); CHKERRQ(ierr); 221864f1c2SLeila Ghaffari ierr = DMSetType(*dm, DMPLEX); CHKERRQ(ierr); 23c32b0260SJed Brown { 24c32b0260SJed Brown PetscBool skip = PETSC_TRUE; 25c32b0260SJed Brown PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, 26c32b0260SJed Brown NULL)); 27c32b0260SJed Brown PetscCall(DMSetMatrixPreallocateSkip(*dm, skip)); 28c32b0260SJed Brown } 294ea65e7bSJed Brown ierr = DMSetMatType(*dm, mat_type); CHKERRQ(ierr); 304ea65e7bSJed Brown ierr = DMSetVecType(*dm, vec_type); CHKERRQ(ierr); 314ea65e7bSJed Brown 321864f1c2SLeila Ghaffari // Set Tensor elements 331864f1c2SLeila Ghaffari ierr = PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"); CHKERRQ(ierr); 343796c488SJed Brown ierr = PetscOptionsSetValue(NULL, "-dm_sparse_localize", "0"); CHKERRQ(ierr); 351864f1c2SLeila Ghaffari // Set CL options 361864f1c2SLeila Ghaffari ierr = DMSetFromOptions(*dm); CHKERRQ(ierr); 3777841947SLeila Ghaffari ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr); 3877841947SLeila Ghaffari PetscFunctionReturn(0); 3977841947SLeila Ghaffari } 4077841947SLeila Ghaffari 4177841947SLeila Ghaffari // Setup DM 4277841947SLeila Ghaffari PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, 43c95f9967SJed Brown SimpleBC bc, Physics phys) { 4477841947SLeila Ghaffari PetscErrorCode ierr; 4577841947SLeila Ghaffari PetscFunctionBeginUser; 4677841947SLeila Ghaffari { 4777841947SLeila Ghaffari // Configure the finite element space and boundary conditions 4877841947SLeila Ghaffari PetscFE fe; 4977841947SLeila Ghaffari PetscInt num_comp_q = 5; 50496f2382SJed Brown DMLabel label; 5177841947SLeila Ghaffari ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q, 5277841947SLeila Ghaffari PETSC_FALSE, degree, PETSC_DECIDE, 5377841947SLeila Ghaffari &fe); CHKERRQ(ierr); 5477841947SLeila Ghaffari ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr); 5577841947SLeila Ghaffari ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr); 5677841947SLeila Ghaffari ierr = DMCreateDS(dm); CHKERRQ(ierr); 57496f2382SJed Brown ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr); 582fe7aee7SLeila Ghaffari // Set wall BCs 592fe7aee7SLeila Ghaffari if (bc->num_wall > 0) { 602fe7aee7SLeila Ghaffari ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 612fe7aee7SLeila Ghaffari bc->num_wall, bc->walls, 0, bc->num_comps, 622fe7aee7SLeila Ghaffari bc->wall_comps, (void(*)(void))problem->bc, 63c95f9967SJed Brown NULL, problem->bc_ctx, NULL); CHKERRQ(ierr); 642fe7aee7SLeila Ghaffari } 652fe7aee7SLeila Ghaffari // Set slip BCs in the x direction 662fe7aee7SLeila Ghaffari if (bc->num_slip[0] > 0) { 672fe7aee7SLeila Ghaffari PetscInt comps[1] = {1}; 682fe7aee7SLeila Ghaffari ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label, 692fe7aee7SLeila Ghaffari bc->num_slip[0], bc->slips[0], 0, 1, comps, 70c95f9967SJed Brown (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr); 712fe7aee7SLeila Ghaffari } 722fe7aee7SLeila Ghaffari // Set slip BCs in the y direction 732fe7aee7SLeila Ghaffari if (bc->num_slip[1] > 0) { 742fe7aee7SLeila Ghaffari PetscInt comps[1] = {2}; 752fe7aee7SLeila Ghaffari ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label, 762fe7aee7SLeila Ghaffari bc->num_slip[1], bc->slips[1], 0, 1, comps, 77c95f9967SJed Brown (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr); 782fe7aee7SLeila Ghaffari } 792fe7aee7SLeila Ghaffari // Set slip BCs in the z direction 802fe7aee7SLeila Ghaffari if (bc->num_slip[2] > 0) { 812fe7aee7SLeila Ghaffari PetscInt comps[1] = {3}; 822fe7aee7SLeila Ghaffari ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label, 832fe7aee7SLeila Ghaffari bc->num_slip[2], bc->slips[2], 0, 1, comps, 84c95f9967SJed Brown (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr); 852fe7aee7SLeila Ghaffari } 86363b60e1SJames Wright { 87363b60e1SJames Wright PetscBool use_strongstg = PETSC_FALSE; 88363b60e1SJames Wright ierr = PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL); 89363b60e1SJames Wright CHKERRQ(ierr); 90363b60e1SJames Wright 91363b60e1SJames Wright if (use_strongstg) { 92192a7459SJames Wright ierr = SetupStrongSTG(dm, bc, problem, phys); CHKERRQ(ierr); 93363b60e1SJames Wright } 94363b60e1SJames Wright } 95363b60e1SJames Wright 9677841947SLeila Ghaffari ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 9777841947SLeila Ghaffari CHKERRQ(ierr); 9877841947SLeila Ghaffari ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 9977841947SLeila Ghaffari } 100dc805cc4SLeila Ghaffari 10177841947SLeila Ghaffari // Empty name for conserved field (because there is only one field) 10277841947SLeila Ghaffari PetscSection section; 10377841947SLeila Ghaffari ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 10477841947SLeila Ghaffari ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr); 105*97baf651SJames Wright switch (phys->state_var) { 106*97baf651SJames Wright case STATEVAR_CONSERVATIVE: 10777841947SLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 0, "Density"); 10877841947SLeila Ghaffari CHKERRQ(ierr); 10977841947SLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 1, "Momentum X"); 11077841947SLeila Ghaffari CHKERRQ(ierr); 11177841947SLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 2, "Momentum Y"); 11277841947SLeila Ghaffari CHKERRQ(ierr); 11377841947SLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 3, "Momentum Z"); 11477841947SLeila Ghaffari CHKERRQ(ierr); 11577841947SLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 4, "Energy Density"); 11677841947SLeila Ghaffari CHKERRQ(ierr); 117*97baf651SJames Wright break; 118*97baf651SJames Wright 119*97baf651SJames Wright case STATEVAR_PRIMITIVE: 120*97baf651SJames Wright ierr = PetscSectionSetComponentName(section, 0, 0, "Pressure"); 121*97baf651SJames Wright CHKERRQ(ierr); 122*97baf651SJames Wright ierr = PetscSectionSetComponentName(section, 0, 1, "Velocity X"); 123*97baf651SJames Wright CHKERRQ(ierr); 124*97baf651SJames Wright ierr = PetscSectionSetComponentName(section, 0, 2, "Velocity Y"); 125*97baf651SJames Wright CHKERRQ(ierr); 126*97baf651SJames Wright ierr = PetscSectionSetComponentName(section, 0, 3, "Velocity Z"); 127*97baf651SJames Wright CHKERRQ(ierr); 128*97baf651SJames Wright ierr = PetscSectionSetComponentName(section, 0, 4, "Temperature"); 129*97baf651SJames Wright CHKERRQ(ierr); 130*97baf651SJames Wright break; 13177841947SLeila Ghaffari } 13277841947SLeila Ghaffari PetscFunctionReturn(0); 13377841947SLeila Ghaffari } 13477841947SLeila Ghaffari 13577841947SLeila Ghaffari // Refine DM for high-order viz 13677841947SLeila Ghaffari PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, 137c95f9967SJed Brown SimpleBC bc, Physics phys) { 13877841947SLeila Ghaffari PetscErrorCode ierr; 13977841947SLeila Ghaffari DM dm_hierarchy[user->app_ctx->viz_refine + 1]; 14077841947SLeila Ghaffari VecType vec_type; 14177841947SLeila Ghaffari PetscFunctionBeginUser; 14277841947SLeila Ghaffari 14377841947SLeila Ghaffari ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr); 14477841947SLeila Ghaffari 14577841947SLeila Ghaffari dm_hierarchy[0] = dm; 14677841947SLeila Ghaffari for (PetscInt i = 0, d = user->app_ctx->degree; 14777841947SLeila Ghaffari i < user->app_ctx->viz_refine; i++) { 14877841947SLeila Ghaffari Mat interp_next; 14977841947SLeila Ghaffari ierr = DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i+1]); 15077841947SLeila Ghaffari CHKERRQ(ierr); 15177841947SLeila Ghaffari ierr = DMClearDS(dm_hierarchy[i+1]); CHKERRQ(ierr); 15277841947SLeila Ghaffari ierr = DMClearFields(dm_hierarchy[i+1]); CHKERRQ(ierr); 15377841947SLeila Ghaffari ierr = DMSetCoarseDM(dm_hierarchy[i+1], dm_hierarchy[i]); CHKERRQ(ierr); 15477841947SLeila Ghaffari d = (d + 1) / 2; 15577841947SLeila Ghaffari if (i + 1 == user->app_ctx->viz_refine) d = 1; 15677841947SLeila Ghaffari ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr); 15777841947SLeila Ghaffari ierr = DMSetVecType(dm_hierarchy[i+1], vec_type); CHKERRQ(ierr); 158c95f9967SJed Brown ierr = SetUpDM(dm_hierarchy[i+1], problem, d, bc, phys); 15977841947SLeila Ghaffari CHKERRQ(ierr); 16077841947SLeila Ghaffari ierr = DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i+1], &interp_next, 16177841947SLeila Ghaffari NULL); CHKERRQ(ierr); 16277841947SLeila Ghaffari if (!i) user->interp_viz = interp_next; 16377841947SLeila Ghaffari else { 16477841947SLeila Ghaffari Mat C; 16577841947SLeila Ghaffari ierr = MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, 16677841947SLeila Ghaffari PETSC_DECIDE, &C); CHKERRQ(ierr); 16777841947SLeila Ghaffari ierr = MatDestroy(&interp_next); CHKERRQ(ierr); 16877841947SLeila Ghaffari ierr = MatDestroy(&user->interp_viz); CHKERRQ(ierr); 16977841947SLeila Ghaffari user->interp_viz = C; 17077841947SLeila Ghaffari } 17177841947SLeila Ghaffari } 17277841947SLeila Ghaffari for (PetscInt i=1; i<user->app_ctx->viz_refine; i++) { 17377841947SLeila Ghaffari ierr = DMDestroy(&dm_hierarchy[i]); CHKERRQ(ierr); 17477841947SLeila Ghaffari } 17577841947SLeila Ghaffari user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine]; 17677841947SLeila Ghaffari 17777841947SLeila Ghaffari PetscFunctionReturn(0); 17877841947SLeila Ghaffari } 179