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" 18a175e481SJames Wright #include "petscmat.h" 1919706a06SJames Wright #include "petscsys.h" 20a175e481SJames Wright #include "petscvec.h" 2151ee423eSJames Wright 2251ee423eSJames Wright PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree, SimpleBC bc) { 23a175e481SJames Wright user->spanstats.num_comp_stats = 22; 24a175e481SJames Wright PetscReal domain_min[3], domain_max[3]; 2551ee423eSJames Wright PetscFunctionBeginUser; 2651ee423eSJames Wright 27a175e481SJames Wright // Get spanwise length 28a175e481SJames Wright PetscCall(DMGetBoundingBox(user->dm, domain_min, domain_max)); 29a175e481SJames Wright user->spanstats.span_width = domain_max[2] - domain_min[1]; 30a175e481SJames 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)); 43a175e481SJames 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, §ion)); 6051ee423eSJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 61a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "Mean Density")); 62a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "Mean Pressure")); 63a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "Mean Pressure Squared")); 64a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "Mean Pressure Velocity X")); 65a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "Mean Pressure Velocity Y")); 66a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 5, "Mean Pressure Velocity Z")); 67a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 6, "Mean Density Temperature")); 68a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 7, "Mean Density Temperature Flux X")); 69a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 8, "Mean Density Temperature Flux Y")); 70a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 9, "Mean Density Temperature Flux Z")); 71a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 10, "Mean Momentum X")); 72a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 11, "Mean Momentum Y")); 73a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 12, "Mean Momentum Z")); 74a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 13, "Mean Momentum Flux XX")); 75a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 14, "Mean Momentum Flux YY")); 76a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 15, "Mean Momentum Flux ZZ")); 77a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 16, "Mean Momentum Flux YZ")); 78a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 17, "Mean Momentum Flux XZ")); 79a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 18, "Mean Momentum Flux XY")); 80a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 19, "Mean Velocity X")); 81a175e481SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 20, "Mean Velocity Y")); 82a175e481SJames 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 104a175e481SJames 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 135a175e481SJames 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 147a175e481SJames 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 175a175e481SJames Wright PetscCall(PetscSFViewFromOptions(statssf, NULL, "-spanstats_sf_view")); 176a175e481SJames 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 183a175e481SJames Wright // Compute mass matrix for statistics projection 184a175e481SJames Wright PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data) { 185a175e481SJames Wright CeedQFunction qf_mass; 186a175e481SJames Wright CeedOperator op_mass; 187a175e481SJames Wright CeedInt num_comp_q, q_data_size; 188a175e481SJames Wright PetscFunctionBeginUser; 189a175e481SJames Wright 190a175e481SJames Wright // CEED Restriction 191a175e481SJames Wright CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_stats, &num_comp_q); 192a175e481SJames Wright CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_qd, &q_data_size); 193a175e481SJames Wright 194a175e481SJames Wright // Create Mass CeedOperator 195a175e481SJames Wright PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass)); 196a175e481SJames Wright CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass); 197a175e481SJames Wright CeedOperatorSetField(op_mass, "q", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE); 198a175e481SJames Wright CeedOperatorSetField(op_mass, "qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, ceed_data->spanstats.q_data); 199a175e481SJames Wright CeedOperatorSetField(op_mass, "v", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE); 200a175e481SJames Wright 201a175e481SJames Wright // Setup KSP for L^2 projection 202a175e481SJames Wright { 203a175e481SJames Wright MatopApplyContext M_ctx; 204a175e481SJames Wright PetscInt l_size, g_size; 205a175e481SJames Wright Mat mat_mass; 206a175e481SJames Wright VecType vec_type; 207a175e481SJames Wright KSP ksp; 208a175e481SJames Wright Vec ones, M_inv; 209a175e481SJames Wright CeedVector x_ceed, y_ceed; 210a175e481SJames Wright 211a175e481SJames Wright PetscCall(DMCreateGlobalVector(user->spanstats.dm, &M_inv)); 212a175e481SJames Wright PetscCall(VecGetLocalSize(M_inv, &l_size)); 213a175e481SJames Wright PetscCall(VecGetSize(M_inv, &g_size)); 214a175e481SJames Wright PetscCall(VecGetType(M_inv, &vec_type)); 215a175e481SJames Wright 216a175e481SJames Wright PetscCall(PetscMalloc1(1, &M_ctx)); 217a175e481SJames Wright PetscCall(MatCreateShell(PETSC_COMM_WORLD, l_size, l_size, g_size, g_size, M_ctx, &mat_mass)); 218a175e481SJames Wright PetscCall(MatShellSetOperation(mat_mass, MATOP_MULT, (void (*)(void))MatMult_Ceed)); 219a175e481SJames Wright PetscCall(MatShellSetOperation(mat_mass, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed)); 220a175e481SJames Wright PetscCall(MatShellSetVecType(mat_mass, vec_type)); 221a175e481SJames Wright 222a175e481SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &x_ceed, NULL); 223a175e481SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &y_ceed, NULL); 224a175e481SJames Wright 225a175e481SJames Wright PetscCall(SetupMatopApplyCtx(PETSC_COMM_WORLD, user->spanstats.dm, user->ceed, op_mass, x_ceed, y_ceed, NULL, M_ctx)); 226a175e481SJames Wright user->spanstats.M_ctx = M_ctx; 227a175e481SJames Wright 228a175e481SJames Wright // Create lumped mass matrix inverse 229a175e481SJames Wright PetscCall(DMGetGlobalVector(user->spanstats.dm, &ones)); 230a175e481SJames Wright PetscCall(VecZeroEntries(M_inv)); 231a175e481SJames Wright PetscCall(VecSet(ones, 1)); 232a175e481SJames Wright PetscCall(MatMult(mat_mass, ones, M_inv)); 233a175e481SJames Wright PetscCall(VecReciprocal(M_inv)); 234a175e481SJames Wright user->spanstats.M_inv = M_inv; 235a175e481SJames Wright PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &ones)); 236a175e481SJames Wright 237a175e481SJames Wright PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp)); 238a175e481SJames Wright PetscCall(KSPSetOptionsPrefix(ksp, "spanstats_")); 239a175e481SJames Wright { 240a175e481SJames Wright PC pc; 241a175e481SJames Wright PetscCall(KSPGetPC(ksp, &pc)); 242a175e481SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 243a175e481SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL)); 244a175e481SJames Wright PetscCall(KSPSetType(ksp, KSPCG)); 245a175e481SJames Wright PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL)); 246a175e481SJames Wright PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 247a175e481SJames Wright } 248a175e481SJames Wright PetscCall(KSPSetOperators(ksp, mat_mass, mat_mass)); 249a175e481SJames Wright PetscCall(KSPSetFromOptions(ksp)); 250a175e481SJames Wright user->spanstats.ksp = ksp; 251a175e481SJames Wright } 252a175e481SJames Wright 253a175e481SJames Wright // Cleanup 254a175e481SJames Wright CeedQFunctionDestroy(&qf_mass); 255a175e481SJames Wright PetscFunctionReturn(0); 256a175e481SJames Wright } 257a175e481SJames Wright 258a175e481SJames Wright // Create CeedOperators and KSP for the statistics collection and processing 259a175e481SJames Wright PetscErrorCode CreateStatisticsOperators(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) { 260a175e481SJames Wright CeedInt num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q; 261a175e481SJames Wright CeedOperator op_setup_sur; 262a175e481SJames Wright PetscFunctionBeginUser; 263a175e481SJames Wright CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q); 264a175e481SJames Wright 265a175e481SJames Wright // Create Operator for statistics collection 266a175e481SJames Wright switch (user->phys->state_var) { 267a175e481SJames Wright case STATEVAR_PRIMITIVE: 268a175e481SJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &ceed_data->spanstats.qf_stats_collect); 269a175e481SJames Wright break; 270a175e481SJames Wright case STATEVAR_CONSERVATIVE: 271a175e481SJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &ceed_data->spanstats.qf_stats_collect); 272a175e481SJames Wright break; 273a175e481SJames Wright } 274a175e481SJames Wright 275a175e481SJames Wright if (user->app_ctx->stats_test) { 276a175e481SJames Wright CeedQFunctionDestroy(&ceed_data->spanstats.qf_stats_collect); 277a175e481SJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionTest, ChildStatsCollectionTest_loc, &ceed_data->spanstats.qf_stats_collect); 278a175e481SJames Wright } 279a175e481SJames Wright 280a175e481SJames Wright CeedQFunctionSetContext(ceed_data->spanstats.qf_stats_collect, problem->apply_vol_ifunction.qfunction_context); 281a175e481SJames Wright CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP); 282a175e481SJames Wright CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE); 283a175e481SJames Wright CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP); 284a175e481SJames Wright CeedQFunctionAddOutput(ceed_data->spanstats.qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE); 285a175e481SJames Wright 286a175e481SJames Wright CeedOperatorCreate(ceed, ceed_data->spanstats.qf_stats_collect, NULL, NULL, &user->spanstats.op_stats_collect); 287a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 288a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 289a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord); 290a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_collect, "v", ceed_data->spanstats.elem_restr_child_colloc, CEED_BASIS_COLLOCATED, 291a175e481SJames Wright CEED_VECTOR_ACTIVE); 292a175e481SJames Wright 293*823a1283SJames Wright // Create Operator for L^2 projection of statistics 294a175e481SJames Wright // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function. 295a175e481SJames Wright // Therefore, an Identity QF is sufficient 296a175e481SJames Wright CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &ceed_data->spanstats.qf_stats_proj); 297a175e481SJames Wright 298a175e481SJames Wright CeedOperatorCreate(ceed, ceed_data->spanstats.qf_stats_proj, NULL, NULL, &user->spanstats.op_stats_proj); 299a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_proj, "input", ceed_data->spanstats.elem_restr_parent_colloc, CEED_BASIS_COLLOCATED, 300a175e481SJames Wright CEED_VECTOR_ACTIVE); 301a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_proj, "output", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, 302a175e481SJames Wright CEED_VECTOR_ACTIVE); 303a175e481SJames Wright 304a175e481SJames Wright // Get q_data for lumped mass matrix formation 305a175e481SJames Wright CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur); 306a175e481SJames Wright CeedOperatorSetField(op_setup_sur, "dx", ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, CEED_VECTOR_ACTIVE); 307a175e481SJames Wright CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->spanstats.basis_x, CEED_VECTOR_NONE); 308a175e481SJames Wright CeedOperatorSetField(op_setup_sur, "surface qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 309a175e481SJames Wright CeedOperatorApply(op_setup_sur, ceed_data->spanstats.x_coord, ceed_data->spanstats.q_data, CEED_REQUEST_IMMEDIATE); 310a175e481SJames Wright 311a175e481SJames Wright CeedOperatorDestroy(&op_setup_sur); 312a175e481SJames Wright PetscFunctionReturn(0); 313a175e481SJames Wright } 314a175e481SJames Wright 315*823a1283SJames Wright PetscErrorCode SetupErrorTesting(Ceed ceed, User user, CeedData ceed_data) { 316*823a1283SJames Wright CeedInt num_comp_stats = user->spanstats.num_comp_stats, num_comp_x; 317*823a1283SJames Wright CeedQFunction qf_error; 318*823a1283SJames Wright CeedOperator op_error; 319*823a1283SJames Wright CeedInt q_data_size; 320*823a1283SJames Wright CeedVector x_ceed, y_ceed; 321*823a1283SJames Wright PetscFunctionBeginUser; 322*823a1283SJames Wright 323*823a1283SJames Wright CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_qd, &q_data_size); 324*823a1283SJames Wright CeedBasisGetNumComponents(ceed_data->spanstats.basis_x, &num_comp_x); 325*823a1283SJames Wright 326*823a1283SJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionTest_Error, ChildStatsCollectionTest_Error_loc, &qf_error); 327*823a1283SJames Wright CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP); 328*823a1283SJames Wright CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE); 329*823a1283SJames Wright CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP); 330*823a1283SJames Wright CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP); 331*823a1283SJames Wright 332*823a1283SJames Wright CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error); 333*823a1283SJames Wright CeedOperatorSetField(op_error, "q", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE); 334*823a1283SJames Wright CeedOperatorSetField(op_error, "qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, ceed_data->spanstats.q_data); 335*823a1283SJames Wright CeedOperatorSetField(op_error, "x", ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, ceed_data->spanstats.x_coord); 336*823a1283SJames Wright CeedOperatorSetField(op_error, "v", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE); 337*823a1283SJames Wright 338*823a1283SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &x_ceed, NULL); 339*823a1283SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &y_ceed, NULL); 340*823a1283SJames Wright 341*823a1283SJames Wright PetscCall(PetscCalloc1(1, &user->spanstats.test_error_ctx)); 342*823a1283SJames Wright PetscCall(SetupMatopApplyCtx(PETSC_COMM_WORLD, user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, user->spanstats.test_error_ctx)); 343*823a1283SJames Wright 344*823a1283SJames Wright PetscFunctionReturn(0); 345*823a1283SJames Wright } 346*823a1283SJames Wright 347a175e481SJames Wright // Setup for statistics collection 34819706a06SJames Wright PetscErrorCode SetupStatsCollection(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) { 34919706a06SJames Wright DM dm = user->spanstats.dm; 35019706a06SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm); 35119706a06SJames Wright CeedInt dim, P, Q, num_comp_x; 35219706a06SJames Wright Vec X_loc; 35319706a06SJames Wright PetscMemType X_loc_memtype; 35419706a06SJames Wright const PetscScalar *X_loc_array; 35519706a06SJames Wright PetscFunctionBeginUser; 35619706a06SJames Wright 35719706a06SJames Wright PetscCall(DMGetDimension(dm, &dim)); 358a175e481SJames Wright CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &Q); 359a175e481SJames Wright CeedBasisGetNumNodes1D(ceed_data->basis_q, &P); 36019706a06SJames Wright 3611737222fSJames Wright PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, problem->q_data_size_sur, &ceed_data->spanstats.elem_restr_parent_stats, 36219706a06SJames Wright &ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.elem_restr_parent_qd)); 36319706a06SJames Wright CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_x, &num_comp_x); 36419706a06SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.x_coord, NULL); 365a175e481SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &user->spanstats.rhs_ceed, NULL); 3661737222fSJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_qd, &ceed_data->spanstats.q_data, NULL); 36719706a06SJames Wright 368a175e481SJames Wright CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS, &ceed_data->spanstats.basis_x); 36919706a06SJames Wright CeedBasisCreateTensorH1Lagrange(ceed, dim, user->spanstats.num_comp_stats, P, Q, CEED_GAUSS, &ceed_data->spanstats.basis_stats); 37019706a06SJames Wright 3711737222fSJames Wright PetscCall(CreateElemRestrColloc(ceed, user->spanstats.num_comp_stats, ceed_data->spanstats.basis_stats, 3721737222fSJames Wright ceed_data->spanstats.elem_restr_parent_stats, &ceed_data->spanstats.elem_restr_parent_colloc, 3731737222fSJames Wright &user->spanstats.parent_stats, NULL)); 3741737222fSJames Wright PetscCall(CreateElemRestrColloc(ceed, user->spanstats.num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q, 3751737222fSJames Wright &ceed_data->spanstats.elem_restr_child_colloc, &user->spanstats.child_stats, NULL)); 376a175e481SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_child_colloc, &user->spanstats.child_inst_stats, NULL); 377a175e481SJames Wright CeedVectorSetValue(user->spanstats.child_stats, 0); 3781737222fSJames Wright 37919706a06SJames Wright // -- Copy DM coordinates into CeedVector 38019706a06SJames Wright { 38119706a06SJames Wright DM cdm; 38219706a06SJames Wright PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 38319706a06SJames Wright if (cdm) { 38419706a06SJames Wright PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 38519706a06SJames Wright } else { 38619706a06SJames Wright PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 38719706a06SJames Wright } 38819706a06SJames Wright } 38919706a06SJames Wright PetscCall(VecScale(X_loc, problem->dm_scale)); 39019706a06SJames Wright PetscCall(VecGetArrayReadAndMemType(X_loc, &X_loc_array, &X_loc_memtype)); 39119706a06SJames Wright CeedVectorSetArray(ceed_data->spanstats.x_coord, MemTypeP2C(X_loc_memtype), CEED_COPY_VALUES, (PetscScalar *)X_loc_array); 39219706a06SJames Wright PetscCall(VecRestoreArrayRead(X_loc, &X_loc_array)); 39319706a06SJames Wright 39419706a06SJames Wright // Create SF for communicating child data back their respective parents 39519706a06SJames Wright PetscCall(PetscSFCreate(comm, &user->spanstats.sf)); 39619706a06SJames Wright PetscCall(CreateStatsSF(ceed, ceed_data, user->dm, user->spanstats.dm, user->spanstats.sf)); 39751ee423eSJames Wright 398a175e481SJames Wright // Create CeedOperators for statistics collection 399a175e481SJames Wright PetscCall(CreateStatisticsOperators(ceed, user, ceed_data, problem)); 400a175e481SJames Wright 401a175e481SJames Wright // Setup KSP and Mat for L^2 projection of statistics 402a175e481SJames Wright PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data)); 403a175e481SJames Wright 404*823a1283SJames Wright if (user->app_ctx->stats_test) { 405*823a1283SJames Wright PetscCall(SetupErrorTesting(ceed, user, ceed_data)); 406*823a1283SJames Wright } 407*823a1283SJames Wright 408a175e481SJames Wright PetscFunctionReturn(0); 409a175e481SJames Wright } 410a175e481SJames Wright 411a175e481SJames Wright // Collect statistics based on the solution Q 412a175e481SJames Wright PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) { 413a175e481SJames Wright PetscMemType q_mem_type; 414a175e481SJames Wright const PetscScalar *q_arr; 415a175e481SJames Wright Vec Q_loc; 416a175e481SJames Wright PetscFunctionBeginUser; 417a175e481SJames Wright 418a175e481SJames Wright PetscCall(DMGetLocalVector(user->dm, &Q_loc)); 419a175e481SJames Wright PetscCall(VecZeroEntries(Q_loc)); 420a175e481SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 421a175e481SJames Wright 422a175e481SJames Wright PetscCall(VecGetArrayReadAndMemType(Q_loc, &q_arr, &q_mem_type)); 423a175e481SJames Wright CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q_arr); 424a175e481SJames Wright 425a175e481SJames Wright CeedOperatorApply(user->spanstats.op_stats_collect, user->q_ceed, user->spanstats.child_inst_stats, CEED_REQUEST_IMMEDIATE); 426a175e481SJames Wright 427a175e481SJames Wright CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL); 428a175e481SJames Wright PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q_arr)); 429a175e481SJames Wright PetscCall(DMRestoreLocalVector(user->dm, &Q_loc)); 430a175e481SJames Wright 431a175e481SJames Wright // Record averaging using left rectangle rule 432a175e481SJames Wright PetscScalar delta_t = solution_time - user->spanstats.prev_time; 433a175e481SJames Wright PetscScalar prev_timeinterval = user->spanstats.prev_time - user->app_ctx->cont_time; 434a175e481SJames Wright CeedVectorScale(user->spanstats.child_stats, prev_timeinterval / (prev_timeinterval + delta_t)); 435a175e481SJames Wright CeedVectorAXPY(user->spanstats.child_stats, delta_t / (prev_timeinterval + delta_t), user->spanstats.child_inst_stats); 436a175e481SJames Wright user->spanstats.prev_time = solution_time; 437a175e481SJames Wright 438a175e481SJames Wright PetscFunctionReturn(0); 439a175e481SJames Wright } 440a175e481SJames Wright 441a175e481SJames Wright // Process the child statistics into parent statistics and project them onto stats 442a175e481SJames Wright PetscErrorCode ProcessStatistics(User user, Vec *stats) { 443a175e481SJames Wright Span_Stats user_stats = user->spanstats; 444a175e481SJames Wright const PetscScalar *child_stats; 445a175e481SJames Wright PetscScalar *parent_stats; 446a175e481SJames Wright MPI_Datatype unit; 447a175e481SJames Wright Vec rhs_loc, rhs; 448a175e481SJames Wright PetscMemType rhs_mem_type; 449a175e481SJames Wright CeedScalar *rhs_arr; 450a175e481SJames Wright CeedMemType ceed_mem_type; 451a175e481SJames Wright PetscFunctionBeginUser; 452a175e481SJames Wright 453a175e481SJames Wright CeedGetPreferredMemType(user->ceed, &ceed_mem_type); 454a175e481SJames Wright CeedVectorSetValue(user_stats.parent_stats, 0); 455a175e481SJames Wright 456a175e481SJames Wright CeedVectorGetArrayRead(user_stats.child_stats, ceed_mem_type, &child_stats); 457a175e481SJames Wright CeedVectorGetArray(user_stats.parent_stats, ceed_mem_type, &parent_stats); 458a175e481SJames Wright 459a175e481SJames Wright if (user_stats.num_comp_stats == 1) unit = MPIU_REAL; 460a175e481SJames Wright else { 461a175e481SJames Wright PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit)); 462a175e481SJames Wright PetscCallMPI(MPI_Type_commit(&unit)); 463a175e481SJames Wright } 464a175e481SJames Wright 465a175e481SJames Wright PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM)); 466a175e481SJames Wright PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM)); 467a175e481SJames Wright 468a175e481SJames Wright CeedVectorRestoreArrayRead(user_stats.child_stats, &child_stats); 469a175e481SJames Wright CeedVectorRestoreArray(user_stats.parent_stats, &parent_stats); 470a175e481SJames Wright PetscCallMPI(MPI_Type_free(&unit)); 471a175e481SJames Wright 472a175e481SJames Wright CeedVectorScale(user_stats.parent_stats, 1 / user_stats.span_width); 473a175e481SJames Wright 474a175e481SJames Wright // L^2 projection with the parent_data 475a175e481SJames Wright PetscCall(DMGetGlobalVector(user_stats.dm, &rhs)); 476a175e481SJames Wright PetscCall(DMGetLocalVector(user_stats.dm, &rhs_loc)); 477a175e481SJames Wright PetscCall(VecZeroEntries(rhs)); 478a175e481SJames Wright PetscCall(VecZeroEntries(rhs_loc)); 479a175e481SJames Wright PetscCall(VecGetArrayWriteAndMemType(rhs_loc, &rhs_arr, &rhs_mem_type)); 480a175e481SJames Wright CeedVectorSetArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), CEED_USE_POINTER, (PetscScalar *)rhs_arr); 481a175e481SJames Wright 482a175e481SJames Wright CeedOperatorApply(user_stats.op_stats_proj, user_stats.parent_stats, user_stats.rhs_ceed, CEED_REQUEST_IMMEDIATE); 483a175e481SJames Wright 484a175e481SJames Wright CeedVectorTakeArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), &rhs_arr); 485a175e481SJames Wright PetscCall(VecRestoreArrayAndMemType(rhs_loc, &rhs_arr)); 486a175e481SJames Wright PetscCall(DMLocalToGlobal(user_stats.dm, rhs_loc, ADD_VALUES, rhs)); 487a175e481SJames Wright 488a175e481SJames Wright PetscCall(VecDuplicate(rhs, stats)); 489a175e481SJames Wright PetscCall(VecPointwiseMult(*stats, rhs, user_stats.M_inv)); 490a175e481SJames Wright 491a175e481SJames Wright PetscCall(KSPSolve(user_stats.ksp, rhs, *stats)); 492a175e481SJames Wright 493a175e481SJames Wright PetscFunctionReturn(0); 494a175e481SJames Wright } 495a175e481SJames Wright 496a175e481SJames Wright // TSMonitor for the statistics collection and processing 497a175e481SJames Wright PetscErrorCode TSMonitor_Statistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) { 498a175e481SJames Wright User user = (User)ctx; 499a175e481SJames Wright Vec stats; 500a175e481SJames Wright PetscFunctionBeginUser; 501a175e481SJames Wright 502a175e481SJames Wright // Do not collect or process on the first step of the run (ie. on the initial condition) 503a175e481SJames Wright if (steps == user->app_ctx->cont_steps) PetscFunctionReturn(0); 504a175e481SJames Wright 505a175e481SJames Wright if (steps % user->app_ctx->stats_collect_interval == 0) { 506a175e481SJames Wright PetscCall(CollectStatistics(user, solution_time, Q)); 507a175e481SJames Wright } 508a175e481SJames Wright 509a175e481SJames Wright if (steps % user->app_ctx->stats_write_interval == 0 && user->app_ctx->stats_write_interval != -1) { 510a175e481SJames Wright PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats)); 511a175e481SJames Wright PetscCall(ProcessStatistics(user, &stats)); 512a175e481SJames Wright PetscCall(VecViewFromOptions(stats, NULL, "-stats_write_view")); 513a175e481SJames Wright PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats)); 514a175e481SJames Wright } 515a175e481SJames Wright PetscFunctionReturn(0); 516a175e481SJames Wright } 517a175e481SJames Wright 518a175e481SJames Wright // Function to be called at the end of a simulation 519a175e481SJames Wright PetscErrorCode StatsCollectFinalCall(User user, PetscReal solution_time, Vec Q) { 520a175e481SJames Wright Vec stats; 521a175e481SJames Wright PetscFunctionBeginUser; 522a175e481SJames Wright 523a175e481SJames Wright PetscCall(CollectStatistics(user, solution_time, Q)); 524a175e481SJames Wright 525a175e481SJames Wright PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats)); 526a175e481SJames Wright PetscCall(ProcessStatistics(user, &stats)); 527a175e481SJames Wright PetscCall(VecViewFromOptions(stats, NULL, "-stats_write_view")); 528a175e481SJames Wright 529*823a1283SJames Wright if (user->app_ctx->stats_test) { 530*823a1283SJames Wright Vec error; 531*823a1283SJames Wright PetscCall(VecDuplicate(stats, &error)); 532*823a1283SJames Wright PetscCall(ApplyLocal_Ceed(stats, error, user->spanstats.test_error_ctx)); 533*823a1283SJames Wright PetscScalar error_sq = 0; 534*823a1283SJames Wright PetscCall(VecSum(error, &error_sq)); 535*823a1283SJames Wright PetscScalar l2_error = sqrt(error_sq); 536*823a1283SJames Wright PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error)); 537*823a1283SJames Wright } 538a175e481SJames Wright PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats)); 539a175e481SJames Wright 54051ee423eSJames Wright PetscFunctionReturn(0); 54151ee423eSJames Wright } 542