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 24*269a910fSJames Wright typedef struct { 25*269a910fSJames Wright CeedElemRestriction elem_restr_parent_x, elem_restr_parent_stats, elem_restr_parent_qd, elem_restr_parent_colloc, elem_restr_child_colloc; 26*269a910fSJames Wright CeedBasis basis_x, basis_stats; 27*269a910fSJames Wright CeedVector x_coord, q_data; 28*269a910fSJames Wright } * SpanStatsSetupData; 29*269a910fSJames Wright 3051ee423eSJames Wright PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree, SimpleBC bc) { 313a4208e6SJames Wright user->spanstats.num_comp_stats = TURB_NUM_COMPONENTS; 32a175e481SJames Wright PetscReal domain_min[3], domain_max[3]; 3378837792SJames Wright PetscFE fe; 3478837792SJames Wright PetscSection section; 354eed8630SJames Wright PetscLogStage stage_stats_setup; 361f595ac1SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)user->dm); 3751ee423eSJames Wright PetscFunctionBeginUser; 3851ee423eSJames Wright 394eed8630SJames Wright PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup)); 404eed8630SJames Wright if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup)); 414eed8630SJames Wright PetscCall(PetscLogStagePush(stage_stats_setup)); 424eed8630SJames Wright 43a175e481SJames Wright // Get spanwise length 44a175e481SJames Wright PetscCall(DMGetBoundingBox(user->dm, domain_min, domain_max)); 45a175e481SJames Wright user->spanstats.span_width = domain_max[2] - domain_min[1]; 46a175e481SJames Wright 47f967ad79SJames Wright { // Get DM from surface 481f595ac1SJames Wright DM parent_distributed_dm; 49c9198418SJames Wright PetscSF isoperiodicface; 5051ee423eSJames Wright DMLabel label; 511f595ac1SJames Wright PetscMPIInt size; 52c9198418SJames Wright 53c9198418SJames Wright PetscCall(DMPlexGetIsoperiodicFaceSF(user->dm, &isoperiodicface)); 54c9198418SJames Wright 55c9198418SJames Wright if (isoperiodicface) { 56c9198418SJames Wright PetscSF inv_isoperiodicface; 57c9198418SJames Wright PetscInt nleaves; 58c9198418SJames Wright const PetscInt *ilocal; 59c9198418SJames Wright 60c9198418SJames Wright PetscCall(PetscSFCreateInverseSF(isoperiodicface, &inv_isoperiodicface)); 61c9198418SJames Wright PetscCall(PetscSFGetGraph(inv_isoperiodicface, NULL, &nleaves, &ilocal, NULL)); 62c9198418SJames Wright PetscCall(DMCreateLabel(user->dm, "Periodic Face")); 63c9198418SJames Wright PetscCall(DMGetLabel(user->dm, "Periodic Face", &label)); 64c9198418SJames Wright for (PetscInt i = 0; i < nleaves; i++) { 65c9198418SJames Wright PetscCall(DMLabelSetValue(label, ilocal[i], 1)); 66c9198418SJames Wright } 67c9198418SJames Wright } else { 6851ee423eSJames Wright PetscCall(DMGetLabel(user->dm, "Face Sets", &label)); 69c9198418SJames Wright } 70c9198418SJames Wright 7151ee423eSJames Wright PetscCall(DMPlexLabelComplete(user->dm, label)); 7251ee423eSJames Wright PetscCall(DMPlexFilter(user->dm, label, 1, &user->spanstats.dm)); 7351ee423eSJames Wright PetscCall(DMProjectCoordinates(user->spanstats.dm, NULL)); // Ensure that a coordinate FE exists 741f595ac1SJames Wright 751f595ac1SJames Wright PetscCall(DMPlexDistribute(user->spanstats.dm, 0, NULL, &parent_distributed_dm)); 761f595ac1SJames Wright PetscCallMPI(MPI_Comm_size(comm, &size)); 771f595ac1SJames Wright if (parent_distributed_dm) { 781f595ac1SJames Wright PetscCall(DMDestroy(&user->spanstats.dm)); 791f595ac1SJames Wright user->spanstats.dm = parent_distributed_dm; 801f595ac1SJames Wright } else if (size > 1) { 811f595ac1SJames Wright PetscCall(PetscPrintf(comm, "WARNING: Turbulent spanwise statistics: parent DM could not be distributed accross %d ranks.\n", size)); 821f595ac1SJames Wright } 8351ee423eSJames Wright } 8451ee423eSJames Wright 8551ee423eSJames Wright PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats")); 86b7d66439SJames Wright PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "turbulence_spanstats_")); 8751ee423eSJames Wright PetscCall(DMSetFromOptions(user->spanstats.dm)); 88a175e481SJames Wright PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view")); // -spanstats_dm_view 89c9198418SJames Wright 9078837792SJames Wright // Create FE space for parent DM 9151ee423eSJames Wright PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim - 1, user->spanstats.num_comp_stats, PETSC_FALSE, degree, PETSC_DECIDE, &fe)); 9251ee423eSJames Wright PetscCall(PetscObjectSetName((PetscObject)fe, "stats")); 9351ee423eSJames Wright PetscCall(DMAddField(user->spanstats.dm, NULL, (PetscObject)fe)); 9451ee423eSJames Wright PetscCall(DMCreateDS(user->spanstats.dm)); 9551ee423eSJames Wright PetscCall(DMPlexSetClosurePermutationTensor(user->spanstats.dm, PETSC_DETERMINE, NULL)); 9651ee423eSJames Wright 9778837792SJames Wright // Create Section for data 9851ee423eSJames Wright PetscCall(DMGetLocalSection(user->spanstats.dm, §ion)); 9951ee423eSJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 1003a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY, "MeanDensity")); 1013a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE, "MeanPressure")); 1023a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_SQUARED, "MeanPressureSquared")); 1033a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_X, "MeanPressureVelocityX")); 1043a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Y, "MeanPressureVelocityY")); 1053a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Z, "MeanPressureVelocityZ")); 1063a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE, "MeanDensityTemperature")); 1073a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_X, "MeanDensityTemperatureFluxX")); 1083a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Y, "MeanDensityTemperatureFluxY")); 1093a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Z, "MeanDensityTemperatureFluxZ")); 1103a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_X, "MeanMomentumX")); 1113a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Y, "MeanMomentumY")); 1123a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Z, "MeanMomentumZ")); 1133a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XX, "MeanMomentumFluxXX")); 1143a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YY, "MeanMomentumFluxYY")); 1153a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_ZZ, "MeanMomentumFluxZZ")); 1163a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YZ, "MeanMomentumFluxYZ")); 1173a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XZ, "MeanMomentumFluxXZ")); 1183a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XY, "MeanMomentumFluxXY")); 1193a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_X, "MeanVelocityX")); 1203a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Y, "MeanVelocityY")); 1213a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Z, "MeanVelocityZ")); 12219706a06SJames Wright 12378837792SJames Wright // Cleanup 12478837792SJames Wright PetscCall(PetscFEDestroy(&fe)); 12578837792SJames Wright 1264eed8630SJames Wright PetscCall(PetscLogStagePop()); 12719706a06SJames Wright PetscFunctionReturn(0); 12819706a06SJames Wright } 12919706a06SJames Wright 1301737222fSJames Wright // Create CeedElemRestriction for collocated data based on associated CeedBasis and CeedElemRestriction 1311737222fSJames Wright // Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction 1321737222fSJames Wright PetscErrorCode CreateElemRestrColloc(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base, 1331737222fSJames Wright CeedElemRestriction *elem_restr_collocated, CeedVector *l_vec, CeedVector *e_vec) { 1341737222fSJames Wright CeedInt num_elem_qpts, loc_num_elem; 1351737222fSJames Wright PetscFunctionBeginUser; 1361737222fSJames Wright 1371737222fSJames Wright CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts); 1381737222fSJames Wright CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem); 1391737222fSJames Wright 1401737222fSJames Wright const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp}; 1411737222fSJames Wright CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides, 1421737222fSJames Wright elem_restr_collocated); 1431737222fSJames Wright CeedElemRestrictionCreateVector(*elem_restr_collocated, l_vec, e_vec); 1441737222fSJames Wright PetscFunctionReturn(0); 1451737222fSJames Wright } 1461737222fSJames Wright 147a175e481SJames Wright // Get coordinates of quadrature points 14819706a06SJames Wright PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, CeedVector *qx_coords, 14919706a06SJames Wright PetscInt *total_nqpnts) { 150*269a910fSJames Wright CeedElemRestriction elem_restr_qx; 15119706a06SJames Wright CeedQFunction qf_quad_coords; 15219706a06SJames Wright CeedOperator op_quad_coords; 15319706a06SJames Wright PetscInt num_comp_x, loc_num_elem, num_elem_qpts; 15419706a06SJames Wright PetscFunctionBeginUser; 15519706a06SJames Wright 15619706a06SJames Wright // Create Element Restriction and CeedVector for quadrature coordinates 15719706a06SJames Wright CeedBasisGetNumQuadraturePoints(basis_x, &num_elem_qpts); 15819706a06SJames Wright CeedElemRestrictionGetNumElements(elem_restr_x, &loc_num_elem); 15919706a06SJames Wright CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x); 16019706a06SJames Wright *total_nqpnts = num_elem_qpts * loc_num_elem; 1611737222fSJames Wright PetscCall(CreateElemRestrColloc(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx, qx_coords, NULL)); 16219706a06SJames Wright 16319706a06SJames Wright // Create QFunction 16419706a06SJames Wright CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords); 16519706a06SJames Wright 16619706a06SJames Wright // Create Operator 16719706a06SJames Wright CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords); 16819706a06SJames Wright CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE); 16919706a06SJames Wright CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 17019706a06SJames Wright 17119706a06SJames Wright CeedOperatorApply(op_quad_coords, x_coords, *qx_coords, CEED_REQUEST_IMMEDIATE); 17219706a06SJames Wright 173*269a910fSJames Wright CeedElemRestrictionDestroy(&elem_restr_qx); 17419706a06SJames Wright CeedQFunctionDestroy(&qf_quad_coords); 17519706a06SJames Wright CeedOperatorDestroy(&op_quad_coords); 17619706a06SJames Wright PetscFunctionReturn(0); 17719706a06SJames Wright } 17819706a06SJames Wright 179*269a910fSJames Wright PetscErrorCode SpanStatsSetupDataCreate(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem, SpanStatsSetupData *stats_data) { 180*269a910fSJames Wright DM dm = user->spanstats.dm; 181*269a910fSJames Wright CeedInt dim, P, Q, num_comp_x, num_comp_stats = user->spanstats.num_comp_stats; 182*269a910fSJames Wright Vec X_loc; 183*269a910fSJames Wright PetscMemType X_loc_memtype; 184*269a910fSJames Wright const PetscScalar *X_loc_array; 185*269a910fSJames Wright PetscFunctionBeginUser; 186*269a910fSJames Wright 187*269a910fSJames Wright PetscCall(PetscNew(stats_data)); 188*269a910fSJames Wright 189*269a910fSJames Wright PetscCall(DMGetDimension(dm, &dim)); 190*269a910fSJames Wright CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &Q); 191*269a910fSJames Wright CeedBasisGetNumNodes1D(ceed_data->basis_q, &P); 192*269a910fSJames Wright 193*269a910fSJames Wright PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, problem->q_data_size_sur, &(*stats_data)->elem_restr_parent_stats, 194*269a910fSJames Wright &(*stats_data)->elem_restr_parent_x, &(*stats_data)->elem_restr_parent_qd)); 195*269a910fSJames Wright CeedElemRestrictionGetNumComponents((*stats_data)->elem_restr_parent_x, &num_comp_x); 196*269a910fSJames Wright CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_x, &(*stats_data)->x_coord, NULL); 197*269a910fSJames Wright CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_qd, &(*stats_data)->q_data, NULL); 198*269a910fSJames Wright 199*269a910fSJames Wright CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS, &(*stats_data)->basis_x); 200*269a910fSJames Wright CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_stats, P, Q, CEED_GAUSS, &(*stats_data)->basis_stats); 201*269a910fSJames Wright 202*269a910fSJames Wright PetscCall(CreateElemRestrColloc(ceed, num_comp_stats, (*stats_data)->basis_stats, (*stats_data)->elem_restr_parent_stats, 203*269a910fSJames Wright &(*stats_data)->elem_restr_parent_colloc, NULL, NULL)); 204*269a910fSJames Wright PetscCall( 205*269a910fSJames Wright CreateElemRestrColloc(ceed, num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q, &(*stats_data)->elem_restr_child_colloc, NULL, NULL)); 206*269a910fSJames Wright 207*269a910fSJames Wright { // -- Copy DM coordinates into CeedVector 208*269a910fSJames Wright DM cdm; 209*269a910fSJames Wright PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 210*269a910fSJames Wright if (cdm) { 211*269a910fSJames Wright PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 212*269a910fSJames Wright } else { 213*269a910fSJames Wright PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 214*269a910fSJames Wright } 215*269a910fSJames Wright } 216*269a910fSJames Wright PetscCall(VecScale(X_loc, problem->dm_scale)); 217*269a910fSJames Wright PetscCall(VecGetArrayReadAndMemType(X_loc, &X_loc_array, &X_loc_memtype)); 218*269a910fSJames Wright CeedVectorSetArray((*stats_data)->x_coord, MemTypeP2C(X_loc_memtype), CEED_COPY_VALUES, (PetscScalar *)X_loc_array); 219*269a910fSJames Wright PetscCall(VecRestoreArrayRead(X_loc, &X_loc_array)); 220*269a910fSJames Wright 221*269a910fSJames Wright PetscFunctionReturn(0); 222*269a910fSJames Wright } 223*269a910fSJames Wright 224*269a910fSJames Wright PetscErrorCode SpanStatsSetupDataDestroy(SpanStatsSetupData data) { 225*269a910fSJames Wright PetscFunctionBeginUser; 226*269a910fSJames Wright 227*269a910fSJames Wright CeedElemRestrictionDestroy(&data->elem_restr_parent_x); 228*269a910fSJames Wright CeedElemRestrictionDestroy(&data->elem_restr_parent_stats); 229*269a910fSJames Wright CeedElemRestrictionDestroy(&data->elem_restr_parent_qd); 230*269a910fSJames Wright CeedElemRestrictionDestroy(&data->elem_restr_parent_colloc); 231*269a910fSJames Wright CeedElemRestrictionDestroy(&data->elem_restr_child_colloc); 232*269a910fSJames Wright 233*269a910fSJames Wright CeedBasisDestroy(&data->basis_x); 234*269a910fSJames Wright CeedBasisDestroy(&data->basis_stats); 235*269a910fSJames Wright 236*269a910fSJames Wright CeedVectorDestroy(&data->x_coord); 237*269a910fSJames Wright CeedVectorDestroy(&data->q_data); 238*269a910fSJames Wright 239*269a910fSJames Wright PetscCall(PetscFree(data)); 240*269a910fSJames Wright PetscFunctionReturn(0); 241*269a910fSJames Wright } 242*269a910fSJames Wright 243a175e481SJames Wright // Create PetscSF for child-to-parent communication 244*269a910fSJames Wright PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, SpanStatsSetupData stats_data, DM parentdm, DM childdm, PetscSF *statssf) { 24519706a06SJames Wright PetscInt child_num_qpnts, parent_num_qpnts, num_comp_x; 24619706a06SJames Wright CeedVector child_qx_coords, parent_qx_coords; 24719706a06SJames Wright PetscReal *child_coords, *parent_coords; 24819706a06SJames Wright PetscFunctionBeginUser; 24919706a06SJames Wright 250*269a910fSJames Wright PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)childdm), statssf)); 251*269a910fSJames Wright 25219706a06SJames Wright // Assume that child and parent have the same number of components 25319706a06SJames Wright CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x); 25419706a06SJames Wright const PetscInt num_comp_sf = num_comp_x - 1; // Number of coord components used in the creation of the SF 25519706a06SJames Wright 25619706a06SJames Wright // Get quad_coords for child DM 257a175e481SJames Wright PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &child_qx_coords, &child_num_qpnts)); 25819706a06SJames Wright 25919706a06SJames Wright // Get quad_coords for parent DM 260*269a910fSJames Wright PetscCall(GetQuadratureCoords(ceed, parentdm, stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord, &parent_qx_coords, 261*269a910fSJames Wright &parent_num_qpnts)); 26219706a06SJames Wright 26319706a06SJames Wright // Remove z component of coordinates for matching 26419706a06SJames Wright { 26519706a06SJames Wright const PetscReal *child_quad_coords, *parent_quad_coords; 26619706a06SJames Wright 26719706a06SJames Wright CeedVectorGetArrayRead(child_qx_coords, CEED_MEM_HOST, &child_quad_coords); 26819706a06SJames Wright CeedVectorGetArrayRead(parent_qx_coords, CEED_MEM_HOST, &parent_quad_coords); 26919706a06SJames Wright 27019706a06SJames Wright PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords)); 27119706a06SJames Wright for (int i = 0; i < child_num_qpnts; i++) { 27219706a06SJames Wright child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x]; 27319706a06SJames Wright child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x]; 27419706a06SJames Wright } 27519706a06SJames Wright for (int i = 0; i < parent_num_qpnts; i++) { 27619706a06SJames Wright parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x]; 27719706a06SJames Wright parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x]; 27819706a06SJames Wright } 27919706a06SJames Wright CeedVectorRestoreArrayRead(child_qx_coords, &child_quad_coords); 28019706a06SJames Wright CeedVectorRestoreArrayRead(parent_qx_coords, &parent_quad_coords); 28119706a06SJames Wright } 28219706a06SJames Wright 283*269a910fSJames Wright PetscCall(PetscSFSetGraphFromCoordinates(*statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords)); 28419706a06SJames Wright 285*269a910fSJames Wright PetscCall(PetscSFViewFromOptions(*statssf, NULL, "-spanstats_sf_view")); 286a175e481SJames Wright 28719706a06SJames Wright PetscCall(PetscFree2(child_coords, parent_coords)); 28819706a06SJames Wright CeedVectorDestroy(&child_qx_coords); 28919706a06SJames Wright CeedVectorDestroy(&parent_qx_coords); 29019706a06SJames Wright PetscFunctionReturn(0); 29119706a06SJames Wright } 29219706a06SJames Wright 293a175e481SJames Wright // Compute mass matrix for statistics projection 294*269a910fSJames Wright PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) { 295*269a910fSJames Wright CeedQFunction qf_mass, qf_stats_proj; 296*269a910fSJames Wright CeedOperator op_mass, op_setup_sur; 297*269a910fSJames Wright CeedInt q_data_size, num_comp_stats = user->spanstats.num_comp_stats; 298f967ad79SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)user->spanstats.dm); 299a175e481SJames Wright PetscFunctionBeginUser; 300a175e481SJames Wright 301*269a910fSJames Wright // Create Operator for L^2 projection of statistics 302*269a910fSJames Wright // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function. 303*269a910fSJames Wright // Therefore, an Identity QF is sufficient 304*269a910fSJames Wright CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &qf_stats_proj); 305*269a910fSJames Wright 306*269a910fSJames Wright CeedOperatorCreate(ceed, qf_stats_proj, NULL, NULL, &user->spanstats.op_stats_proj); 307*269a910fSJames Wright CeedOperatorSetField(user->spanstats.op_stats_proj, "input", stats_data->elem_restr_parent_colloc, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 308*269a910fSJames Wright CeedOperatorSetField(user->spanstats.op_stats_proj, "output", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE); 309*269a910fSJames Wright 310*269a910fSJames Wright // Get q_data for mass matrix operator 311*269a910fSJames Wright CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur); 312*269a910fSJames Wright CeedOperatorSetField(op_setup_sur, "dx", stats_data->elem_restr_parent_x, stats_data->basis_x, CEED_VECTOR_ACTIVE); 313*269a910fSJames Wright CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, stats_data->basis_x, CEED_VECTOR_NONE); 314*269a910fSJames Wright CeedOperatorSetField(op_setup_sur, "surface qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 315*269a910fSJames Wright CeedOperatorApply(op_setup_sur, stats_data->x_coord, stats_data->q_data, CEED_REQUEST_IMMEDIATE); 316*269a910fSJames Wright 317a175e481SJames Wright // CEED Restriction 318*269a910fSJames Wright CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size); 319a175e481SJames Wright 320a175e481SJames Wright // Create Mass CeedOperator 321*269a910fSJames Wright PetscCall(CreateMassQFunction(ceed, num_comp_stats, q_data_size, &qf_mass)); 322a175e481SJames Wright CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass); 323*269a910fSJames Wright CeedOperatorSetField(op_mass, "u", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE); 324*269a910fSJames Wright CeedOperatorSetField(op_mass, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_COLLOCATED, stats_data->q_data); 325*269a910fSJames Wright CeedOperatorSetField(op_mass, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE); 326a175e481SJames Wright 327f967ad79SJames Wright { // Setup KSP for L^2 projection 328a175e481SJames Wright MatopApplyContext M_ctx; 329a175e481SJames Wright PetscInt l_size, g_size; 330a175e481SJames Wright Mat mat_mass; 331a175e481SJames Wright VecType vec_type; 332a175e481SJames Wright KSP ksp; 333a175e481SJames Wright Vec ones, M_inv; 334a175e481SJames Wright CeedVector x_ceed, y_ceed; 335a175e481SJames Wright 336a175e481SJames Wright PetscCall(DMCreateGlobalVector(user->spanstats.dm, &M_inv)); 337a175e481SJames Wright PetscCall(VecGetLocalSize(M_inv, &l_size)); 338a175e481SJames Wright PetscCall(VecGetSize(M_inv, &g_size)); 339a175e481SJames Wright PetscCall(VecGetType(M_inv, &vec_type)); 340a175e481SJames Wright 341*269a910fSJames Wright CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &x_ceed, NULL); 342*269a910fSJames Wright CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &y_ceed, NULL); 3436665b873SJames Wright PetscCall(MatopApplyContextCreate(user->spanstats.dm, user->ceed, op_mass, x_ceed, y_ceed, NULL, &M_ctx)); 3446665b873SJames Wright CeedVectorDestroy(&x_ceed); 3456665b873SJames Wright CeedVectorDestroy(&y_ceed); 3466665b873SJames Wright 347f967ad79SJames Wright PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, M_ctx, &mat_mass)); 3486665b873SJames Wright PetscCall(MatShellSetContextDestroy(mat_mass, (PetscErrorCode(*)(void *))MatopApplyContextDestroy)); 349a175e481SJames Wright PetscCall(MatShellSetOperation(mat_mass, MATOP_MULT, (void (*)(void))MatMult_Ceed)); 350a175e481SJames Wright PetscCall(MatShellSetOperation(mat_mass, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed)); 351a175e481SJames Wright PetscCall(MatShellSetVecType(mat_mass, vec_type)); 352a175e481SJames Wright 353a175e481SJames Wright // Create lumped mass matrix inverse 354a175e481SJames Wright PetscCall(DMGetGlobalVector(user->spanstats.dm, &ones)); 355a175e481SJames Wright PetscCall(VecZeroEntries(M_inv)); 356a175e481SJames Wright PetscCall(VecSet(ones, 1)); 357a175e481SJames Wright PetscCall(MatMult(mat_mass, ones, M_inv)); 358a175e481SJames Wright PetscCall(VecReciprocal(M_inv)); 359a175e481SJames Wright user->spanstats.M_inv = M_inv; 360a175e481SJames Wright PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &ones)); 361a175e481SJames Wright 362f967ad79SJames Wright PetscCall(KSPCreate(comm, &ksp)); 363b7d66439SJames Wright PetscCall(KSPSetOptionsPrefix(ksp, "turbulence_spanstats_")); 364a175e481SJames Wright { 365a175e481SJames Wright PC pc; 366a175e481SJames Wright PetscCall(KSPGetPC(ksp, &pc)); 367a175e481SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 368a175e481SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL)); 369a175e481SJames Wright PetscCall(KSPSetType(ksp, KSPCG)); 370a175e481SJames Wright PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL)); 371a175e481SJames Wright PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 372a175e481SJames Wright } 373a175e481SJames Wright PetscCall(KSPSetOperators(ksp, mat_mass, mat_mass)); 374a175e481SJames Wright PetscCall(KSPSetFromOptions(ksp)); 375a175e481SJames Wright user->spanstats.ksp = ksp; 376a175e481SJames Wright } 377a175e481SJames Wright 378a175e481SJames Wright // Cleanup 379a175e481SJames Wright CeedQFunctionDestroy(&qf_mass); 380*269a910fSJames Wright CeedQFunctionDestroy(&qf_stats_proj); 3816665b873SJames Wright CeedOperatorDestroy(&op_mass); 382*269a910fSJames Wright CeedOperatorDestroy(&op_setup_sur); 383a175e481SJames Wright PetscFunctionReturn(0); 384a175e481SJames Wright } 385a175e481SJames Wright 386*269a910fSJames Wright // Create CeedOperator for statistics collection 387*269a910fSJames Wright PetscErrorCode CreateStatisticCollectionOperator(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data, ProblemData *problem) { 388a175e481SJames Wright CeedInt num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q; 389495b9769SJames Wright Turbulence_SpanStatsContext collect_ctx; 390495b9769SJames Wright NewtonianIdealGasContext newtonian_ig_ctx; 391495b9769SJames Wright CeedQFunctionContext collect_context; 392*269a910fSJames Wright CeedQFunction qf_stats_collect; 393a175e481SJames Wright PetscFunctionBeginUser; 394a175e481SJames Wright CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q); 395a175e481SJames Wright 396a175e481SJames Wright // Create Operator for statistics collection 397a175e481SJames Wright switch (user->phys->state_var) { 398a175e481SJames Wright case STATEVAR_PRIMITIVE: 399*269a910fSJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &qf_stats_collect); 400a175e481SJames Wright break; 401a175e481SJames Wright case STATEVAR_CONSERVATIVE: 402*269a910fSJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &qf_stats_collect); 403a175e481SJames Wright break; 404*269a910fSJames Wright default: 405*269a910fSJames Wright SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No statisics collection available for chosen state variable"); 406a175e481SJames Wright } 407a175e481SJames Wright 408b7d66439SJames Wright if (user->spanstats.do_mms_test) { 409*269a910fSJames Wright CeedQFunctionDestroy(&qf_stats_collect); 410*269a910fSJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest, ChildStatsCollectionMMSTest_loc, &qf_stats_collect); 411a175e481SJames Wright } 412a175e481SJames Wright 413*269a910fSJames Wright { // Setup Collection Context 414495b9769SJames Wright PetscCall(PetscNew(&collect_ctx)); 415495b9769SJames Wright CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx); 416495b9769SJames Wright collect_ctx->gas = *newtonian_ig_ctx; 417495b9769SJames Wright 418495b9769SJames Wright CeedQFunctionContextCreate(user->ceed, &collect_context); 419495b9769SJames Wright CeedQFunctionContextSetData(collect_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx); 420495b9769SJames Wright CeedQFunctionContextSetDataDestroy(collect_context, CEED_MEM_HOST, FreeContextPetsc); 421495b9769SJames Wright 422495b9769SJames Wright CeedQFunctionContextRegisterDouble(collect_context, "solution time", offsetof(struct Turbulence_SpanStatsContext_, solution_time), 1, 423495b9769SJames Wright "Current solution time"); 424495b9769SJames Wright CeedQFunctionContextRegisterDouble(collect_context, "previous time", offsetof(struct Turbulence_SpanStatsContext_, previous_time), 1, 425495b9769SJames Wright "Previous time statistics collection was done"); 426495b9769SJames Wright 427495b9769SJames Wright CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx); 428495b9769SJames Wright } 429495b9769SJames Wright 430*269a910fSJames Wright CeedQFunctionSetContext(qf_stats_collect, collect_context); 431495b9769SJames Wright CeedQFunctionContextDestroy(&collect_context); 432*269a910fSJames Wright CeedQFunctionAddInput(qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP); 433*269a910fSJames Wright CeedQFunctionAddInput(qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE); 434*269a910fSJames Wright CeedQFunctionAddInput(qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP); 435*269a910fSJames Wright CeedQFunctionAddOutput(qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE); 436a175e481SJames Wright 437*269a910fSJames Wright CeedOperatorCreate(ceed, qf_stats_collect, NULL, NULL, &user->spanstats.op_stats_collect); 438a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 439a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 440a175e481SJames Wright CeedOperatorSetField(user->spanstats.op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord); 441*269a910fSJames Wright CeedOperatorSetField(user->spanstats.op_stats_collect, "v", stats_data->elem_restr_child_colloc, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 442a175e481SJames Wright 443495b9769SJames Wright CeedOperatorContextGetFieldLabel(user->spanstats.op_stats_collect, "solution time", &user->spanstats.solution_time_label); 444495b9769SJames Wright CeedOperatorContextGetFieldLabel(user->spanstats.op_stats_collect, "previous time", &user->spanstats.previous_time_label); 445495b9769SJames Wright 446*269a910fSJames Wright CeedQFunctionDestroy(&qf_stats_collect); 447a175e481SJames Wright PetscFunctionReturn(0); 448a175e481SJames Wright } 449a175e481SJames Wright 450b7d66439SJames Wright // Creates operator for calculating error of method of manufactured solution (MMS) test 451*269a910fSJames Wright PetscErrorCode SetupMMSErrorChecking(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) { 452*269a910fSJames Wright CeedInt num_comp_stats = user->spanstats.num_comp_stats, num_comp_x, q_data_size; 453823a1283SJames Wright CeedQFunction qf_error; 454823a1283SJames Wright CeedOperator op_error; 455823a1283SJames Wright CeedVector x_ceed, y_ceed; 456823a1283SJames Wright PetscFunctionBeginUser; 457823a1283SJames Wright 458*269a910fSJames Wright CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size); 459*269a910fSJames Wright CeedBasisGetNumComponents(stats_data->basis_x, &num_comp_x); 460823a1283SJames Wright 461b7d66439SJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest_Error, ChildStatsCollectionMMSTest_Error_loc, &qf_error); 462823a1283SJames Wright CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP); 463823a1283SJames Wright CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE); 464823a1283SJames Wright CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP); 465823a1283SJames Wright CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP); 466823a1283SJames Wright 467823a1283SJames Wright CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error); 468*269a910fSJames Wright CeedOperatorSetField(op_error, "q", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE); 469*269a910fSJames Wright CeedOperatorSetField(op_error, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_COLLOCATED, stats_data->q_data); 470*269a910fSJames Wright CeedOperatorSetField(op_error, "x", stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord); 471*269a910fSJames Wright CeedOperatorSetField(op_error, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE); 472823a1283SJames Wright 473*269a910fSJames Wright CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &x_ceed, NULL); 474*269a910fSJames Wright CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &y_ceed, NULL); 4756665b873SJames Wright PetscCall(MatopApplyContextCreate(user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, &user->spanstats.mms_error_ctx)); 476823a1283SJames Wright 4776665b873SJames Wright CeedOperatorDestroy(&op_error); 478*269a910fSJames Wright CeedQFunctionDestroy(&qf_error); 4796665b873SJames Wright CeedVectorDestroy(&x_ceed); 4806665b873SJames Wright CeedVectorDestroy(&y_ceed); 481823a1283SJames Wright PetscFunctionReturn(0); 482823a1283SJames Wright } 483823a1283SJames Wright 484a175e481SJames Wright // Setup for statistics collection 48519706a06SJames Wright PetscErrorCode SetupStatsCollection(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) { 486*269a910fSJames Wright SpanStatsSetupData stats_data; 4874eed8630SJames Wright PetscLogStage stage_stats_setup; 48819706a06SJames Wright PetscFunctionBeginUser; 4894eed8630SJames Wright PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup)); 4904eed8630SJames Wright if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup)); 4914eed8630SJames Wright PetscCall(PetscLogStagePush(stage_stats_setup)); 4924eed8630SJames Wright 493*269a910fSJames Wright // Create necessary CeedObjects for setting up statistics 494*269a910fSJames Wright PetscCall(SpanStatsSetupDataCreate(ceed, user, ceed_data, problem, &stats_data)); 495*269a910fSJames Wright CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &user->spanstats.rhs_ceed, NULL); 496*269a910fSJames Wright CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_colloc, &user->spanstats.parent_stats, NULL); 497*269a910fSJames Wright CeedElemRestrictionCreateVector(stats_data->elem_restr_child_colloc, &user->spanstats.child_stats, NULL); 498a175e481SJames Wright CeedVectorSetValue(user->spanstats.child_stats, 0); 4991737222fSJames Wright 50019706a06SJames Wright // Create SF for communicating child data back their respective parents 501*269a910fSJames Wright PetscCall(CreateStatsSF(ceed, ceed_data, stats_data, user->dm, user->spanstats.dm, &user->spanstats.sf)); 50251ee423eSJames Wright 503a175e481SJames Wright // Create CeedOperators for statistics collection 504*269a910fSJames Wright PetscCall(CreateStatisticCollectionOperator(ceed, user, ceed_data, stats_data, problem)); 505a175e481SJames Wright 506a175e481SJames Wright // Setup KSP and Mat for L^2 projection of statistics 507*269a910fSJames Wright PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data, stats_data)); 508a175e481SJames Wright 509b7d66439SJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_turbulence_spanstats_mms", &user->spanstats.do_mms_test, NULL)); 510b7d66439SJames Wright if (user->spanstats.do_mms_test) { 511*269a910fSJames Wright PetscCall(SetupMMSErrorChecking(ceed, user, ceed_data, stats_data)); 512823a1283SJames Wright } 513823a1283SJames Wright 514f967ad79SJames Wright { // Setup stats viewer with prefix 5158ed52730SJames Wright PetscViewerType viewer_type; 516b7d66439SJames Wright PetscCall(PetscViewerGetType(user->app_ctx->turb_spanstats_viewer, &viewer_type)); 517b7d66439SJames Wright PetscCall(PetscOptionsSetValue(NULL, "-ts_monitor_turbulence_spanstats_viewer_type", viewer_type)); 5188ed52730SJames Wright 519b7d66439SJames Wright PetscCall(PetscViewerSetOptionsPrefix(user->app_ctx->turb_spanstats_viewer, "ts_monitor_turbulence_spanstats_")); 520b7d66439SJames Wright PetscCall(PetscViewerSetFromOptions(user->app_ctx->turb_spanstats_viewer)); 5218ed52730SJames Wright } 5224eed8630SJames Wright 523*269a910fSJames Wright PetscCall(SpanStatsSetupDataDestroy(stats_data)); 5244eed8630SJames Wright PetscCall(PetscLogStagePop()); 525a175e481SJames Wright PetscFunctionReturn(0); 526a175e481SJames Wright } 527a175e481SJames Wright 528a175e481SJames Wright // Collect statistics based on the solution Q 529a175e481SJames Wright PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) { 530a175e481SJames Wright PetscMemType q_mem_type; 531a175e481SJames Wright const PetscScalar *q_arr; 532a175e481SJames Wright PetscFunctionBeginUser; 533a175e481SJames Wright 5346dcea3beSJames Wright PetscLogStage stage_stats_collect; 5356dcea3beSJames Wright PetscCall(PetscLogStageGetId("Stats Collect", &stage_stats_collect)); 5366dcea3beSJames Wright if (stage_stats_collect == -1) PetscCall(PetscLogStageRegister("Stats Collect", &stage_stats_collect)); 5376dcea3beSJames Wright PetscCall(PetscLogStagePush(stage_stats_collect)); 5386dcea3beSJames Wright 539a0b9a424SJames Wright PetscCall(UpdateBoundaryValues(user, user->Q_loc, solution_time)); 540495b9769SJames Wright CeedOperatorContextSetDouble(user->spanstats.op_stats_collect, user->spanstats.solution_time_label, &solution_time); 541a0b9a424SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, user->Q_loc)); 542a0b9a424SJames Wright PetscCall(VecGetArrayReadAndMemType(user->Q_loc, &q_arr, &q_mem_type)); 543a175e481SJames Wright CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q_arr); 544a175e481SJames Wright 545495b9769SJames Wright CeedOperatorApplyAdd(user->spanstats.op_stats_collect, user->q_ceed, user->spanstats.child_stats, CEED_REQUEST_IMMEDIATE); 546a175e481SJames Wright 547a175e481SJames Wright CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL); 548a0b9a424SJames Wright PetscCall(VecRestoreArrayReadAndMemType(user->Q_loc, &q_arr)); 549a175e481SJames Wright 550495b9769SJames Wright CeedOperatorContextSetDouble(user->spanstats.op_stats_collect, user->spanstats.previous_time_label, &solution_time); 551a175e481SJames Wright 5526dcea3beSJames Wright PetscCall(PetscLogStagePop()); 553a175e481SJames Wright PetscFunctionReturn(0); 554a175e481SJames Wright } 555a175e481SJames Wright 556a175e481SJames Wright // Process the child statistics into parent statistics and project them onto stats 55778bbfb6fSJed Brown PetscErrorCode ProcessStatistics(User user, Vec stats) { 558a175e481SJames Wright Span_Stats user_stats = user->spanstats; 559a175e481SJames Wright const PetscScalar *child_stats; 560a175e481SJames Wright PetscScalar *parent_stats; 561a175e481SJames Wright MPI_Datatype unit; 562a175e481SJames Wright Vec rhs_loc, rhs; 563a175e481SJames Wright PetscMemType rhs_mem_type; 564a175e481SJames Wright CeedScalar *rhs_arr; 565a175e481SJames Wright CeedMemType ceed_mem_type; 566a175e481SJames Wright PetscFunctionBeginUser; 567a175e481SJames Wright 5686dcea3beSJames Wright PetscLogStage stage_stats_process; 5696dcea3beSJames Wright PetscCall(PetscLogStageGetId("Stats Process", &stage_stats_process)); 5706dcea3beSJames Wright if (stage_stats_process == -1) PetscCall(PetscLogStageRegister("Stats Process", &stage_stats_process)); 5716dcea3beSJames Wright PetscCall(PetscLogStagePush(stage_stats_process)); 5726dcea3beSJames Wright 573a175e481SJames Wright CeedGetPreferredMemType(user->ceed, &ceed_mem_type); 574a175e481SJames Wright CeedVectorSetValue(user_stats.parent_stats, 0); 575a175e481SJames Wright 576a175e481SJames Wright CeedVectorGetArrayRead(user_stats.child_stats, ceed_mem_type, &child_stats); 577a175e481SJames Wright CeedVectorGetArray(user_stats.parent_stats, ceed_mem_type, &parent_stats); 578a175e481SJames Wright 579a175e481SJames Wright if (user_stats.num_comp_stats == 1) unit = MPIU_REAL; 580a175e481SJames Wright else { 581a175e481SJames Wright PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit)); 582a175e481SJames Wright PetscCallMPI(MPI_Type_commit(&unit)); 583a175e481SJames Wright } 584a175e481SJames Wright 585a175e481SJames Wright PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM)); 586a175e481SJames Wright PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM)); 587a175e481SJames Wright 588a175e481SJames Wright CeedVectorRestoreArrayRead(user_stats.child_stats, &child_stats); 589a175e481SJames Wright CeedVectorRestoreArray(user_stats.parent_stats, &parent_stats); 590a175e481SJames Wright PetscCallMPI(MPI_Type_free(&unit)); 591a175e481SJames Wright 59278bbfb6fSJed Brown PetscReal solution_time; 59378bbfb6fSJed Brown PetscCall(DMGetOutputSequenceNumber(user_stats.dm, NULL, &solution_time)); 59478bbfb6fSJed Brown PetscReal summing_duration = solution_time - user->app_ctx->cont_time; 59578bbfb6fSJed Brown CeedVectorScale(user_stats.parent_stats, 1 / (summing_duration * user_stats.span_width)); 596a175e481SJames Wright 597a175e481SJames Wright // L^2 projection with the parent_data 598a175e481SJames Wright PetscCall(DMGetLocalVector(user_stats.dm, &rhs_loc)); 599a175e481SJames Wright PetscCall(VecZeroEntries(rhs_loc)); 600a175e481SJames Wright PetscCall(VecGetArrayWriteAndMemType(rhs_loc, &rhs_arr, &rhs_mem_type)); 601a175e481SJames Wright CeedVectorSetArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), CEED_USE_POINTER, (PetscScalar *)rhs_arr); 602a175e481SJames Wright 603a175e481SJames Wright CeedOperatorApply(user_stats.op_stats_proj, user_stats.parent_stats, user_stats.rhs_ceed, CEED_REQUEST_IMMEDIATE); 604a175e481SJames Wright 605a175e481SJames Wright CeedVectorTakeArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), &rhs_arr); 606a175e481SJames Wright PetscCall(VecRestoreArrayAndMemType(rhs_loc, &rhs_arr)); 60778bbfb6fSJed Brown 60878bbfb6fSJed Brown PetscCall(DMGetGlobalVector(user_stats.dm, &rhs)); 60978bbfb6fSJed Brown PetscCall(VecZeroEntries(rhs)); 610a175e481SJames Wright PetscCall(DMLocalToGlobal(user_stats.dm, rhs_loc, ADD_VALUES, rhs)); 61178bbfb6fSJed Brown PetscCall(DMRestoreLocalVector(user_stats.dm, &rhs_loc)); 612a175e481SJames Wright 61378bbfb6fSJed Brown PetscCall(KSPSolve(user_stats.ksp, rhs, stats)); 614a175e481SJames Wright 61578bbfb6fSJed Brown PetscCall(DMRestoreGlobalVector(user_stats.dm, &rhs)); 6166dcea3beSJames Wright PetscCall(PetscLogStagePop()); 617a175e481SJames Wright PetscFunctionReturn(0); 618a175e481SJames Wright } 619a175e481SJames Wright 620a175e481SJames Wright // TSMonitor for the statistics collection and processing 621a175e481SJames Wright PetscErrorCode TSMonitor_Statistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) { 622a175e481SJames Wright User user = (User)ctx; 623a175e481SJames Wright Vec stats; 62478bbfb6fSJed Brown TSConvergedReason reason; 625b7d66439SJames Wright PetscInt collect_interval = user->app_ctx->turb_spanstats_collect_interval, viewer_interval = user->app_ctx->turb_spanstats_viewer_interval; 626a175e481SJames Wright PetscFunctionBeginUser; 62778bbfb6fSJed Brown PetscCall(TSGetConvergedReason(ts, &reason)); 628a175e481SJames Wright // Do not collect or process on the first step of the run (ie. on the initial condition) 62978bbfb6fSJed Brown if (steps == user->app_ctx->cont_steps && reason == TS_CONVERGED_ITERATING) PetscFunctionReturn(0); 630a175e481SJames Wright 63127240365SJames Wright PetscBool run_processing_and_viewer = (steps % viewer_interval == 0 && viewer_interval != -1) || reason != TS_CONVERGED_ITERATING; 632a175e481SJames Wright 63327240365SJames Wright if (steps % collect_interval == 0 || run_processing_and_viewer) { 63427240365SJames Wright PetscCall(CollectStatistics(user, solution_time, Q)); 63527240365SJames Wright 63627240365SJames Wright if (run_processing_and_viewer) { 6378ed52730SJames Wright PetscCall(DMSetOutputSequenceNumber(user->spanstats.dm, steps, solution_time)); 63878bbfb6fSJed Brown PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats)); 63978bbfb6fSJed Brown PetscCall(ProcessStatistics(user, stats)); 6408fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_NONE) { 641b7d66439SJames Wright PetscCall(PetscViewerPushFormat(user->app_ctx->turb_spanstats_viewer, user->app_ctx->turb_spanstats_viewer_format)); 642b7d66439SJames Wright PetscCall(VecView(stats, user->app_ctx->turb_spanstats_viewer)); 643b7d66439SJames Wright PetscCall(PetscViewerPopFormat(user->app_ctx->turb_spanstats_viewer)); 6448fb33541SJames Wright } 6458fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_TURB_SPANSTATS && reason != TS_CONVERGED_ITERATING) { 6468fb33541SJames Wright PetscCall(RegressionTests_NS(user->app_ctx, stats)); 6478fb33541SJames Wright } 648b7d66439SJames Wright if (user->spanstats.do_mms_test && reason != TS_CONVERGED_ITERATING) { 649823a1283SJames Wright Vec error; 650823a1283SJames Wright PetscCall(VecDuplicate(stats, &error)); 651b7d66439SJames Wright PetscCall(ApplyLocal_Ceed(stats, error, user->spanstats.mms_error_ctx)); 652823a1283SJames Wright PetscScalar error_sq = 0; 653823a1283SJames Wright PetscCall(VecSum(error, &error_sq)); 654823a1283SJames Wright PetscScalar l2_error = sqrt(error_sq); 655823a1283SJames Wright PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error)); 656823a1283SJames Wright } 657a175e481SJames Wright PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats)); 658522ee345SJames Wright } 65927240365SJames Wright } 66051ee423eSJames Wright PetscFunctionReturn(0); 66151ee423eSJames Wright } 662fd170fd0SJames Wright 6636665b873SJames Wright PetscErrorCode DestroyStats(User user, CeedData ceed_data) { 664fd170fd0SJames Wright PetscFunctionBeginUser; 665fd170fd0SJames Wright 666fd170fd0SJames Wright // -- CeedVectors 667fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.child_stats); 668fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.parent_stats); 669fd170fd0SJames Wright CeedVectorDestroy(&user->spanstats.rhs_ceed); 670fd170fd0SJames Wright 671fd170fd0SJames Wright // -- CeedOperators 672fd170fd0SJames Wright CeedOperatorDestroy(&user->spanstats.op_stats_collect); 673fd170fd0SJames Wright CeedOperatorDestroy(&user->spanstats.op_stats_proj); 6746665b873SJames Wright PetscCall(MatopApplyContextDestroy(user->spanstats.mms_error_ctx)); 675fd170fd0SJames Wright 676fd170fd0SJames Wright // -- Vec 677fd170fd0SJames Wright PetscCall(VecDestroy(&user->spanstats.M_inv)); 678fd170fd0SJames Wright 679fd170fd0SJames Wright // -- KSP 680fd170fd0SJames Wright PetscCall(KSPDestroy(&user->spanstats.ksp)); 681fd170fd0SJames Wright 682fd170fd0SJames Wright // -- SF 683fd170fd0SJames Wright PetscCall(PetscSFDestroy(&user->spanstats.sf)); 684fd170fd0SJames Wright 685fd170fd0SJames Wright // -- DM 686fd170fd0SJames Wright PetscCall(DMDestroy(&user->spanstats.dm)); 687fd170fd0SJames Wright 688fd170fd0SJames Wright PetscFunctionReturn(0); 689fd170fd0SJames Wright } 690