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 8a175e481SJames Wright /// Functions for setting up and performing statistics collection 9a175e481SJames Wright 10a175e481SJames Wright #include "../qfunctions/turb_spanstats.h" 1151ee423eSJames Wright 1219706a06SJames Wright #include <petscsf.h> 1319706a06SJames Wright 14a175e481SJames Wright #include "../include/matops.h" 1551ee423eSJames Wright #include "../navierstokes.h" 16a175e481SJames Wright #include "ceed/ceed.h" 17a175e481SJames Wright #include "petscerror.h" 186dcea3beSJames Wright #include "petsclog.h" 19a175e481SJames Wright #include "petscmat.h" 2019706a06SJames Wright #include "petscsys.h" 21a175e481SJames Wright #include "petscvec.h" 2251ee423eSJames Wright 2351ee423eSJames Wright PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree, SimpleBC bc) { 24a175e481SJames Wright user->spanstats.num_comp_stats = 22; 25a175e481SJames Wright PetscReal domain_min[3], domain_max[3]; 2651ee423eSJames Wright PetscFunctionBeginUser; 2751ee423eSJames Wright 28a175e481SJames Wright // Get spanwise length 29a175e481SJames Wright PetscCall(DMGetBoundingBox(user->dm, domain_min, domain_max)); 30a175e481SJames Wright user->spanstats.span_width = domain_max[2] - domain_min[1]; 31a175e481SJames Wright 3251ee423eSJames Wright // Get DM from surface 3351ee423eSJames Wright { 3451ee423eSJames Wright DMLabel label; 3551ee423eSJames Wright PetscCall(DMGetLabel(user->dm, "Face Sets", &label)); 3651ee423eSJames Wright PetscCall(DMPlexLabelComplete(user->dm, label)); 3751ee423eSJames Wright PetscCall(DMPlexFilter(user->dm, label, 1, &user->spanstats.dm)); 3851ee423eSJames Wright PetscCall(DMProjectCoordinates(user->spanstats.dm, NULL)); // Ensure that a coordinate FE exists 3951ee423eSJames Wright } 4051ee423eSJames Wright 4151ee423eSJames Wright PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats")); 4251ee423eSJames Wright PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "spanstats_")); 4351ee423eSJames Wright PetscCall(DMSetFromOptions(user->spanstats.dm)); 44a175e481SJames Wright PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view")); // -spanstats_dm_view 4551ee423eSJames Wright { 4651ee423eSJames Wright PetscFE fe; 4751ee423eSJames Wright DMLabel label; 4851ee423eSJames Wright 4951ee423eSJames Wright PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim - 1, user->spanstats.num_comp_stats, PETSC_FALSE, degree, PETSC_DECIDE, &fe)); 5051ee423eSJames Wright PetscCall(PetscObjectSetName((PetscObject)fe, "stats")); 5151ee423eSJames Wright PetscCall(DMAddField(user->spanstats.dm, NULL, (PetscObject)fe)); 5251ee423eSJames Wright PetscCall(DMCreateDS(user->spanstats.dm)); 5351ee423eSJames Wright PetscCall(DMGetLabel(user->spanstats.dm, "Face Sets", &label)); 5451ee423eSJames Wright 5551ee423eSJames Wright PetscCall(DMPlexSetClosurePermutationTensor(user->spanstats.dm, PETSC_DETERMINE, NULL)); 5651ee423eSJames Wright PetscCall(PetscFEDestroy(&fe)); 5751ee423eSJames Wright } 5851ee423eSJames Wright 5951ee423eSJames Wright PetscSection section; 6051ee423eSJames Wright PetscCall(DMGetLocalSection(user->spanstats.dm, §ion)); 6151ee423eSJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 62a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "Mean Density")); 63a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "Mean Pressure")); 64a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "Mean Pressure Squared")); 65a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "Mean Pressure Velocity X")); 66a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "Mean Pressure Velocity Y")); 67a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 5, "Mean Pressure Velocity Z")); 68a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 6, "Mean Density Temperature")); 69a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 7, "Mean Density Temperature Flux X")); 70a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 8, "Mean Density Temperature Flux Y")); 71a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 9, "Mean Density Temperature Flux Z")); 72a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 10, "Mean Momentum X")); 73a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 11, "Mean Momentum Y")); 74a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 12, "Mean Momentum Z")); 75a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 13, "Mean Momentum Flux XX")); 76a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 14, "Mean Momentum Flux YY")); 77a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 15, "Mean Momentum Flux ZZ")); 78a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 16, "Mean Momentum Flux YZ")); 79a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 17, "Mean Momentum Flux XZ")); 80a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 18, "Mean Momentum Flux XY")); 81a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 19, "Mean Velocity X")); 82a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 20, "Mean Velocity Y")); 83a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 21, "Mean Velocity Z")); 8419706a06SJames Wright 8519706a06SJames Wright PetscFunctionReturn(0); 8619706a06SJames Wright } 8719706a06SJames Wright 881737222fSJames Wright // Create CeedElemRestriction for collocated data based on associated CeedBasis and CeedElemRestriction 891737222fSJames Wright // Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction 901737222fSJames Wright PetscErrorCode CreateElemRestrColloc(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base, 911737222fSJames Wright CeedElemRestriction *elem_restr_collocated, CeedVector *l_vec, CeedVector *e_vec) { 921737222fSJames Wright CeedInt num_elem_qpts, loc_num_elem; 931737222fSJames Wright PetscFunctionBeginUser; 941737222fSJames Wright 951737222fSJames Wright CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts); 961737222fSJames Wright CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem); 971737222fSJames Wright 981737222fSJames Wright const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp}; 991737222fSJames Wright CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides, 1001737222fSJames Wright elem_restr_collocated); 1011737222fSJames Wright CeedElemRestrictionCreateVector(*elem_restr_collocated, l_vec, e_vec); 1021737222fSJames Wright PetscFunctionReturn(0); 1031737222fSJames Wright } 1041737222fSJames Wright 105a175e481SJames Wright // Get coordinates of quadrature points 10619706a06SJames Wright PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, CeedVector *qx_coords, 10719706a06SJames Wright PetscInt *total_nqpnts) { 10819706a06SJames Wright CeedQFunction qf_quad_coords; 10919706a06SJames Wright CeedOperator op_quad_coords; 11019706a06SJames Wright PetscInt num_comp_x, loc_num_elem, num_elem_qpts; 11119706a06SJames Wright CeedElemRestriction elem_restr_qx; 11219706a06SJames Wright PetscFunctionBeginUser; 11319706a06SJames Wright 11419706a06SJames Wright // Create Element Restriction and CeedVector for quadrature coordinates 11519706a06SJames Wright CeedBasisGetNumQuadraturePoints(basis_x, &num_elem_qpts); 11619706a06SJames Wright CeedElemRestrictionGetNumElements(elem_restr_x, &loc_num_elem); 11719706a06SJames Wright CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x); 11819706a06SJames Wright *total_nqpnts = num_elem_qpts * loc_num_elem; 1191737222fSJames Wright PetscCall(CreateElemRestrColloc(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx, qx_coords, NULL)); 12019706a06SJames Wright 12119706a06SJames Wright // Create QFunction 12219706a06SJames Wright CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords); 12319706a06SJames Wright 12419706a06SJames Wright // Create Operator 12519706a06SJames Wright CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords); 12619706a06SJames Wright CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE); 12719706a06SJames Wright CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 12819706a06SJames Wright 12919706a06SJames Wright CeedOperatorApply(op_quad_coords, x_coords, *qx_coords, CEED_REQUEST_IMMEDIATE); 13019706a06SJames Wright 13119706a06SJames Wright CeedQFunctionDestroy(&qf_quad_coords); 13219706a06SJames Wright CeedOperatorDestroy(&op_quad_coords); 13319706a06SJames Wright PetscFunctionReturn(0); 13419706a06SJames Wright } 13519706a06SJames Wright 136a175e481SJames Wright // Create PetscSF for child-to-parent communication 13719706a06SJames Wright PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, DM parentdm, DM childdm, PetscSF statssf) { 13819706a06SJames Wright PetscInt child_num_qpnts, parent_num_qpnts, num_comp_x; 13919706a06SJames Wright CeedVector child_qx_coords, parent_qx_coords; 14019706a06SJames Wright PetscReal *child_coords, *parent_coords; 14119706a06SJames Wright PetscFunctionBeginUser; 14219706a06SJames Wright 14319706a06SJames Wright // Assume that child and parent have the same number of components 14419706a06SJames Wright CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x); 14519706a06SJames Wright const PetscInt num_comp_sf = num_comp_x - 1; // Number of coord components used in the creation of the SF 14619706a06SJames Wright 14719706a06SJames Wright // Get quad_coords for child DM 148a175e481SJames Wright PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &child_qx_coords, &child_num_qpnts)); 14919706a06SJames Wright 15019706a06SJames Wright // Get quad_coords for parent DM 15119706a06SJames Wright PetscCall(GetQuadratureCoords(ceed, parentdm, ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, ceed_data->spanstats.x_coord, 15219706a06SJames Wright &parent_qx_coords, &parent_num_qpnts)); 15319706a06SJames Wright 15419706a06SJames Wright // Remove z component of coordinates for matching 15519706a06SJames Wright { 15619706a06SJames Wright const PetscReal *child_quad_coords, *parent_quad_coords; 15719706a06SJames Wright 15819706a06SJames Wright CeedVectorGetArrayRead(child_qx_coords, CEED_MEM_HOST, &child_quad_coords); 15919706a06SJames Wright CeedVectorGetArrayRead(parent_qx_coords, CEED_MEM_HOST, &parent_quad_coords); 16019706a06SJames Wright 16119706a06SJames Wright PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords)); 16219706a06SJames Wright for (int i = 0; i < child_num_qpnts; i++) { 16319706a06SJames Wright child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x]; 16419706a06SJames Wright child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x]; 16519706a06SJames Wright } 16619706a06SJames Wright for (int i = 0; i < parent_num_qpnts; i++) { 16719706a06SJames Wright parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x]; 16819706a06SJames Wright parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x]; 16919706a06SJames Wright } 17019706a06SJames Wright CeedVectorRestoreArrayRead(child_qx_coords, &child_quad_coords); 17119706a06SJames Wright CeedVectorRestoreArrayRead(parent_qx_coords, &parent_quad_coords); 17219706a06SJames Wright } 17319706a06SJames Wright 17419706a06SJames Wright PetscCall(PetscSFSetGraphFromCoordinates(statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords)); 17519706a06SJames Wright 176a175e481SJames Wright PetscCall(PetscSFViewFromOptions(statssf, NULL, "-spanstats_sf_view")); 177a175e481SJames Wright 17819706a06SJames Wright PetscCall(PetscFree2(child_coords, parent_coords)); 17919706a06SJames Wright CeedVectorDestroy(&child_qx_coords); 18019706a06SJames Wright CeedVectorDestroy(&parent_qx_coords); 18119706a06SJames Wright PetscFunctionReturn(0); 18219706a06SJames Wright } 18319706a06SJames Wright 184a175e481SJames Wright // Compute mass matrix for statistics projection 185a175e481SJames Wright PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data) { 186a175e481SJames Wright CeedQFunction qf_mass; 187a175e481SJames Wright CeedOperator op_mass; 188a175e481SJames Wright CeedInt num_comp_q, q_data_size; 189a175e481SJames Wright PetscFunctionBeginUser; 190a175e481SJames Wright 191a175e481SJames Wright // CEED Restriction 192a175e481SJames Wright CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_stats, &num_comp_q); 193a175e481SJames Wright CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_qd, &q_data_size); 194a175e481SJames Wright 195a175e481SJames Wright // Create Mass CeedOperator 196a175e481SJames Wright PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass)); 197a175e481SJames Wright CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass); 198a175e481SJames Wright CeedOperatorSetField(op_mass, "q", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE); 199a175e481SJames Wright CeedOperatorSetField(op_mass, "qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, ceed_data->spanstats.q_data); 200a175e481SJames Wright CeedOperatorSetField(op_mass, "v", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE); 201a175e481SJames Wright 202a175e481SJames Wright // Setup KSP for L^2 projection 203a175e481SJames Wright { 204a175e481SJames Wright MatopApplyContext M_ctx; 205a175e481SJames Wright PetscInt l_size, g_size; 206a175e481SJames Wright Mat mat_mass; 207a175e481SJames Wright VecType vec_type; 208a175e481SJames Wright KSP ksp; 209a175e481SJames Wright Vec ones, M_inv; 210a175e481SJames Wright CeedVector x_ceed, y_ceed; 211a175e481SJames Wright 212a175e481SJames Wright PetscCall(DMCreateGlobalVector(user->spanstats.dm, &M_inv)); 213a175e481SJames Wright PetscCall(VecGetLocalSize(M_inv, &l_size)); 214a175e481SJames Wright PetscCall(VecGetSize(M_inv, &g_size)); 215a175e481SJames Wright PetscCall(VecGetType(M_inv, &vec_type)); 216a175e481SJames Wright 217a175e481SJames Wright PetscCall(PetscMalloc1(1, &M_ctx)); 218a175e481SJames Wright PetscCall(MatCreateShell(PETSC_COMM_WORLD, l_size, l_size, g_size, g_size, M_ctx, &mat_mass)); 219a175e481SJames Wright PetscCall(MatShellSetOperation(mat_mass, MATOP_MULT, (void (*)(void))MatMult_Ceed)); 220a175e481SJames Wright PetscCall(MatShellSetOperation(mat_mass, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed)); 221a175e481SJames Wright PetscCall(MatShellSetVecType(mat_mass, vec_type)); 222a175e481SJames Wright 223a175e481SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &x_ceed, NULL); 224a175e481SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &y_ceed, NULL); 225a175e481SJames Wright 226a175e481SJames Wright PetscCall(SetupMatopApplyCtx(PETSC_COMM_WORLD, user->spanstats.dm, user->ceed, op_mass, x_ceed, y_ceed, NULL, M_ctx)); 227a175e481SJames Wright user->spanstats.M_ctx = M_ctx; 228a175e481SJames Wright 229a175e481SJames Wright // Create lumped mass matrix inverse 230a175e481SJames Wright PetscCall(DMGetGlobalVector(user->spanstats.dm, &ones)); 231a175e481SJames Wright PetscCall(VecZeroEntries(M_inv)); 232a175e481SJames Wright PetscCall(VecSet(ones, 1)); 233a175e481SJames Wright PetscCall(MatMult(mat_mass, ones, M_inv)); 234a175e481SJames Wright PetscCall(VecReciprocal(M_inv)); 235a175e481SJames Wright user->spanstats.M_inv = M_inv; 236a175e481SJames Wright PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &ones)); 237a175e481SJames Wright 238a175e481SJames Wright PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp)); 239a175e481SJames Wright PetscCall(KSPSetOptionsPrefix(ksp, "spanstats_")); 240a175e481SJames Wright { 241a175e481SJames Wright PC pc; 242a175e481SJames Wright PetscCall(KSPGetPC(ksp, &pc)); 243a175e481SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 244a175e481SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL)); 245a175e481SJames Wright PetscCall(KSPSetType(ksp, KSPCG)); 246a175e481SJames Wright PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL)); 247a175e481SJames Wright PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 248a175e481SJames Wright } 249a175e481SJames Wright PetscCall(KSPSetOperators(ksp, mat_mass, mat_mass)); 250a175e481SJames Wright PetscCall(KSPSetFromOptions(ksp)); 251a175e481SJames Wright user->spanstats.ksp = ksp; 252a175e481SJames Wright } 253a175e481SJames Wright 254a175e481SJames Wright // Cleanup 255a175e481SJames Wright CeedQFunctionDestroy(&qf_mass); 256a175e481SJames Wright PetscFunctionReturn(0); 257a175e481SJames Wright } 258a175e481SJames Wright 259a175e481SJames Wright // Create CeedOperators and KSP for the statistics collection and processing 260a175e481SJames Wright PetscErrorCode CreateStatisticsOperators(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) { 261a175e481SJames Wright CeedInt num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q; 262a175e481SJames Wright CeedOperator op_setup_sur; 263a175e481SJames Wright PetscFunctionBeginUser; 264a175e481SJames Wright CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q); 265a175e481SJames Wright 266a175e481SJames Wright // Create Operator for statistics collection 267a175e481SJames Wright switch (user->phys->state_var) { 268a175e481SJames Wright case STATEVAR_PRIMITIVE: 269a175e481SJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &ceed_data->spanstats.qf_stats_collect); 270a175e481SJames Wright break; 271a175e481SJames Wright case STATEVAR_CONSERVATIVE: 272a175e481SJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &ceed_data->spanstats.qf_stats_collect); 273a175e481SJames Wright break; 274a175e481SJames Wright } 275a175e481SJames Wright 276a175e481SJames Wright if (user->app_ctx->stats_test) { 277a175e481SJames Wright CeedQFunctionDestroy(&ceed_data->spanstats.qf_stats_collect); 278a175e481SJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionTest, ChildStatsCollectionTest_loc, &ceed_data->spanstats.qf_stats_collect); 279a175e481SJames Wright } 280a175e481SJames Wright 281a175e481SJames Wright CeedQFunctionSetContext(ceed_data->spanstats.qf_stats_collect, problem->apply_vol_ifunction.qfunction_context); 282a175e481SJames Wright CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP); 283a175e481SJames Wright CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE); 284a175e481SJames Wright CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP); 285a175e481SJames Wright CeedQFunctionAddOutput(ceed_data->spanstats.qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE); 286a175e481SJames Wright 287a175e481SJames Wright CeedOperatorCreate(ceed, ceed_data->spanstats.qf_stats_collect, NULL, NULL, &user->spanstats.op_stats_collect); 288a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 289a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 290a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord); 291a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_collect, "v", ceed_data->spanstats.elem_restr_child_colloc, CEED_BASIS_COLLOCATED, 292a175e481SJames Wright CEED_VECTOR_ACTIVE); 293a175e481SJames Wright 294823a1283SJames Wright // Create Operator for L^2 projection of statistics 295a175e481SJames Wright // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function. 296a175e481SJames Wright // Therefore, an Identity QF is sufficient 297a175e481SJames Wright CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &ceed_data->spanstats.qf_stats_proj); 298a175e481SJames Wright 299a175e481SJames Wright CeedOperatorCreate(ceed, ceed_data->spanstats.qf_stats_proj, NULL, NULL, &user->spanstats.op_stats_proj); 300a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_proj, "input", ceed_data->spanstats.elem_restr_parent_colloc, CEED_BASIS_COLLOCATED, 301a175e481SJames Wright CEED_VECTOR_ACTIVE); 302a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_proj, "output", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, 303a175e481SJames Wright CEED_VECTOR_ACTIVE); 304a175e481SJames Wright 305a175e481SJames Wright // Get q_data for lumped mass matrix formation 306a175e481SJames Wright CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur); 307a175e481SJames Wright CeedOperatorSetField(op_setup_sur, "dx", ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, CEED_VECTOR_ACTIVE); 308a175e481SJames Wright CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->spanstats.basis_x, CEED_VECTOR_NONE); 309a175e481SJames Wright CeedOperatorSetField(op_setup_sur, "surface qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 310a175e481SJames Wright CeedOperatorApply(op_setup_sur, ceed_data->spanstats.x_coord, ceed_data->spanstats.q_data, CEED_REQUEST_IMMEDIATE); 311a175e481SJames Wright 312a175e481SJames Wright CeedOperatorDestroy(&op_setup_sur); 313a175e481SJames Wright PetscFunctionReturn(0); 314a175e481SJames Wright } 315a175e481SJames Wright 316823a1283SJames Wright PetscErrorCode SetupErrorTesting(Ceed ceed, User user, CeedData ceed_data) { 317823a1283SJames Wright CeedInt num_comp_stats = user->spanstats.num_comp_stats, num_comp_x; 318823a1283SJames Wright CeedQFunction qf_error; 319823a1283SJames Wright CeedOperator op_error; 320823a1283SJames Wright CeedInt q_data_size; 321823a1283SJames Wright CeedVector x_ceed, y_ceed; 322823a1283SJames Wright PetscFunctionBeginUser; 323823a1283SJames Wright 324823a1283SJames Wright CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_qd, &q_data_size); 325823a1283SJames Wright CeedBasisGetNumComponents(ceed_data->spanstats.basis_x, &num_comp_x); 326823a1283SJames Wright 327823a1283SJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionTest_Error, ChildStatsCollectionTest_Error_loc, &qf_error); 328823a1283SJames Wright CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP); 329823a1283SJames Wright CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE); 330823a1283SJames Wright CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP); 331823a1283SJames Wright CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP); 332823a1283SJames Wright 333823a1283SJames Wright CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error); 334823a1283SJames Wright CeedOperatorSetField(op_error, "q", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE); 335823a1283SJames Wright CeedOperatorSetField(op_error, "qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, ceed_data->spanstats.q_data); 336823a1283SJames Wright CeedOperatorSetField(op_error, "x", ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, ceed_data->spanstats.x_coord); 337823a1283SJames Wright CeedOperatorSetField(op_error, "v", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE); 338823a1283SJames Wright 339823a1283SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &x_ceed, NULL); 340823a1283SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &y_ceed, NULL); 341823a1283SJames Wright 342823a1283SJames Wright PetscCall(PetscCalloc1(1, &user->spanstats.test_error_ctx)); 343823a1283SJames Wright PetscCall(SetupMatopApplyCtx(PETSC_COMM_WORLD, user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, user->spanstats.test_error_ctx)); 344823a1283SJames Wright 345823a1283SJames Wright PetscFunctionReturn(0); 346823a1283SJames Wright } 347823a1283SJames Wright 348a175e481SJames Wright // Setup for statistics collection 34919706a06SJames Wright PetscErrorCode SetupStatsCollection(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) { 35019706a06SJames Wright DM dm = user->spanstats.dm; 35119706a06SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm); 35219706a06SJames Wright CeedInt dim, P, Q, num_comp_x; 35319706a06SJames Wright Vec X_loc; 35419706a06SJames Wright PetscMemType X_loc_memtype; 35519706a06SJames Wright const PetscScalar *X_loc_array; 35619706a06SJames Wright PetscFunctionBeginUser; 35719706a06SJames Wright 35819706a06SJames Wright PetscCall(DMGetDimension(dm, &dim)); 359a175e481SJames Wright CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &Q); 360a175e481SJames Wright CeedBasisGetNumNodes1D(ceed_data->basis_q, &P); 36119706a06SJames Wright 3621737222fSJames Wright PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, problem->q_data_size_sur, &ceed_data->spanstats.elem_restr_parent_stats, 36319706a06SJames Wright &ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.elem_restr_parent_qd)); 36419706a06SJames Wright CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_x, &num_comp_x); 36519706a06SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.x_coord, NULL); 366a175e481SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &user->spanstats.rhs_ceed, NULL); 3671737222fSJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_qd, &ceed_data->spanstats.q_data, NULL); 36819706a06SJames Wright 369a175e481SJames Wright CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS, &ceed_data->spanstats.basis_x); 37019706a06SJames Wright CeedBasisCreateTensorH1Lagrange(ceed, dim, user->spanstats.num_comp_stats, P, Q, CEED_GAUSS, &ceed_data->spanstats.basis_stats); 37119706a06SJames Wright 3721737222fSJames Wright PetscCall(CreateElemRestrColloc(ceed, user->spanstats.num_comp_stats, ceed_data->spanstats.basis_stats, 3731737222fSJames Wright ceed_data->spanstats.elem_restr_parent_stats, &ceed_data->spanstats.elem_restr_parent_colloc, 3741737222fSJames Wright &user->spanstats.parent_stats, NULL)); 3751737222fSJames Wright PetscCall(CreateElemRestrColloc(ceed, user->spanstats.num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q, 3761737222fSJames Wright &ceed_data->spanstats.elem_restr_child_colloc, &user->spanstats.child_stats, NULL)); 377a175e481SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_child_colloc, &user->spanstats.child_inst_stats, NULL); 378a175e481SJames Wright CeedVectorSetValue(user->spanstats.child_stats, 0); 3791737222fSJames Wright 38019706a06SJames Wright // -- Copy DM coordinates into CeedVector 38119706a06SJames Wright { 38219706a06SJames Wright DM cdm; 38319706a06SJames Wright PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 38419706a06SJames Wright if (cdm) { 38519706a06SJames Wright PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 38619706a06SJames Wright } else { 38719706a06SJames Wright PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 38819706a06SJames Wright } 38919706a06SJames Wright } 39019706a06SJames Wright PetscCall(VecScale(X_loc, problem->dm_scale)); 39119706a06SJames Wright PetscCall(VecGetArrayReadAndMemType(X_loc, &X_loc_array, &X_loc_memtype)); 39219706a06SJames Wright CeedVectorSetArray(ceed_data->spanstats.x_coord, MemTypeP2C(X_loc_memtype), CEED_COPY_VALUES, (PetscScalar *)X_loc_array); 39319706a06SJames Wright PetscCall(VecRestoreArrayRead(X_loc, &X_loc_array)); 39419706a06SJames Wright 39519706a06SJames Wright // Create SF for communicating child data back their respective parents 39619706a06SJames Wright PetscCall(PetscSFCreate(comm, &user->spanstats.sf)); 39719706a06SJames Wright PetscCall(CreateStatsSF(ceed, ceed_data, user->dm, user->spanstats.dm, user->spanstats.sf)); 39851ee423eSJames Wright 399a175e481SJames Wright // Create CeedOperators for statistics collection 400a175e481SJames Wright PetscCall(CreateStatisticsOperators(ceed, user, ceed_data, problem)); 401a175e481SJames Wright 402a175e481SJames Wright // Setup KSP and Mat for L^2 projection of statistics 403a175e481SJames Wright PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data)); 404a175e481SJames Wright 405823a1283SJames Wright if (user->app_ctx->stats_test) { 406823a1283SJames Wright PetscCall(SetupErrorTesting(ceed, user, ceed_data)); 407823a1283SJames Wright } 408823a1283SJames Wright 409a175e481SJames Wright PetscFunctionReturn(0); 410a175e481SJames Wright } 411a175e481SJames Wright 412a175e481SJames Wright // Collect statistics based on the solution Q 413a175e481SJames Wright PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) { 414a175e481SJames Wright PetscMemType q_mem_type; 415a175e481SJames Wright const PetscScalar *q_arr; 416a175e481SJames Wright Vec Q_loc; 417a175e481SJames Wright PetscFunctionBeginUser; 418a175e481SJames Wright 4196dcea3beSJames Wright PetscLogStage stage_stats_collect; 4206dcea3beSJames Wright PetscCall(PetscLogStageGetId("Stats Collect", &stage_stats_collect)); 4216dcea3beSJames Wright if (stage_stats_collect == -1) PetscCall(PetscLogStageRegister("Stats Collect", &stage_stats_collect)); 4226dcea3beSJames Wright PetscCall(PetscLogStagePush(stage_stats_collect)); 4236dcea3beSJames Wright 424a175e481SJames Wright PetscCall(DMGetLocalVector(user->dm, &Q_loc)); 425a175e481SJames Wright PetscCall(VecZeroEntries(Q_loc)); 426a175e481SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 427a175e481SJames Wright 428a175e481SJames Wright PetscCall(VecGetArrayReadAndMemType(Q_loc, &q_arr, &q_mem_type)); 429a175e481SJames Wright CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q_arr); 430a175e481SJames Wright 431a175e481SJames Wright CeedOperatorApply(user->spanstats.op_stats_collect, user->q_ceed, user->spanstats.child_inst_stats, CEED_REQUEST_IMMEDIATE); 432a175e481SJames Wright 433a175e481SJames Wright CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL); 434a175e481SJames Wright PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q_arr)); 435a175e481SJames Wright PetscCall(DMRestoreLocalVector(user->dm, &Q_loc)); 436a175e481SJames Wright 437a175e481SJames Wright // Record averaging using left rectangle rule 438a175e481SJames Wright PetscScalar delta_t = solution_time - user->spanstats.prev_time; 439a175e481SJames Wright PetscScalar prev_timeinterval = user->spanstats.prev_time - user->app_ctx->cont_time; 440a175e481SJames Wright CeedVectorScale(user->spanstats.child_stats, prev_timeinterval / (prev_timeinterval + delta_t)); 441a175e481SJames Wright CeedVectorAXPY(user->spanstats.child_stats, delta_t / (prev_timeinterval + delta_t), user->spanstats.child_inst_stats); 442a175e481SJames Wright user->spanstats.prev_time = solution_time; 443a175e481SJames Wright 4446dcea3beSJames Wright PetscCall(PetscLogStagePop()); 445a175e481SJames Wright PetscFunctionReturn(0); 446a175e481SJames Wright } 447a175e481SJames Wright 448a175e481SJames Wright // Process the child statistics into parent statistics and project them onto stats 449a175e481SJames Wright PetscErrorCode ProcessStatistics(User user, Vec *stats) { 450a175e481SJames Wright Span_Stats user_stats = user->spanstats; 451a175e481SJames Wright const PetscScalar *child_stats; 452a175e481SJames Wright PetscScalar *parent_stats; 453a175e481SJames Wright MPI_Datatype unit; 454a175e481SJames Wright Vec rhs_loc, rhs; 455a175e481SJames Wright PetscMemType rhs_mem_type; 456a175e481SJames Wright CeedScalar *rhs_arr; 457a175e481SJames Wright CeedMemType ceed_mem_type; 458a175e481SJames Wright PetscFunctionBeginUser; 459a175e481SJames Wright 4606dcea3beSJames Wright PetscLogStage stage_stats_process; 4616dcea3beSJames Wright PetscCall(PetscLogStageGetId("Stats Process", &stage_stats_process)); 4626dcea3beSJames Wright if (stage_stats_process == -1) PetscCall(PetscLogStageRegister("Stats Process", &stage_stats_process)); 4636dcea3beSJames Wright PetscCall(PetscLogStagePush(stage_stats_process)); 4646dcea3beSJames Wright 465a175e481SJames Wright CeedGetPreferredMemType(user->ceed, &ceed_mem_type); 466a175e481SJames Wright CeedVectorSetValue(user_stats.parent_stats, 0); 467a175e481SJames Wright 468a175e481SJames Wright CeedVectorGetArrayRead(user_stats.child_stats, ceed_mem_type, &child_stats); 469a175e481SJames Wright CeedVectorGetArray(user_stats.parent_stats, ceed_mem_type, &parent_stats); 470a175e481SJames Wright 471a175e481SJames Wright if (user_stats.num_comp_stats == 1) unit = MPIU_REAL; 472a175e481SJames Wright else { 473a175e481SJames Wright PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit)); 474a175e481SJames Wright PetscCallMPI(MPI_Type_commit(&unit)); 475a175e481SJames Wright } 476a175e481SJames Wright 477a175e481SJames Wright PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM)); 478a175e481SJames Wright PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM)); 479a175e481SJames Wright 480a175e481SJames Wright CeedVectorRestoreArrayRead(user_stats.child_stats, &child_stats); 481a175e481SJames Wright CeedVectorRestoreArray(user_stats.parent_stats, &parent_stats); 482a175e481SJames Wright PetscCallMPI(MPI_Type_free(&unit)); 483a175e481SJames Wright 484a175e481SJames Wright CeedVectorScale(user_stats.parent_stats, 1 / user_stats.span_width); 485a175e481SJames Wright 486a175e481SJames Wright // L^2 projection with the parent_data 487a175e481SJames Wright PetscCall(DMGetGlobalVector(user_stats.dm, &rhs)); 488a175e481SJames Wright PetscCall(DMGetLocalVector(user_stats.dm, &rhs_loc)); 489a175e481SJames Wright PetscCall(VecZeroEntries(rhs)); 490a175e481SJames Wright PetscCall(VecZeroEntries(rhs_loc)); 491a175e481SJames Wright PetscCall(VecGetArrayWriteAndMemType(rhs_loc, &rhs_arr, &rhs_mem_type)); 492a175e481SJames Wright CeedVectorSetArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), CEED_USE_POINTER, (PetscScalar *)rhs_arr); 493a175e481SJames Wright 494a175e481SJames Wright CeedOperatorApply(user_stats.op_stats_proj, user_stats.parent_stats, user_stats.rhs_ceed, CEED_REQUEST_IMMEDIATE); 495a175e481SJames Wright 496a175e481SJames Wright CeedVectorTakeArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), &rhs_arr); 497a175e481SJames Wright PetscCall(VecRestoreArrayAndMemType(rhs_loc, &rhs_arr)); 498a175e481SJames Wright PetscCall(DMLocalToGlobal(user_stats.dm, rhs_loc, ADD_VALUES, rhs)); 499a175e481SJames Wright 500a175e481SJames Wright PetscCall(VecDuplicate(rhs, stats)); 501a175e481SJames Wright PetscCall(VecPointwiseMult(*stats, rhs, user_stats.M_inv)); 502a175e481SJames Wright 503a175e481SJames Wright PetscCall(KSPSolve(user_stats.ksp, rhs, *stats)); 504a175e481SJames Wright 5056dcea3beSJames Wright PetscCall(PetscLogStagePop()); 506a175e481SJames Wright PetscFunctionReturn(0); 507a175e481SJames Wright } 508a175e481SJames Wright 509a175e481SJames Wright // TSMonitor for the statistics collection and processing 510a175e481SJames Wright PetscErrorCode TSMonitor_Statistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) { 511a175e481SJames Wright User user = (User)ctx; 512a175e481SJames Wright Vec stats; 513a175e481SJames Wright PetscFunctionBeginUser; 514a175e481SJames Wright 515a175e481SJames Wright // Do not collect or process on the first step of the run (ie. on the initial condition) 516a175e481SJames Wright if (steps == user->app_ctx->cont_steps) PetscFunctionReturn(0); 517a175e481SJames Wright 518a175e481SJames Wright if (steps % user->app_ctx->stats_collect_interval == 0) { 519a175e481SJames Wright PetscCall(CollectStatistics(user, solution_time, Q)); 520a175e481SJames Wright } 521a175e481SJames Wright 522a175e481SJames Wright if (steps % user->app_ctx->stats_write_interval == 0 && user->app_ctx->stats_write_interval != -1) { 523a175e481SJames Wright PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats)); 524a175e481SJames Wright PetscCall(ProcessStatistics(user, &stats)); 525a175e481SJames Wright PetscCall(VecViewFromOptions(stats, NULL, "-stats_write_view")); 526a175e481SJames Wright PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats)); 527a175e481SJames Wright } 528a175e481SJames Wright PetscFunctionReturn(0); 529a175e481SJames Wright } 530a175e481SJames Wright // Function to be called at the end of a simulation 531a175e481SJames Wright PetscErrorCode StatsCollectFinalCall(User user, PetscReal solution_time, Vec Q) { 532a175e481SJames Wright Vec stats; 533a175e481SJames Wright PetscFunctionBeginUser; 534a175e481SJames Wright 535a175e481SJames Wright PetscCall(CollectStatistics(user, solution_time, Q)); 536a175e481SJames Wright 537a175e481SJames Wright PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats)); 538a175e481SJames Wright PetscCall(ProcessStatistics(user, &stats)); 539a175e481SJames Wright PetscCall(VecViewFromOptions(stats, NULL, "-stats_write_view")); 540a175e481SJames Wright 541823a1283SJames Wright if (user->app_ctx->stats_test) { 542823a1283SJames Wright Vec error; 543823a1283SJames Wright PetscCall(VecDuplicate(stats, &error)); 544823a1283SJames Wright PetscCall(ApplyLocal_Ceed(stats, error, user->spanstats.test_error_ctx)); 545823a1283SJames Wright PetscScalar error_sq = 0; 546823a1283SJames Wright PetscCall(VecSum(error, &error_sq)); 547823a1283SJames Wright PetscScalar l2_error = sqrt(error_sq); 548823a1283SJames Wright PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error)); 549823a1283SJames Wright } 550a175e481SJames Wright PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats)); 551a175e481SJames Wright 55251ee423eSJames Wright PetscFunctionReturn(0); 55351ee423eSJames Wright } 554*fd170fd0SJames Wright 555*fd170fd0SJames Wright PetscErrorCode CleanupStats(User user, CeedData ceed_data) { 556*fd170fd0SJames Wright PetscFunctionBeginUser; 557*fd170fd0SJames Wright 558*fd170fd0SJames Wright // -- CeedVectors 559*fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.child_stats); 560*fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.child_inst_stats); 561*fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.parent_stats); 562*fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.rhs_ceed); 563*fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.x_ceed); 564*fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.y_ceed); 565*fd170fd0SJames Wright CeedVectorDestroy(&ceed_data->spanstats.x_coord); 566*fd170fd0SJames Wright CeedVectorDestroy(&ceed_data->spanstats.q_data); 567*fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.M_ctx->x_ceed); 568*fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.M_ctx->y_ceed); 569*fd170fd0SJames Wright if (user->app_ctx->stats_test) { 570*fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.test_error_ctx->x_ceed); 571*fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.test_error_ctx->y_ceed); 572*fd170fd0SJames Wright } 573*fd170fd0SJames Wright 574*fd170fd0SJames Wright // -- QFunctions 575*fd170fd0SJames Wright CeedQFunctionDestroy(&ceed_data->spanstats.qf_stats_collect); 576*fd170fd0SJames Wright CeedQFunctionDestroy(&ceed_data->spanstats.qf_stats_proj); 577*fd170fd0SJames Wright 578*fd170fd0SJames Wright // -- CeedBasis 579*fd170fd0SJames Wright CeedBasisDestroy(&ceed_data->spanstats.basis_stats); 580*fd170fd0SJames Wright CeedBasisDestroy(&ceed_data->spanstats.basis_x); 581*fd170fd0SJames Wright 582*fd170fd0SJames Wright // -- CeedElemRestriction 583*fd170fd0SJames Wright CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_parent_stats); 584*fd170fd0SJames Wright CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_parent_qd); 585*fd170fd0SJames Wright CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_parent_x); 586*fd170fd0SJames Wright CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_parent_colloc); 587*fd170fd0SJames Wright CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_child_colloc); 588*fd170fd0SJames Wright 589*fd170fd0SJames Wright // -- CeedOperators 590*fd170fd0SJames Wright CeedOperatorDestroy(&user->spanstats.op_stats_collect); 591*fd170fd0SJames Wright CeedOperatorDestroy(&user->spanstats.op_stats_proj); 592*fd170fd0SJames Wright CeedOperatorDestroy(&user->spanstats.M_ctx->op); 593*fd170fd0SJames Wright if (user->app_ctx->stats_test) CeedOperatorDestroy(&user->spanstats.test_error_ctx->op); 594*fd170fd0SJames Wright 595*fd170fd0SJames Wright // -- Vec 596*fd170fd0SJames Wright PetscCall(VecDestroy(&user->spanstats.M_inv)); 597*fd170fd0SJames Wright 598*fd170fd0SJames Wright // -- KSP 599*fd170fd0SJames Wright PetscCall(KSPDestroy(&user->spanstats.ksp)); 600*fd170fd0SJames Wright 601*fd170fd0SJames Wright // -- SF 602*fd170fd0SJames Wright PetscCall(PetscSFDestroy(&user->spanstats.sf)); 603*fd170fd0SJames Wright 604*fd170fd0SJames Wright // -- DM 605*fd170fd0SJames Wright PetscCall(DMDestroy(&user->spanstats.dm)); 606*fd170fd0SJames Wright 607*fd170fd0SJames Wright PetscFunctionReturn(0); 608*fd170fd0SJames Wright } 609