1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3a515125bSLeila Ghaffari // 4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5a515125bSLeila Ghaffari // 6727da7e7SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7a515125bSLeila Ghaffari 8a515125bSLeila Ghaffari /// @file 9a515125bSLeila Ghaffari /// Miscellaneous utility functions 10a515125bSLeila Ghaffari 11e419654dSJeremy L Thompson #include <ceed.h> 12e419654dSJeremy L Thompson #include <petscdm.h> 13e419654dSJeremy L Thompson #include <petscts.h> 14e419654dSJeremy L Thompson 15a515125bSLeila Ghaffari #include "../navierstokes.h" 169f59f36eSJames Wright #include "../qfunctions/mass.h" 17a515125bSLeila Ghaffari 182b916ea7SJeremy L Thompson PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) { 19a515125bSLeila Ghaffari PetscFunctionBeginUser; 20a515125bSLeila Ghaffari 21a515125bSLeila Ghaffari // --------------------------------------------------------------------------- 22b7f03f12SJed Brown // Update time for evaluation 23a515125bSLeila Ghaffari // --------------------------------------------------------------------------- 248f18bb8bSJames Wright if (user->phys->ics_time_label) CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time); 25a515125bSLeila Ghaffari 26a515125bSLeila Ghaffari // --------------------------------------------------------------------------- 27a515125bSLeila Ghaffari // ICs 28a515125bSLeila Ghaffari // --------------------------------------------------------------------------- 29a515125bSLeila Ghaffari // -- CEED Restriction 30a515125bSLeila Ghaffari CeedVector q0_ceed; 31a515125bSLeila Ghaffari CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &q0_ceed, NULL); 32a515125bSLeila Ghaffari 33a515125bSLeila Ghaffari // -- Place PETSc vector in CEED vector 348f18bb8bSJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx)); 35a515125bSLeila Ghaffari 36a515125bSLeila Ghaffari // --------------------------------------------------------------------------- 37a515125bSLeila Ghaffari // Fix multiplicity for output of ICs 38a515125bSLeila Ghaffari // --------------------------------------------------------------------------- 39a515125bSLeila Ghaffari // -- CEED Restriction 40a515125bSLeila Ghaffari CeedVector mult_vec; 41a515125bSLeila Ghaffari CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL); 42a515125bSLeila Ghaffari 43a515125bSLeila Ghaffari // -- Place PETSc vector in CEED vector 44a515125bSLeila Ghaffari PetscMemType m_mem_type; 45a515125bSLeila Ghaffari Vec multiplicity_loc; 462b916ea7SJeremy L Thompson PetscCall(DMGetLocalVector(dm, &multiplicity_loc)); 47fd969b44SJames Wright PetscCall(VecP2C(multiplicity_loc, &m_mem_type, mult_vec)); 48a515125bSLeila Ghaffari 49a515125bSLeila Ghaffari // -- Get multiplicity 50a515125bSLeila Ghaffari CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec); 51a515125bSLeila Ghaffari 52a515125bSLeila Ghaffari // -- Restore vectors 53fd969b44SJames Wright PetscCall(VecC2P(mult_vec, m_mem_type, multiplicity_loc)); 54a515125bSLeila Ghaffari 55a515125bSLeila Ghaffari // -- Local-to-Global 56a515125bSLeila Ghaffari Vec multiplicity; 572b916ea7SJeremy L Thompson PetscCall(DMGetGlobalVector(dm, &multiplicity)); 582b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(multiplicity)); 592b916ea7SJeremy L Thompson PetscCall(DMLocalToGlobal(dm, multiplicity_loc, ADD_VALUES, multiplicity)); 60a515125bSLeila Ghaffari 61a515125bSLeila Ghaffari // -- Fix multiplicity 622b916ea7SJeremy L Thompson PetscCall(VecPointwiseDivide(Q, Q, multiplicity)); 632b916ea7SJeremy L Thompson PetscCall(VecPointwiseDivide(Q_loc, Q_loc, multiplicity_loc)); 64a515125bSLeila Ghaffari 65a515125bSLeila Ghaffari // -- Restore vectors 662b916ea7SJeremy L Thompson PetscCall(DMRestoreLocalVector(dm, &multiplicity_loc)); 672b916ea7SJeremy L Thompson PetscCall(DMRestoreGlobalVector(dm, &multiplicity)); 68a515125bSLeila Ghaffari 69a515125bSLeila Ghaffari // Cleanup 70a515125bSLeila Ghaffari CeedVectorDestroy(&mult_vec); 71a515125bSLeila Ghaffari CeedVectorDestroy(&q0_ceed); 72a515125bSLeila Ghaffari 73d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 74a515125bSLeila Ghaffari } 75a515125bSLeila Ghaffari 762b916ea7SJeremy L Thompson PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 772b916ea7SJeremy L Thompson Vec grad_FVM) { 789d437337SJames Wright Vec Qbc, boundary_mask; 79a515125bSLeila Ghaffari PetscFunctionBegin; 80a515125bSLeila Ghaffari 812eb7bf1fSJames Wright // Mask (zero) Strong BC entries 829d437337SJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 839d437337SJames Wright PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask)); 849d437337SJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 859d437337SJames Wright 862b916ea7SJeremy L Thompson PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 872b916ea7SJeremy L Thompson PetscCall(VecAXPY(Q_loc, 1., Qbc)); 882b916ea7SJeremy L Thompson PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 89a515125bSLeila Ghaffari 90d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 91a515125bSLeila Ghaffari } 92a515125bSLeila Ghaffari 93e7754af5SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number 94e7754af5SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) { 95*e1233009SJames Wright PetscInt file_step_number; 96*e1233009SJames Wright PetscInt32 token; 97e7754af5SKenneth E. Jansen PetscReal file_time; 98e7754af5SKenneth E. Jansen PetscFunctionBeginUser; 99e7754af5SKenneth E. Jansen 100e7754af5SKenneth E. Jansen // Attempt 101*e1233009SJames Wright PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32)); 102*e1233009SJames Wright if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 || 103*e1233009SJames Wright token == FLUIDS_FILE_TOKEN) { // New style format; we're reading a file with step number and time in the header 104e7754af5SKenneth E. Jansen PetscCall(PetscViewerBinaryRead(viewer, &file_step_number, 1, NULL, PETSC_INT)); 105e7754af5SKenneth E. Jansen PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL)); 106e7754af5SKenneth E. Jansen if (time) *time = file_time; 107e7754af5SKenneth E. Jansen if (step_number) *step_number = file_step_number; 108e7754af5SKenneth E. Jansen } else if (token == VEC_FILE_CLASSID) { // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ] 109e7754af5SKenneth E. Jansen PetscInt length, N; 110e7754af5SKenneth E. Jansen PetscCall(PetscViewerBinaryRead(viewer, &length, 1, NULL, PETSC_INT)); 111e7754af5SKenneth E. Jansen PetscCall(VecGetSize(Q, &N)); 112e7754af5SKenneth E. Jansen PetscCheck(length == N, comm, PETSC_ERR_ARG_INCOMP, "File Vec has length %" PetscInt_FMT " but DM has global Vec size %" PetscInt_FMT, length, N); 113e7754af5SKenneth E. Jansen PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE)); 114e7754af5SKenneth E. Jansen } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file"); 115e7754af5SKenneth E. Jansen 116e7754af5SKenneth E. Jansen // Load Q from existent solution 117e7754af5SKenneth E. Jansen PetscCall(VecLoad(Q, viewer)); 118e7754af5SKenneth E. Jansen 119d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 120e7754af5SKenneth E. Jansen } 121e7754af5SKenneth E. Jansen 122a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI 123a515125bSLeila Ghaffari PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q) { 124a515125bSLeila Ghaffari Vec Qref; 125a515125bSLeila Ghaffari PetscViewer viewer; 126a515125bSLeila Ghaffari PetscReal error, Qrefnorm; 127e7754af5SKenneth E. Jansen MPI_Comm comm = PetscObjectComm((PetscObject)Q); 128a515125bSLeila Ghaffari PetscFunctionBegin; 129a515125bSLeila Ghaffari 130a515125bSLeila Ghaffari // Read reference file 1312b916ea7SJeremy L Thompson PetscCall(VecDuplicate(Q, &Qref)); 132e7754af5SKenneth E. Jansen PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer)); 133e7754af5SKenneth E. Jansen PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL)); 134a515125bSLeila Ghaffari 135a515125bSLeila Ghaffari // Compute error with respect to reference solution 1362b916ea7SJeremy L Thompson PetscCall(VecAXPY(Q, -1.0, Qref)); 1372b916ea7SJeremy L Thompson PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm)); 1382b916ea7SJeremy L Thompson PetscCall(VecScale(Q, 1. / Qrefnorm)); 1392b916ea7SJeremy L Thompson PetscCall(VecNorm(Q, NORM_MAX, &error)); 140a515125bSLeila Ghaffari 141a515125bSLeila Ghaffari // Check error 142a515125bSLeila Ghaffari if (error > app_ctx->test_tol) { 1432b916ea7SJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error)); 144a515125bSLeila Ghaffari } 145a515125bSLeila Ghaffari 146a515125bSLeila Ghaffari // Cleanup 1472b916ea7SJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 1482b916ea7SJeremy L Thompson PetscCall(VecDestroy(&Qref)); 149a515125bSLeila Ghaffari 150d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 151a515125bSLeila Ghaffari } 152a515125bSLeila Ghaffari 153a515125bSLeila Ghaffari // Get error for problems with exact solutions 1542b916ea7SJeremy L Thompson PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) { 155a515125bSLeila Ghaffari PetscInt loc_nodes; 156a515125bSLeila Ghaffari Vec Q_exact, Q_exact_loc; 157a515125bSLeila Ghaffari PetscReal rel_error, norm_error, norm_exact; 158a515125bSLeila Ghaffari PetscFunctionBegin; 159a515125bSLeila Ghaffari 160a515125bSLeila Ghaffari // Get exact solution at final time 1612b916ea7SJeremy L Thompson PetscCall(DMCreateGlobalVector(dm, &Q_exact)); 1622b916ea7SJeremy L Thompson PetscCall(DMGetLocalVector(dm, &Q_exact_loc)); 1632b916ea7SJeremy L Thompson PetscCall(VecGetSize(Q_exact_loc, &loc_nodes)); 1642b916ea7SJeremy L Thompson PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time)); 165a515125bSLeila Ghaffari 166a515125bSLeila Ghaffari // Get |exact solution - obtained solution| 1672b916ea7SJeremy L Thompson PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact)); 1682b916ea7SJeremy L Thompson PetscCall(VecAXPY(Q, -1.0, Q_exact)); 1692b916ea7SJeremy L Thompson PetscCall(VecNorm(Q, NORM_1, &norm_error)); 170a515125bSLeila Ghaffari 171a515125bSLeila Ghaffari // Compute relative error 172a515125bSLeila Ghaffari rel_error = norm_error / norm_exact; 173a515125bSLeila Ghaffari 174a515125bSLeila Ghaffari // Output relative error 1752b916ea7SJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error)); 176a515125bSLeila Ghaffari // Cleanup 1772b916ea7SJeremy L Thompson PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc)); 1782b916ea7SJeremy L Thompson PetscCall(VecDestroy(&Q_exact)); 179a515125bSLeila Ghaffari 180d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 181a515125bSLeila Ghaffari } 182a515125bSLeila Ghaffari 183a515125bSLeila Ghaffari // Post-processing 1842b916ea7SJeremy L Thompson PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) { 185a515125bSLeila Ghaffari PetscInt steps; 186f0784ed3SJed Brown TSConvergedReason reason; 187a515125bSLeila Ghaffari PetscFunctionBegin; 188a515125bSLeila Ghaffari 189a515125bSLeila Ghaffari // Print relative error 1900e1e9333SJames Wright if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) { 1912b916ea7SJeremy L Thompson PetscCall(GetError_NS(ceed_data, dm, user, Q, final_time)); 192a515125bSLeila Ghaffari } 193a515125bSLeila Ghaffari 194a515125bSLeila Ghaffari // Print final time and number of steps 1952b916ea7SJeremy L Thompson PetscCall(TSGetStepNumber(ts, &steps)); 196f0784ed3SJed Brown PetscCall(TSGetConvergedReason(ts, &reason)); 1970e1e9333SJames Wright if (user->app_ctx->test_type == TESTTYPE_NONE) { 198f0784ed3SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason], 199f0784ed3SJed Brown steps, (double)final_time)); 200a515125bSLeila Ghaffari } 201a515125bSLeila Ghaffari 202a515125bSLeila Ghaffari // Output numerical values from command line 2032b916ea7SJeremy L Thompson PetscCall(VecViewFromOptions(Q, NULL, "-vec_view")); 204a515125bSLeila Ghaffari 205a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI 2060e1e9333SJames Wright if (user->app_ctx->test_type == TESTTYPE_SOLVER) { 2072b916ea7SJeremy L Thompson PetscCall(RegressionTests_NS(user->app_ctx, Q)); 208a515125bSLeila Ghaffari } 209d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 210a515125bSLeila Ghaffari } 211a515125bSLeila Ghaffari 212*e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN = 0xceedf00; // for backwards compatibility 213*e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32; 214*e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64; 2159293eaa1SJed Brown 216a515125bSLeila Ghaffari // Gather initial Q values in case of continuation of simulation 217a515125bSLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) { 218a515125bSLeila Ghaffari PetscViewer viewer; 2192b916ea7SJeremy L Thompson 220a515125bSLeila Ghaffari PetscFunctionBegin; 221a515125bSLeila Ghaffari 2222b916ea7SJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer)); 223e7754af5SKenneth E. Jansen PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps)); 2242b916ea7SJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 225a515125bSLeila Ghaffari 226d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 227a515125bSLeila Ghaffari } 228a515125bSLeila Ghaffari 229a515125bSLeila Ghaffari // Record boundary values from initial condition 230a515125bSLeila Ghaffari PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc) { 2319d437337SJames Wright Vec Qbc, boundary_mask; 232a515125bSLeila Ghaffari PetscFunctionBegin; 233a515125bSLeila Ghaffari 2342b916ea7SJeremy L Thompson PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 2352b916ea7SJeremy L Thompson PetscCall(VecCopy(Q_loc, Qbc)); 2362b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(Q_loc)); 2372b916ea7SJeremy L Thompson PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc)); 2382b916ea7SJeremy L Thompson PetscCall(VecAXPY(Qbc, -1., Q_loc)); 2392b916ea7SJeremy L Thompson PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 2402b916ea7SJeremy L Thompson PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS)); 241a515125bSLeila Ghaffari 2429d437337SJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 2439d437337SJames Wright PetscCall(DMGetGlobalVector(dm, &Q)); 2449d437337SJames Wright PetscCall(VecZeroEntries(boundary_mask)); 2459d437337SJames Wright PetscCall(VecSet(Q, 1.0)); 2469d437337SJames Wright PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask)); 2479d437337SJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 2489d437337SJames Wright 249d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 250a515125bSLeila Ghaffari } 25115a3537eSJed Brown 25215a3537eSJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes 25315a3537eSJed Brown int FreeContextPetsc(void *data) { 2542b916ea7SJeremy L Thompson if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed"); 25515a3537eSJed Brown return CEED_ERROR_SUCCESS; 25615a3537eSJed Brown } 2579f59f36eSJames Wright 2589f59f36eSJames Wright // Return mass qfunction specification for number of components N 2599f59f36eSJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) { 2609f59f36eSJames Wright PetscFunctionBeginUser; 2619f59f36eSJames Wright 2629f59f36eSJames Wright switch (N) { 2639f59f36eSJames Wright case 1: 2646f539f71SJames Wright CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf); 2659f59f36eSJames Wright break; 2669f59f36eSJames Wright case 5: 2676f539f71SJames Wright CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf); 2689f59f36eSJames Wright break; 269c38c977aSJames Wright case 7: 2706f539f71SJames Wright CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf); 271c38c977aSJames Wright break; 2729f59f36eSJames Wright case 9: 2736f539f71SJames Wright CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf); 2749f59f36eSJames Wright break; 2759f59f36eSJames Wright case 22: 2766f539f71SJames Wright CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf); 2779f59f36eSJames Wright break; 2789f59f36eSJames Wright default: 2796f539f71SJames Wright SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N); 2809f59f36eSJames Wright } 2819f59f36eSJames Wright 2829f59f36eSJames Wright CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP); 2839f59f36eSJames Wright CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE); 2849f59f36eSJames Wright CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP); 285d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2869f59f36eSJames Wright } 287e5e81594SJames Wright 288e5e81594SJames Wright /* @brief L^2 Projection of a source FEM function to a target FEM space 289e5e81594SJames Wright * 290e5e81594SJames Wright * To solve system using a lumped mass matrix, pass a KSP object with ksp_type=preonly, pc_type=jacobi, pc_jacobi_type=rowsum. 291e5e81594SJames Wright * 292e5e81594SJames Wright * @param[in] source_vec Global Vec of the source FEM function. NULL indicates using rhs_matop_ctx->X_loc 293e5e81594SJames Wright * @param[out] target_vec Global Vec of the target (result) FEM function. NULL indicates using rhs_matop_ctx->Y_loc 294e5e81594SJames Wright * @param[in] rhs_matop_ctx MatopApplyContext for performing the RHS evaluation 295e5e81594SJames Wright * @param[in] ksp KSP for solving the consistent projection problem 296e5e81594SJames Wright */ 297f7325489SJames Wright PetscErrorCode ComputeL2Projection(Vec source_vec, Vec target_vec, OperatorApplyContext rhs_matop_ctx, KSP ksp) { 298e5e81594SJames Wright PetscFunctionBeginUser; 299e5e81594SJames Wright 300f7325489SJames Wright PetscCall(ApplyCeedOperatorGlobalToGlobal(source_vec, target_vec, rhs_matop_ctx)); 301e5e81594SJames Wright PetscCall(KSPSolve(ksp, target_vec, target_vec)); 302e5e81594SJames Wright 303d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 304e5e81594SJames Wright } 305676080b4SJames Wright 306457a5831SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) { 307457a5831SJames Wright PetscFunctionBeginUser; 308d949ddfcSJames Wright if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS); 309457a5831SJames Wright 310457a5831SJames Wright PetscCall(DMDestroy(&context->dm)); 311457a5831SJames Wright PetscCall(KSPDestroy(&context->ksp)); 312457a5831SJames Wright 313457a5831SJames Wright PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx)); 314457a5831SJames Wright 315457a5831SJames Wright PetscCall(PetscFree(context)); 316457a5831SJames Wright 317d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 318457a5831SJames Wright } 319457a5831SJames Wright 320676080b4SJames Wright /* 321676080b4SJames Wright * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer 322676080b4SJames Wright * 323676080b4SJames Wright * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`. 324676080b4SJames Wright * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example). 325676080b4SJames Wright * 326676080b4SJames Wright * Assumes that the first line of the file has the number of rows and columns as the only two entries, separated by a single space. 327676080b4SJames Wright * 328676080b4SJames Wright * @param[in] comm MPI_Comm for the program 329676080b4SJames Wright * @param[in] path Path to the file 330676080b4SJames Wright * @param[in] char_array_len Length of the character array that should contain each line 331676080b4SJames Wright * @param[out] dims Dimensions of the file, taken from the first line of the file 332676080b4SJames Wright * @param[out] fp File pointer to the opened file 333676080b4SJames Wright */ 334676080b4SJames Wright PetscErrorCode PHASTADatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2], 335676080b4SJames Wright FILE **fp) { 336defe8520SJames Wright int ndims; 337676080b4SJames Wright char line[char_array_len]; 338676080b4SJames Wright char **array; 339676080b4SJames Wright 340676080b4SJames Wright PetscFunctionBeginUser; 341676080b4SJames Wright PetscCall(PetscFOpen(comm, path, "r", fp)); 342676080b4SJames Wright PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line)); 343676080b4SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 344defe8520SJames Wright PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path); 345676080b4SJames Wright 346676080b4SJames Wright for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]); 347676080b4SJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array)); 348676080b4SJames Wright 349d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 350676080b4SJames Wright } 351676080b4SJames Wright 352676080b4SJames Wright /* 353676080b4SJames Wright * @brief Get the number of rows for the PHASTA file at path. 354676080b4SJames Wright * 355676080b4SJames Wright * Assumes that the first line of the file has the number of rows and columns as the only two entries, separated by a single space. 356676080b4SJames Wright * 357676080b4SJames Wright * @param[in] comm MPI_Comm for the program 358676080b4SJames Wright * @param[in] path Path to the file 359676080b4SJames Wright * @param[out] nrows Number of rows 360676080b4SJames Wright */ 361676080b4SJames Wright PetscErrorCode PHASTADatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) { 362676080b4SJames Wright const PetscInt char_array_len = 512; 363676080b4SJames Wright PetscInt dims[2]; 364676080b4SJames Wright FILE *fp; 365676080b4SJames Wright 366676080b4SJames Wright PetscFunctionBeginUser; 367676080b4SJames Wright PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp)); 368676080b4SJames Wright *nrows = dims[0]; 369676080b4SJames Wright PetscCall(PetscFClose(comm, fp)); 370676080b4SJames Wright 371d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 372676080b4SJames Wright } 37362b7942eSJames Wright 37462b7942eSJames Wright PetscErrorCode PHASTADatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) { 375defe8520SJames Wright PetscInt dims[2]; 376defe8520SJames Wright int ndims; 37762b7942eSJames Wright FILE *fp; 37862b7942eSJames Wright const PetscInt char_array_len = 512; 37962b7942eSJames Wright char line[char_array_len]; 38062b7942eSJames Wright char **row_array; 38162b7942eSJames Wright PetscFunctionBeginUser; 38262b7942eSJames Wright 38362b7942eSJames Wright PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp)); 38462b7942eSJames Wright 38562b7942eSJames Wright for (PetscInt i = 0; i < dims[0]; i++) { 38662b7942eSJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 38762b7942eSJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array)); 3885d28dccaSJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 389defe8520SJames Wright "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 39062b7942eSJames Wright 39162b7942eSJames Wright for (PetscInt j = 0; j < dims[1]; j++) { 39262b7942eSJames Wright array[i * dims[1] + j] = (PetscReal)atof(row_array[j]); 39362b7942eSJames Wright } 39462b7942eSJames Wright } 39562b7942eSJames Wright 39662b7942eSJames Wright PetscCall(PetscFClose(comm, fp)); 39762b7942eSJames Wright 398d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 39962b7942eSJames Wright } 4007eedc94cSJames Wright 4017eedc94cSJames Wright PetscLogEvent FLUIDS_CeedOperatorApply; 4027eedc94cSJames Wright PetscLogEvent FLUIDS_CeedOperatorAssemble; 4037eedc94cSJames Wright PetscLogEvent FLUIDS_CeedOperatorAssembleDiagonal; 4047eedc94cSJames Wright PetscLogEvent FLUIDS_CeedOperatorAssemblePointBlockDiagonal; 4057eedc94cSJames Wright static PetscClassId libCEED_classid; 4067eedc94cSJames Wright 4077eedc94cSJames Wright PetscErrorCode RegisterLogEvents() { 4087eedc94cSJames Wright PetscFunctionBeginUser; 4097eedc94cSJames Wright PetscCall(PetscClassIdRegister("libCEED", &libCEED_classid)); 4107eedc94cSJames Wright PetscCall(PetscLogEventRegister("CeedOpApply", libCEED_classid, &FLUIDS_CeedOperatorApply)); 4117eedc94cSJames Wright PetscCall(PetscLogEventRegister("CeedOpAsm", libCEED_classid, &FLUIDS_CeedOperatorAssemble)); 4127eedc94cSJames Wright PetscCall(PetscLogEventRegister("CeedOpAsmD", libCEED_classid, &FLUIDS_CeedOperatorAssembleDiagonal)); 4137eedc94cSJames Wright PetscCall(PetscLogEventRegister("CeedOpAsmPBD", libCEED_classid, &FLUIDS_CeedOperatorAssemblePointBlockDiagonal)); 414d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4157eedc94cSJames Wright } 416defe8520SJames Wright 417defe8520SJames Wright /** 418defe8520SJames Wright @brief Translate array of CeedInt to PetscInt. 419defe8520SJames Wright If the types differ, `array_ceed` is freed with `free()` and `array_petsc` is allocated with `malloc()`. 420defe8520SJames Wright Caller is responsible for freeing `array_petsc` with `free()`. 421defe8520SJames Wright 422defe8520SJames Wright @param[in] num_entries Number of array entries 423defe8520SJames Wright @param[in,out] array_ceed Array of CeedInts 424defe8520SJames Wright @param[out] array_petsc Array of PetscInts 425defe8520SJames Wright **/ 426defe8520SJames Wright PetscErrorCode IntArrayC2P(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc) { 427defe8520SJames Wright CeedInt int_c = 0; 428defe8520SJames Wright PetscInt int_p = 0; 429defe8520SJames Wright 430defe8520SJames Wright PetscFunctionBeginUser; 431defe8520SJames Wright if (sizeof(int_c) == sizeof(int_p)) { 432defe8520SJames Wright *array_petsc = (PetscInt *)*array_ceed; 433defe8520SJames Wright } else { 434defe8520SJames Wright *array_petsc = malloc(num_entries * sizeof(PetscInt)); 435defe8520SJames Wright for (PetscInt i = 0; i < num_entries; i++) (*array_petsc)[i] = (*array_ceed)[i]; 436defe8520SJames Wright free(*array_ceed); 437defe8520SJames Wright } 438defe8520SJames Wright *array_ceed = NULL; 439defe8520SJames Wright 440defe8520SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 441defe8520SJames Wright } 442defe8520SJames Wright 443defe8520SJames Wright /** 444defe8520SJames Wright @brief Translate array of PetscInt to CeedInt. 445defe8520SJames Wright If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`. 446defe8520SJames Wright Caller is responsible for freeing `array_ceed` with `PetscFree()`. 447defe8520SJames Wright 448defe8520SJames Wright @param[in] num_entries Number of array entries 449defe8520SJames Wright @param[in,out] array_petsc Array of PetscInts 450defe8520SJames Wright @param[out] array_ceed Array of CeedInts 451defe8520SJames Wright **/ 452defe8520SJames Wright PetscErrorCode IntArrayP2C(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed) { 453defe8520SJames Wright CeedInt int_c = 0; 454defe8520SJames Wright PetscInt int_p = 0; 455defe8520SJames Wright 456defe8520SJames Wright PetscFunctionBeginUser; 457defe8520SJames Wright if (sizeof(int_c) == sizeof(int_p)) { 458defe8520SJames Wright *array_ceed = (CeedInt *)*array_petsc; 459defe8520SJames Wright } else { 460defe8520SJames Wright PetscCall(PetscMalloc1(num_entries, array_ceed)); 461defe8520SJames Wright for (PetscInt i = 0; i < num_entries; i++) (*array_ceed)[i] = (*array_petsc)[i]; 462defe8520SJames Wright PetscCall(PetscFree(*array_petsc)); 463defe8520SJames Wright } 464defe8520SJames Wright *array_petsc = NULL; 465defe8520SJames Wright 466defe8520SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 467defe8520SJames Wright } 468