13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 377841947SLeila Ghaffari // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 577841947SLeila Ghaffari // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 777841947SLeila Ghaffari 877841947SLeila Ghaffari /// @file 977841947SLeila Ghaffari /// Miscellaneous utility functions 1077841947SLeila Ghaffari 1149aac155SJeremy L Thompson #include <ceed.h> 1249aac155SJeremy L Thompson #include <petscdm.h> 1349aac155SJeremy L Thompson #include <petscts.h> 1449aac155SJeremy L Thompson 1577841947SLeila Ghaffari #include "../navierstokes.h" 16ef080ff9SJames Wright #include "../qfunctions/mass.h" 1777841947SLeila Ghaffari 182b730f8bSJeremy L Thompson PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) { 1977841947SLeila Ghaffari PetscFunctionBeginUser; 2077841947SLeila Ghaffari 2177841947SLeila Ghaffari // --------------------------------------------------------------------------- 22a0add3c9SJed Brown // Update time for evaluation 2377841947SLeila Ghaffari // --------------------------------------------------------------------------- 245263e9c6SJames Wright if (user->phys->ics_time_label) CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time); 2577841947SLeila Ghaffari 2677841947SLeila Ghaffari // --------------------------------------------------------------------------- 2777841947SLeila Ghaffari // ICs 2877841947SLeila Ghaffari // --------------------------------------------------------------------------- 2977841947SLeila Ghaffari // -- CEED Restriction 3077841947SLeila Ghaffari CeedVector q0_ceed; 3177841947SLeila Ghaffari CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &q0_ceed, NULL); 3277841947SLeila Ghaffari 3377841947SLeila Ghaffari // -- Place PETSc vector in CEED vector 345263e9c6SJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx)); 3577841947SLeila Ghaffari 3677841947SLeila Ghaffari // --------------------------------------------------------------------------- 3777841947SLeila Ghaffari // Fix multiplicity for output of ICs 3877841947SLeila Ghaffari // --------------------------------------------------------------------------- 3977841947SLeila Ghaffari // -- CEED Restriction 4077841947SLeila Ghaffari CeedVector mult_vec; 4177841947SLeila Ghaffari CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL); 4277841947SLeila Ghaffari 4377841947SLeila Ghaffari // -- Place PETSc vector in CEED vector 4477841947SLeila Ghaffari PetscMemType m_mem_type; 4577841947SLeila Ghaffari Vec multiplicity_loc; 462b730f8bSJeremy L Thompson PetscCall(DMGetLocalVector(dm, &multiplicity_loc)); 47c798d55aSJames Wright PetscCall(VecP2C(multiplicity_loc, &m_mem_type, mult_vec)); 4877841947SLeila Ghaffari 4977841947SLeila Ghaffari // -- Get multiplicity 5077841947SLeila Ghaffari CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec); 5177841947SLeila Ghaffari 5277841947SLeila Ghaffari // -- Restore vectors 53c798d55aSJames Wright PetscCall(VecC2P(mult_vec, m_mem_type, multiplicity_loc)); 5477841947SLeila Ghaffari 5577841947SLeila Ghaffari // -- Local-to-Global 5677841947SLeila Ghaffari Vec multiplicity; 572b730f8bSJeremy L Thompson PetscCall(DMGetGlobalVector(dm, &multiplicity)); 582b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(multiplicity)); 592b730f8bSJeremy L Thompson PetscCall(DMLocalToGlobal(dm, multiplicity_loc, ADD_VALUES, multiplicity)); 6077841947SLeila Ghaffari 6177841947SLeila Ghaffari // -- Fix multiplicity 622b730f8bSJeremy L Thompson PetscCall(VecPointwiseDivide(Q, Q, multiplicity)); 632b730f8bSJeremy L Thompson PetscCall(VecPointwiseDivide(Q_loc, Q_loc, multiplicity_loc)); 6477841947SLeila Ghaffari 6577841947SLeila Ghaffari // -- Restore vectors 662b730f8bSJeremy L Thompson PetscCall(DMRestoreLocalVector(dm, &multiplicity_loc)); 672b730f8bSJeremy L Thompson PetscCall(DMRestoreGlobalVector(dm, &multiplicity)); 6877841947SLeila Ghaffari 6977841947SLeila Ghaffari // Cleanup 7077841947SLeila Ghaffari CeedVectorDestroy(&mult_vec); 7177841947SLeila Ghaffari CeedVectorDestroy(&q0_ceed); 7277841947SLeila Ghaffari 73ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 7477841947SLeila Ghaffari } 7577841947SLeila Ghaffari 762b730f8bSJeremy L Thompson PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 772b730f8bSJeremy L Thompson Vec grad_FVM) { 785571c6fdSJames Wright Vec Qbc, boundary_mask; 7977841947SLeila Ghaffari PetscFunctionBegin; 8077841947SLeila Ghaffari 813722cd23SJames Wright // Mask (zero) Strong BC entries 825571c6fdSJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 835571c6fdSJames Wright PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask)); 845571c6fdSJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 855571c6fdSJames Wright 862b730f8bSJeremy L Thompson PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 872b730f8bSJeremy L Thompson PetscCall(VecAXPY(Q_loc, 1., Qbc)); 882b730f8bSJeremy L Thompson PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 8977841947SLeila Ghaffari 90ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 9177841947SLeila Ghaffari } 9277841947SLeila Ghaffari 93530ad8c4SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number 94530ad8c4SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) { 95*0de6a49fSJames Wright PetscInt file_step_number; 96*0de6a49fSJames Wright PetscInt32 token; 97530ad8c4SKenneth E. Jansen PetscReal file_time; 98530ad8c4SKenneth E. Jansen PetscFunctionBeginUser; 99530ad8c4SKenneth E. Jansen 100530ad8c4SKenneth E. Jansen // Attempt 101*0de6a49fSJames Wright PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32)); 102*0de6a49fSJames Wright if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 || 103*0de6a49fSJames Wright token == FLUIDS_FILE_TOKEN) { // New style format; we're reading a file with step number and time in the header 104530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinaryRead(viewer, &file_step_number, 1, NULL, PETSC_INT)); 105530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL)); 106530ad8c4SKenneth E. Jansen if (time) *time = file_time; 107530ad8c4SKenneth E. Jansen if (step_number) *step_number = file_step_number; 108530ad8c4SKenneth E. Jansen } else if (token == VEC_FILE_CLASSID) { // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ] 109530ad8c4SKenneth E. Jansen PetscInt length, N; 110530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinaryRead(viewer, &length, 1, NULL, PETSC_INT)); 111530ad8c4SKenneth E. Jansen PetscCall(VecGetSize(Q, &N)); 112530ad8c4SKenneth 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); 113530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE)); 114530ad8c4SKenneth E. Jansen } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file"); 115530ad8c4SKenneth E. Jansen 116530ad8c4SKenneth E. Jansen // Load Q from existent solution 117530ad8c4SKenneth E. Jansen PetscCall(VecLoad(Q, viewer)); 118530ad8c4SKenneth E. Jansen 119ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 120530ad8c4SKenneth E. Jansen } 121530ad8c4SKenneth E. Jansen 12277841947SLeila Ghaffari // Compare reference solution values with current test run for CI 12377841947SLeila Ghaffari PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q) { 12477841947SLeila Ghaffari Vec Qref; 12577841947SLeila Ghaffari PetscViewer viewer; 12677841947SLeila Ghaffari PetscReal error, Qrefnorm; 127530ad8c4SKenneth E. Jansen MPI_Comm comm = PetscObjectComm((PetscObject)Q); 12877841947SLeila Ghaffari PetscFunctionBegin; 12977841947SLeila Ghaffari 13077841947SLeila Ghaffari // Read reference file 1312b730f8bSJeremy L Thompson PetscCall(VecDuplicate(Q, &Qref)); 132530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer)); 133530ad8c4SKenneth E. Jansen PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL)); 13477841947SLeila Ghaffari 13577841947SLeila Ghaffari // Compute error with respect to reference solution 1362b730f8bSJeremy L Thompson PetscCall(VecAXPY(Q, -1.0, Qref)); 1372b730f8bSJeremy L Thompson PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm)); 1382b730f8bSJeremy L Thompson PetscCall(VecScale(Q, 1. / Qrefnorm)); 1392b730f8bSJeremy L Thompson PetscCall(VecNorm(Q, NORM_MAX, &error)); 14077841947SLeila Ghaffari 14177841947SLeila Ghaffari // Check error 14277841947SLeila Ghaffari if (error > app_ctx->test_tol) { 1432b730f8bSJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error)); 14477841947SLeila Ghaffari } 14577841947SLeila Ghaffari 14677841947SLeila Ghaffari // Cleanup 1472b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 1482b730f8bSJeremy L Thompson PetscCall(VecDestroy(&Qref)); 14977841947SLeila Ghaffari 150ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 15177841947SLeila Ghaffari } 15277841947SLeila Ghaffari 15377841947SLeila Ghaffari // Get error for problems with exact solutions 1542b730f8bSJeremy L Thompson PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) { 15577841947SLeila Ghaffari PetscInt loc_nodes; 15677841947SLeila Ghaffari Vec Q_exact, Q_exact_loc; 15777841947SLeila Ghaffari PetscReal rel_error, norm_error, norm_exact; 15877841947SLeila Ghaffari PetscFunctionBegin; 15977841947SLeila Ghaffari 16077841947SLeila Ghaffari // Get exact solution at final time 1612b730f8bSJeremy L Thompson PetscCall(DMCreateGlobalVector(dm, &Q_exact)); 1622b730f8bSJeremy L Thompson PetscCall(DMGetLocalVector(dm, &Q_exact_loc)); 1632b730f8bSJeremy L Thompson PetscCall(VecGetSize(Q_exact_loc, &loc_nodes)); 1642b730f8bSJeremy L Thompson PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time)); 16577841947SLeila Ghaffari 16677841947SLeila Ghaffari // Get |exact solution - obtained solution| 1672b730f8bSJeremy L Thompson PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact)); 1682b730f8bSJeremy L Thompson PetscCall(VecAXPY(Q, -1.0, Q_exact)); 1692b730f8bSJeremy L Thompson PetscCall(VecNorm(Q, NORM_1, &norm_error)); 17077841947SLeila Ghaffari 17177841947SLeila Ghaffari // Compute relative error 17277841947SLeila Ghaffari rel_error = norm_error / norm_exact; 17377841947SLeila Ghaffari 17477841947SLeila Ghaffari // Output relative error 1752b730f8bSJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error)); 17677841947SLeila Ghaffari // Cleanup 1772b730f8bSJeremy L Thompson PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc)); 1782b730f8bSJeremy L Thompson PetscCall(VecDestroy(&Q_exact)); 17977841947SLeila Ghaffari 180ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 18177841947SLeila Ghaffari } 18277841947SLeila Ghaffari 18377841947SLeila Ghaffari // Post-processing 1842b730f8bSJeremy L Thompson PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) { 18577841947SLeila Ghaffari PetscInt steps; 186cf7a0454SJed Brown TSConvergedReason reason; 18777841947SLeila Ghaffari PetscFunctionBegin; 18877841947SLeila Ghaffari 18977841947SLeila Ghaffari // Print relative error 1908fb33541SJames Wright if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) { 1912b730f8bSJeremy L Thompson PetscCall(GetError_NS(ceed_data, dm, user, Q, final_time)); 19277841947SLeila Ghaffari } 19377841947SLeila Ghaffari 19477841947SLeila Ghaffari // Print final time and number of steps 1952b730f8bSJeremy L Thompson PetscCall(TSGetStepNumber(ts, &steps)); 196cf7a0454SJed Brown PetscCall(TSGetConvergedReason(ts, &reason)); 1978fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_NONE) { 198cf7a0454SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason], 199cf7a0454SJed Brown steps, (double)final_time)); 20077841947SLeila Ghaffari } 20177841947SLeila Ghaffari 20277841947SLeila Ghaffari // Output numerical values from command line 2032b730f8bSJeremy L Thompson PetscCall(VecViewFromOptions(Q, NULL, "-vec_view")); 20477841947SLeila Ghaffari 20577841947SLeila Ghaffari // Compare reference solution values with current test run for CI 2068fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_SOLVER) { 2072b730f8bSJeremy L Thompson PetscCall(RegressionTests_NS(user->app_ctx, Q)); 20877841947SLeila Ghaffari } 209ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 21077841947SLeila Ghaffari } 21177841947SLeila Ghaffari 212*0de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN = 0xceedf00; // for backwards compatibility 213*0de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32; 214*0de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64; 2154de8550aSJed Brown 21677841947SLeila Ghaffari // Gather initial Q values in case of continuation of simulation 21777841947SLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) { 21877841947SLeila Ghaffari PetscViewer viewer; 2192b730f8bSJeremy L Thompson 22077841947SLeila Ghaffari PetscFunctionBegin; 22177841947SLeila Ghaffari 2222b730f8bSJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer)); 223530ad8c4SKenneth E. Jansen PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps)); 2242b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 22577841947SLeila Ghaffari 226ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 22777841947SLeila Ghaffari } 22877841947SLeila Ghaffari 22977841947SLeila Ghaffari // Record boundary values from initial condition 23077841947SLeila Ghaffari PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc) { 2315571c6fdSJames Wright Vec Qbc, boundary_mask; 23277841947SLeila Ghaffari PetscFunctionBegin; 23377841947SLeila Ghaffari 2342b730f8bSJeremy L Thompson PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 2352b730f8bSJeremy L Thompson PetscCall(VecCopy(Q_loc, Qbc)); 2362b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(Q_loc)); 2372b730f8bSJeremy L Thompson PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc)); 2382b730f8bSJeremy L Thompson PetscCall(VecAXPY(Qbc, -1., Q_loc)); 2392b730f8bSJeremy L Thompson PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 2402b730f8bSJeremy L Thompson PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS)); 24177841947SLeila Ghaffari 2425571c6fdSJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 2435571c6fdSJames Wright PetscCall(DMGetGlobalVector(dm, &Q)); 2445571c6fdSJames Wright PetscCall(VecZeroEntries(boundary_mask)); 2455571c6fdSJames Wright PetscCall(VecSet(Q, 1.0)); 2465571c6fdSJames Wright PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask)); 2475571c6fdSJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 2485571c6fdSJames Wright 249ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 25077841947SLeila Ghaffari } 251841e4c73SJed Brown 252841e4c73SJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes 253841e4c73SJed Brown int FreeContextPetsc(void *data) { 2542b730f8bSJeremy L Thompson if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed"); 255841e4c73SJed Brown return CEED_ERROR_SUCCESS; 256841e4c73SJed Brown } 257ef080ff9SJames Wright 258ef080ff9SJames Wright // Return mass qfunction specification for number of components N 259ef080ff9SJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) { 260ef080ff9SJames Wright PetscFunctionBeginUser; 261ef080ff9SJames Wright 262ef080ff9SJames Wright switch (N) { 263ef080ff9SJames Wright case 1: 26483ae4962SJames Wright CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf); 265ef080ff9SJames Wright break; 266ef080ff9SJames Wright case 5: 26783ae4962SJames Wright CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf); 268ef080ff9SJames Wright break; 269052409adSJames Wright case 7: 27083ae4962SJames Wright CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf); 271052409adSJames Wright break; 272ef080ff9SJames Wright case 9: 27383ae4962SJames Wright CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf); 274ef080ff9SJames Wright break; 275ef080ff9SJames Wright case 22: 27683ae4962SJames Wright CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf); 277ef080ff9SJames Wright break; 278ef080ff9SJames Wright default: 27983ae4962SJames Wright SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N); 280ef080ff9SJames Wright } 281ef080ff9SJames Wright 282ef080ff9SJames Wright CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP); 283ef080ff9SJames Wright CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE); 284ef080ff9SJames Wright CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP); 285ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 286ef080ff9SJames Wright } 2870df12fefSJames Wright 2880df12fefSJames Wright /* @brief L^2 Projection of a source FEM function to a target FEM space 2890df12fefSJames Wright * 2900df12fefSJames Wright * To solve system using a lumped mass matrix, pass a KSP object with ksp_type=preonly, pc_type=jacobi, pc_jacobi_type=rowsum. 2910df12fefSJames Wright * 2920df12fefSJames Wright * @param[in] source_vec Global Vec of the source FEM function. NULL indicates using rhs_matop_ctx->X_loc 2930df12fefSJames Wright * @param[out] target_vec Global Vec of the target (result) FEM function. NULL indicates using rhs_matop_ctx->Y_loc 2940df12fefSJames Wright * @param[in] rhs_matop_ctx MatopApplyContext for performing the RHS evaluation 2950df12fefSJames Wright * @param[in] ksp KSP for solving the consistent projection problem 2960df12fefSJames Wright */ 2974021610dSJames Wright PetscErrorCode ComputeL2Projection(Vec source_vec, Vec target_vec, OperatorApplyContext rhs_matop_ctx, KSP ksp) { 2980df12fefSJames Wright PetscFunctionBeginUser; 2990df12fefSJames Wright 3004021610dSJames Wright PetscCall(ApplyCeedOperatorGlobalToGlobal(source_vec, target_vec, rhs_matop_ctx)); 3010df12fefSJames Wright PetscCall(KSPSolve(ksp, target_vec, target_vec)); 3020df12fefSJames Wright 303ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3040df12fefSJames Wright } 3058b892a05SJames Wright 306999ff5c7SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) { 307999ff5c7SJames Wright PetscFunctionBeginUser; 308ee4ca9cbSJames Wright if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS); 309999ff5c7SJames Wright 310999ff5c7SJames Wright PetscCall(DMDestroy(&context->dm)); 311999ff5c7SJames Wright PetscCall(KSPDestroy(&context->ksp)); 312999ff5c7SJames Wright 313999ff5c7SJames Wright PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx)); 314999ff5c7SJames Wright 315999ff5c7SJames Wright PetscCall(PetscFree(context)); 316999ff5c7SJames Wright 317ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 318999ff5c7SJames Wright } 319999ff5c7SJames Wright 3208b892a05SJames Wright /* 3218b892a05SJames Wright * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer 3228b892a05SJames Wright * 3238b892a05SJames Wright * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`. 3248b892a05SJames Wright * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example). 3258b892a05SJames Wright * 3268b892a05SJames 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. 3278b892a05SJames Wright * 3288b892a05SJames Wright * @param[in] comm MPI_Comm for the program 3298b892a05SJames Wright * @param[in] path Path to the file 3308b892a05SJames Wright * @param[in] char_array_len Length of the character array that should contain each line 3318b892a05SJames Wright * @param[out] dims Dimensions of the file, taken from the first line of the file 3328b892a05SJames Wright * @param[out] fp File pointer to the opened file 3338b892a05SJames Wright */ 3348b892a05SJames Wright PetscErrorCode PHASTADatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2], 3358b892a05SJames Wright FILE **fp) { 336457e73b2SJames Wright int ndims; 3378b892a05SJames Wright char line[char_array_len]; 3388b892a05SJames Wright char **array; 3398b892a05SJames Wright 3408b892a05SJames Wright PetscFunctionBeginUser; 3418b892a05SJames Wright PetscCall(PetscFOpen(comm, path, "r", fp)); 3428b892a05SJames Wright PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line)); 3438b892a05SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 344457e73b2SJames Wright PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path); 3458b892a05SJames Wright 3468b892a05SJames Wright for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]); 3478b892a05SJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array)); 3488b892a05SJames Wright 349ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3508b892a05SJames Wright } 3518b892a05SJames Wright 3528b892a05SJames Wright /* 3538b892a05SJames Wright * @brief Get the number of rows for the PHASTA file at path. 3548b892a05SJames Wright * 3558b892a05SJames 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. 3568b892a05SJames Wright * 3578b892a05SJames Wright * @param[in] comm MPI_Comm for the program 3588b892a05SJames Wright * @param[in] path Path to the file 3598b892a05SJames Wright * @param[out] nrows Number of rows 3608b892a05SJames Wright */ 3618b892a05SJames Wright PetscErrorCode PHASTADatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) { 3628b892a05SJames Wright const PetscInt char_array_len = 512; 3638b892a05SJames Wright PetscInt dims[2]; 3648b892a05SJames Wright FILE *fp; 3658b892a05SJames Wright 3668b892a05SJames Wright PetscFunctionBeginUser; 3678b892a05SJames Wright PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp)); 3688b892a05SJames Wright *nrows = dims[0]; 3698b892a05SJames Wright PetscCall(PetscFClose(comm, fp)); 3708b892a05SJames Wright 371ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3728b892a05SJames Wright } 3734cc9442bSJames Wright 3744cc9442bSJames Wright PetscErrorCode PHASTADatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) { 375457e73b2SJames Wright PetscInt dims[2]; 376457e73b2SJames Wright int ndims; 3774cc9442bSJames Wright FILE *fp; 3784cc9442bSJames Wright const PetscInt char_array_len = 512; 3794cc9442bSJames Wright char line[char_array_len]; 3804cc9442bSJames Wright char **row_array; 3814cc9442bSJames Wright PetscFunctionBeginUser; 3824cc9442bSJames Wright 3834cc9442bSJames Wright PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp)); 3844cc9442bSJames Wright 3854cc9442bSJames Wright for (PetscInt i = 0; i < dims[0]; i++) { 3864cc9442bSJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 3874cc9442bSJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array)); 3880e654f56SJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 389457e73b2SJames Wright "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 3904cc9442bSJames Wright 3914cc9442bSJames Wright for (PetscInt j = 0; j < dims[1]; j++) { 3924cc9442bSJames Wright array[i * dims[1] + j] = (PetscReal)atof(row_array[j]); 3934cc9442bSJames Wright } 3944cc9442bSJames Wright } 3954cc9442bSJames Wright 3964cc9442bSJames Wright PetscCall(PetscFClose(comm, fp)); 3974cc9442bSJames Wright 398ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3994cc9442bSJames Wright } 40075d1798cSJames Wright 40175d1798cSJames Wright PetscLogEvent FLUIDS_CeedOperatorApply; 40275d1798cSJames Wright PetscLogEvent FLUIDS_CeedOperatorAssemble; 40375d1798cSJames Wright PetscLogEvent FLUIDS_CeedOperatorAssembleDiagonal; 40475d1798cSJames Wright PetscLogEvent FLUIDS_CeedOperatorAssemblePointBlockDiagonal; 40575d1798cSJames Wright static PetscClassId libCEED_classid; 40675d1798cSJames Wright 40775d1798cSJames Wright PetscErrorCode RegisterLogEvents() { 40875d1798cSJames Wright PetscFunctionBeginUser; 40975d1798cSJames Wright PetscCall(PetscClassIdRegister("libCEED", &libCEED_classid)); 41075d1798cSJames Wright PetscCall(PetscLogEventRegister("CeedOpApply", libCEED_classid, &FLUIDS_CeedOperatorApply)); 41175d1798cSJames Wright PetscCall(PetscLogEventRegister("CeedOpAsm", libCEED_classid, &FLUIDS_CeedOperatorAssemble)); 41275d1798cSJames Wright PetscCall(PetscLogEventRegister("CeedOpAsmD", libCEED_classid, &FLUIDS_CeedOperatorAssembleDiagonal)); 41375d1798cSJames Wright PetscCall(PetscLogEventRegister("CeedOpAsmPBD", libCEED_classid, &FLUIDS_CeedOperatorAssemblePointBlockDiagonal)); 414ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 41575d1798cSJames Wright } 416457e73b2SJames Wright 417457e73b2SJames Wright /** 418457e73b2SJames Wright @brief Translate array of CeedInt to PetscInt. 419457e73b2SJames Wright If the types differ, `array_ceed` is freed with `free()` and `array_petsc` is allocated with `malloc()`. 420457e73b2SJames Wright Caller is responsible for freeing `array_petsc` with `free()`. 421457e73b2SJames Wright 422457e73b2SJames Wright @param[in] num_entries Number of array entries 423457e73b2SJames Wright @param[in,out] array_ceed Array of CeedInts 424457e73b2SJames Wright @param[out] array_petsc Array of PetscInts 425457e73b2SJames Wright **/ 426457e73b2SJames Wright PetscErrorCode IntArrayC2P(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc) { 427457e73b2SJames Wright CeedInt int_c = 0; 428457e73b2SJames Wright PetscInt int_p = 0; 429457e73b2SJames Wright 430457e73b2SJames Wright PetscFunctionBeginUser; 431457e73b2SJames Wright if (sizeof(int_c) == sizeof(int_p)) { 432457e73b2SJames Wright *array_petsc = (PetscInt *)*array_ceed; 433457e73b2SJames Wright } else { 434457e73b2SJames Wright *array_petsc = malloc(num_entries * sizeof(PetscInt)); 435457e73b2SJames Wright for (PetscInt i = 0; i < num_entries; i++) (*array_petsc)[i] = (*array_ceed)[i]; 436457e73b2SJames Wright free(*array_ceed); 437457e73b2SJames Wright } 438457e73b2SJames Wright *array_ceed = NULL; 439457e73b2SJames Wright 440457e73b2SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 441457e73b2SJames Wright } 442457e73b2SJames Wright 443457e73b2SJames Wright /** 444457e73b2SJames Wright @brief Translate array of PetscInt to CeedInt. 445457e73b2SJames Wright If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`. 446457e73b2SJames Wright Caller is responsible for freeing `array_ceed` with `PetscFree()`. 447457e73b2SJames Wright 448457e73b2SJames Wright @param[in] num_entries Number of array entries 449457e73b2SJames Wright @param[in,out] array_petsc Array of PetscInts 450457e73b2SJames Wright @param[out] array_ceed Array of CeedInts 451457e73b2SJames Wright **/ 452457e73b2SJames Wright PetscErrorCode IntArrayP2C(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed) { 453457e73b2SJames Wright CeedInt int_c = 0; 454457e73b2SJames Wright PetscInt int_p = 0; 455457e73b2SJames Wright 456457e73b2SJames Wright PetscFunctionBeginUser; 457457e73b2SJames Wright if (sizeof(int_c) == sizeof(int_p)) { 458457e73b2SJames Wright *array_ceed = (CeedInt *)*array_petsc; 459457e73b2SJames Wright } else { 460457e73b2SJames Wright PetscCall(PetscMalloc1(num_entries, array_ceed)); 461457e73b2SJames Wright for (PetscInt i = 0; i < num_entries; i++) (*array_ceed)[i] = (*array_petsc)[i]; 462457e73b2SJames Wright PetscCall(PetscFree(*array_petsc)); 463457e73b2SJames Wright } 464457e73b2SJames Wright *array_petsc = NULL; 465457e73b2SJames Wright 466457e73b2SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 467457e73b2SJames Wright } 468