xref: /libCEED/examples/fluids/src/setupdm.c (revision 363b60e10ca170e7545568eddd43e64f5967d744)
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"
12*363b60e1SJames Wright #include "../problems/stg_shur14.h"
1377841947SLeila Ghaffari 
141864f1c2SLeila Ghaffari // Create mesh
151864f1c2SLeila Ghaffari PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, DM *dm) {
1677841947SLeila Ghaffari   PetscErrorCode   ierr;
1777841947SLeila Ghaffari   PetscFunctionBeginUser;
181864f1c2SLeila Ghaffari   // Create DMPLEX
191864f1c2SLeila Ghaffari   ierr = DMCreate(comm, dm); CHKERRQ(ierr);
201864f1c2SLeila Ghaffari   ierr = DMSetType(*dm, DMPLEX); CHKERRQ(ierr);
211864f1c2SLeila Ghaffari   // Set Tensor elements
221864f1c2SLeila Ghaffari   ierr = PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"); CHKERRQ(ierr);
231864f1c2SLeila Ghaffari   // Set CL options
241864f1c2SLeila Ghaffari   ierr = DMSetFromOptions(*dm); CHKERRQ(ierr);
2577841947SLeila Ghaffari   ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr);
2677841947SLeila Ghaffari   PetscFunctionReturn(0);
2777841947SLeila Ghaffari }
2877841947SLeila Ghaffari 
2977841947SLeila Ghaffari // Setup DM
3077841947SLeila Ghaffari PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree,
31c95f9967SJed Brown                        SimpleBC bc, Physics phys) {
3277841947SLeila Ghaffari   PetscErrorCode ierr;
3377841947SLeila Ghaffari   PetscFunctionBeginUser;
3477841947SLeila Ghaffari   {
3577841947SLeila Ghaffari     // Configure the finite element space and boundary conditions
3677841947SLeila Ghaffari     PetscFE  fe;
3777841947SLeila Ghaffari     PetscInt num_comp_q = 5;
38496f2382SJed Brown     DMLabel label;
3977841947SLeila Ghaffari     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q,
4077841947SLeila Ghaffari                                  PETSC_FALSE, degree, PETSC_DECIDE,
4177841947SLeila Ghaffari                                  &fe); CHKERRQ(ierr);
4277841947SLeila Ghaffari     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
4377841947SLeila Ghaffari     ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr);
4477841947SLeila Ghaffari     ierr = DMCreateDS(dm); CHKERRQ(ierr);
45496f2382SJed Brown     ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr);
462fe7aee7SLeila Ghaffari     // Set wall BCs
472fe7aee7SLeila Ghaffari     if (bc->num_wall > 0) {
482fe7aee7SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label,
492fe7aee7SLeila Ghaffari                            bc->num_wall, bc->walls, 0, bc->num_comps,
502fe7aee7SLeila Ghaffari                            bc->wall_comps, (void(*)(void))problem->bc,
51c95f9967SJed Brown                            NULL, problem->bc_ctx, NULL);  CHKERRQ(ierr);
522fe7aee7SLeila Ghaffari     }
532fe7aee7SLeila Ghaffari     // Set slip BCs in the x direction
542fe7aee7SLeila Ghaffari     if (bc->num_slip[0] > 0) {
552fe7aee7SLeila Ghaffari       PetscInt comps[1] = {1};
562fe7aee7SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label,
572fe7aee7SLeila Ghaffari                            bc->num_slip[0], bc->slips[0], 0, 1, comps,
58c95f9967SJed Brown                            (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr);
592fe7aee7SLeila Ghaffari     }
602fe7aee7SLeila Ghaffari     // Set slip BCs in the y direction
612fe7aee7SLeila Ghaffari     if (bc->num_slip[1] > 0) {
622fe7aee7SLeila Ghaffari       PetscInt comps[1] = {2};
632fe7aee7SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label,
642fe7aee7SLeila Ghaffari                            bc->num_slip[1], bc->slips[1], 0, 1, comps,
65c95f9967SJed Brown                            (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr);
662fe7aee7SLeila Ghaffari     }
672fe7aee7SLeila Ghaffari     // Set slip BCs in the z direction
682fe7aee7SLeila Ghaffari     if (bc->num_slip[2] > 0) {
692fe7aee7SLeila Ghaffari       PetscInt comps[1] = {3};
702fe7aee7SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label,
712fe7aee7SLeila Ghaffari                            bc->num_slip[2], bc->slips[2], 0, 1, comps,
72c95f9967SJed Brown                            (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr);
732fe7aee7SLeila Ghaffari     }
74*363b60e1SJames Wright     {
75*363b60e1SJames Wright       PetscBool use_strongstg = PETSC_FALSE;
76*363b60e1SJames Wright       ierr = PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL);
77*363b60e1SJames Wright       CHKERRQ(ierr);
78*363b60e1SJames Wright       STGShur14Context stg_ctx;
79*363b60e1SJames Wright 
80*363b60e1SJames Wright       if (use_strongstg) {
81*363b60e1SJames Wright         CeedQFunctionContextGetData(problem->apply_inflow.qfunction_context,
82*363b60e1SJames Wright                                     CEED_MEM_HOST, &stg_ctx);
83*363b60e1SJames Wright         ierr = SetupStrongSTG(dm, bc, problem, stg_ctx); CHKERRQ(ierr);
84*363b60e1SJames Wright         CeedQFunctionContextRestoreData(problem->apply_inflow.qfunction_context,
85*363b60e1SJames Wright                                         &stg_ctx);
86*363b60e1SJames Wright       }
87*363b60e1SJames Wright     }
88*363b60e1SJames Wright 
8977841947SLeila Ghaffari     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
9077841947SLeila Ghaffari     CHKERRQ(ierr);
9177841947SLeila Ghaffari     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
9277841947SLeila Ghaffari   }
9377841947SLeila Ghaffari   {
9477841947SLeila Ghaffari     // Empty name for conserved field (because there is only one field)
9577841947SLeila Ghaffari     PetscSection section;
9677841947SLeila Ghaffari     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
9777841947SLeila Ghaffari     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
9877841947SLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
9977841947SLeila Ghaffari     CHKERRQ(ierr);
10077841947SLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 1, "Momentum X");
10177841947SLeila Ghaffari     CHKERRQ(ierr);
10277841947SLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 2, "Momentum Y");
10377841947SLeila Ghaffari     CHKERRQ(ierr);
10477841947SLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 3, "Momentum Z");
10577841947SLeila Ghaffari     CHKERRQ(ierr);
10677841947SLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 4, "Energy Density");
10777841947SLeila Ghaffari     CHKERRQ(ierr);
10877841947SLeila Ghaffari   }
10977841947SLeila Ghaffari   PetscFunctionReturn(0);
11077841947SLeila Ghaffari }
11177841947SLeila Ghaffari 
11277841947SLeila Ghaffari // Refine DM for high-order viz
11377841947SLeila Ghaffari PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem,
114c95f9967SJed Brown                            SimpleBC bc, Physics phys) {
11577841947SLeila Ghaffari   PetscErrorCode ierr;
11677841947SLeila Ghaffari   DM             dm_hierarchy[user->app_ctx->viz_refine + 1];
11777841947SLeila Ghaffari   VecType        vec_type;
11877841947SLeila Ghaffari   PetscFunctionBeginUser;
11977841947SLeila Ghaffari 
12077841947SLeila Ghaffari   ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
12177841947SLeila Ghaffari 
12277841947SLeila Ghaffari   dm_hierarchy[0] = dm;
12377841947SLeila Ghaffari   for (PetscInt i = 0, d = user->app_ctx->degree;
12477841947SLeila Ghaffari        i < user->app_ctx->viz_refine; i++) {
12577841947SLeila Ghaffari     Mat interp_next;
12677841947SLeila Ghaffari     ierr = DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i+1]);
12777841947SLeila Ghaffari     CHKERRQ(ierr);
12877841947SLeila Ghaffari     ierr = DMClearDS(dm_hierarchy[i+1]); CHKERRQ(ierr);
12977841947SLeila Ghaffari     ierr = DMClearFields(dm_hierarchy[i+1]); CHKERRQ(ierr);
13077841947SLeila Ghaffari     ierr = DMSetCoarseDM(dm_hierarchy[i+1], dm_hierarchy[i]); CHKERRQ(ierr);
13177841947SLeila Ghaffari     d = (d + 1) / 2;
13277841947SLeila Ghaffari     if (i + 1 == user->app_ctx->viz_refine) d = 1;
13377841947SLeila Ghaffari     ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr);
13477841947SLeila Ghaffari     ierr = DMSetVecType(dm_hierarchy[i+1], vec_type); CHKERRQ(ierr);
135c95f9967SJed Brown     ierr = SetUpDM(dm_hierarchy[i+1], problem, d, bc, phys);
13677841947SLeila Ghaffari     CHKERRQ(ierr);
13777841947SLeila Ghaffari     ierr = DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i+1], &interp_next,
13877841947SLeila Ghaffari                                  NULL); CHKERRQ(ierr);
13977841947SLeila Ghaffari     if (!i) user->interp_viz = interp_next;
14077841947SLeila Ghaffari     else {
14177841947SLeila Ghaffari       Mat C;
14277841947SLeila Ghaffari       ierr = MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX,
14377841947SLeila Ghaffari                         PETSC_DECIDE, &C); CHKERRQ(ierr);
14477841947SLeila Ghaffari       ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
14577841947SLeila Ghaffari       ierr = MatDestroy(&user->interp_viz); CHKERRQ(ierr);
14677841947SLeila Ghaffari       user->interp_viz = C;
14777841947SLeila Ghaffari     }
14877841947SLeila Ghaffari   }
14977841947SLeila Ghaffari   for (PetscInt i=1; i<user->app_ctx->viz_refine; i++) {
15077841947SLeila Ghaffari     ierr = DMDestroy(&dm_hierarchy[i]); CHKERRQ(ierr);
15177841947SLeila Ghaffari   }
15277841947SLeila Ghaffari   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
15377841947SLeila Ghaffari 
15477841947SLeila Ghaffari   PetscFunctionReturn(0);
15577841947SLeila Ghaffari }
156