xref: /libCEED/examples/fluids/src/setupdm.c (revision 778419476db2fdf3952b1b9a5faed3f766f76223)
1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3 // reserved. See files LICENSE and NOTICE for details.
4 //
5 // This file is part of CEED, a collection of benchmarks, miniapps, software
6 // libraries and APIs for efficient high-order finite element and spectral
7 // element discretizations for exascale applications. For more information and
8 // source code availability see http://github.com/ceed.
9 //
10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11 // a collaborative effort of two U.S. Department of Energy organizations (Office
12 // of Science and the National Nuclear Security Administration) responsible for
13 // the planning and preparation of a capable exascale ecosystem, including
14 // software, applications, hardware, advanced system engineering and early
15 // testbed platforms, in support of the nation's exascale computing imperative.
16 
17 /// @file
18 /// Setup DM for Navier-Stokes example using PETSc
19 
20 #include "../navierstokes.h"
21 
22 // Read mesh and distribute DM in parallel
23 PetscErrorCode CreateDistributedDM(MPI_Comm comm, ProblemData *problem,
24                                    SetupContext setup_ctx, DM *dm) {
25   DM               dist_mesh = NULL;
26   PetscPartitioner part;
27   PetscInt         dim = problem->dim;
28   const PetscReal  scale[3] = {setup_ctx->lx, setup_ctx->ly, setup_ctx->lz};
29   PetscErrorCode   ierr;
30   PetscFunctionBeginUser;
31 
32   ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale,
33                              NULL, PETSC_TRUE, dm); CHKERRQ(ierr);
34 
35   // Distribute DM in parallel
36   ierr = DMPlexGetPartitioner(*dm, &part); CHKERRQ(ierr);
37   ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
38   ierr = DMPlexDistribute(*dm, 0, NULL, &dist_mesh); CHKERRQ(ierr);
39   if (dist_mesh) {
40     ierr = DMDestroy(dm); CHKERRQ(ierr);
41     *dm  = dist_mesh;
42   }
43   ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr);
44 
45   PetscFunctionReturn(0);
46 }
47 
48 // Setup DM
49 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree,
50                        SimpleBC bc, Physics phys, void *setup_ctx) {
51   PetscErrorCode ierr;
52   PetscFunctionBeginUser;
53   {
54     // Configure the finite element space and boundary conditions
55     PetscFE  fe;
56     PetscInt num_comp_q = 5;
57     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q,
58                                  PETSC_FALSE, degree, PETSC_DECIDE,
59                                  &fe); CHKERRQ(ierr);
60     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
61     ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr);
62     ierr = DMCreateDS(dm); CHKERRQ(ierr);
63     ierr = problem->bc_func(dm, bc, phys, setup_ctx);
64     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
65     CHKERRQ(ierr);
66     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
67   }
68   {
69     // Empty name for conserved field (because there is only one field)
70     PetscSection section;
71     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
72     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
73     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
74     CHKERRQ(ierr);
75     ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX");
76     CHKERRQ(ierr);
77     ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY");
78     CHKERRQ(ierr);
79     ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ");
80     CHKERRQ(ierr);
81     ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity");
82     CHKERRQ(ierr);
83   }
84   PetscFunctionReturn(0);
85 }
86 
87 // Refine DM for high-order viz
88 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem,
89                            SimpleBC bc, Physics phys, void *setup_ctx) {
90   PetscErrorCode ierr;
91   DM             dm_hierarchy[user->app_ctx->viz_refine + 1];
92   VecType        vec_type;
93   PetscFunctionBeginUser;
94 
95   ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
96 
97   dm_hierarchy[0] = dm;
98   for (PetscInt i = 0, d = user->app_ctx->degree;
99        i < user->app_ctx->viz_refine; i++) {
100     Mat interp_next;
101     ierr = DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i+1]);
102     CHKERRQ(ierr);
103     ierr = DMClearDS(dm_hierarchy[i+1]); CHKERRQ(ierr);
104     ierr = DMClearFields(dm_hierarchy[i+1]); CHKERRQ(ierr);
105     ierr = DMSetCoarseDM(dm_hierarchy[i+1], dm_hierarchy[i]); CHKERRQ(ierr);
106     d = (d + 1) / 2;
107     if (i + 1 == user->app_ctx->viz_refine) d = 1;
108     ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr);
109     ierr = DMSetVecType(dm_hierarchy[i+1], vec_type); CHKERRQ(ierr);
110     ierr = SetUpDM(dm_hierarchy[i+1], problem, d, bc, phys, setup_ctx);
111     CHKERRQ(ierr);
112     ierr = DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i+1], &interp_next,
113                                  NULL); CHKERRQ(ierr);
114     if (!i) user->interp_viz = interp_next;
115     else {
116       Mat C;
117       ierr = MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX,
118                         PETSC_DECIDE, &C); CHKERRQ(ierr);
119       ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
120       ierr = MatDestroy(&user->interp_viz); CHKERRQ(ierr);
121       user->interp_viz = C;
122     }
123   }
124   for (PetscInt i=1; i<user->app_ctx->viz_refine; i++) {
125     ierr = DMDestroy(&dm_hierarchy[i]); CHKERRQ(ierr);
126   }
127   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
128 
129   PetscFunctionReturn(0);
130 }
131