1*3d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2*3d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 377841947SLeila Ghaffari // 4*3d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 577841947SLeila Ghaffari // 6*3d8e8822SJeremy 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" 1277841947SLeila Ghaffari 131864f1c2SLeila Ghaffari // Create mesh 141864f1c2SLeila Ghaffari PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, DM *dm) { 1577841947SLeila Ghaffari PetscErrorCode ierr; 1677841947SLeila Ghaffari PetscFunctionBeginUser; 171864f1c2SLeila Ghaffari // Create DMPLEX 181864f1c2SLeila Ghaffari ierr = DMCreate(comm, dm); CHKERRQ(ierr); 191864f1c2SLeila Ghaffari ierr = DMSetType(*dm, DMPLEX); CHKERRQ(ierr); 201864f1c2SLeila Ghaffari // Set Tensor elements 211864f1c2SLeila Ghaffari ierr = PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"); CHKERRQ(ierr); 221864f1c2SLeila Ghaffari // Set CL options 231864f1c2SLeila Ghaffari ierr = DMSetFromOptions(*dm); CHKERRQ(ierr); 2477841947SLeila Ghaffari ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr); 2577841947SLeila Ghaffari PetscFunctionReturn(0); 2677841947SLeila Ghaffari } 2777841947SLeila Ghaffari 2877841947SLeila Ghaffari // Setup DM 2977841947SLeila Ghaffari PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, 3077841947SLeila Ghaffari SimpleBC bc, Physics phys, void *setup_ctx) { 3177841947SLeila Ghaffari PetscErrorCode ierr; 3277841947SLeila Ghaffari PetscFunctionBeginUser; 3377841947SLeila Ghaffari { 3477841947SLeila Ghaffari // Configure the finite element space and boundary conditions 3577841947SLeila Ghaffari PetscFE fe; 3677841947SLeila Ghaffari PetscInt num_comp_q = 5; 37496f2382SJed Brown DMLabel label; 3877841947SLeila Ghaffari ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q, 3977841947SLeila Ghaffari PETSC_FALSE, degree, PETSC_DECIDE, 4077841947SLeila Ghaffari &fe); CHKERRQ(ierr); 4177841947SLeila Ghaffari ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr); 4277841947SLeila Ghaffari ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr); 4377841947SLeila Ghaffari ierr = DMCreateDS(dm); CHKERRQ(ierr); 447ed3e4cdSJeremy L Thompson { 457ed3e4cdSJeremy L Thompson /* create FE field for coordinates */ 467ed3e4cdSJeremy L Thompson PetscFE fe_coords; 477ed3e4cdSJeremy L Thompson PetscInt num_comp_coord; 487ed3e4cdSJeremy L Thompson ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr); 497ed3e4cdSJeremy L Thompson ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_coord, 507ed3e4cdSJeremy L Thompson PETSC_FALSE, 1, 1, &fe_coords); CHKERRQ(ierr); 517ed3e4cdSJeremy L Thompson ierr = DMProjectCoordinates(dm, fe_coords); CHKERRQ(ierr); 527ed3e4cdSJeremy L Thompson ierr = PetscFEDestroy(&fe_coords); CHKERRQ(ierr); 537ed3e4cdSJeremy L Thompson } 54496f2382SJed Brown ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr); 552fe7aee7SLeila Ghaffari // Set wall BCs 562fe7aee7SLeila Ghaffari if (bc->num_wall > 0) { 572fe7aee7SLeila Ghaffari ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 582fe7aee7SLeila Ghaffari bc->num_wall, bc->walls, 0, bc->num_comps, 592fe7aee7SLeila Ghaffari bc->wall_comps, (void(*)(void))problem->bc, 602fe7aee7SLeila Ghaffari NULL, setup_ctx, NULL); CHKERRQ(ierr); 612fe7aee7SLeila Ghaffari } 622fe7aee7SLeila Ghaffari // Set slip BCs in the x direction 632fe7aee7SLeila Ghaffari if (bc->num_slip[0] > 0) { 642fe7aee7SLeila Ghaffari PetscInt comps[1] = {1}; 652fe7aee7SLeila Ghaffari ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label, 662fe7aee7SLeila Ghaffari bc->num_slip[0], bc->slips[0], 0, 1, comps, 672fe7aee7SLeila Ghaffari (void(*)(void))NULL, NULL, setup_ctx, NULL); CHKERRQ(ierr); 682fe7aee7SLeila Ghaffari } 692fe7aee7SLeila Ghaffari // Set slip BCs in the y direction 702fe7aee7SLeila Ghaffari if (bc->num_slip[1] > 0) { 712fe7aee7SLeila Ghaffari PetscInt comps[1] = {2}; 722fe7aee7SLeila Ghaffari ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label, 732fe7aee7SLeila Ghaffari bc->num_slip[1], bc->slips[1], 0, 1, comps, 742fe7aee7SLeila Ghaffari (void(*)(void))NULL, NULL, setup_ctx, NULL); CHKERRQ(ierr); 752fe7aee7SLeila Ghaffari } 762fe7aee7SLeila Ghaffari // Set slip BCs in the z direction 772fe7aee7SLeila Ghaffari if (bc->num_slip[2] > 0) { 782fe7aee7SLeila Ghaffari PetscInt comps[1] = {3}; 792fe7aee7SLeila Ghaffari ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label, 802fe7aee7SLeila Ghaffari bc->num_slip[2], bc->slips[2], 0, 1, comps, 812fe7aee7SLeila Ghaffari (void(*)(void))NULL, NULL, setup_ctx, NULL); CHKERRQ(ierr); 822fe7aee7SLeila Ghaffari } 8377841947SLeila Ghaffari ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 8477841947SLeila Ghaffari CHKERRQ(ierr); 8577841947SLeila Ghaffari ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 8677841947SLeila Ghaffari } 8777841947SLeila Ghaffari { 8877841947SLeila Ghaffari // Empty name for conserved field (because there is only one field) 8977841947SLeila Ghaffari PetscSection section; 9077841947SLeila Ghaffari ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 9177841947SLeila Ghaffari ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr); 9277841947SLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 0, "Density"); 9377841947SLeila Ghaffari CHKERRQ(ierr); 9477841947SLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 1, "Momentum X"); 9577841947SLeila Ghaffari CHKERRQ(ierr); 9677841947SLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 2, "Momentum Y"); 9777841947SLeila Ghaffari CHKERRQ(ierr); 9877841947SLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 3, "Momentum Z"); 9977841947SLeila Ghaffari CHKERRQ(ierr); 10077841947SLeila Ghaffari ierr = PetscSectionSetComponentName(section, 0, 4, "Energy Density"); 10177841947SLeila Ghaffari CHKERRQ(ierr); 10277841947SLeila Ghaffari } 10377841947SLeila Ghaffari PetscFunctionReturn(0); 10477841947SLeila Ghaffari } 10577841947SLeila Ghaffari 10677841947SLeila Ghaffari // Refine DM for high-order viz 10777841947SLeila Ghaffari PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, 10877841947SLeila Ghaffari SimpleBC bc, Physics phys, void *setup_ctx) { 10977841947SLeila Ghaffari PetscErrorCode ierr; 11077841947SLeila Ghaffari DM dm_hierarchy[user->app_ctx->viz_refine + 1]; 11177841947SLeila Ghaffari VecType vec_type; 11277841947SLeila Ghaffari PetscFunctionBeginUser; 11377841947SLeila Ghaffari 11477841947SLeila Ghaffari ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr); 11577841947SLeila Ghaffari 11677841947SLeila Ghaffari dm_hierarchy[0] = dm; 11777841947SLeila Ghaffari for (PetscInt i = 0, d = user->app_ctx->degree; 11877841947SLeila Ghaffari i < user->app_ctx->viz_refine; i++) { 11977841947SLeila Ghaffari Mat interp_next; 12077841947SLeila Ghaffari ierr = DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i+1]); 12177841947SLeila Ghaffari CHKERRQ(ierr); 12277841947SLeila Ghaffari ierr = DMClearDS(dm_hierarchy[i+1]); CHKERRQ(ierr); 12377841947SLeila Ghaffari ierr = DMClearFields(dm_hierarchy[i+1]); CHKERRQ(ierr); 12477841947SLeila Ghaffari ierr = DMSetCoarseDM(dm_hierarchy[i+1], dm_hierarchy[i]); CHKERRQ(ierr); 12577841947SLeila Ghaffari d = (d + 1) / 2; 12677841947SLeila Ghaffari if (i + 1 == user->app_ctx->viz_refine) d = 1; 12777841947SLeila Ghaffari ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr); 12877841947SLeila Ghaffari ierr = DMSetVecType(dm_hierarchy[i+1], vec_type); CHKERRQ(ierr); 12977841947SLeila Ghaffari ierr = SetUpDM(dm_hierarchy[i+1], problem, d, bc, phys, setup_ctx); 13077841947SLeila Ghaffari CHKERRQ(ierr); 13177841947SLeila Ghaffari ierr = DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i+1], &interp_next, 13277841947SLeila Ghaffari NULL); CHKERRQ(ierr); 13377841947SLeila Ghaffari if (!i) user->interp_viz = interp_next; 13477841947SLeila Ghaffari else { 13577841947SLeila Ghaffari Mat C; 13677841947SLeila Ghaffari ierr = MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, 13777841947SLeila Ghaffari PETSC_DECIDE, &C); CHKERRQ(ierr); 13877841947SLeila Ghaffari ierr = MatDestroy(&interp_next); CHKERRQ(ierr); 13977841947SLeila Ghaffari ierr = MatDestroy(&user->interp_viz); CHKERRQ(ierr); 14077841947SLeila Ghaffari user->interp_viz = C; 14177841947SLeila Ghaffari } 14277841947SLeila Ghaffari } 14377841947SLeila Ghaffari for (PetscInt i=1; i<user->app_ctx->viz_refine; i++) { 14477841947SLeila Ghaffari ierr = DMDestroy(&dm_hierarchy[i]); CHKERRQ(ierr); 14577841947SLeila Ghaffari } 14677841947SLeila Ghaffari user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine]; 14777841947SLeila Ghaffari 14877841947SLeila Ghaffari PetscFunctionReturn(0); 14977841947SLeila Ghaffari } 150