xref: /libCEED/examples/fluids/src/turb_spanstats.c (revision 19706a06bf84b1c610357895ba1ad37c35991712)
151ee423eSJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
251ee423eSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
351ee423eSJames Wright //
451ee423eSJames Wright // SPDX-License-Identifier: BSD-2-Clause
551ee423eSJames Wright //
651ee423eSJames Wright // This file is part of CEED:  http://github.com/ceed
751ee423eSJames Wright 
851ee423eSJames Wright /// @file
951ee423eSJames Wright /// Utility functions for setting up statistics collection
1051ee423eSJames Wright 
11*19706a06SJames Wright #include <petscsf.h>
12*19706a06SJames Wright 
1351ee423eSJames Wright #include "../navierstokes.h"
14*19706a06SJames Wright #include "petscsys.h"
1551ee423eSJames Wright 
1651ee423eSJames Wright PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree, SimpleBC bc) {
1751ee423eSJames Wright   user->spanstats.num_comp_stats = 1;
1851ee423eSJames Wright   PetscFunctionBeginUser;
1951ee423eSJames Wright 
2051ee423eSJames Wright   // Get DM from surface
2151ee423eSJames Wright   {
2251ee423eSJames Wright     DMLabel label;
2351ee423eSJames Wright     PetscCall(DMGetLabel(user->dm, "Face Sets", &label));
2451ee423eSJames Wright     PetscCall(DMPlexLabelComplete(user->dm, label));
2551ee423eSJames Wright     PetscCall(DMPlexFilter(user->dm, label, 1, &user->spanstats.dm));
2651ee423eSJames Wright     PetscCall(DMProjectCoordinates(user->spanstats.dm, NULL));  // Ensure that a coordinate FE exists
2751ee423eSJames Wright   }
2851ee423eSJames Wright 
2951ee423eSJames Wright   PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats"));
3051ee423eSJames Wright   PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "spanstats_"));
3151ee423eSJames Wright   PetscCall(PetscOptionsSetValue(NULL, "-spanstats_dm_sparse_localize", "0"));  // [Jed] Not relevant because not periodic in this direction
3251ee423eSJames Wright 
3351ee423eSJames Wright   PetscCall(DMSetFromOptions(user->spanstats.dm));
3451ee423eSJames Wright   PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view"));  // -spanstats_dm_view (option includes prefix)
3551ee423eSJames Wright   {
3651ee423eSJames Wright     PetscFE fe;
3751ee423eSJames Wright     DMLabel label;
3851ee423eSJames Wright 
3951ee423eSJames Wright     PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim - 1, user->spanstats.num_comp_stats, PETSC_FALSE, degree, PETSC_DECIDE, &fe));
4051ee423eSJames Wright     PetscCall(PetscObjectSetName((PetscObject)fe, "stats"));
4151ee423eSJames Wright     PetscCall(DMAddField(user->spanstats.dm, NULL, (PetscObject)fe));
4251ee423eSJames Wright     PetscCall(DMCreateDS(user->spanstats.dm));
4351ee423eSJames Wright     PetscCall(DMGetLabel(user->spanstats.dm, "Face Sets", &label));
4451ee423eSJames Wright 
4551ee423eSJames Wright     // // Set wall BCs
4651ee423eSJames Wright     // if (bc->num_wall > 0) {
4751ee423eSJames Wright     //   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, bc->num_wall, bc->walls, 0, bc->num_comps, bc->wall_comps,
4851ee423eSJames Wright     //                           (void (*)(void))problem->bc, NULL, problem->bc_ctx, NULL));
4951ee423eSJames Wright     // }
5051ee423eSJames Wright     // // Set slip BCs in the x direction
5151ee423eSJames Wright     // if (bc->num_slip[0] > 0) {
5251ee423eSJames Wright     //   PetscInt comps[1] = {1};
5351ee423eSJames Wright     //   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label, bc->num_slip[0], bc->slips[0], 0, 1, comps, (void (*)(void))NULL, NULL,
5451ee423eSJames Wright     //                           problem->bc_ctx, NULL));
5551ee423eSJames Wright     // }
5651ee423eSJames Wright     // // Set slip BCs in the y direction
5751ee423eSJames Wright     // if (bc->num_slip[1] > 0) {
5851ee423eSJames Wright     //   PetscInt comps[1] = {2};
5951ee423eSJames Wright     //   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label, bc->num_slip[1], bc->slips[1], 0, 1, comps, (void (*)(void))NULL, NULL,
6051ee423eSJames Wright     //                           problem->bc_ctx, NULL));
6151ee423eSJames Wright     // }
6251ee423eSJames Wright     // // Set slip BCs in the z direction
6351ee423eSJames Wright     // if (bc->num_slip[2] > 0) {
6451ee423eSJames Wright     //   PetscInt comps[1] = {3};
6551ee423eSJames Wright     //   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label, bc->num_slip[2], bc->slips[2], 0, 1, comps, (void (*)(void))NULL, NULL,
6651ee423eSJames Wright     //                           problem->bc_ctx, NULL));
6751ee423eSJames Wright     // }
6851ee423eSJames Wright 
6951ee423eSJames Wright     PetscCall(DMPlexSetClosurePermutationTensor(user->spanstats.dm, PETSC_DETERMINE, NULL));
7051ee423eSJames Wright     PetscCall(PetscFEDestroy(&fe));
7151ee423eSJames Wright   }
7251ee423eSJames Wright 
7351ee423eSJames Wright   PetscSection section;
7451ee423eSJames Wright   PetscCall(DMGetLocalSection(user->spanstats.dm, &section));
7551ee423eSJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
7651ee423eSJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "Test"));
7751ee423eSJames Wright   // PetscCall(PetscSectionSetComponentName(section, 0, 0, "Mean Velocity Products XX"));
7851ee423eSJames Wright   // PetscCall(PetscSectionSetComponentName(section, 0, 1, "Mean Velocity Products YY"));
7951ee423eSJames Wright   // PetscCall(PetscSectionSetComponentName(section, 0, 2, "Mean Velocity Products ZZ"));
8051ee423eSJames Wright   // PetscCall(PetscSectionSetComponentName(section, 0, 3, "Mean Velocity Products YZ"));
8151ee423eSJames Wright   // PetscCall(PetscSectionSetComponentName(section, 0, 4, "Mean Velocity Products XZ"));
8251ee423eSJames Wright   // PetscCall(PetscSectionSetComponentName(section, 0, 5, "Mean Velocity Products XY"));
8351ee423eSJames Wright 
84*19706a06SJames Wright   Vec test;
85*19706a06SJames Wright   PetscCall(DMCreateLocalVector(user->spanstats.dm, &test));
86*19706a06SJames Wright   PetscCall(VecZeroEntries(test));
87*19706a06SJames Wright   PetscCall(VecViewFromOptions(test, NULL, "-test_view"));
88*19706a06SJames Wright 
89*19706a06SJames Wright   PetscFunctionReturn(0);
90*19706a06SJames Wright }
91*19706a06SJames Wright 
92*19706a06SJames Wright PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, CeedVector *qx_coords,
93*19706a06SJames Wright                                    PetscInt *total_nqpnts) {
94*19706a06SJames Wright   CeedQFunction       qf_quad_coords;
95*19706a06SJames Wright   CeedOperator        op_quad_coords;
96*19706a06SJames Wright   PetscInt            num_comp_x, loc_num_elem, num_elem_qpts;
97*19706a06SJames Wright   CeedElemRestriction elem_restr_qx;
98*19706a06SJames Wright   PetscFunctionBeginUser;
99*19706a06SJames Wright 
100*19706a06SJames Wright   // Create Element Restriction and CeedVector for quadrature coordinates
101*19706a06SJames Wright   CeedBasisGetNumQuadraturePoints(basis_x, &num_elem_qpts);
102*19706a06SJames Wright   CeedElemRestrictionGetNumElements(elem_restr_x, &loc_num_elem);
103*19706a06SJames Wright   CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x);
104*19706a06SJames Wright   *total_nqpnts           = num_elem_qpts * loc_num_elem;
105*19706a06SJames Wright   const CeedInt strides[] = {num_comp_x, 1, num_elem_qpts * num_comp_x};
106*19706a06SJames Wright   CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp_x, num_comp_x * loc_num_elem * num_elem_qpts, strides, &elem_restr_qx);
107*19706a06SJames Wright   CeedElemRestrictionCreateVector(elem_restr_qx, qx_coords, NULL);
108*19706a06SJames Wright 
109*19706a06SJames Wright   // Create QFunction
110*19706a06SJames Wright   CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords);
111*19706a06SJames Wright 
112*19706a06SJames Wright   // Create Operator
113*19706a06SJames Wright   CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords);
114*19706a06SJames Wright   CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
115*19706a06SJames Wright   CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
116*19706a06SJames Wright 
117*19706a06SJames Wright   CeedOperatorApply(op_quad_coords, x_coords, *qx_coords, CEED_REQUEST_IMMEDIATE);
118*19706a06SJames Wright 
119*19706a06SJames Wright   CeedQFunctionDestroy(&qf_quad_coords);
120*19706a06SJames Wright   CeedOperatorDestroy(&op_quad_coords);
121*19706a06SJames Wright   PetscFunctionReturn(0);
122*19706a06SJames Wright }
123*19706a06SJames Wright 
124*19706a06SJames Wright PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, DM parentdm, DM childdm, PetscSF statssf) {
125*19706a06SJames Wright   PetscInt   child_num_qpnts, parent_num_qpnts, num_comp_x;
126*19706a06SJames Wright   CeedVector child_qx_coords, parent_qx_coords;
127*19706a06SJames Wright   PetscReal *child_coords, *parent_coords;
128*19706a06SJames Wright   PetscFunctionBeginUser;
129*19706a06SJames Wright 
130*19706a06SJames Wright   // Assume that child and parent have the same number of components
131*19706a06SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x);
132*19706a06SJames Wright   const PetscInt num_comp_sf = num_comp_x - 1;  // Number of coord components used in the creation of the SF
133*19706a06SJames Wright 
134*19706a06SJames Wright   // Get quad_coords for child DM
135*19706a06SJames Wright   PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_xc, ceed_data->x_coord, &child_qx_coords, &child_num_qpnts));
136*19706a06SJames Wright 
137*19706a06SJames Wright   // Get quad_coords for parent DM
138*19706a06SJames Wright   PetscCall(GetQuadratureCoords(ceed, parentdm, ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, ceed_data->spanstats.x_coord,
139*19706a06SJames Wright                                 &parent_qx_coords, &parent_num_qpnts));
140*19706a06SJames Wright 
141*19706a06SJames Wright   // Remove z component of coordinates for matching
142*19706a06SJames Wright   {
143*19706a06SJames Wright     const PetscReal *child_quad_coords, *parent_quad_coords;
144*19706a06SJames Wright 
145*19706a06SJames Wright     CeedVectorGetArrayRead(child_qx_coords, CEED_MEM_HOST, &child_quad_coords);
146*19706a06SJames Wright     CeedVectorGetArrayRead(parent_qx_coords, CEED_MEM_HOST, &parent_quad_coords);
147*19706a06SJames Wright 
148*19706a06SJames Wright     PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords));
149*19706a06SJames Wright     for (int i = 0; i < child_num_qpnts; i++) {
150*19706a06SJames Wright       child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x];
151*19706a06SJames Wright       child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x];
152*19706a06SJames Wright     }
153*19706a06SJames Wright     for (int i = 0; i < parent_num_qpnts; i++) {
154*19706a06SJames Wright       parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x];
155*19706a06SJames Wright       parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x];
156*19706a06SJames Wright     }
157*19706a06SJames Wright     CeedVectorRestoreArrayRead(child_qx_coords, &child_quad_coords);
158*19706a06SJames Wright     CeedVectorRestoreArrayRead(parent_qx_coords, &parent_quad_coords);
159*19706a06SJames Wright   }
160*19706a06SJames Wright 
161*19706a06SJames Wright   // Only check the first two components of the coordinates
162*19706a06SJames Wright   PetscCall(PetscSFSetGraphFromCoordinates(statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords));
163*19706a06SJames Wright 
164*19706a06SJames Wright   PetscCall(PetscFree2(child_coords, parent_coords));
165*19706a06SJames Wright   CeedVectorDestroy(&ceed_data->spanstats.x_coord);
166*19706a06SJames Wright   CeedVectorDestroy(&child_qx_coords);
167*19706a06SJames Wright   CeedVectorDestroy(&parent_qx_coords);
168*19706a06SJames Wright   PetscFunctionReturn(0);
169*19706a06SJames Wright }
170*19706a06SJames Wright 
171*19706a06SJames Wright PetscErrorCode SetupStatsCollection(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
172*19706a06SJames Wright   DM                 dm   = user->spanstats.dm;
173*19706a06SJames Wright   MPI_Comm           comm = PetscObjectComm((PetscObject)dm);
174*19706a06SJames Wright   CeedInt            dim, P, Q, num_comp_x;
175*19706a06SJames Wright   Vec                X_loc;
176*19706a06SJames Wright   PetscMemType       X_loc_memtype;
177*19706a06SJames Wright   const PetscScalar *X_loc_array;
178*19706a06SJames Wright   PetscFunctionBeginUser;
179*19706a06SJames Wright 
180*19706a06SJames Wright   PetscCall(DMGetDimension(dm, &dim));
181*19706a06SJames Wright   CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_x, &Q);
182*19706a06SJames Wright   CeedBasisGetNumNodes1D(ceed_data->basis_x, &P);
183*19706a06SJames Wright 
184*19706a06SJames Wright   // TODO: Possibly need to create a elem_restr_qcollocated for the global domain as well
185*19706a06SJames Wright   PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, user->spanstats.num_comp_stats, &ceed_data->spanstats.elem_restr_parent_stats,
186*19706a06SJames Wright                                     &ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.elem_restr_parent_qd));
187*19706a06SJames Wright   CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_x, &num_comp_x);
188*19706a06SJames Wright   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.x_coord, NULL);
189*19706a06SJames Wright 
190*19706a06SJames Wright   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS_LOBATTO, &ceed_data->spanstats.basis_x);
191*19706a06SJames Wright   CeedBasisCreateTensorH1Lagrange(ceed, dim, user->spanstats.num_comp_stats, P, Q, CEED_GAUSS, &ceed_data->spanstats.basis_stats);
192*19706a06SJames Wright 
193*19706a06SJames Wright   // -- Copy DM coordinates into CeedVector
194*19706a06SJames Wright   {
195*19706a06SJames Wright     DM cdm;
196*19706a06SJames Wright     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
197*19706a06SJames Wright     if (cdm) {
198*19706a06SJames Wright       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
199*19706a06SJames Wright     } else {
200*19706a06SJames Wright       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
201*19706a06SJames Wright     }
202*19706a06SJames Wright   }
203*19706a06SJames Wright   PetscCall(VecScale(X_loc, problem->dm_scale));
204*19706a06SJames Wright   PetscCall(VecGetArrayReadAndMemType(X_loc, &X_loc_array, &X_loc_memtype));
205*19706a06SJames Wright   CeedVectorSetArray(ceed_data->spanstats.x_coord, MemTypeP2C(X_loc_memtype), CEED_COPY_VALUES, (PetscScalar *)X_loc_array);
206*19706a06SJames Wright   PetscCall(VecRestoreArrayRead(X_loc, &X_loc_array));
207*19706a06SJames Wright 
208*19706a06SJames Wright   // Create SF for communicating child data back their respective parents
209*19706a06SJames Wright   PetscCall(PetscSFCreate(comm, &user->spanstats.sf));
210*19706a06SJames Wright   PetscCall(CreateStatsSF(ceed, ceed_data, user->dm, user->spanstats.dm, user->spanstats.sf));
21151ee423eSJames Wright 
21251ee423eSJames Wright   PetscFunctionReturn(0);
21351ee423eSJames Wright }
214