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 851ee423eSJames Wright /// @file 951ee423eSJames Wright /// Utility functions for setting up statistics collection 1051ee423eSJames Wright 1119706a06SJames Wright #include <petscsf.h> 1219706a06SJames Wright 1351ee423eSJames Wright #include "../navierstokes.h" 1419706a06SJames Wright #include "petscsys.h" 1551ee423eSJames Wright 1651ee423eSJames Wright PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree, SimpleBC bc) { 1751ee423eSJames Wright user->spanstats.num_comp_stats = 1; 1851ee423eSJames Wright PetscFunctionBeginUser; 1951ee423eSJames Wright 2051ee423eSJames Wright // Get DM from surface 2151ee423eSJames Wright { 2251ee423eSJames Wright DMLabel label; 2351ee423eSJames Wright PetscCall(DMGetLabel(user->dm, "Face Sets", &label)); 2451ee423eSJames Wright PetscCall(DMPlexLabelComplete(user->dm, label)); 2551ee423eSJames Wright PetscCall(DMPlexFilter(user->dm, label, 1, &user->spanstats.dm)); 2651ee423eSJames Wright PetscCall(DMProjectCoordinates(user->spanstats.dm, NULL)); // Ensure that a coordinate FE exists 2751ee423eSJames Wright } 2851ee423eSJames Wright 2951ee423eSJames Wright PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats")); 3051ee423eSJames Wright PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "spanstats_")); 3151ee423eSJames Wright PetscCall(PetscOptionsSetValue(NULL, "-spanstats_dm_sparse_localize", "0")); // [Jed] Not relevant because not periodic in this direction 3251ee423eSJames Wright 3351ee423eSJames Wright PetscCall(DMSetFromOptions(user->spanstats.dm)); 3451ee423eSJames Wright PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view")); // -spanstats_dm_view (option includes prefix) 3551ee423eSJames Wright { 3651ee423eSJames Wright PetscFE fe; 3751ee423eSJames Wright DMLabel label; 3851ee423eSJames Wright 3951ee423eSJames Wright PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim - 1, user->spanstats.num_comp_stats, PETSC_FALSE, degree, PETSC_DECIDE, &fe)); 4051ee423eSJames Wright PetscCall(PetscObjectSetName((PetscObject)fe, "stats")); 4151ee423eSJames Wright PetscCall(DMAddField(user->spanstats.dm, NULL, (PetscObject)fe)); 4251ee423eSJames Wright PetscCall(DMCreateDS(user->spanstats.dm)); 4351ee423eSJames Wright PetscCall(DMGetLabel(user->spanstats.dm, "Face Sets", &label)); 4451ee423eSJames Wright 4551ee423eSJames Wright // // Set wall BCs 4651ee423eSJames Wright // if (bc->num_wall > 0) { 4751ee423eSJames Wright // PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, bc->num_wall, bc->walls, 0, bc->num_comps, bc->wall_comps, 4851ee423eSJames Wright // (void (*)(void))problem->bc, NULL, problem->bc_ctx, NULL)); 4951ee423eSJames Wright // } 5051ee423eSJames Wright // // Set slip BCs in the x direction 5151ee423eSJames Wright // if (bc->num_slip[0] > 0) { 5251ee423eSJames Wright // PetscInt comps[1] = {1}; 5351ee423eSJames Wright // PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label, bc->num_slip[0], bc->slips[0], 0, 1, comps, (void (*)(void))NULL, NULL, 5451ee423eSJames Wright // problem->bc_ctx, NULL)); 5551ee423eSJames Wright // } 5651ee423eSJames Wright // // Set slip BCs in the y direction 5751ee423eSJames Wright // if (bc->num_slip[1] > 0) { 5851ee423eSJames Wright // PetscInt comps[1] = {2}; 5951ee423eSJames Wright // PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label, bc->num_slip[1], bc->slips[1], 0, 1, comps, (void (*)(void))NULL, NULL, 6051ee423eSJames Wright // problem->bc_ctx, NULL)); 6151ee423eSJames Wright // } 6251ee423eSJames Wright // // Set slip BCs in the z direction 6351ee423eSJames Wright // if (bc->num_slip[2] > 0) { 6451ee423eSJames Wright // PetscInt comps[1] = {3}; 6551ee423eSJames Wright // PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label, bc->num_slip[2], bc->slips[2], 0, 1, comps, (void (*)(void))NULL, NULL, 6651ee423eSJames Wright // problem->bc_ctx, NULL)); 6751ee423eSJames Wright // } 6851ee423eSJames Wright 6951ee423eSJames Wright PetscCall(DMPlexSetClosurePermutationTensor(user->spanstats.dm, PETSC_DETERMINE, NULL)); 7051ee423eSJames Wright PetscCall(PetscFEDestroy(&fe)); 7151ee423eSJames Wright } 7251ee423eSJames Wright 7351ee423eSJames Wright PetscSection section; 7451ee423eSJames Wright PetscCall(DMGetLocalSection(user->spanstats.dm, §ion)); 7551ee423eSJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 7651ee423eSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "Test")); 7751ee423eSJames Wright // PetscCall(PetscSectionSetComponentName(section, 0, 0, "Mean Velocity Products XX")); 7851ee423eSJames Wright // PetscCall(PetscSectionSetComponentName(section, 0, 1, "Mean Velocity Products YY")); 7951ee423eSJames Wright // PetscCall(PetscSectionSetComponentName(section, 0, 2, "Mean Velocity Products ZZ")); 8051ee423eSJames Wright // PetscCall(PetscSectionSetComponentName(section, 0, 3, "Mean Velocity Products YZ")); 8151ee423eSJames Wright // PetscCall(PetscSectionSetComponentName(section, 0, 4, "Mean Velocity Products XZ")); 8251ee423eSJames Wright // PetscCall(PetscSectionSetComponentName(section, 0, 5, "Mean Velocity Products XY")); 8351ee423eSJames Wright 8419706a06SJames Wright Vec test; 8519706a06SJames Wright PetscCall(DMCreateLocalVector(user->spanstats.dm, &test)); 8619706a06SJames Wright PetscCall(VecZeroEntries(test)); 8719706a06SJames Wright PetscCall(VecViewFromOptions(test, NULL, "-test_view")); 8819706a06SJames Wright 8919706a06SJames Wright PetscFunctionReturn(0); 9019706a06SJames Wright } 9119706a06SJames Wright 92*1737222fSJames Wright // Create CeedElemRestriction for collocated data based on associated CeedBasis and CeedElemRestriction 93*1737222fSJames Wright // Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction 94*1737222fSJames Wright PetscErrorCode CreateElemRestrColloc(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base, 95*1737222fSJames Wright CeedElemRestriction *elem_restr_collocated, CeedVector *l_vec, CeedVector *e_vec) { 96*1737222fSJames Wright CeedInt num_elem_qpts, loc_num_elem; 97*1737222fSJames Wright PetscFunctionBeginUser; 98*1737222fSJames Wright 99*1737222fSJames Wright CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts); 100*1737222fSJames Wright CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem); 101*1737222fSJames Wright 102*1737222fSJames Wright const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp}; 103*1737222fSJames Wright CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides, 104*1737222fSJames Wright elem_restr_collocated); 105*1737222fSJames Wright CeedElemRestrictionCreateVector(*elem_restr_collocated, l_vec, e_vec); 106*1737222fSJames Wright PetscFunctionReturn(0); 107*1737222fSJames Wright } 108*1737222fSJames Wright 10919706a06SJames Wright PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, CeedVector *qx_coords, 11019706a06SJames Wright PetscInt *total_nqpnts) { 11119706a06SJames Wright CeedQFunction qf_quad_coords; 11219706a06SJames Wright CeedOperator op_quad_coords; 11319706a06SJames Wright PetscInt num_comp_x, loc_num_elem, num_elem_qpts; 11419706a06SJames Wright CeedElemRestriction elem_restr_qx; 11519706a06SJames Wright PetscFunctionBeginUser; 11619706a06SJames Wright 11719706a06SJames Wright // Create Element Restriction and CeedVector for quadrature coordinates 11819706a06SJames Wright CeedBasisGetNumQuadraturePoints(basis_x, &num_elem_qpts); 11919706a06SJames Wright CeedElemRestrictionGetNumElements(elem_restr_x, &loc_num_elem); 12019706a06SJames Wright CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x); 12119706a06SJames Wright *total_nqpnts = num_elem_qpts * loc_num_elem; 122*1737222fSJames Wright PetscCall(CreateElemRestrColloc(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx, qx_coords, NULL)); 12319706a06SJames Wright 12419706a06SJames Wright // Create QFunction 12519706a06SJames Wright CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords); 12619706a06SJames Wright 12719706a06SJames Wright // Create Operator 12819706a06SJames Wright CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords); 12919706a06SJames Wright CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE); 13019706a06SJames Wright CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 13119706a06SJames Wright 13219706a06SJames Wright CeedOperatorApply(op_quad_coords, x_coords, *qx_coords, CEED_REQUEST_IMMEDIATE); 13319706a06SJames Wright 13419706a06SJames Wright CeedQFunctionDestroy(&qf_quad_coords); 13519706a06SJames Wright CeedOperatorDestroy(&op_quad_coords); 13619706a06SJames Wright PetscFunctionReturn(0); 13719706a06SJames Wright } 13819706a06SJames Wright 13919706a06SJames Wright PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, DM parentdm, DM childdm, PetscSF statssf) { 14019706a06SJames Wright PetscInt child_num_qpnts, parent_num_qpnts, num_comp_x; 14119706a06SJames Wright CeedVector child_qx_coords, parent_qx_coords; 14219706a06SJames Wright PetscReal *child_coords, *parent_coords; 14319706a06SJames Wright PetscFunctionBeginUser; 14419706a06SJames Wright 14519706a06SJames Wright // Assume that child and parent have the same number of components 14619706a06SJames Wright CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x); 14719706a06SJames Wright const PetscInt num_comp_sf = num_comp_x - 1; // Number of coord components used in the creation of the SF 14819706a06SJames Wright 14919706a06SJames Wright // Get quad_coords for child DM 15019706a06SJames Wright PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_xc, ceed_data->x_coord, &child_qx_coords, &child_num_qpnts)); 15119706a06SJames Wright 15219706a06SJames Wright // Get quad_coords for parent DM 15319706a06SJames Wright PetscCall(GetQuadratureCoords(ceed, parentdm, ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, ceed_data->spanstats.x_coord, 15419706a06SJames Wright &parent_qx_coords, &parent_num_qpnts)); 15519706a06SJames Wright 15619706a06SJames Wright // Remove z component of coordinates for matching 15719706a06SJames Wright { 15819706a06SJames Wright const PetscReal *child_quad_coords, *parent_quad_coords; 15919706a06SJames Wright 16019706a06SJames Wright CeedVectorGetArrayRead(child_qx_coords, CEED_MEM_HOST, &child_quad_coords); 16119706a06SJames Wright CeedVectorGetArrayRead(parent_qx_coords, CEED_MEM_HOST, &parent_quad_coords); 16219706a06SJames Wright 16319706a06SJames Wright PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords)); 16419706a06SJames Wright for (int i = 0; i < child_num_qpnts; i++) { 16519706a06SJames Wright child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x]; 16619706a06SJames Wright child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x]; 16719706a06SJames Wright } 16819706a06SJames Wright for (int i = 0; i < parent_num_qpnts; i++) { 16919706a06SJames Wright parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x]; 17019706a06SJames Wright parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x]; 17119706a06SJames Wright } 17219706a06SJames Wright CeedVectorRestoreArrayRead(child_qx_coords, &child_quad_coords); 17319706a06SJames Wright CeedVectorRestoreArrayRead(parent_qx_coords, &parent_quad_coords); 17419706a06SJames Wright } 17519706a06SJames Wright 17619706a06SJames Wright // Only check the first two components of the coordinates 17719706a06SJames Wright PetscCall(PetscSFSetGraphFromCoordinates(statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords)); 17819706a06SJames Wright 17919706a06SJames Wright PetscCall(PetscFree2(child_coords, parent_coords)); 18019706a06SJames Wright CeedVectorDestroy(&ceed_data->spanstats.x_coord); 18119706a06SJames Wright CeedVectorDestroy(&child_qx_coords); 18219706a06SJames Wright CeedVectorDestroy(&parent_qx_coords); 18319706a06SJames Wright PetscFunctionReturn(0); 18419706a06SJames Wright } 18519706a06SJames Wright 18619706a06SJames Wright PetscErrorCode SetupStatsCollection(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) { 18719706a06SJames Wright DM dm = user->spanstats.dm; 18819706a06SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm); 18919706a06SJames Wright CeedInt dim, P, Q, num_comp_x; 19019706a06SJames Wright Vec X_loc; 19119706a06SJames Wright PetscMemType X_loc_memtype; 19219706a06SJames Wright const PetscScalar *X_loc_array; 19319706a06SJames Wright PetscFunctionBeginUser; 19419706a06SJames Wright 19519706a06SJames Wright PetscCall(DMGetDimension(dm, &dim)); 19619706a06SJames Wright CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_x, &Q); 19719706a06SJames Wright CeedBasisGetNumNodes1D(ceed_data->basis_x, &P); 19819706a06SJames Wright 199*1737222fSJames Wright PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, problem->q_data_size_sur, &ceed_data->spanstats.elem_restr_parent_stats, 20019706a06SJames Wright &ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.elem_restr_parent_qd)); 20119706a06SJames Wright CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_x, &num_comp_x); 20219706a06SJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.x_coord, NULL); 203*1737222fSJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &user->spanstats.stats_ceed, NULL); 204*1737222fSJames Wright CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_qd, &ceed_data->spanstats.q_data, NULL); 20519706a06SJames Wright 20619706a06SJames Wright CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS_LOBATTO, &ceed_data->spanstats.basis_x); 20719706a06SJames Wright CeedBasisCreateTensorH1Lagrange(ceed, dim, user->spanstats.num_comp_stats, P, Q, CEED_GAUSS, &ceed_data->spanstats.basis_stats); 20819706a06SJames Wright 209*1737222fSJames Wright PetscCall(CreateElemRestrColloc(ceed, user->spanstats.num_comp_stats, ceed_data->spanstats.basis_stats, 210*1737222fSJames Wright ceed_data->spanstats.elem_restr_parent_stats, &ceed_data->spanstats.elem_restr_parent_colloc, 211*1737222fSJames Wright &user->spanstats.parent_stats, NULL)); 212*1737222fSJames Wright PetscCall(CreateElemRestrColloc(ceed, user->spanstats.num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q, 213*1737222fSJames Wright &ceed_data->spanstats.elem_restr_child_colloc, &user->spanstats.child_stats, NULL)); 214*1737222fSJames Wright 21519706a06SJames Wright // -- Copy DM coordinates into CeedVector 21619706a06SJames Wright { 21719706a06SJames Wright DM cdm; 21819706a06SJames Wright PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 21919706a06SJames Wright if (cdm) { 22019706a06SJames Wright PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 22119706a06SJames Wright } else { 22219706a06SJames Wright PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 22319706a06SJames Wright } 22419706a06SJames Wright } 22519706a06SJames Wright PetscCall(VecScale(X_loc, problem->dm_scale)); 22619706a06SJames Wright PetscCall(VecGetArrayReadAndMemType(X_loc, &X_loc_array, &X_loc_memtype)); 22719706a06SJames Wright CeedVectorSetArray(ceed_data->spanstats.x_coord, MemTypeP2C(X_loc_memtype), CEED_COPY_VALUES, (PetscScalar *)X_loc_array); 22819706a06SJames Wright PetscCall(VecRestoreArrayRead(X_loc, &X_loc_array)); 22919706a06SJames Wright 23019706a06SJames Wright // Create SF for communicating child data back their respective parents 23119706a06SJames Wright PetscCall(PetscSFCreate(comm, &user->spanstats.sf)); 23219706a06SJames Wright PetscCall(CreateStatsSF(ceed, ceed_data, user->dm, user->spanstats.dm, user->spanstats.sf)); 23351ee423eSJames Wright 23451ee423eSJames Wright PetscFunctionReturn(0); 23551ee423eSJames Wright } 236