xref: /honee/src/setupdm.c (revision 05a512bda999ca8e4cad9319b5f88d1ff577e38b)
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 
22*05a512bdSLeila Ghaffari // Create mesh
23*05a512bdSLeila Ghaffari PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, DM *dm) {
24a515125bSLeila Ghaffari   PetscErrorCode   ierr;
25a515125bSLeila Ghaffari   PetscFunctionBeginUser;
26*05a512bdSLeila Ghaffari   // Create DMPLEX
27*05a512bdSLeila Ghaffari   ierr = DMCreate(comm, dm); CHKERRQ(ierr);
28*05a512bdSLeila Ghaffari   ierr = DMSetType(*dm, DMPLEX); CHKERRQ(ierr);
29*05a512bdSLeila Ghaffari   // Set Tensor elements
30*05a512bdSLeila Ghaffari   ierr = PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"); CHKERRQ(ierr);
31*05a512bdSLeila Ghaffari   // Set CL options
32*05a512bdSLeila 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;
46a515125bSLeila Ghaffari     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q,
47a515125bSLeila Ghaffari                                  PETSC_FALSE, degree, PETSC_DECIDE,
48a515125bSLeila Ghaffari                                  &fe); CHKERRQ(ierr);
49a515125bSLeila Ghaffari     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
50a515125bSLeila Ghaffari     ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr);
51a515125bSLeila Ghaffari     ierr = DMCreateDS(dm); CHKERRQ(ierr);
526b1ccf21SJeremy L Thompson     {
536b1ccf21SJeremy L Thompson       /* create FE field for coordinates */
546b1ccf21SJeremy L Thompson       PetscFE fe_coords;
556b1ccf21SJeremy L Thompson       PetscInt num_comp_coord;
566b1ccf21SJeremy L Thompson       ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr);
576b1ccf21SJeremy L Thompson       ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_coord,
586b1ccf21SJeremy L Thompson                                    PETSC_FALSE, 1, 1, &fe_coords); CHKERRQ(ierr);
596b1ccf21SJeremy L Thompson       ierr = DMProjectCoordinates(dm, fe_coords); CHKERRQ(ierr);
606b1ccf21SJeremy L Thompson       ierr = PetscFEDestroy(&fe_coords); CHKERRQ(ierr);
616b1ccf21SJeremy L Thompson     }
62a515125bSLeila Ghaffari     ierr = problem->bc_func(dm, bc, phys, setup_ctx);
63a515125bSLeila Ghaffari     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
64a515125bSLeila Ghaffari     CHKERRQ(ierr);
65a515125bSLeila Ghaffari     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
66a515125bSLeila Ghaffari   }
67a515125bSLeila Ghaffari   {
68a515125bSLeila Ghaffari     // Empty name for conserved field (because there is only one field)
69a515125bSLeila Ghaffari     PetscSection section;
70a515125bSLeila Ghaffari     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
71a515125bSLeila Ghaffari     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
72a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
73a515125bSLeila Ghaffari     CHKERRQ(ierr);
74a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 1, "Momentum X");
75a515125bSLeila Ghaffari     CHKERRQ(ierr);
76a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 2, "Momentum Y");
77a515125bSLeila Ghaffari     CHKERRQ(ierr);
78a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 3, "Momentum Z");
79a515125bSLeila Ghaffari     CHKERRQ(ierr);
80a515125bSLeila Ghaffari     ierr = PetscSectionSetComponentName(section, 0, 4, "Energy Density");
81a515125bSLeila Ghaffari     CHKERRQ(ierr);
82a515125bSLeila Ghaffari   }
83a515125bSLeila Ghaffari   PetscFunctionReturn(0);
84a515125bSLeila Ghaffari }
85a515125bSLeila Ghaffari 
86a515125bSLeila Ghaffari // Refine DM for high-order viz
87a515125bSLeila Ghaffari PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem,
88a515125bSLeila Ghaffari                            SimpleBC bc, Physics phys, void *setup_ctx) {
89a515125bSLeila Ghaffari   PetscErrorCode ierr;
90a515125bSLeila Ghaffari   DM             dm_hierarchy[user->app_ctx->viz_refine + 1];
91a515125bSLeila Ghaffari   VecType        vec_type;
92a515125bSLeila Ghaffari   PetscFunctionBeginUser;
93a515125bSLeila Ghaffari 
94a515125bSLeila Ghaffari   ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
95a515125bSLeila Ghaffari 
96a515125bSLeila Ghaffari   dm_hierarchy[0] = dm;
97a515125bSLeila Ghaffari   for (PetscInt i = 0, d = user->app_ctx->degree;
98a515125bSLeila Ghaffari        i < user->app_ctx->viz_refine; i++) {
99a515125bSLeila Ghaffari     Mat interp_next;
100a515125bSLeila Ghaffari     ierr = DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i+1]);
101a515125bSLeila Ghaffari     CHKERRQ(ierr);
102a515125bSLeila Ghaffari     ierr = DMClearDS(dm_hierarchy[i+1]); CHKERRQ(ierr);
103a515125bSLeila Ghaffari     ierr = DMClearFields(dm_hierarchy[i+1]); CHKERRQ(ierr);
104a515125bSLeila Ghaffari     ierr = DMSetCoarseDM(dm_hierarchy[i+1], dm_hierarchy[i]); CHKERRQ(ierr);
105a515125bSLeila Ghaffari     d = (d + 1) / 2;
106a515125bSLeila Ghaffari     if (i + 1 == user->app_ctx->viz_refine) d = 1;
107a515125bSLeila Ghaffari     ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr);
108a515125bSLeila Ghaffari     ierr = DMSetVecType(dm_hierarchy[i+1], vec_type); CHKERRQ(ierr);
109a515125bSLeila Ghaffari     ierr = SetUpDM(dm_hierarchy[i+1], problem, d, bc, phys, setup_ctx);
110a515125bSLeila Ghaffari     CHKERRQ(ierr);
111a515125bSLeila Ghaffari     ierr = DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i+1], &interp_next,
112a515125bSLeila Ghaffari                                  NULL); CHKERRQ(ierr);
113a515125bSLeila Ghaffari     if (!i) user->interp_viz = interp_next;
114a515125bSLeila Ghaffari     else {
115a515125bSLeila Ghaffari       Mat C;
116a515125bSLeila Ghaffari       ierr = MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX,
117a515125bSLeila Ghaffari                         PETSC_DECIDE, &C); CHKERRQ(ierr);
118a515125bSLeila Ghaffari       ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
119a515125bSLeila Ghaffari       ierr = MatDestroy(&user->interp_viz); CHKERRQ(ierr);
120a515125bSLeila Ghaffari       user->interp_viz = C;
121a515125bSLeila Ghaffari     }
122a515125bSLeila Ghaffari   }
123a515125bSLeila Ghaffari   for (PetscInt i=1; i<user->app_ctx->viz_refine; i++) {
124a515125bSLeila Ghaffari     ierr = DMDestroy(&dm_hierarchy[i]); CHKERRQ(ierr);
125a515125bSLeila Ghaffari   }
126a515125bSLeila Ghaffari   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
127a515125bSLeila Ghaffari 
128a515125bSLeila Ghaffari   PetscFunctionReturn(0);
129a515125bSLeila Ghaffari }
130