xref: /honee/src/setupdm.c (revision 00359db47665a79ecb0241f6ccbf886b649022df)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3a515125bSLeila Ghaffari 
4a515125bSLeila Ghaffari /// @file
5ea615d4cSJames Wright /// Setup DM for HONEE
6a515125bSLeila Ghaffari 
7e419654dSJeremy L Thompson #include <ceed.h>
8e419654dSJeremy L Thompson #include <petscdmplex.h>
967263decSKenneth E. Jansen #include <petscds.h>
10e419654dSJeremy L Thompson 
11149fb536SJames Wright #include <navierstokes.h>
12866f23abSJames Wright #include "../problems/stg_shur14.h"
13a515125bSLeila Ghaffari 
1405a512bdSLeila Ghaffari // Create mesh
CreateDM(Honee honee,ProblemData problem,MatType mat_type,VecType vec_type,DM * dm)15*4d9179f2SJames Wright PetscErrorCode CreateDM(Honee honee, ProblemData problem, MatType mat_type, VecType vec_type, DM *dm) {
16*4d9179f2SJames Wright   MPI_Comm  comm = honee->comm;
17eb9b4fe1SJames Wright   PetscBool isCGNS;
18eb9b4fe1SJames Wright   char      filename[PETSC_MAX_PATH_LEN] = "";
19eb9b4fe1SJames Wright 
20a515125bSLeila Ghaffari   PetscFunctionBeginUser;
21eb9b4fe1SJames Wright   PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_filename", filename, sizeof(filename), NULL));
22eb9b4fe1SJames Wright   PetscCall(HoneeCheckFilenameExtension(comm, filename, ".cgns", &isCGNS));
23eb9b4fe1SJames Wright 
24eb9b4fe1SJames Wright   if (isCGNS) {
25eb9b4fe1SJames Wright     // Must create from file to keep the sfNatural field in tact
26eb9b4fe1SJames Wright     PetscCall(DMPlexCreateFromFile(comm, filename, "HONEE", PETSC_TRUE, dm));
27eb9b4fe1SJames Wright     PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_filename", ""));
28eb9b4fe1SJames Wright   } else {
292b916ea7SJeremy L Thompson     PetscCall(DMCreate(comm, dm));
302b916ea7SJeremy L Thompson     PetscCall(DMSetType(*dm, DMPLEX));
31eb9b4fe1SJames Wright   }
32eb9b4fe1SJames Wright 
33b150a244SJed Brown   {
34b150a244SJed Brown     PetscBool skip = PETSC_TRUE;
352b916ea7SJeremy L Thompson     PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL));
36b150a244SJed Brown     PetscCall(DMSetMatrixPreallocateSkip(*dm, skip));
37b150a244SJed Brown   }
382b916ea7SJeremy L Thompson   PetscCall(DMSetMatType(*dm, mat_type));
392b916ea7SJeremy L Thompson   PetscCall(DMSetVecType(*dm, vec_type));
4059572072SJames Wright   PetscCall(HoneeOptionsSetValueDefault(NULL, "-dm_plex_simplex", "0"));
41e747eef9SJames Wright   PetscCall(HoneeOptionsSetValueDefault(NULL, "-dm_localize", "0"));  // Delay localization until DMSetupByOrderEnd_FEM
42517ae6cdSJames Wright   {  // Use Mat graph generation for partitioning. See https://gitlab.com/petsc/petsc/-/merge_requests/8336
43517ae6cdSJames Wright     PetscBool   distribute = PETSC_FALSE;
44517ae6cdSJames Wright     PetscMPIInt num_ranks;
45517ae6cdSJames Wright     PetscCall(DMPlexDistributeGetDefault(*dm, &distribute));
46517ae6cdSJames Wright     PetscCallMPI(MPI_Comm_size(comm, &num_ranks));
47517ae6cdSJames Wright     if (distribute && num_ranks > 1) PetscCall(HoneeOptionsSetValueDefault(NULL, "-dm_plex_csr_alg", "mat"));
48517ae6cdSJames Wright   }
49e747eef9SJames Wright   PetscCall(DMSetSparseLocalize(*dm, PETSC_FALSE));
50e747eef9SJames Wright   PetscCall(DMSetBlockingType(*dm, DM_BLOCKING_FIELD_NODE));
51e7934809SJames Wright   PetscCall(DMPlexSetPartitionBalance(*dm, PETSC_TRUE));
52ebfabadfSJames Wright 
5305a512bdSLeila Ghaffari   // Set CL options
542b916ea7SJeremy L Thompson   PetscCall(DMSetFromOptions(*dm));
552b916ea7SJeremy L Thompson   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
56eb9b4fe1SJames Wright   {  // Force isoperiodic point SF to be created to update sfNatural.
57eb9b4fe1SJames Wright     // Needs to be done before removing the field corresponding to sfNatural, but after IsoperiodicFaces have been set
58c4a0f6c7SJames Wright     PetscInt num_fields;
59c4a0f6c7SJames Wright     PetscCall(DMGetNumFields(*dm, &num_fields));
60c4a0f6c7SJames Wright     if (num_fields) {
61eb9b4fe1SJames Wright       PetscSection dummy;
62eb9b4fe1SJames Wright       PetscCall(DMGetGlobalSection(*dm, &dummy));
63eb9b4fe1SJames Wright     }
64c4a0f6c7SJames Wright   }
65*4d9179f2SJames Wright 
66*4d9179f2SJames Wright   PetscCall(DMPlexSetScale(*dm, PETSC_UNIT_LENGTH, honee->units->meter));
67*4d9179f2SJames Wright   PetscCall(DMPlexSetScale(*dm, PETSC_UNIT_MASS, honee->units->kilogram));
68*4d9179f2SJames Wright   PetscCall(DMPlexSetScale(*dm, PETSC_UNIT_TIME, honee->units->second));
69*4d9179f2SJames Wright   PetscCall(DMPlexSetScale(*dm, PETSC_UNIT_TEMPERATURE, honee->units->Kelvin));
70d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
71a515125bSLeila Ghaffari }
72a515125bSLeila Ghaffari 
73a515125bSLeila Ghaffari // Setup DM
SetUpDM(DM dm,ProblemData problem,PetscInt degree,PetscInt q_extra,Physics phys)74d3c60affSJames Wright PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, Physics phys) {
75a515125bSLeila Ghaffari   PetscInt num_comp_q = 5;
76da4ca0cfSJames Wright 
77eb9b4fe1SJames Wright   PetscFunctionBeginUser;
78eb9b4fe1SJames Wright   PetscCall(DMClearFields(dm));
79eb9b4fe1SJames Wright   PetscCall(DMSetLocalSection(dm, NULL));
80b2e5b5b3SJames Wright   PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, PETSC_DECIDE, q_extra, 1, &num_comp_q, dm));
81da4ca0cfSJames Wright 
82d6cac220SJames Wright   {
83d6cac220SJames Wright     PetscBool use_strongstg = PETSC_FALSE;
84d6cac220SJames Wright     PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
85d6cac220SJames Wright     if (use_strongstg) PetscCall(SetupStrongStg(dm, problem, phys));
86d6cac220SJames Wright   }
87d6cac220SJames Wright 
88da4ca0cfSJames Wright   {  // Add strong boundary conditions to DM
89047c2946SJed Brown     DMLabel label;
902b916ea7SJeremy L Thompson     PetscCall(DMGetLabel(dm, "Face Sets", &label));
915892ddc8SJames Wright     PetscCheck(label, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DM does not have 'Face Sets' label");
92c5e9980aSAdeleke O. Bankole     PetscCall(DMPlexLabelComplete(dm, label));
93487a3b6eSJames Wright 
94487a3b6eSJames Wright     for (PetscInt i = 0; i < problem->num_bc_defs; i++) {
95487a3b6eSJames Wright       BCDefinition    bc_def = problem->bc_defs[i];
96487a3b6eSJames Wright       PetscInt        num_essential_comps, num_label_values;
97487a3b6eSJames Wright       const PetscInt *essential_comps, *label_values;
98487a3b6eSJames Wright       const char     *name;
99487a3b6eSJames Wright 
10009240e0aSJames Wright       PetscCall(BCDefinitionSetDM(problem->bc_defs[i], dm));
101487a3b6eSJames Wright       PetscCall(BCDefinitionGetEssential(bc_def, &num_essential_comps, &essential_comps));
102487a3b6eSJames Wright       if (essential_comps > 0) {
103487a3b6eSJames Wright         PetscCall(BCDefinitionGetInfo(bc_def, &name, &num_label_values, &label_values));
104487a3b6eSJames Wright         PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, name, label, num_label_values, label_values, 0, num_essential_comps, essential_comps, NULL, NULL,
105487a3b6eSJames Wright                                 NULL, NULL));
106002797a3SLeila Ghaffari       }
107002797a3SLeila Ghaffari     }
10867263decSKenneth E. Jansen   }
10967263decSKenneth E. Jansen 
110da4ca0cfSJames Wright   PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm));
111cbe60e31SLeila Ghaffari 
112a515125bSLeila Ghaffari   // Empty name for conserved field (because there is only one field)
113a515125bSLeila Ghaffari   PetscSection section;
1142b916ea7SJeremy L Thompson   PetscCall(DMGetLocalSection(dm, &section));
1152b916ea7SJeremy L Thompson   PetscCall(PetscSectionSetFieldName(section, 0, ""));
116f27c2204SJames Wright   for (PetscInt i = 0; i < problem->num_components; i++) {
117f27c2204SJames Wright     PetscCall(PetscSectionSetComponentName(section, 0, i, problem->component_names[i]));
118a515125bSLeila Ghaffari   }
119d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
120a515125bSLeila Ghaffari }
121a515125bSLeila Ghaffari 
122a515125bSLeila Ghaffari // Refine DM for high-order viz
VizRefineDM(DM dm,Honee honee,ProblemData problem,Physics phys)123d3c60affSJames Wright PetscErrorCode VizRefineDM(DM dm, Honee honee, ProblemData problem, Physics phys) {
1240c373b74SJames Wright   DM      dm_hierarchy[honee->app_ctx->viz_refine + 1];
125a515125bSLeila Ghaffari   VecType vec_type;
126a515125bSLeila Ghaffari 
12706f41313SJames Wright   PetscFunctionBeginUser;
1282b916ea7SJeremy L Thompson   PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
129a515125bSLeila Ghaffari 
130a515125bSLeila Ghaffari   dm_hierarchy[0] = dm;
1310c373b74SJames Wright   for (PetscInt i = 0, d = honee->app_ctx->degree; i < honee->app_ctx->viz_refine; i++) {
132a515125bSLeila Ghaffari     Mat interp_next;
1332b916ea7SJeremy L Thompson     PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1]));
1342b916ea7SJeremy L Thompson     PetscCall(DMClearDS(dm_hierarchy[i + 1]));
1352b916ea7SJeremy L Thompson     PetscCall(DMClearFields(dm_hierarchy[i + 1]));
1362b916ea7SJeremy L Thompson     PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i]));
137a515125bSLeila Ghaffari     d                = (d + 1) / 2;
1380c373b74SJames Wright     PetscInt q_order = d + honee->app_ctx->q_extra;
1390c373b74SJames Wright     if (i + 1 == honee->app_ctx->viz_refine) d = 1;
1402b916ea7SJeremy L Thompson     PetscCall(DMGetVecType(dm, &vec_type));
1412b916ea7SJeremy L Thompson     PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type));
142d3c60affSJames Wright     PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, phys));
1432b916ea7SJeremy L Thompson     PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL));
1440c373b74SJames Wright     if (!i) honee->interp_viz = interp_next;
145a515125bSLeila Ghaffari     else {
146a515125bSLeila Ghaffari       Mat C;
1470c373b74SJames Wright       PetscCall(MatMatMult(interp_next, honee->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C));
1482b916ea7SJeremy L Thompson       PetscCall(MatDestroy(&interp_next));
1490c373b74SJames Wright       PetscCall(MatDestroy(&honee->interp_viz));
1500c373b74SJames Wright       honee->interp_viz = C;
151a515125bSLeila Ghaffari     }
152a515125bSLeila Ghaffari   }
1530c373b74SJames Wright   for (PetscInt i = 1; i < honee->app_ctx->viz_refine; i++) {
1542b916ea7SJeremy L Thompson     PetscCall(DMDestroy(&dm_hierarchy[i]));
155a515125bSLeila Ghaffari   }
1560c373b74SJames Wright   honee->dm_viz = dm_hierarchy[honee->app_ctx->viz_refine];
157d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
158a515125bSLeila Ghaffari }
159