xref: /libCEED/examples/fluids/src/turb_spanstats.c (revision a175e4813653f41fdd11086020b9aeb03072d2d4)
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 /// @file
8*a175e481SJames Wright /// Functions for setting up and performing statistics collection
9*a175e481SJames Wright 
10*a175e481SJames Wright #include "../qfunctions/turb_spanstats.h"
1151ee423eSJames Wright 
1219706a06SJames Wright #include <petscsf.h>
1319706a06SJames Wright 
14*a175e481SJames Wright #include "../include/matops.h"
1551ee423eSJames Wright #include "../navierstokes.h"
16*a175e481SJames Wright #include "ceed/ceed.h"
17*a175e481SJames Wright #include "petscerror.h"
18*a175e481SJames Wright #include "petscmat.h"
1919706a06SJames Wright #include "petscsys.h"
20*a175e481SJames Wright #include "petscvec.h"
2151ee423eSJames Wright 
2251ee423eSJames Wright PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree, SimpleBC bc) {
23*a175e481SJames Wright   user->spanstats.num_comp_stats = 22;
24*a175e481SJames Wright   PetscReal domain_min[3], domain_max[3];
2551ee423eSJames Wright   PetscFunctionBeginUser;
2651ee423eSJames Wright 
27*a175e481SJames Wright   // Get spanwise length
28*a175e481SJames Wright   PetscCall(DMGetBoundingBox(user->dm, domain_min, domain_max));
29*a175e481SJames Wright   user->spanstats.span_width = domain_max[2] - domain_min[1];
30*a175e481SJames Wright 
3151ee423eSJames Wright   // Get DM from surface
3251ee423eSJames Wright   {
3351ee423eSJames Wright     DMLabel label;
3451ee423eSJames Wright     PetscCall(DMGetLabel(user->dm, "Face Sets", &label));
3551ee423eSJames Wright     PetscCall(DMPlexLabelComplete(user->dm, label));
3651ee423eSJames Wright     PetscCall(DMPlexFilter(user->dm, label, 1, &user->spanstats.dm));
3751ee423eSJames Wright     PetscCall(DMProjectCoordinates(user->spanstats.dm, NULL));  // Ensure that a coordinate FE exists
3851ee423eSJames Wright   }
3951ee423eSJames Wright 
4051ee423eSJames Wright   PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats"));
4151ee423eSJames Wright   PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "spanstats_"));
4251ee423eSJames Wright   PetscCall(DMSetFromOptions(user->spanstats.dm));
43*a175e481SJames Wright   PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view"));  // -spanstats_dm_view
4451ee423eSJames Wright   {
4551ee423eSJames Wright     PetscFE fe;
4651ee423eSJames Wright     DMLabel label;
4751ee423eSJames Wright 
4851ee423eSJames Wright     PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim - 1, user->spanstats.num_comp_stats, PETSC_FALSE, degree, PETSC_DECIDE, &fe));
4951ee423eSJames Wright     PetscCall(PetscObjectSetName((PetscObject)fe, "stats"));
5051ee423eSJames Wright     PetscCall(DMAddField(user->spanstats.dm, NULL, (PetscObject)fe));
5151ee423eSJames Wright     PetscCall(DMCreateDS(user->spanstats.dm));
5251ee423eSJames Wright     PetscCall(DMGetLabel(user->spanstats.dm, "Face Sets", &label));
5351ee423eSJames Wright 
5451ee423eSJames Wright     PetscCall(DMPlexSetClosurePermutationTensor(user->spanstats.dm, PETSC_DETERMINE, NULL));
5551ee423eSJames Wright     PetscCall(PetscFEDestroy(&fe));
5651ee423eSJames Wright   }
5751ee423eSJames Wright 
5851ee423eSJames Wright   PetscSection section;
5951ee423eSJames Wright   PetscCall(DMGetLocalSection(user->spanstats.dm, &section));
6051ee423eSJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
61*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "Mean Density"));
62*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "Mean Pressure"));
63*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "Mean Pressure Squared"));
64*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "Mean Pressure Velocity X"));
65*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "Mean Pressure Velocity Y"));
66*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "Mean Pressure Velocity Z"));
67*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 6, "Mean Density Temperature"));
68*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 7, "Mean Density Temperature Flux X"));
69*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 8, "Mean Density Temperature Flux Y"));
70*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 9, "Mean Density Temperature Flux Z"));
71*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 10, "Mean Momentum X"));
72*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 11, "Mean Momentum Y"));
73*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 12, "Mean Momentum Z"));
74*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 13, "Mean Momentum Flux XX"));
75*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 14, "Mean Momentum Flux YY"));
76*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 15, "Mean Momentum Flux ZZ"));
77*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 16, "Mean Momentum Flux YZ"));
78*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 17, "Mean Momentum Flux XZ"));
79*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 18, "Mean Momentum Flux XY"));
80*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 19, "Mean Velocity X"));
81*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 20, "Mean Velocity Y"));
82*a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 21, "Mean Velocity Z"));
8319706a06SJames Wright 
8419706a06SJames Wright   PetscFunctionReturn(0);
8519706a06SJames Wright }
8619706a06SJames Wright 
871737222fSJames Wright // Create CeedElemRestriction for collocated data based on associated CeedBasis and CeedElemRestriction
881737222fSJames Wright // Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction
891737222fSJames Wright PetscErrorCode CreateElemRestrColloc(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base,
901737222fSJames Wright                                      CeedElemRestriction *elem_restr_collocated, CeedVector *l_vec, CeedVector *e_vec) {
911737222fSJames Wright   CeedInt num_elem_qpts, loc_num_elem;
921737222fSJames Wright   PetscFunctionBeginUser;
931737222fSJames Wright 
941737222fSJames Wright   CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts);
951737222fSJames Wright   CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem);
961737222fSJames Wright 
971737222fSJames Wright   const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp};
981737222fSJames Wright   CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides,
991737222fSJames Wright                                    elem_restr_collocated);
1001737222fSJames Wright   CeedElemRestrictionCreateVector(*elem_restr_collocated, l_vec, e_vec);
1011737222fSJames Wright   PetscFunctionReturn(0);
1021737222fSJames Wright }
1031737222fSJames Wright 
104*a175e481SJames Wright // Get coordinates of quadrature points
10519706a06SJames Wright PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, CeedVector *qx_coords,
10619706a06SJames Wright                                    PetscInt *total_nqpnts) {
10719706a06SJames Wright   CeedQFunction       qf_quad_coords;
10819706a06SJames Wright   CeedOperator        op_quad_coords;
10919706a06SJames Wright   PetscInt            num_comp_x, loc_num_elem, num_elem_qpts;
11019706a06SJames Wright   CeedElemRestriction elem_restr_qx;
11119706a06SJames Wright   PetscFunctionBeginUser;
11219706a06SJames Wright 
11319706a06SJames Wright   // Create Element Restriction and CeedVector for quadrature coordinates
11419706a06SJames Wright   CeedBasisGetNumQuadraturePoints(basis_x, &num_elem_qpts);
11519706a06SJames Wright   CeedElemRestrictionGetNumElements(elem_restr_x, &loc_num_elem);
11619706a06SJames Wright   CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x);
11719706a06SJames Wright   *total_nqpnts = num_elem_qpts * loc_num_elem;
1181737222fSJames Wright   PetscCall(CreateElemRestrColloc(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx, qx_coords, NULL));
11919706a06SJames Wright 
12019706a06SJames Wright   // Create QFunction
12119706a06SJames Wright   CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords);
12219706a06SJames Wright 
12319706a06SJames Wright   // Create Operator
12419706a06SJames Wright   CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords);
12519706a06SJames Wright   CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
12619706a06SJames Wright   CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
12719706a06SJames Wright 
12819706a06SJames Wright   CeedOperatorApply(op_quad_coords, x_coords, *qx_coords, CEED_REQUEST_IMMEDIATE);
12919706a06SJames Wright 
13019706a06SJames Wright   CeedQFunctionDestroy(&qf_quad_coords);
13119706a06SJames Wright   CeedOperatorDestroy(&op_quad_coords);
13219706a06SJames Wright   PetscFunctionReturn(0);
13319706a06SJames Wright }
13419706a06SJames Wright 
135*a175e481SJames Wright // Create PetscSF for child-to-parent communication
13619706a06SJames Wright PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, DM parentdm, DM childdm, PetscSF statssf) {
13719706a06SJames Wright   PetscInt   child_num_qpnts, parent_num_qpnts, num_comp_x;
13819706a06SJames Wright   CeedVector child_qx_coords, parent_qx_coords;
13919706a06SJames Wright   PetscReal *child_coords, *parent_coords;
14019706a06SJames Wright   PetscFunctionBeginUser;
14119706a06SJames Wright 
14219706a06SJames Wright   // Assume that child and parent have the same number of components
14319706a06SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x);
14419706a06SJames Wright   const PetscInt num_comp_sf = num_comp_x - 1;  // Number of coord components used in the creation of the SF
14519706a06SJames Wright 
14619706a06SJames Wright   // Get quad_coords for child DM
147*a175e481SJames Wright   PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &child_qx_coords, &child_num_qpnts));
14819706a06SJames Wright 
14919706a06SJames Wright   // Get quad_coords for parent DM
15019706a06SJames Wright   PetscCall(GetQuadratureCoords(ceed, parentdm, ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, ceed_data->spanstats.x_coord,
15119706a06SJames Wright                                 &parent_qx_coords, &parent_num_qpnts));
15219706a06SJames Wright 
15319706a06SJames Wright   // Remove z component of coordinates for matching
15419706a06SJames Wright   {
15519706a06SJames Wright     const PetscReal *child_quad_coords, *parent_quad_coords;
15619706a06SJames Wright 
15719706a06SJames Wright     CeedVectorGetArrayRead(child_qx_coords, CEED_MEM_HOST, &child_quad_coords);
15819706a06SJames Wright     CeedVectorGetArrayRead(parent_qx_coords, CEED_MEM_HOST, &parent_quad_coords);
15919706a06SJames Wright 
16019706a06SJames Wright     PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords));
16119706a06SJames Wright     for (int i = 0; i < child_num_qpnts; i++) {
16219706a06SJames Wright       child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x];
16319706a06SJames Wright       child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x];
16419706a06SJames Wright     }
16519706a06SJames Wright     for (int i = 0; i < parent_num_qpnts; i++) {
16619706a06SJames Wright       parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x];
16719706a06SJames Wright       parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x];
16819706a06SJames Wright     }
16919706a06SJames Wright     CeedVectorRestoreArrayRead(child_qx_coords, &child_quad_coords);
17019706a06SJames Wright     CeedVectorRestoreArrayRead(parent_qx_coords, &parent_quad_coords);
17119706a06SJames Wright   }
17219706a06SJames Wright 
17319706a06SJames Wright   PetscCall(PetscSFSetGraphFromCoordinates(statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords));
17419706a06SJames Wright 
175*a175e481SJames Wright   PetscCall(PetscSFViewFromOptions(statssf, NULL, "-spanstats_sf_view"));
176*a175e481SJames Wright 
17719706a06SJames Wright   PetscCall(PetscFree2(child_coords, parent_coords));
17819706a06SJames Wright   CeedVectorDestroy(&child_qx_coords);
17919706a06SJames Wright   CeedVectorDestroy(&parent_qx_coords);
18019706a06SJames Wright   PetscFunctionReturn(0);
18119706a06SJames Wright }
18219706a06SJames Wright 
183*a175e481SJames Wright // Compute mass matrix for statistics projection
184*a175e481SJames Wright PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data) {
185*a175e481SJames Wright   CeedQFunction qf_mass;
186*a175e481SJames Wright   CeedOperator  op_mass;
187*a175e481SJames Wright   CeedInt       num_comp_q, q_data_size;
188*a175e481SJames Wright   PetscFunctionBeginUser;
189*a175e481SJames Wright 
190*a175e481SJames Wright   // CEED Restriction
191*a175e481SJames Wright   CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_stats, &num_comp_q);
192*a175e481SJames Wright   CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_qd, &q_data_size);
193*a175e481SJames Wright 
194*a175e481SJames Wright   // Create Mass CeedOperator
195*a175e481SJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
196*a175e481SJames Wright   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
197*a175e481SJames Wright   CeedOperatorSetField(op_mass, "q", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE);
198*a175e481SJames Wright   CeedOperatorSetField(op_mass, "qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, ceed_data->spanstats.q_data);
199*a175e481SJames Wright   CeedOperatorSetField(op_mass, "v", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE);
200*a175e481SJames Wright 
201*a175e481SJames Wright   // Setup KSP for L^2 projection
202*a175e481SJames Wright   {
203*a175e481SJames Wright     MatopApplyContext M_ctx;
204*a175e481SJames Wright     PetscInt          l_size, g_size;
205*a175e481SJames Wright     Mat               mat_mass;
206*a175e481SJames Wright     VecType           vec_type;
207*a175e481SJames Wright     KSP               ksp;
208*a175e481SJames Wright     Vec               ones, M_inv;
209*a175e481SJames Wright     CeedVector        x_ceed, y_ceed;
210*a175e481SJames Wright 
211*a175e481SJames Wright     PetscCall(DMCreateGlobalVector(user->spanstats.dm, &M_inv));
212*a175e481SJames Wright     PetscCall(VecGetLocalSize(M_inv, &l_size));
213*a175e481SJames Wright     PetscCall(VecGetSize(M_inv, &g_size));
214*a175e481SJames Wright     PetscCall(VecGetType(M_inv, &vec_type));
215*a175e481SJames Wright 
216*a175e481SJames Wright     PetscCall(PetscMalloc1(1, &M_ctx));
217*a175e481SJames Wright     PetscCall(MatCreateShell(PETSC_COMM_WORLD, l_size, l_size, g_size, g_size, M_ctx, &mat_mass));
218*a175e481SJames Wright     PetscCall(MatShellSetOperation(mat_mass, MATOP_MULT, (void (*)(void))MatMult_Ceed));
219*a175e481SJames Wright     PetscCall(MatShellSetOperation(mat_mass, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed));
220*a175e481SJames Wright     PetscCall(MatShellSetVecType(mat_mass, vec_type));
221*a175e481SJames Wright 
222*a175e481SJames Wright     CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &x_ceed, NULL);
223*a175e481SJames Wright     CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &y_ceed, NULL);
224*a175e481SJames Wright 
225*a175e481SJames Wright     PetscCall(SetupMatopApplyCtx(PETSC_COMM_WORLD, user->spanstats.dm, user->ceed, op_mass, x_ceed, y_ceed, NULL, M_ctx));
226*a175e481SJames Wright     user->spanstats.M_ctx = M_ctx;
227*a175e481SJames Wright 
228*a175e481SJames Wright     // Create lumped mass matrix inverse
229*a175e481SJames Wright     PetscCall(DMGetGlobalVector(user->spanstats.dm, &ones));
230*a175e481SJames Wright     PetscCall(VecZeroEntries(M_inv));
231*a175e481SJames Wright     PetscCall(VecSet(ones, 1));
232*a175e481SJames Wright     PetscCall(MatMult(mat_mass, ones, M_inv));
233*a175e481SJames Wright     PetscCall(VecReciprocal(M_inv));
234*a175e481SJames Wright     user->spanstats.M_inv = M_inv;
235*a175e481SJames Wright     PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &ones));
236*a175e481SJames Wright 
237*a175e481SJames Wright     PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
238*a175e481SJames Wright     PetscCall(KSPSetOptionsPrefix(ksp, "spanstats_"));
239*a175e481SJames Wright     {
240*a175e481SJames Wright       PC pc;
241*a175e481SJames Wright       PetscCall(KSPGetPC(ksp, &pc));
242*a175e481SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
243*a175e481SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
244*a175e481SJames Wright       PetscCall(KSPSetType(ksp, KSPCG));
245*a175e481SJames Wright       PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
246*a175e481SJames Wright       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
247*a175e481SJames Wright     }
248*a175e481SJames Wright     PetscCall(KSPSetOperators(ksp, mat_mass, mat_mass));
249*a175e481SJames Wright     PetscCall(KSPSetFromOptions(ksp));
250*a175e481SJames Wright     user->spanstats.ksp = ksp;
251*a175e481SJames Wright   }
252*a175e481SJames Wright 
253*a175e481SJames Wright   // Cleanup
254*a175e481SJames Wright   CeedQFunctionDestroy(&qf_mass);
255*a175e481SJames Wright   PetscFunctionReturn(0);
256*a175e481SJames Wright }
257*a175e481SJames Wright 
258*a175e481SJames Wright // Create CeedOperators and KSP for the statistics collection and processing
259*a175e481SJames Wright PetscErrorCode CreateStatisticsOperators(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
260*a175e481SJames Wright   CeedInt      num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q;
261*a175e481SJames Wright   CeedOperator op_setup_sur;
262*a175e481SJames Wright   PetscFunctionBeginUser;
263*a175e481SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q);
264*a175e481SJames Wright 
265*a175e481SJames Wright   // Create Operator for statistics collection
266*a175e481SJames Wright   switch (user->phys->state_var) {
267*a175e481SJames Wright     case STATEVAR_PRIMITIVE:
268*a175e481SJames Wright       CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &ceed_data->spanstats.qf_stats_collect);
269*a175e481SJames Wright       break;
270*a175e481SJames Wright     case STATEVAR_CONSERVATIVE:
271*a175e481SJames Wright       CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &ceed_data->spanstats.qf_stats_collect);
272*a175e481SJames Wright       break;
273*a175e481SJames Wright   }
274*a175e481SJames Wright 
275*a175e481SJames Wright   if (user->app_ctx->stats_test) {
276*a175e481SJames Wright     CeedQFunctionDestroy(&ceed_data->spanstats.qf_stats_collect);
277*a175e481SJames Wright     CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionTest, ChildStatsCollectionTest_loc, &ceed_data->spanstats.qf_stats_collect);
278*a175e481SJames Wright   }
279*a175e481SJames Wright 
280*a175e481SJames Wright   CeedQFunctionSetContext(ceed_data->spanstats.qf_stats_collect, problem->apply_vol_ifunction.qfunction_context);
281*a175e481SJames Wright   CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP);
282*a175e481SJames Wright   CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE);
283*a175e481SJames Wright   CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP);
284*a175e481SJames Wright   CeedQFunctionAddOutput(ceed_data->spanstats.qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE);
285*a175e481SJames Wright 
286*a175e481SJames Wright   CeedOperatorCreate(ceed, ceed_data->spanstats.qf_stats_collect, NULL, NULL, &user->spanstats.op_stats_collect);
287*a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
288*a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
289*a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
290*a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_collect, "v", ceed_data->spanstats.elem_restr_child_colloc, CEED_BASIS_COLLOCATED,
291*a175e481SJames Wright                        CEED_VECTOR_ACTIVE);
292*a175e481SJames Wright 
293*a175e481SJames Wright   // Create Operator for statistics projection / processing
294*a175e481SJames Wright   // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function.
295*a175e481SJames Wright   // Therefore, an Identity QF is sufficient
296*a175e481SJames Wright   CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &ceed_data->spanstats.qf_stats_proj);
297*a175e481SJames Wright 
298*a175e481SJames Wright   CeedOperatorCreate(ceed, ceed_data->spanstats.qf_stats_proj, NULL, NULL, &user->spanstats.op_stats_proj);
299*a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_proj, "input", ceed_data->spanstats.elem_restr_parent_colloc, CEED_BASIS_COLLOCATED,
300*a175e481SJames Wright                        CEED_VECTOR_ACTIVE);
301*a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_proj, "output", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats,
302*a175e481SJames Wright                        CEED_VECTOR_ACTIVE);
303*a175e481SJames Wright 
304*a175e481SJames Wright   // Get q_data for lumped mass matrix formation
305*a175e481SJames Wright   CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur);
306*a175e481SJames Wright   CeedOperatorSetField(op_setup_sur, "dx", ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, CEED_VECTOR_ACTIVE);
307*a175e481SJames Wright   CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->spanstats.basis_x, CEED_VECTOR_NONE);
308*a175e481SJames Wright   CeedOperatorSetField(op_setup_sur, "surface qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
309*a175e481SJames Wright   CeedOperatorApply(op_setup_sur, ceed_data->spanstats.x_coord, ceed_data->spanstats.q_data, CEED_REQUEST_IMMEDIATE);
310*a175e481SJames Wright 
311*a175e481SJames Wright   CeedOperatorDestroy(&op_setup_sur);
312*a175e481SJames Wright   PetscFunctionReturn(0);
313*a175e481SJames Wright }
314*a175e481SJames Wright 
315*a175e481SJames Wright // Setup for statistics collection
31619706a06SJames Wright PetscErrorCode SetupStatsCollection(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
31719706a06SJames Wright   DM                 dm   = user->spanstats.dm;
31819706a06SJames Wright   MPI_Comm           comm = PetscObjectComm((PetscObject)dm);
31919706a06SJames Wright   CeedInt            dim, P, Q, num_comp_x;
32019706a06SJames Wright   Vec                X_loc;
32119706a06SJames Wright   PetscMemType       X_loc_memtype;
32219706a06SJames Wright   const PetscScalar *X_loc_array;
32319706a06SJames Wright   PetscFunctionBeginUser;
32419706a06SJames Wright 
32519706a06SJames Wright   PetscCall(DMGetDimension(dm, &dim));
326*a175e481SJames Wright   CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &Q);
327*a175e481SJames Wright   CeedBasisGetNumNodes1D(ceed_data->basis_q, &P);
32819706a06SJames Wright 
3291737222fSJames Wright   PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, problem->q_data_size_sur, &ceed_data->spanstats.elem_restr_parent_stats,
33019706a06SJames Wright                                     &ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.elem_restr_parent_qd));
33119706a06SJames Wright   CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_x, &num_comp_x);
33219706a06SJames Wright   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.x_coord, NULL);
333*a175e481SJames Wright   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &user->spanstats.rhs_ceed, NULL);
3341737222fSJames Wright   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_qd, &ceed_data->spanstats.q_data, NULL);
33519706a06SJames Wright 
336*a175e481SJames Wright   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS, &ceed_data->spanstats.basis_x);
33719706a06SJames Wright   CeedBasisCreateTensorH1Lagrange(ceed, dim, user->spanstats.num_comp_stats, P, Q, CEED_GAUSS, &ceed_data->spanstats.basis_stats);
33819706a06SJames Wright 
3391737222fSJames Wright   PetscCall(CreateElemRestrColloc(ceed, user->spanstats.num_comp_stats, ceed_data->spanstats.basis_stats,
3401737222fSJames Wright                                   ceed_data->spanstats.elem_restr_parent_stats, &ceed_data->spanstats.elem_restr_parent_colloc,
3411737222fSJames Wright                                   &user->spanstats.parent_stats, NULL));
3421737222fSJames Wright   PetscCall(CreateElemRestrColloc(ceed, user->spanstats.num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q,
3431737222fSJames Wright                                   &ceed_data->spanstats.elem_restr_child_colloc, &user->spanstats.child_stats, NULL));
344*a175e481SJames Wright   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_child_colloc, &user->spanstats.child_inst_stats, NULL);
345*a175e481SJames Wright   CeedVectorSetValue(user->spanstats.child_stats, 0);
3461737222fSJames Wright 
34719706a06SJames Wright   // -- Copy DM coordinates into CeedVector
34819706a06SJames Wright   {
34919706a06SJames Wright     DM cdm;
35019706a06SJames Wright     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
35119706a06SJames Wright     if (cdm) {
35219706a06SJames Wright       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
35319706a06SJames Wright     } else {
35419706a06SJames Wright       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
35519706a06SJames Wright     }
35619706a06SJames Wright   }
35719706a06SJames Wright   PetscCall(VecScale(X_loc, problem->dm_scale));
35819706a06SJames Wright   PetscCall(VecGetArrayReadAndMemType(X_loc, &X_loc_array, &X_loc_memtype));
35919706a06SJames Wright   CeedVectorSetArray(ceed_data->spanstats.x_coord, MemTypeP2C(X_loc_memtype), CEED_COPY_VALUES, (PetscScalar *)X_loc_array);
36019706a06SJames Wright   PetscCall(VecRestoreArrayRead(X_loc, &X_loc_array));
36119706a06SJames Wright 
36219706a06SJames Wright   // Create SF for communicating child data back their respective parents
36319706a06SJames Wright   PetscCall(PetscSFCreate(comm, &user->spanstats.sf));
36419706a06SJames Wright   PetscCall(CreateStatsSF(ceed, ceed_data, user->dm, user->spanstats.dm, user->spanstats.sf));
36551ee423eSJames Wright 
366*a175e481SJames Wright   // Create CeedOperators for statistics collection
367*a175e481SJames Wright   PetscCall(CreateStatisticsOperators(ceed, user, ceed_data, problem));
368*a175e481SJames Wright 
369*a175e481SJames Wright   // Setup KSP and Mat for L^2 projection of statistics
370*a175e481SJames Wright   PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data));
371*a175e481SJames Wright 
372*a175e481SJames Wright   PetscFunctionReturn(0);
373*a175e481SJames Wright }
374*a175e481SJames Wright 
375*a175e481SJames Wright // Collect statistics based on the solution Q
376*a175e481SJames Wright PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) {
377*a175e481SJames Wright   PetscMemType       q_mem_type;
378*a175e481SJames Wright   const PetscScalar *q_arr;
379*a175e481SJames Wright   Vec                Q_loc;
380*a175e481SJames Wright   PetscFunctionBeginUser;
381*a175e481SJames Wright 
382*a175e481SJames Wright   PetscCall(DMGetLocalVector(user->dm, &Q_loc));
383*a175e481SJames Wright   PetscCall(VecZeroEntries(Q_loc));
384*a175e481SJames Wright   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
385*a175e481SJames Wright 
386*a175e481SJames Wright   PetscCall(VecGetArrayReadAndMemType(Q_loc, &q_arr, &q_mem_type));
387*a175e481SJames Wright   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q_arr);
388*a175e481SJames Wright 
389*a175e481SJames Wright   CeedOperatorApply(user->spanstats.op_stats_collect, user->q_ceed, user->spanstats.child_inst_stats, CEED_REQUEST_IMMEDIATE);
390*a175e481SJames Wright 
391*a175e481SJames Wright   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
392*a175e481SJames Wright   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q_arr));
393*a175e481SJames Wright   PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
394*a175e481SJames Wright 
395*a175e481SJames Wright   // Record averaging using left rectangle rule
396*a175e481SJames Wright   PetscScalar delta_t           = solution_time - user->spanstats.prev_time;
397*a175e481SJames Wright   PetscScalar prev_timeinterval = user->spanstats.prev_time - user->app_ctx->cont_time;
398*a175e481SJames Wright   CeedVectorScale(user->spanstats.child_stats, prev_timeinterval / (prev_timeinterval + delta_t));
399*a175e481SJames Wright   CeedVectorAXPY(user->spanstats.child_stats, delta_t / (prev_timeinterval + delta_t), user->spanstats.child_inst_stats);
400*a175e481SJames Wright   user->spanstats.prev_time = solution_time;
401*a175e481SJames Wright 
402*a175e481SJames Wright   PetscFunctionReturn(0);
403*a175e481SJames Wright }
404*a175e481SJames Wright 
405*a175e481SJames Wright // Process the child statistics into parent statistics and project them onto stats
406*a175e481SJames Wright PetscErrorCode ProcessStatistics(User user, Vec *stats) {
407*a175e481SJames Wright   Span_Stats         user_stats = user->spanstats;
408*a175e481SJames Wright   const PetscScalar *child_stats;
409*a175e481SJames Wright   PetscScalar       *parent_stats;
410*a175e481SJames Wright   MPI_Datatype       unit;
411*a175e481SJames Wright   Vec                rhs_loc, rhs;
412*a175e481SJames Wright   PetscMemType       rhs_mem_type;
413*a175e481SJames Wright   CeedScalar        *rhs_arr;
414*a175e481SJames Wright   CeedMemType        ceed_mem_type;
415*a175e481SJames Wright   PetscFunctionBeginUser;
416*a175e481SJames Wright 
417*a175e481SJames Wright   CeedGetPreferredMemType(user->ceed, &ceed_mem_type);
418*a175e481SJames Wright   CeedVectorSetValue(user_stats.parent_stats, 0);
419*a175e481SJames Wright 
420*a175e481SJames Wright   CeedVectorGetArrayRead(user_stats.child_stats, ceed_mem_type, &child_stats);
421*a175e481SJames Wright   CeedVectorGetArray(user_stats.parent_stats, ceed_mem_type, &parent_stats);
422*a175e481SJames Wright 
423*a175e481SJames Wright   if (user_stats.num_comp_stats == 1) unit = MPIU_REAL;
424*a175e481SJames Wright   else {
425*a175e481SJames Wright     PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit));
426*a175e481SJames Wright     PetscCallMPI(MPI_Type_commit(&unit));
427*a175e481SJames Wright   }
428*a175e481SJames Wright 
429*a175e481SJames Wright   PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
430*a175e481SJames Wright   PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
431*a175e481SJames Wright 
432*a175e481SJames Wright   CeedVectorRestoreArrayRead(user_stats.child_stats, &child_stats);
433*a175e481SJames Wright   CeedVectorRestoreArray(user_stats.parent_stats, &parent_stats);
434*a175e481SJames Wright   PetscCallMPI(MPI_Type_free(&unit));
435*a175e481SJames Wright 
436*a175e481SJames Wright   CeedVectorScale(user_stats.parent_stats, 1 / user_stats.span_width);
437*a175e481SJames Wright 
438*a175e481SJames Wright   // L^2 projection with the parent_data
439*a175e481SJames Wright   PetscCall(DMGetGlobalVector(user_stats.dm, &rhs));
440*a175e481SJames Wright   PetscCall(DMGetLocalVector(user_stats.dm, &rhs_loc));
441*a175e481SJames Wright   PetscCall(VecZeroEntries(rhs));
442*a175e481SJames Wright   PetscCall(VecZeroEntries(rhs_loc));
443*a175e481SJames Wright   PetscCall(VecGetArrayWriteAndMemType(rhs_loc, &rhs_arr, &rhs_mem_type));
444*a175e481SJames Wright   CeedVectorSetArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), CEED_USE_POINTER, (PetscScalar *)rhs_arr);
445*a175e481SJames Wright 
446*a175e481SJames Wright   CeedOperatorApply(user_stats.op_stats_proj, user_stats.parent_stats, user_stats.rhs_ceed, CEED_REQUEST_IMMEDIATE);
447*a175e481SJames Wright 
448*a175e481SJames Wright   CeedVectorTakeArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), &rhs_arr);
449*a175e481SJames Wright   PetscCall(VecRestoreArrayAndMemType(rhs_loc, &rhs_arr));
450*a175e481SJames Wright   PetscCall(DMLocalToGlobal(user_stats.dm, rhs_loc, ADD_VALUES, rhs));
451*a175e481SJames Wright 
452*a175e481SJames Wright   PetscCall(VecDuplicate(rhs, stats));
453*a175e481SJames Wright   PetscCall(VecPointwiseMult(*stats, rhs, user_stats.M_inv));
454*a175e481SJames Wright 
455*a175e481SJames Wright   PetscCall(KSPSolve(user_stats.ksp, rhs, *stats));
456*a175e481SJames Wright 
457*a175e481SJames Wright   PetscFunctionReturn(0);
458*a175e481SJames Wright }
459*a175e481SJames Wright 
460*a175e481SJames Wright // TSMonitor for the statistics collection and processing
461*a175e481SJames Wright PetscErrorCode TSMonitor_Statistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) {
462*a175e481SJames Wright   User user = (User)ctx;
463*a175e481SJames Wright   Vec  stats;
464*a175e481SJames Wright   PetscFunctionBeginUser;
465*a175e481SJames Wright 
466*a175e481SJames Wright   // Do not collect or process on the first step of the run (ie. on the initial condition)
467*a175e481SJames Wright   if (steps == user->app_ctx->cont_steps) PetscFunctionReturn(0);
468*a175e481SJames Wright 
469*a175e481SJames Wright   if (steps % user->app_ctx->stats_collect_interval == 0) {
470*a175e481SJames Wright     PetscCall(CollectStatistics(user, solution_time, Q));
471*a175e481SJames Wright   }
472*a175e481SJames Wright 
473*a175e481SJames Wright   if (steps % user->app_ctx->stats_write_interval == 0 && user->app_ctx->stats_write_interval != -1) {
474*a175e481SJames Wright     PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats));
475*a175e481SJames Wright     PetscCall(ProcessStatistics(user, &stats));
476*a175e481SJames Wright     PetscCall(VecViewFromOptions(stats, NULL, "-stats_write_view"));
477*a175e481SJames Wright     PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats));
478*a175e481SJames Wright   }
479*a175e481SJames Wright   PetscFunctionReturn(0);
480*a175e481SJames Wright }
481*a175e481SJames Wright 
482*a175e481SJames Wright // Function to be called at the end of a simulation
483*a175e481SJames Wright PetscErrorCode StatsCollectFinalCall(User user, PetscReal solution_time, Vec Q) {
484*a175e481SJames Wright   Vec stats;
485*a175e481SJames Wright   PetscFunctionBeginUser;
486*a175e481SJames Wright 
487*a175e481SJames Wright   PetscCall(CollectStatistics(user, solution_time, Q));
488*a175e481SJames Wright 
489*a175e481SJames Wright   PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats));
490*a175e481SJames Wright   PetscCall(ProcessStatistics(user, &stats));
491*a175e481SJames Wright   PetscCall(VecViewFromOptions(stats, NULL, "-stats_write_view"));
492*a175e481SJames Wright 
493*a175e481SJames Wright   PetscScalar *stats_arr;
494*a175e481SJames Wright   PetscCall(VecGetArray(stats, &stats_arr));
495*a175e481SJames Wright   PetscCall(VecRestoreArray(stats, &stats_arr));
496*a175e481SJames Wright   PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats));
497*a175e481SJames Wright 
49851ee423eSJames Wright   PetscFunctionReturn(0);
49951ee423eSJames Wright }
500