xref: /honee/src/setupdm.c (revision 047c2946ef37c736435bf19d58b41537a935a06a)
1a515125bSLeila Ghaffari // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2a515125bSLeila Ghaffari // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3a515125bSLeila Ghaffari // reserved. See files LICENSE and NOTICE for details.
4a515125bSLeila Ghaffari //
5a515125bSLeila Ghaffari // This file is part of CEED, a collection of benchmarks, miniapps, software
6a515125bSLeila Ghaffari // libraries and APIs for efficient high-order finite element and spectral
7a515125bSLeila Ghaffari // element discretizations for exascale applications. For more information and
8a515125bSLeila Ghaffari // source code availability see http://github.com/ceed.
9a515125bSLeila Ghaffari //
10a515125bSLeila Ghaffari // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11a515125bSLeila Ghaffari // a collaborative effort of two U.S. Department of Energy organizations (Office
12a515125bSLeila Ghaffari // of Science and the National Nuclear Security Administration) responsible for
13a515125bSLeila Ghaffari // the planning and preparation of a capable exascale ecosystem, including
14a515125bSLeila Ghaffari // software, applications, hardware, advanced system engineering and early
15a515125bSLeila Ghaffari // testbed platforms, in support of the nation's exascale computing imperative.
16a515125bSLeila Ghaffari 
17a515125bSLeila Ghaffari /// @file
18a515125bSLeila Ghaffari /// Setup DM for Navier-Stokes example using PETSc
19a515125bSLeila Ghaffari 
20a515125bSLeila Ghaffari #include "../navierstokes.h"
21a515125bSLeila Ghaffari 
2205a512bdSLeila Ghaffari // Create mesh
2305a512bdSLeila Ghaffari PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, DM *dm) {
24a515125bSLeila Ghaffari   PetscErrorCode   ierr;
25a515125bSLeila Ghaffari   PetscFunctionBeginUser;
2605a512bdSLeila Ghaffari   // Create DMPLEX
2705a512bdSLeila Ghaffari   ierr = DMCreate(comm, dm); CHKERRQ(ierr);
2805a512bdSLeila Ghaffari   ierr = DMSetType(*dm, DMPLEX); CHKERRQ(ierr);
2905a512bdSLeila Ghaffari   // Set Tensor elements
3005a512bdSLeila Ghaffari   ierr = PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"); CHKERRQ(ierr);
3105a512bdSLeila Ghaffari   // Set CL options
3205a512bdSLeila Ghaffari   ierr = DMSetFromOptions(*dm); CHKERRQ(ierr);
33a515125bSLeila Ghaffari   ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr);
34a515125bSLeila Ghaffari   PetscFunctionReturn(0);
35a515125bSLeila Ghaffari }
36a515125bSLeila Ghaffari 
37a515125bSLeila Ghaffari // Setup DM
38a515125bSLeila Ghaffari PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree,
39a515125bSLeila Ghaffari                        SimpleBC bc, Physics phys, void *setup_ctx) {
40a515125bSLeila Ghaffari   PetscErrorCode ierr;
41a515125bSLeila Ghaffari   PetscFunctionBeginUser;
42a515125bSLeila Ghaffari   {
43a515125bSLeila Ghaffari     // Configure the finite element space and boundary conditions
44a515125bSLeila Ghaffari     PetscFE  fe;
45a515125bSLeila Ghaffari     PetscInt num_comp_q = 5;
46*047c2946SJed Brown     DMLabel label;
47a515125bSLeila Ghaffari     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q,
48a515125bSLeila Ghaffari                                  PETSC_FALSE, degree, PETSC_DECIDE,
49a515125bSLeila Ghaffari                                  &fe); CHKERRQ(ierr);
50a515125bSLeila Ghaffari     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
51a515125bSLeila Ghaffari     ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr);
52a515125bSLeila Ghaffari     ierr = DMCreateDS(dm); CHKERRQ(ierr);
536b1ccf21SJeremy L Thompson     {
546b1ccf21SJeremy L Thompson       /* create FE field for coordinates */
556b1ccf21SJeremy L Thompson       PetscFE fe_coords;
566b1ccf21SJeremy L Thompson       PetscInt num_comp_coord;
576b1ccf21SJeremy L Thompson       ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr);
586b1ccf21SJeremy L Thompson       ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_coord,
596b1ccf21SJeremy L Thompson                                    PETSC_FALSE, 1, 1, &fe_coords); CHKERRQ(ierr);
606b1ccf21SJeremy L Thompson       ierr = DMProjectCoordinates(dm, fe_coords); CHKERRQ(ierr);
616b1ccf21SJeremy L Thompson       ierr = PetscFEDestroy(&fe_coords); CHKERRQ(ierr);
626b1ccf21SJeremy L Thompson     }
63*047c2946SJed Brown     ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr);
64002797a3SLeila Ghaffari     // Set wall BCs
65002797a3SLeila Ghaffari     if (bc->num_wall > 0) {
66002797a3SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label,
67002797a3SLeila Ghaffari                            bc->num_wall, bc->walls, 0, bc->num_comps,
68002797a3SLeila Ghaffari                            bc->wall_comps, (void(*)(void))problem->bc,
69002797a3SLeila Ghaffari                            NULL, setup_ctx, NULL);  CHKERRQ(ierr);
70002797a3SLeila Ghaffari     }
71002797a3SLeila Ghaffari     // Set slip BCs in the x direction
72002797a3SLeila Ghaffari     if (bc->num_slip[0] > 0) {
73002797a3SLeila Ghaffari       PetscInt comps[1] = {1};
74002797a3SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label,
75002797a3SLeila Ghaffari                            bc->num_slip[0], bc->slips[0], 0, 1, comps,
76002797a3SLeila Ghaffari                            (void(*)(void))NULL, NULL, setup_ctx, NULL); CHKERRQ(ierr);
77002797a3SLeila Ghaffari     }
78002797a3SLeila Ghaffari     // Set slip BCs in the y direction
79002797a3SLeila Ghaffari     if (bc->num_slip[1] > 0) {
80002797a3SLeila Ghaffari       PetscInt comps[1] = {2};
81002797a3SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label,
82002797a3SLeila Ghaffari                            bc->num_slip[1], bc->slips[1], 0, 1, comps,
83002797a3SLeila Ghaffari                            (void(*)(void))NULL, NULL, setup_ctx, NULL); CHKERRQ(ierr);
84002797a3SLeila Ghaffari     }
85002797a3SLeila Ghaffari     // Set slip BCs in the z direction
86002797a3SLeila Ghaffari     if (bc->num_slip[2] > 0) {
87002797a3SLeila Ghaffari       PetscInt comps[1] = {3};
88002797a3SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label,
89002797a3SLeila Ghaffari                            bc->num_slip[2], bc->slips[2], 0, 1, comps,
90002797a3SLeila Ghaffari                            (void(*)(void))NULL, NULL, setup_ctx, NULL); CHKERRQ(ierr);
91002797a3SLeila Ghaffari     }
92a515125bSLeila Ghaffari     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
93a515125bSLeila Ghaffari     CHKERRQ(ierr);
94a515125bSLeila Ghaffari     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
95a515125bSLeila Ghaffari   }
96a515125bSLeila Ghaffari   {
97a515125bSLeila Ghaffari     // Empty name for conserved field (because there is only one field)
98a515125bSLeila Ghaffari     PetscSection section;
99a515125bSLeila Ghaffari     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
100a515125bSLeila Ghaffari     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
101a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
102a515125bSLeila Ghaffari     CHKERRQ(ierr);
103a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 1, "Momentum X");
104a515125bSLeila Ghaffari     CHKERRQ(ierr);
105a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 2, "Momentum Y");
106a515125bSLeila Ghaffari     CHKERRQ(ierr);
107a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 3, "Momentum Z");
108a515125bSLeila Ghaffari     CHKERRQ(ierr);
109a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 4, "Energy Density");
110a515125bSLeila Ghaffari     CHKERRQ(ierr);
111a515125bSLeila Ghaffari   }
112a515125bSLeila Ghaffari   PetscFunctionReturn(0);
113a515125bSLeila Ghaffari }
114a515125bSLeila Ghaffari 
115a515125bSLeila Ghaffari // Refine DM for high-order viz
116a515125bSLeila Ghaffari PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem,
117a515125bSLeila Ghaffari                            SimpleBC bc, Physics phys, void *setup_ctx) {
118a515125bSLeila Ghaffari   PetscErrorCode ierr;
119a515125bSLeila Ghaffari   DM             dm_hierarchy[user->app_ctx->viz_refine + 1];
120a515125bSLeila Ghaffari   VecType        vec_type;
121a515125bSLeila Ghaffari   PetscFunctionBeginUser;
122a515125bSLeila Ghaffari 
123a515125bSLeila Ghaffari   ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
124a515125bSLeila Ghaffari 
125a515125bSLeila Ghaffari   dm_hierarchy[0] = dm;
126a515125bSLeila Ghaffari   for (PetscInt i = 0, d = user->app_ctx->degree;
127a515125bSLeila Ghaffari        i < user->app_ctx->viz_refine; i++) {
128a515125bSLeila Ghaffari     Mat interp_next;
129a515125bSLeila Ghaffari     ierr = DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i+1]);
130a515125bSLeila Ghaffari     CHKERRQ(ierr);
131a515125bSLeila Ghaffari     ierr = DMClearDS(dm_hierarchy[i+1]); CHKERRQ(ierr);
132a515125bSLeila Ghaffari     ierr = DMClearFields(dm_hierarchy[i+1]); CHKERRQ(ierr);
133a515125bSLeila Ghaffari     ierr = DMSetCoarseDM(dm_hierarchy[i+1], dm_hierarchy[i]); CHKERRQ(ierr);
134a515125bSLeila Ghaffari     d = (d + 1) / 2;
135a515125bSLeila Ghaffari     if (i + 1 == user->app_ctx->viz_refine) d = 1;
136a515125bSLeila Ghaffari     ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr);
137a515125bSLeila Ghaffari     ierr = DMSetVecType(dm_hierarchy[i+1], vec_type); CHKERRQ(ierr);
138a515125bSLeila Ghaffari     ierr = SetUpDM(dm_hierarchy[i+1], problem, d, bc, phys, setup_ctx);
139a515125bSLeila Ghaffari     CHKERRQ(ierr);
140a515125bSLeila Ghaffari     ierr = DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i+1], &interp_next,
141a515125bSLeila Ghaffari                                  NULL); CHKERRQ(ierr);
142a515125bSLeila Ghaffari     if (!i) user->interp_viz = interp_next;
143a515125bSLeila Ghaffari     else {
144a515125bSLeila Ghaffari       Mat C;
145a515125bSLeila Ghaffari       ierr = MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX,
146a515125bSLeila Ghaffari                         PETSC_DECIDE, &C); CHKERRQ(ierr);
147a515125bSLeila Ghaffari       ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
148a515125bSLeila Ghaffari       ierr = MatDestroy(&user->interp_viz); CHKERRQ(ierr);
149a515125bSLeila Ghaffari       user->interp_viz = C;
150a515125bSLeila Ghaffari     }
151a515125bSLeila Ghaffari   }
152a515125bSLeila Ghaffari   for (PetscInt i=1; i<user->app_ctx->viz_refine; i++) {
153a515125bSLeila Ghaffari     ierr = DMDestroy(&dm_hierarchy[i]); CHKERRQ(ierr);
154a515125bSLeila Ghaffari   }
155a515125bSLeila Ghaffari   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
156a515125bSLeila Ghaffari 
157a515125bSLeila Ghaffari   PetscFunctionReturn(0);
158a515125bSLeila Ghaffari }
159