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" 228ed52730SJames Wright #include "petscviewer.h" 2351ee423eSJames Wright 2451ee423eSJames Wright PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree, SimpleBC bc) { 253a4208e6SJames Wright user->spanstats.num_comp_stats = TURB_NUM_COMPONENTS; 26a175e481SJames Wright PetscReal domain_min[3], domain_max[3]; 2778837792SJames Wright PetscFE fe; 2878837792SJames Wright PetscSection section; 294eed8630SJames Wright PetscLogStage stage_stats_setup; 301f595ac1SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)user->dm); 3151ee423eSJames Wright PetscFunctionBeginUser; 3251ee423eSJames Wright 334eed8630SJames Wright PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup)); 344eed8630SJames Wright if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup)); 354eed8630SJames Wright PetscCall(PetscLogStagePush(stage_stats_setup)); 364eed8630SJames Wright 37a175e481SJames Wright // Get spanwise length 38a175e481SJames Wright PetscCall(DMGetBoundingBox(user->dm, domain_min, domain_max)); 39a175e481SJames Wright user->spanstats.span_width = domain_max[2] - domain_min[1]; 40a175e481SJames Wright 41f967ad79SJames Wright { // Get DM from surface 421f595ac1SJames Wright DM parent_distributed_dm; 43c9198418SJames Wright PetscSF isoperiodicface; 4451ee423eSJames Wright DMLabel label; 451f595ac1SJames Wright PetscMPIInt size; 46c9198418SJames Wright 47c9198418SJames Wright PetscCall(DMPlexGetIsoperiodicFaceSF(user->dm, &isoperiodicface)); 48c9198418SJames Wright 49c9198418SJames Wright if (isoperiodicface) { 50c9198418SJames Wright PetscSF inv_isoperiodicface; 51c9198418SJames Wright PetscInt nleaves; 52c9198418SJames Wright const PetscInt *ilocal; 53c9198418SJames Wright 54c9198418SJames Wright PetscCall(PetscSFCreateInverseSF(isoperiodicface, &inv_isoperiodicface)); 55c9198418SJames Wright PetscCall(PetscSFGetGraph(inv_isoperiodicface, NULL, &nleaves, &ilocal, NULL)); 56c9198418SJames Wright PetscCall(DMCreateLabel(user->dm, "Periodic Face")); 57c9198418SJames Wright PetscCall(DMGetLabel(user->dm, "Periodic Face", &label)); 58c9198418SJames Wright for (PetscInt i = 0; i < nleaves; i++) { 59c9198418SJames Wright PetscCall(DMLabelSetValue(label, ilocal[i], 1)); 60c9198418SJames Wright } 61c9198418SJames Wright } else { 6251ee423eSJames Wright PetscCall(DMGetLabel(user->dm, "Face Sets", &label)); 63c9198418SJames Wright } 64c9198418SJames Wright 6551ee423eSJames Wright PetscCall(DMPlexLabelComplete(user->dm, label)); 6651ee423eSJames Wright PetscCall(DMPlexFilter(user->dm, label, 1, &user->spanstats.dm)); 6751ee423eSJames Wright PetscCall(DMProjectCoordinates(user->spanstats.dm, NULL)); // Ensure that a coordinate FE exists 681f595ac1SJames Wright 691f595ac1SJames Wright PetscCall(DMPlexDistribute(user->spanstats.dm, 0, NULL, &parent_distributed_dm)); 701f595ac1SJames Wright PetscCallMPI(MPI_Comm_size(comm, &size)); 711f595ac1SJames Wright if (parent_distributed_dm) { 721f595ac1SJames Wright PetscCall(DMDestroy(&user->spanstats.dm)); 731f595ac1SJames Wright user->spanstats.dm = parent_distributed_dm; 741f595ac1SJames Wright } else if (size > 1) { 751f595ac1SJames Wright PetscCall(PetscPrintf(comm, "WARNING: Turbulent spanwise statistics: parent DM could not be distributed accross %d ranks.\n", size)); 761f595ac1SJames Wright } 7751ee423eSJames Wright } 7851ee423eSJames Wright 7951ee423eSJames Wright PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats")); 80b7d66439SJames Wright PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "turbulence_spanstats_")); 8151ee423eSJames Wright PetscCall(DMSetFromOptions(user->spanstats.dm)); 82a175e481SJames Wright PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view")); // -spanstats_dm_view 83c9198418SJames Wright 8478837792SJames Wright // Create FE space for parent DM 8551ee423eSJames Wright PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim - 1, user->spanstats.num_comp_stats, PETSC_FALSE, degree, PETSC_DECIDE, &fe)); 8651ee423eSJames Wright PetscCall(PetscObjectSetName((PetscObject)fe, "stats")); 8751ee423eSJames Wright PetscCall(DMAddField(user->spanstats.dm, NULL, (PetscObject)fe)); 8851ee423eSJames Wright PetscCall(DMCreateDS(user->spanstats.dm)); 8951ee423eSJames Wright PetscCall(DMPlexSetClosurePermutationTensor(user->spanstats.dm, PETSC_DETERMINE, NULL)); 9051ee423eSJames Wright 9178837792SJames Wright // Create Section for data 9251ee423eSJames Wright PetscCall(DMGetLocalSection(user->spanstats.dm, §ion)); 9351ee423eSJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 943a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY, "MeanDensity")); 953a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE, "MeanPressure")); 963a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_SQUARED, "MeanPressureSquared")); 973a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_X, "MeanPressureVelocityX")); 983a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Y, "MeanPressureVelocityY")); 993a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Z, "MeanPressureVelocityZ")); 1003a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE, "MeanDensityTemperature")); 1013a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_X, "MeanDensityTemperatureFluxX")); 1023a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Y, "MeanDensityTemperatureFluxY")); 1033a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Z, "MeanDensityTemperatureFluxZ")); 1043a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_X, "MeanMomentumX")); 1053a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Y, "MeanMomentumY")); 1063a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Z, "MeanMomentumZ")); 1073a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XX, "MeanMomentumFluxXX")); 1083a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YY, "MeanMomentumFluxYY")); 1093a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_ZZ, "MeanMomentumFluxZZ")); 1103a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YZ, "MeanMomentumFluxYZ")); 1113a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XZ, "MeanMomentumFluxXZ")); 1123a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XY, "MeanMomentumFluxXY")); 1133a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_X, "MeanVelocityX")); 1143a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Y, "MeanVelocityY")); 1153a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Z, "MeanVelocityZ")); 11619706a06SJames Wright 11778837792SJames Wright // Cleanup 11878837792SJames Wright PetscCall(PetscFEDestroy(&fe)); 11978837792SJames Wright 1204eed8630SJames Wright PetscCall(PetscLogStagePop()); 12119706a06SJames Wright PetscFunctionReturn(0); 12219706a06SJames Wright } 12319706a06SJames Wright 1241737222fSJames Wright // Create CeedElemRestriction for collocated data based on associated CeedBasis and CeedElemRestriction 1251737222fSJames Wright // Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction 1261737222fSJames Wright PetscErrorCode CreateElemRestrColloc(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base, 1271737222fSJames Wright CeedElemRestriction *elem_restr_collocated, CeedVector *l_vec, CeedVector *e_vec) { 1281737222fSJames Wright CeedInt num_elem_qpts, loc_num_elem; 1291737222fSJames Wright PetscFunctionBeginUser; 1301737222fSJames Wright 1311737222fSJames Wright CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts); 1321737222fSJames Wright CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem); 1331737222fSJames Wright 1341737222fSJames Wright const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp}; 1351737222fSJames Wright CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides, 1361737222fSJames Wright elem_restr_collocated); 1371737222fSJames Wright CeedElemRestrictionCreateVector(*elem_restr_collocated, l_vec, e_vec); 1381737222fSJames Wright PetscFunctionReturn(0); 1391737222fSJames Wright } 1401737222fSJames Wright 141a175e481SJames Wright // Get coordinates of quadrature points 14219706a06SJames Wright PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, CeedVector *qx_coords, 14319706a06SJames Wright PetscInt *total_nqpnts) { 14419706a06SJames Wright CeedQFunction qf_quad_coords; 14519706a06SJames Wright CeedOperator op_quad_coords; 14619706a06SJames Wright PetscInt num_comp_x, loc_num_elem, num_elem_qpts; 14719706a06SJames Wright CeedElemRestriction elem_restr_qx; 14819706a06SJames Wright PetscFunctionBeginUser; 14919706a06SJames Wright 15019706a06SJames Wright // Create Element Restriction and CeedVector for quadrature coordinates 15119706a06SJames Wright CeedBasisGetNumQuadraturePoints(basis_x, &num_elem_qpts); 15219706a06SJames Wright CeedElemRestrictionGetNumElements(elem_restr_x, &loc_num_elem); 15319706a06SJames Wright CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x); 15419706a06SJames Wright *total_nqpnts = num_elem_qpts * loc_num_elem; 1551737222fSJames Wright PetscCall(CreateElemRestrColloc(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx, qx_coords, NULL)); 15619706a06SJames Wright 15719706a06SJames Wright // Create QFunction 15819706a06SJames Wright CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords); 15919706a06SJames Wright 16019706a06SJames Wright // Create Operator 16119706a06SJames Wright CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords); 16219706a06SJames Wright CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE); 16319706a06SJames Wright CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 16419706a06SJames Wright 16519706a06SJames Wright CeedOperatorApply(op_quad_coords, x_coords, *qx_coords, CEED_REQUEST_IMMEDIATE); 16619706a06SJames Wright 16719706a06SJames Wright CeedQFunctionDestroy(&qf_quad_coords); 16819706a06SJames Wright CeedOperatorDestroy(&op_quad_coords); 16919706a06SJames Wright PetscFunctionReturn(0); 17019706a06SJames Wright } 17119706a06SJames Wright 172a175e481SJames Wright // Create PetscSF for child-to-parent communication 17319706a06SJames Wright PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, DM parentdm, DM childdm, PetscSF statssf) { 17419706a06SJames Wright PetscInt child_num_qpnts, parent_num_qpnts, num_comp_x; 17519706a06SJames Wright CeedVector child_qx_coords, parent_qx_coords; 17619706a06SJames Wright PetscReal *child_coords, *parent_coords; 17719706a06SJames Wright PetscFunctionBeginUser; 17819706a06SJames Wright 17919706a06SJames Wright // Assume that child and parent have the same number of components 18019706a06SJames Wright CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x); 18119706a06SJames Wright const PetscInt num_comp_sf = num_comp_x - 1; // Number of coord components used in the creation of the SF 18219706a06SJames Wright 18319706a06SJames Wright // Get quad_coords for child DM 184a175e481SJames Wright PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &child_qx_coords, &child_num_qpnts)); 18519706a06SJames Wright 18619706a06SJames Wright // Get quad_coords for parent DM 18719706a06SJames Wright PetscCall(GetQuadratureCoords(ceed, parentdm, ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, ceed_data->spanstats.x_coord, 18819706a06SJames Wright &parent_qx_coords, &parent_num_qpnts)); 18919706a06SJames Wright 19019706a06SJames Wright // Remove z component of coordinates for matching 19119706a06SJames Wright { 19219706a06SJames Wright const PetscReal *child_quad_coords, *parent_quad_coords; 19319706a06SJames Wright 19419706a06SJames Wright CeedVectorGetArrayRead(child_qx_coords, CEED_MEM_HOST, &child_quad_coords); 19519706a06SJames Wright CeedVectorGetArrayRead(parent_qx_coords, CEED_MEM_HOST, &parent_quad_coords); 19619706a06SJames Wright 19719706a06SJames Wright PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords)); 19819706a06SJames Wright for (int i = 0; i < child_num_qpnts; i++) { 19919706a06SJames Wright child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x]; 20019706a06SJames Wright child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x]; 20119706a06SJames Wright } 20219706a06SJames Wright for (int i = 0; i < parent_num_qpnts; i++) { 20319706a06SJames Wright parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x]; 20419706a06SJames Wright parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x]; 20519706a06SJames Wright } 20619706a06SJames Wright CeedVectorRestoreArrayRead(child_qx_coords, &child_quad_coords); 20719706a06SJames Wright CeedVectorRestoreArrayRead(parent_qx_coords, &parent_quad_coords); 20819706a06SJames Wright } 20919706a06SJames Wright 21019706a06SJames Wright PetscCall(PetscSFSetGraphFromCoordinates(statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords)); 21119706a06SJames Wright 212a175e481SJames Wright PetscCall(PetscSFViewFromOptions(statssf, NULL, "-spanstats_sf_view")); 213a175e481SJames Wright 21419706a06SJames Wright PetscCall(PetscFree2(child_coords, parent_coords)); 21519706a06SJames Wright CeedVectorDestroy(&child_qx_coords); 21619706a06SJames Wright CeedVectorDestroy(&parent_qx_coords); 21719706a06SJames Wright PetscFunctionReturn(0); 21819706a06SJames Wright } 21919706a06SJames Wright 220a175e481SJames Wright // Compute mass matrix for statistics projection 221a175e481SJames Wright PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data) { 222a175e481SJames Wright CeedQFunction qf_mass; 223a175e481SJames Wright CeedOperator op_mass; 224a175e481SJames Wright CeedInt num_comp_q, q_data_size; 225f967ad79SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)user->spanstats.dm); 226a175e481SJames Wright PetscFunctionBeginUser; 227a175e481SJames Wright 228a175e481SJames Wright // CEED Restriction 229a175e481SJames Wright CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_stats, &num_comp_q); 230a175e481SJames Wright CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_qd, &q_data_size); 231a175e481SJames Wright 232a175e481SJames Wright // Create Mass CeedOperator 233a175e481SJames Wright PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass)); 234a175e481SJames Wright CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass); 235a175e481SJames Wright CeedOperatorSetField(op_mass, "q", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE); 236a175e481SJames Wright CeedOperatorSetField(op_mass, "qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, ceed_data->spanstats.q_data); 237a175e481SJames Wright CeedOperatorSetField(op_mass, "v", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE); 238a175e481SJames Wright 239f967ad79SJames Wright { // Setup KSP for L^2 projection 240a175e481SJames Wright MatopApplyContext M_ctx; 241a175e481SJames Wright PetscInt l_size, g_size; 242a175e481SJames Wright Mat mat_mass; 243a175e481SJames Wright VecType vec_type; 244a175e481SJames Wright KSP ksp; 245a175e481SJames Wright Vec ones, M_inv; 246a175e481SJames Wright CeedVector x_ceed, y_ceed; 247a175e481SJames Wright 248a175e481SJames Wright PetscCall(DMCreateGlobalVector(user->spanstats.dm, &M_inv)); 249a175e481SJames Wright PetscCall(VecGetLocalSize(M_inv, &l_size)); 250a175e481SJames Wright PetscCall(VecGetSize(M_inv, &g_size)); 251a175e481SJames Wright PetscCall(VecGetType(M_inv, &vec_type)); 252a175e481SJames Wright 253*6665b873SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &x_ceed, NULL); 254*6665b873SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &y_ceed, NULL); 255*6665b873SJames Wright PetscCall(MatopApplyContextCreate(user->spanstats.dm, user->ceed, op_mass, x_ceed, y_ceed, NULL, &M_ctx)); 256*6665b873SJames Wright CeedVectorDestroy(&x_ceed); 257*6665b873SJames Wright CeedVectorDestroy(&y_ceed); 258*6665b873SJames Wright 259f967ad79SJames Wright PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, M_ctx, &mat_mass)); 260*6665b873SJames Wright PetscCall(MatShellSetContextDestroy(mat_mass, (PetscErrorCode(*)(void *))MatopApplyContextDestroy)); 261a175e481SJames Wright PetscCall(MatShellSetOperation(mat_mass, MATOP_MULT, (void (*)(void))MatMult_Ceed)); 262a175e481SJames Wright PetscCall(MatShellSetOperation(mat_mass, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed)); 263a175e481SJames Wright PetscCall(MatShellSetVecType(mat_mass, vec_type)); 264a175e481SJames Wright 265a175e481SJames Wright // Create lumped mass matrix inverse 266a175e481SJames Wright PetscCall(DMGetGlobalVector(user->spanstats.dm, &ones)); 267a175e481SJames Wright PetscCall(VecZeroEntries(M_inv)); 268a175e481SJames Wright PetscCall(VecSet(ones, 1)); 269a175e481SJames Wright PetscCall(MatMult(mat_mass, ones, M_inv)); 270a175e481SJames Wright PetscCall(VecReciprocal(M_inv)); 271a175e481SJames Wright user->spanstats.M_inv = M_inv; 272a175e481SJames Wright PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &ones)); 273a175e481SJames Wright 274f967ad79SJames Wright PetscCall(KSPCreate(comm, &ksp)); 275b7d66439SJames Wright PetscCall(KSPSetOptionsPrefix(ksp, "turbulence_spanstats_")); 276a175e481SJames Wright { 277a175e481SJames Wright PC pc; 278a175e481SJames Wright PetscCall(KSPGetPC(ksp, &pc)); 279a175e481SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 280a175e481SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL)); 281a175e481SJames Wright PetscCall(KSPSetType(ksp, KSPCG)); 282a175e481SJames Wright PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL)); 283a175e481SJames Wright PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 284a175e481SJames Wright } 285a175e481SJames Wright PetscCall(KSPSetOperators(ksp, mat_mass, mat_mass)); 286a175e481SJames Wright PetscCall(KSPSetFromOptions(ksp)); 287a175e481SJames Wright user->spanstats.ksp = ksp; 288a175e481SJames Wright } 289a175e481SJames Wright 290a175e481SJames Wright // Cleanup 291a175e481SJames Wright CeedQFunctionDestroy(&qf_mass); 292*6665b873SJames Wright CeedOperatorDestroy(&op_mass); 293a175e481SJames Wright PetscFunctionReturn(0); 294a175e481SJames Wright } 295a175e481SJames Wright 296a175e481SJames Wright // Create CeedOperators and KSP for the statistics collection and processing 297a175e481SJames Wright PetscErrorCode CreateStatisticsOperators(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) { 298a175e481SJames Wright CeedInt num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q; 299a175e481SJames Wright CeedOperator op_setup_sur; 300495b9769SJames Wright Turbulence_SpanStatsContext collect_ctx; 301495b9769SJames Wright NewtonianIdealGasContext newtonian_ig_ctx; 302495b9769SJames Wright CeedQFunctionContext collect_context; 303a175e481SJames Wright PetscFunctionBeginUser; 304a175e481SJames Wright CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q); 305a175e481SJames Wright 306a175e481SJames Wright // Create Operator for statistics collection 307a175e481SJames Wright switch (user->phys->state_var) { 308a175e481SJames Wright case STATEVAR_PRIMITIVE: 309a175e481SJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &ceed_data->spanstats.qf_stats_collect); 310a175e481SJames Wright break; 311a175e481SJames Wright case STATEVAR_CONSERVATIVE: 312a175e481SJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &ceed_data->spanstats.qf_stats_collect); 313a175e481SJames Wright break; 314a175e481SJames Wright } 315a175e481SJames Wright 316b7d66439SJames Wright if (user->spanstats.do_mms_test) { 317a175e481SJames Wright CeedQFunctionDestroy(&ceed_data->spanstats.qf_stats_collect); 318b7d66439SJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest, ChildStatsCollectionMMSTest_loc, &ceed_data->spanstats.qf_stats_collect); 319a175e481SJames Wright } 320a175e481SJames Wright 321495b9769SJames Wright // Setup Collection Context 322495b9769SJames Wright { 323495b9769SJames Wright PetscCall(PetscNew(&collect_ctx)); 324495b9769SJames Wright CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx); 325495b9769SJames Wright collect_ctx->gas = *newtonian_ig_ctx; 326495b9769SJames Wright 327495b9769SJames Wright CeedQFunctionContextCreate(user->ceed, &collect_context); 328495b9769SJames Wright CeedQFunctionContextSetData(collect_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx); 329495b9769SJames Wright CeedQFunctionContextSetDataDestroy(collect_context, CEED_MEM_HOST, FreeContextPetsc); 330495b9769SJames Wright 331495b9769SJames Wright CeedQFunctionContextRegisterDouble(collect_context, "solution time", offsetof(struct Turbulence_SpanStatsContext_, solution_time), 1, 332495b9769SJames Wright "Current solution time"); 333495b9769SJames Wright CeedQFunctionContextRegisterDouble(collect_context, "previous time", offsetof(struct Turbulence_SpanStatsContext_, previous_time), 1, 334495b9769SJames Wright "Previous time statistics collection was done"); 335495b9769SJames Wright 336495b9769SJames Wright CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx); 337495b9769SJames Wright } 338495b9769SJames Wright 339495b9769SJames Wright CeedQFunctionSetContext(ceed_data->spanstats.qf_stats_collect, collect_context); 340495b9769SJames Wright CeedQFunctionContextDestroy(&collect_context); 341a175e481SJames Wright CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP); 342a175e481SJames Wright CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE); 343a175e481SJames Wright CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP); 344a175e481SJames Wright CeedQFunctionAddOutput(ceed_data->spanstats.qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE); 345a175e481SJames Wright 346a175e481SJames Wright CeedOperatorCreate(ceed, ceed_data->spanstats.qf_stats_collect, NULL, NULL, &user->spanstats.op_stats_collect); 347a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 348a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 349a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord); 350a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_collect, "v", ceed_data->spanstats.elem_restr_child_colloc, CEED_BASIS_COLLOCATED, 351a175e481SJames Wright CEED_VECTOR_ACTIVE); 352a175e481SJames Wright 353495b9769SJames Wright CeedOperatorContextGetFieldLabel(user->spanstats.op_stats_collect, "solution time", &user->spanstats.solution_time_label); 354495b9769SJames Wright CeedOperatorContextGetFieldLabel(user->spanstats.op_stats_collect, "previous time", &user->spanstats.previous_time_label); 355495b9769SJames Wright 356823a1283SJames Wright // Create Operator for L^2 projection of statistics 357a175e481SJames Wright // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function. 358a175e481SJames Wright // Therefore, an Identity QF is sufficient 359a175e481SJames Wright CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &ceed_data->spanstats.qf_stats_proj); 360a175e481SJames Wright 361a175e481SJames Wright CeedOperatorCreate(ceed, ceed_data->spanstats.qf_stats_proj, NULL, NULL, &user->spanstats.op_stats_proj); 362a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_proj, "input", ceed_data->spanstats.elem_restr_parent_colloc, CEED_BASIS_COLLOCATED, 363a175e481SJames Wright CEED_VECTOR_ACTIVE); 364a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_proj, "output", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, 365a175e481SJames Wright CEED_VECTOR_ACTIVE); 366a175e481SJames Wright 367a175e481SJames Wright // Get q_data for lumped mass matrix formation 368a175e481SJames Wright CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur); 369a175e481SJames Wright CeedOperatorSetField(op_setup_sur, "dx", ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, CEED_VECTOR_ACTIVE); 370a175e481SJames Wright CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->spanstats.basis_x, CEED_VECTOR_NONE); 371a175e481SJames Wright CeedOperatorSetField(op_setup_sur, "surface qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 372a175e481SJames Wright CeedOperatorApply(op_setup_sur, ceed_data->spanstats.x_coord, ceed_data->spanstats.q_data, CEED_REQUEST_IMMEDIATE); 373a175e481SJames Wright 374a175e481SJames Wright CeedOperatorDestroy(&op_setup_sur); 375a175e481SJames Wright PetscFunctionReturn(0); 376a175e481SJames Wright } 377a175e481SJames Wright 378b7d66439SJames Wright // Creates operator for calculating error of method of manufactured solution (MMS) test 379b7d66439SJames Wright PetscErrorCode SetupMMSErrorChecking(Ceed ceed, User user, CeedData ceed_data) { 380823a1283SJames Wright CeedInt num_comp_stats = user->spanstats.num_comp_stats, num_comp_x; 381823a1283SJames Wright CeedQFunction qf_error; 382823a1283SJames Wright CeedOperator op_error; 383823a1283SJames Wright CeedInt q_data_size; 384823a1283SJames Wright CeedVector x_ceed, y_ceed; 385823a1283SJames Wright PetscFunctionBeginUser; 386823a1283SJames Wright 387823a1283SJames Wright CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_qd, &q_data_size); 388823a1283SJames Wright CeedBasisGetNumComponents(ceed_data->spanstats.basis_x, &num_comp_x); 389823a1283SJames Wright 390b7d66439SJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest_Error, ChildStatsCollectionMMSTest_Error_loc, &qf_error); 391823a1283SJames Wright CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP); 392823a1283SJames Wright CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE); 393823a1283SJames Wright CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP); 394823a1283SJames Wright CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP); 395823a1283SJames Wright 396823a1283SJames Wright CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error); 397823a1283SJames Wright CeedOperatorSetField(op_error, "q", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE); 398823a1283SJames Wright CeedOperatorSetField(op_error, "qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, ceed_data->spanstats.q_data); 399823a1283SJames Wright CeedOperatorSetField(op_error, "x", ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, ceed_data->spanstats.x_coord); 400823a1283SJames Wright CeedOperatorSetField(op_error, "v", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE); 401823a1283SJames Wright 402823a1283SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &x_ceed, NULL); 403823a1283SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &y_ceed, NULL); 404*6665b873SJames Wright PetscCall(MatopApplyContextCreate(user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, &user->spanstats.mms_error_ctx)); 405823a1283SJames Wright 406*6665b873SJames Wright CeedOperatorDestroy(&op_error); 407*6665b873SJames Wright CeedVectorDestroy(&x_ceed); 408*6665b873SJames Wright CeedVectorDestroy(&y_ceed); 409823a1283SJames Wright PetscFunctionReturn(0); 410823a1283SJames Wright } 411823a1283SJames Wright 412a175e481SJames Wright // Setup for statistics collection 41319706a06SJames Wright PetscErrorCode SetupStatsCollection(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) { 41419706a06SJames Wright DM dm = user->spanstats.dm; 41519706a06SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm); 41619706a06SJames Wright CeedInt dim, P, Q, num_comp_x; 41719706a06SJames Wright Vec X_loc; 41819706a06SJames Wright PetscMemType X_loc_memtype; 41919706a06SJames Wright const PetscScalar *X_loc_array; 4204eed8630SJames Wright PetscLogStage stage_stats_setup; 42119706a06SJames Wright PetscFunctionBeginUser; 42219706a06SJames Wright 4234eed8630SJames Wright PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup)); 4244eed8630SJames Wright if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup)); 4254eed8630SJames Wright PetscCall(PetscLogStagePush(stage_stats_setup)); 4264eed8630SJames Wright 42719706a06SJames Wright PetscCall(DMGetDimension(dm, &dim)); 428a175e481SJames Wright CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &Q); 429a175e481SJames Wright CeedBasisGetNumNodes1D(ceed_data->basis_q, &P); 43019706a06SJames Wright 4311737222fSJames Wright PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, problem->q_data_size_sur, &ceed_data->spanstats.elem_restr_parent_stats, 43219706a06SJames Wright &ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.elem_restr_parent_qd)); 43319706a06SJames Wright CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_x, &num_comp_x); 43419706a06SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.x_coord, NULL); 435a175e481SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &user->spanstats.rhs_ceed, NULL); 4361737222fSJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_qd, &ceed_data->spanstats.q_data, NULL); 43719706a06SJames Wright 438a175e481SJames Wright CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS, &ceed_data->spanstats.basis_x); 43919706a06SJames Wright CeedBasisCreateTensorH1Lagrange(ceed, dim, user->spanstats.num_comp_stats, P, Q, CEED_GAUSS, &ceed_data->spanstats.basis_stats); 44019706a06SJames Wright 4411737222fSJames Wright PetscCall(CreateElemRestrColloc(ceed, user->spanstats.num_comp_stats, ceed_data->spanstats.basis_stats, 4421737222fSJames Wright ceed_data->spanstats.elem_restr_parent_stats, &ceed_data->spanstats.elem_restr_parent_colloc, 4431737222fSJames Wright &user->spanstats.parent_stats, NULL)); 4441737222fSJames Wright PetscCall(CreateElemRestrColloc(ceed, user->spanstats.num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q, 4451737222fSJames Wright &ceed_data->spanstats.elem_restr_child_colloc, &user->spanstats.child_stats, NULL)); 446a175e481SJames Wright CeedVectorSetValue(user->spanstats.child_stats, 0); 4471737222fSJames Wright 448f967ad79SJames Wright { // -- Copy DM coordinates into CeedVector 44919706a06SJames Wright DM cdm; 45019706a06SJames Wright PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 45119706a06SJames Wright if (cdm) { 45219706a06SJames Wright PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 45319706a06SJames Wright } else { 45419706a06SJames Wright PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 45519706a06SJames Wright } 45619706a06SJames Wright } 45719706a06SJames Wright PetscCall(VecScale(X_loc, problem->dm_scale)); 45819706a06SJames Wright PetscCall(VecGetArrayReadAndMemType(X_loc, &X_loc_array, &X_loc_memtype)); 45919706a06SJames Wright CeedVectorSetArray(ceed_data->spanstats.x_coord, MemTypeP2C(X_loc_memtype), CEED_COPY_VALUES, (PetscScalar *)X_loc_array); 46019706a06SJames Wright PetscCall(VecRestoreArrayRead(X_loc, &X_loc_array)); 46119706a06SJames Wright 46219706a06SJames Wright // Create SF for communicating child data back their respective parents 46319706a06SJames Wright PetscCall(PetscSFCreate(comm, &user->spanstats.sf)); 46419706a06SJames Wright PetscCall(CreateStatsSF(ceed, ceed_data, user->dm, user->spanstats.dm, user->spanstats.sf)); 46551ee423eSJames Wright 466a175e481SJames Wright // Create CeedOperators for statistics collection 467a175e481SJames Wright PetscCall(CreateStatisticsOperators(ceed, user, ceed_data, problem)); 468a175e481SJames Wright 469a175e481SJames Wright // Setup KSP and Mat for L^2 projection of statistics 470a175e481SJames Wright PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data)); 471a175e481SJames Wright 472b7d66439SJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_turbulence_spanstats_mms", &user->spanstats.do_mms_test, NULL)); 473b7d66439SJames Wright if (user->spanstats.do_mms_test) { 474b7d66439SJames Wright PetscCall(SetupMMSErrorChecking(ceed, user, ceed_data)); 475823a1283SJames Wright } 476823a1283SJames Wright 477f967ad79SJames Wright { // Setup stats viewer with prefix 4788ed52730SJames Wright PetscViewerType viewer_type; 479b7d66439SJames Wright PetscCall(PetscViewerGetType(user->app_ctx->turb_spanstats_viewer, &viewer_type)); 480b7d66439SJames Wright PetscCall(PetscOptionsSetValue(NULL, "-ts_monitor_turbulence_spanstats_viewer_type", viewer_type)); 4818ed52730SJames Wright 482b7d66439SJames Wright PetscCall(PetscViewerSetOptionsPrefix(user->app_ctx->turb_spanstats_viewer, "ts_monitor_turbulence_spanstats_")); 483b7d66439SJames Wright PetscCall(PetscViewerSetFromOptions(user->app_ctx->turb_spanstats_viewer)); 4848ed52730SJames Wright } 4854eed8630SJames Wright 4864eed8630SJames Wright PetscCall(PetscLogStagePop()); 487a175e481SJames Wright PetscFunctionReturn(0); 488a175e481SJames Wright } 489a175e481SJames Wright 490a175e481SJames Wright // Collect statistics based on the solution Q 491a175e481SJames Wright PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) { 492a175e481SJames Wright PetscMemType q_mem_type; 493a175e481SJames Wright const PetscScalar *q_arr; 494a175e481SJames Wright PetscFunctionBeginUser; 495a175e481SJames Wright 4966dcea3beSJames Wright PetscLogStage stage_stats_collect; 4976dcea3beSJames Wright PetscCall(PetscLogStageGetId("Stats Collect", &stage_stats_collect)); 4986dcea3beSJames Wright if (stage_stats_collect == -1) PetscCall(PetscLogStageRegister("Stats Collect", &stage_stats_collect)); 4996dcea3beSJames Wright PetscCall(PetscLogStagePush(stage_stats_collect)); 5006dcea3beSJames Wright 501a0b9a424SJames Wright PetscCall(UpdateBoundaryValues(user, user->Q_loc, solution_time)); 502495b9769SJames Wright CeedOperatorContextSetDouble(user->spanstats.op_stats_collect, user->spanstats.solution_time_label, &solution_time); 503a0b9a424SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, user->Q_loc)); 504a0b9a424SJames Wright PetscCall(VecGetArrayReadAndMemType(user->Q_loc, &q_arr, &q_mem_type)); 505a175e481SJames Wright CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q_arr); 506a175e481SJames Wright 507495b9769SJames Wright CeedOperatorApplyAdd(user->spanstats.op_stats_collect, user->q_ceed, user->spanstats.child_stats, CEED_REQUEST_IMMEDIATE); 508a175e481SJames Wright 509a175e481SJames Wright CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL); 510a0b9a424SJames Wright PetscCall(VecRestoreArrayReadAndMemType(user->Q_loc, &q_arr)); 511a175e481SJames Wright 512495b9769SJames Wright CeedOperatorContextSetDouble(user->spanstats.op_stats_collect, user->spanstats.previous_time_label, &solution_time); 513a175e481SJames Wright 5146dcea3beSJames Wright PetscCall(PetscLogStagePop()); 515a175e481SJames Wright PetscFunctionReturn(0); 516a175e481SJames Wright } 517a175e481SJames Wright 518a175e481SJames Wright // Process the child statistics into parent statistics and project them onto stats 51978bbfb6fSJed Brown PetscErrorCode ProcessStatistics(User user, Vec stats) { 520a175e481SJames Wright Span_Stats user_stats = user->spanstats; 521a175e481SJames Wright const PetscScalar *child_stats; 522a175e481SJames Wright PetscScalar *parent_stats; 523a175e481SJames Wright MPI_Datatype unit; 524a175e481SJames Wright Vec rhs_loc, rhs; 525a175e481SJames Wright PetscMemType rhs_mem_type; 526a175e481SJames Wright CeedScalar *rhs_arr; 527a175e481SJames Wright CeedMemType ceed_mem_type; 528a175e481SJames Wright PetscFunctionBeginUser; 529a175e481SJames Wright 5306dcea3beSJames Wright PetscLogStage stage_stats_process; 5316dcea3beSJames Wright PetscCall(PetscLogStageGetId("Stats Process", &stage_stats_process)); 5326dcea3beSJames Wright if (stage_stats_process == -1) PetscCall(PetscLogStageRegister("Stats Process", &stage_stats_process)); 5336dcea3beSJames Wright PetscCall(PetscLogStagePush(stage_stats_process)); 5346dcea3beSJames Wright 535a175e481SJames Wright CeedGetPreferredMemType(user->ceed, &ceed_mem_type); 536a175e481SJames Wright CeedVectorSetValue(user_stats.parent_stats, 0); 537a175e481SJames Wright 538a175e481SJames Wright CeedVectorGetArrayRead(user_stats.child_stats, ceed_mem_type, &child_stats); 539a175e481SJames Wright CeedVectorGetArray(user_stats.parent_stats, ceed_mem_type, &parent_stats); 540a175e481SJames Wright 541a175e481SJames Wright if (user_stats.num_comp_stats == 1) unit = MPIU_REAL; 542a175e481SJames Wright else { 543a175e481SJames Wright PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit)); 544a175e481SJames Wright PetscCallMPI(MPI_Type_commit(&unit)); 545a175e481SJames Wright } 546a175e481SJames Wright 547a175e481SJames Wright PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM)); 548a175e481SJames Wright PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM)); 549a175e481SJames Wright 550a175e481SJames Wright CeedVectorRestoreArrayRead(user_stats.child_stats, &child_stats); 551a175e481SJames Wright CeedVectorRestoreArray(user_stats.parent_stats, &parent_stats); 552a175e481SJames Wright PetscCallMPI(MPI_Type_free(&unit)); 553a175e481SJames Wright 55478bbfb6fSJed Brown PetscReal solution_time; 55578bbfb6fSJed Brown PetscCall(DMGetOutputSequenceNumber(user_stats.dm, NULL, &solution_time)); 55678bbfb6fSJed Brown PetscReal summing_duration = solution_time - user->app_ctx->cont_time; 55778bbfb6fSJed Brown CeedVectorScale(user_stats.parent_stats, 1 / (summing_duration * user_stats.span_width)); 558a175e481SJames Wright 559a175e481SJames Wright // L^2 projection with the parent_data 560a175e481SJames Wright PetscCall(DMGetLocalVector(user_stats.dm, &rhs_loc)); 561a175e481SJames Wright PetscCall(VecZeroEntries(rhs_loc)); 562a175e481SJames Wright PetscCall(VecGetArrayWriteAndMemType(rhs_loc, &rhs_arr, &rhs_mem_type)); 563a175e481SJames Wright CeedVectorSetArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), CEED_USE_POINTER, (PetscScalar *)rhs_arr); 564a175e481SJames Wright 565a175e481SJames Wright CeedOperatorApply(user_stats.op_stats_proj, user_stats.parent_stats, user_stats.rhs_ceed, CEED_REQUEST_IMMEDIATE); 566a175e481SJames Wright 567a175e481SJames Wright CeedVectorTakeArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), &rhs_arr); 568a175e481SJames Wright PetscCall(VecRestoreArrayAndMemType(rhs_loc, &rhs_arr)); 56978bbfb6fSJed Brown 57078bbfb6fSJed Brown PetscCall(DMGetGlobalVector(user_stats.dm, &rhs)); 57178bbfb6fSJed Brown PetscCall(VecZeroEntries(rhs)); 572a175e481SJames Wright PetscCall(DMLocalToGlobal(user_stats.dm, rhs_loc, ADD_VALUES, rhs)); 57378bbfb6fSJed Brown PetscCall(DMRestoreLocalVector(user_stats.dm, &rhs_loc)); 574a175e481SJames Wright 57578bbfb6fSJed Brown PetscCall(KSPSolve(user_stats.ksp, rhs, stats)); 576a175e481SJames Wright 57778bbfb6fSJed Brown PetscCall(DMRestoreGlobalVector(user_stats.dm, &rhs)); 5786dcea3beSJames Wright PetscCall(PetscLogStagePop()); 579a175e481SJames Wright PetscFunctionReturn(0); 580a175e481SJames Wright } 581a175e481SJames Wright 582a175e481SJames Wright // TSMonitor for the statistics collection and processing 583a175e481SJames Wright PetscErrorCode TSMonitor_Statistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) { 584a175e481SJames Wright User user = (User)ctx; 585a175e481SJames Wright Vec stats; 58678bbfb6fSJed Brown TSConvergedReason reason; 587b7d66439SJames Wright PetscInt collect_interval = user->app_ctx->turb_spanstats_collect_interval, viewer_interval = user->app_ctx->turb_spanstats_viewer_interval; 588a175e481SJames Wright PetscFunctionBeginUser; 58978bbfb6fSJed Brown PetscCall(TSGetConvergedReason(ts, &reason)); 590a175e481SJames Wright // Do not collect or process on the first step of the run (ie. on the initial condition) 59178bbfb6fSJed Brown if (steps == user->app_ctx->cont_steps && reason == TS_CONVERGED_ITERATING) PetscFunctionReturn(0); 592a175e481SJames Wright 59327240365SJames Wright PetscBool run_processing_and_viewer = (steps % viewer_interval == 0 && viewer_interval != -1) || reason != TS_CONVERGED_ITERATING; 594a175e481SJames Wright 59527240365SJames Wright if (steps % collect_interval == 0 || run_processing_and_viewer) { 59627240365SJames Wright PetscCall(CollectStatistics(user, solution_time, Q)); 59727240365SJames Wright 59827240365SJames Wright if (run_processing_and_viewer) { 5998ed52730SJames Wright PetscCall(DMSetOutputSequenceNumber(user->spanstats.dm, steps, solution_time)); 60078bbfb6fSJed Brown PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats)); 60178bbfb6fSJed Brown PetscCall(ProcessStatistics(user, stats)); 6028fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_NONE) { 603b7d66439SJames Wright PetscCall(PetscViewerPushFormat(user->app_ctx->turb_spanstats_viewer, user->app_ctx->turb_spanstats_viewer_format)); 604b7d66439SJames Wright PetscCall(VecView(stats, user->app_ctx->turb_spanstats_viewer)); 605b7d66439SJames Wright PetscCall(PetscViewerPopFormat(user->app_ctx->turb_spanstats_viewer)); 6068fb33541SJames Wright } 6078fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_TURB_SPANSTATS && reason != TS_CONVERGED_ITERATING) { 6088fb33541SJames Wright PetscCall(RegressionTests_NS(user->app_ctx, stats)); 6098fb33541SJames Wright } 610b7d66439SJames Wright if (user->spanstats.do_mms_test && reason != TS_CONVERGED_ITERATING) { 611823a1283SJames Wright Vec error; 612823a1283SJames Wright PetscCall(VecDuplicate(stats, &error)); 613b7d66439SJames Wright PetscCall(ApplyLocal_Ceed(stats, error, user->spanstats.mms_error_ctx)); 614823a1283SJames Wright PetscScalar error_sq = 0; 615823a1283SJames Wright PetscCall(VecSum(error, &error_sq)); 616823a1283SJames Wright PetscScalar l2_error = sqrt(error_sq); 617823a1283SJames Wright PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error)); 618823a1283SJames Wright } 619a175e481SJames Wright PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats)); 620522ee345SJames Wright } 62127240365SJames Wright } 62251ee423eSJames Wright PetscFunctionReturn(0); 62351ee423eSJames Wright } 624fd170fd0SJames Wright 625*6665b873SJames Wright PetscErrorCode DestroyStats(User user, CeedData ceed_data) { 626fd170fd0SJames Wright PetscFunctionBeginUser; 627fd170fd0SJames Wright 628fd170fd0SJames Wright // -- CeedVectors 629fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.child_stats); 630fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.parent_stats); 631fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.rhs_ceed); 632fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.x_ceed); 633fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.y_ceed); 634fd170fd0SJames Wright CeedVectorDestroy(&ceed_data->spanstats.x_coord); 635fd170fd0SJames Wright CeedVectorDestroy(&ceed_data->spanstats.q_data); 636fd170fd0SJames Wright 637fd170fd0SJames Wright // -- QFunctions 638fd170fd0SJames Wright CeedQFunctionDestroy(&ceed_data->spanstats.qf_stats_collect); 639fd170fd0SJames Wright CeedQFunctionDestroy(&ceed_data->spanstats.qf_stats_proj); 640fd170fd0SJames Wright 641fd170fd0SJames Wright // -- CeedBasis 642fd170fd0SJames Wright CeedBasisDestroy(&ceed_data->spanstats.basis_stats); 643fd170fd0SJames Wright CeedBasisDestroy(&ceed_data->spanstats.basis_x); 644fd170fd0SJames Wright 645fd170fd0SJames Wright // -- CeedElemRestriction 646fd170fd0SJames Wright CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_parent_stats); 647fd170fd0SJames Wright CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_parent_qd); 648fd170fd0SJames Wright CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_parent_x); 649fd170fd0SJames Wright CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_parent_colloc); 650fd170fd0SJames Wright CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_child_colloc); 651fd170fd0SJames Wright 652fd170fd0SJames Wright // -- CeedOperators 653fd170fd0SJames Wright CeedOperatorDestroy(&user->spanstats.op_stats_collect); 654fd170fd0SJames Wright CeedOperatorDestroy(&user->spanstats.op_stats_proj); 655*6665b873SJames Wright PetscCall(MatopApplyContextDestroy(user->spanstats.mms_error_ctx)); 656fd170fd0SJames Wright 657fd170fd0SJames Wright // -- Vec 658fd170fd0SJames Wright PetscCall(VecDestroy(&user->spanstats.M_inv)); 659fd170fd0SJames Wright 660fd170fd0SJames Wright // -- KSP 661fd170fd0SJames Wright PetscCall(KSPDestroy(&user->spanstats.ksp)); 662fd170fd0SJames Wright 663fd170fd0SJames Wright // -- SF 664fd170fd0SJames Wright PetscCall(PetscSFDestroy(&user->spanstats.sf)); 665fd170fd0SJames Wright 666fd170fd0SJames Wright // -- DM 667fd170fd0SJames Wright PetscCall(DMDestroy(&user->spanstats.dm)); 668fd170fd0SJames Wright 669fd170fd0SJames Wright PetscFunctionReturn(0); 670fd170fd0SJames Wright } 671