15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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> 132e31f396SJames Wright #include <petscsf.h> 1449aac155SJeremy L Thompson #include <petscts.h> 1549aac155SJeremy L Thompson 1677841947SLeila Ghaffari #include "../navierstokes.h" 17ef080ff9SJames Wright #include "../qfunctions/mass.h" 1877841947SLeila Ghaffari 192b730f8bSJeremy L Thompson PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) { 20a424bcd0SJames Wright Ceed ceed = user->ceed; 21f6ce2b0aSJames Wright CeedVector mult_vec; 22f6ce2b0aSJames Wright PetscMemType m_mem_type; 23f6ce2b0aSJames Wright Vec Multiplicity, Multiplicity_loc; 24f6ce2b0aSJames Wright 2577841947SLeila Ghaffari PetscFunctionBeginUser; 26a424bcd0SJames Wright if (user->phys->ics_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time)); 275263e9c6SJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx)); 2877841947SLeila Ghaffari 29a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL)); 3077841947SLeila Ghaffari 3177841947SLeila Ghaffari // -- Get multiplicity 32f6ce2b0aSJames Wright PetscCall(DMGetLocalVector(dm, &Multiplicity_loc)); 33d0593705SJames Wright PetscCall(VecPetscToCeed(Multiplicity_loc, &m_mem_type, mult_vec)); 34a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec)); 35d0593705SJames Wright PetscCall(VecCeedToPetsc(mult_vec, m_mem_type, Multiplicity_loc)); 3677841947SLeila Ghaffari 37f6ce2b0aSJames Wright PetscCall(DMGetGlobalVector(dm, &Multiplicity)); 38f6ce2b0aSJames Wright PetscCall(VecZeroEntries(Multiplicity)); 39f6ce2b0aSJames Wright PetscCall(DMLocalToGlobal(dm, Multiplicity_loc, ADD_VALUES, Multiplicity)); 4077841947SLeila Ghaffari 4177841947SLeila Ghaffari // -- Fix multiplicity 42f6ce2b0aSJames Wright PetscCall(VecPointwiseDivide(Q, Q, Multiplicity)); 43f6ce2b0aSJames Wright PetscCall(VecPointwiseDivide(Q_loc, Q_loc, Multiplicity_loc)); 4477841947SLeila Ghaffari 45f6ce2b0aSJames Wright PetscCall(DMRestoreLocalVector(dm, &Multiplicity_loc)); 46f6ce2b0aSJames Wright PetscCall(DMRestoreGlobalVector(dm, &Multiplicity)); 47a424bcd0SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&mult_vec)); 48ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4977841947SLeila Ghaffari } 5077841947SLeila Ghaffari 513c9e7ad1SJames Wright // Record boundary values from initial condition 523c9e7ad1SJames Wright PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc) { 533c9e7ad1SJames Wright PetscFunctionBeginUser; 54251425b7SJames Wright { // Capture initial condition values in Qbc 55251425b7SJames Wright Vec Qbc; 56251425b7SJames Wright 573c9e7ad1SJames Wright PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 583c9e7ad1SJames Wright PetscCall(VecCopy(Q_loc, Qbc)); 593c9e7ad1SJames Wright PetscCall(VecZeroEntries(Q_loc)); 603c9e7ad1SJames Wright PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc)); 613c9e7ad1SJames Wright PetscCall(VecAXPY(Qbc, -1., Q_loc)); 623c9e7ad1SJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 63251425b7SJames Wright } 643c9e7ad1SJames Wright PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_FromICs)); 653c9e7ad1SJames Wright 66251425b7SJames Wright { // Set boundary mask to zero out essential BCs 67251425b7SJames Wright Vec boundary_mask, ones; 68251425b7SJames Wright 693c9e7ad1SJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 70251425b7SJames Wright PetscCall(DMGetGlobalVector(dm, &ones)); 713c9e7ad1SJames Wright PetscCall(VecZeroEntries(boundary_mask)); 72251425b7SJames Wright PetscCall(VecSet(ones, 1.0)); 73251425b7SJames Wright PetscCall(DMGlobalToLocal(dm, ones, INSERT_VALUES, boundary_mask)); 743c9e7ad1SJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 75251425b7SJames Wright PetscCall(DMRestoreGlobalVector(dm, &ones)); 76251425b7SJames Wright } 773c9e7ad1SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 783c9e7ad1SJames Wright } 793c9e7ad1SJames Wright 803c9e7ad1SJames Wright PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 812b730f8bSJeremy L Thompson Vec grad_FVM) { 825571c6fdSJames Wright Vec Qbc, boundary_mask; 8377841947SLeila Ghaffari 84f17d818dSJames Wright PetscFunctionBeginUser; 853722cd23SJames Wright // Mask (zero) Strong BC entries 865571c6fdSJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 875571c6fdSJames Wright PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask)); 885571c6fdSJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 895571c6fdSJames Wright 902b730f8bSJeremy L Thompson PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 912b730f8bSJeremy L Thompson PetscCall(VecAXPY(Q_loc, 1., Qbc)); 922b730f8bSJeremy L Thompson PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 93ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 9477841947SLeila Ghaffari } 9577841947SLeila Ghaffari 967ceb6423SJed Brown static PetscErrorCode BinaryReadIntoInt(PetscViewer viewer, PetscInt *out, PetscDataType file_type) { 977ceb6423SJed Brown PetscFunctionBeginUser; 987ceb6423SJed Brown if (file_type == PETSC_INT32) { 997ceb6423SJed Brown PetscInt32 val; 1007ceb6423SJed Brown PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT32)); 1017ceb6423SJed Brown *out = val; 1027ceb6423SJed Brown } else if (file_type == PETSC_INT64) { 1037ceb6423SJed Brown PetscInt64 val; 1047ceb6423SJed Brown PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT64)); 1057ceb6423SJed Brown *out = val; 1067ceb6423SJed Brown } else { 1077ceb6423SJed Brown PetscCall(PetscViewerBinaryRead(viewer, out, 1, NULL, PETSC_INT)); 1087ceb6423SJed Brown } 1097ceb6423SJed Brown PetscFunctionReturn(PETSC_SUCCESS); 1107ceb6423SJed Brown } 1117ceb6423SJed Brown 112530ad8c4SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number 113530ad8c4SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) { 1140de6a49fSJames Wright PetscInt file_step_number; 1150de6a49fSJames Wright PetscInt32 token; 116530ad8c4SKenneth E. Jansen PetscReal file_time; 1177ceb6423SJed Brown PetscDataType file_type = PETSC_INT32; 118530ad8c4SKenneth E. Jansen 119f17d818dSJames Wright PetscFunctionBeginUser; 1200de6a49fSJames Wright PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32)); 1210de6a49fSJames Wright if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 || 1220de6a49fSJames Wright token == FLUIDS_FILE_TOKEN) { // New style format; we're reading a file with step number and time in the header 1237ceb6423SJed Brown if (token == FLUIDS_FILE_TOKEN_32) file_type = PETSC_INT32; 1247ceb6423SJed Brown else if (token == FLUIDS_FILE_TOKEN_64) file_type = PETSC_INT64; 1257ceb6423SJed Brown PetscCall(BinaryReadIntoInt(viewer, &file_step_number, file_type)); 126530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL)); 127530ad8c4SKenneth E. Jansen if (time) *time = file_time; 128530ad8c4SKenneth E. Jansen if (step_number) *step_number = file_step_number; 129530ad8c4SKenneth E. Jansen } else if (token == VEC_FILE_CLASSID) { // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ] 130530ad8c4SKenneth E. Jansen PetscInt length, N; 1317ceb6423SJed Brown PetscCall(BinaryReadIntoInt(viewer, &length, file_type)); 132530ad8c4SKenneth E. Jansen PetscCall(VecGetSize(Q, &N)); 133530ad8c4SKenneth 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); 134530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE)); 135530ad8c4SKenneth E. Jansen } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file"); 136530ad8c4SKenneth E. Jansen 137530ad8c4SKenneth E. Jansen PetscCall(VecLoad(Q, viewer)); 138ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 139530ad8c4SKenneth E. Jansen } 140530ad8c4SKenneth E. Jansen 14177841947SLeila Ghaffari // Compare reference solution values with current test run for CI 1423c9e7ad1SJames Wright PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) { 14377841947SLeila Ghaffari Vec Qref; 14477841947SLeila Ghaffari PetscViewer viewer; 14577841947SLeila Ghaffari PetscReal error, Qrefnorm; 146530ad8c4SKenneth E. Jansen MPI_Comm comm = PetscObjectComm((PetscObject)Q); 14777841947SLeila Ghaffari 148f17d818dSJames Wright PetscFunctionBeginUser; 14977841947SLeila Ghaffari // Read reference file 1502b730f8bSJeremy L Thompson PetscCall(VecDuplicate(Q, &Qref)); 151013a5551SJames Wright PetscCheck(strcmp(app_ctx->test_file_path, "") != 0, comm, PETSC_ERR_FILE_READ, "File for regression test not given"); 152530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer)); 153530ad8c4SKenneth E. Jansen PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL)); 15477841947SLeila Ghaffari 15577841947SLeila Ghaffari // Compute error with respect to reference solution 1562b730f8bSJeremy L Thompson PetscCall(VecAXPY(Q, -1.0, Qref)); 1572b730f8bSJeremy L Thompson PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm)); 1582b730f8bSJeremy L Thompson PetscCall(VecScale(Q, 1. / Qrefnorm)); 1592b730f8bSJeremy L Thompson PetscCall(VecNorm(Q, NORM_MAX, &error)); 16077841947SLeila Ghaffari 16177841947SLeila Ghaffari // Check error 16277841947SLeila Ghaffari if (error > app_ctx->test_tol) { 1632b730f8bSJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error)); 16477841947SLeila Ghaffari } 16577841947SLeila Ghaffari 16677841947SLeila Ghaffari // Cleanup 1672b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 1682b730f8bSJeremy L Thompson PetscCall(VecDestroy(&Qref)); 169ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 17077841947SLeila Ghaffari } 17177841947SLeila Ghaffari 17277841947SLeila Ghaffari // Get error for problems with exact solutions 1733c9e7ad1SJames Wright PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) { 17477841947SLeila Ghaffari PetscInt loc_nodes; 17577841947SLeila Ghaffari Vec Q_exact, Q_exact_loc; 17677841947SLeila Ghaffari PetscReal rel_error, norm_error, norm_exact; 17777841947SLeila Ghaffari 178f17d818dSJames Wright PetscFunctionBeginUser; 17977841947SLeila Ghaffari // Get exact solution at final time 180f6ce2b0aSJames Wright PetscCall(DMGetGlobalVector(dm, &Q_exact)); 1812b730f8bSJeremy L Thompson PetscCall(DMGetLocalVector(dm, &Q_exact_loc)); 1822b730f8bSJeremy L Thompson PetscCall(VecGetSize(Q_exact_loc, &loc_nodes)); 1832b730f8bSJeremy L Thompson PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time)); 18477841947SLeila Ghaffari 18577841947SLeila Ghaffari // Get |exact solution - obtained solution| 1862b730f8bSJeremy L Thompson PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact)); 1872b730f8bSJeremy L Thompson PetscCall(VecAXPY(Q, -1.0, Q_exact)); 1882b730f8bSJeremy L Thompson PetscCall(VecNorm(Q, NORM_1, &norm_error)); 18977841947SLeila Ghaffari 19077841947SLeila Ghaffari rel_error = norm_error / norm_exact; 1912b730f8bSJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error)); 1922b730f8bSJeremy L Thompson PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc)); 193f6ce2b0aSJames Wright PetscCall(DMRestoreGlobalVector(dm, &Q_exact)); 194ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 19577841947SLeila Ghaffari } 19677841947SLeila Ghaffari 19777841947SLeila Ghaffari // Post-processing 198731c13d7SJames Wright PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData problem, User user, Vec Q, PetscScalar final_time) { 19977841947SLeila Ghaffari PetscInt steps; 200cf7a0454SJed Brown TSConvergedReason reason; 20177841947SLeila Ghaffari 202f17d818dSJames Wright PetscFunctionBeginUser; 20377841947SLeila Ghaffari // Print relative error 204de327db4SJames Wright if (problem->compute_exact_solution_error && user->app_ctx->test_type == TESTTYPE_NONE) { 2053c9e7ad1SJames Wright PetscCall(PrintError(ceed_data, dm, user, Q, final_time)); 20677841947SLeila Ghaffari } 20777841947SLeila Ghaffari 20877841947SLeila Ghaffari // Print final time and number of steps 2092b730f8bSJeremy L Thompson PetscCall(TSGetStepNumber(ts, &steps)); 210cf7a0454SJed Brown PetscCall(TSGetConvergedReason(ts, &reason)); 2118fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_NONE) { 212cf7a0454SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason], 213cf7a0454SJed Brown steps, (double)final_time)); 21477841947SLeila Ghaffari } 21577841947SLeila Ghaffari 21677841947SLeila Ghaffari // Output numerical values from command line 2172b730f8bSJeremy L Thompson PetscCall(VecViewFromOptions(Q, NULL, "-vec_view")); 21877841947SLeila Ghaffari 21977841947SLeila Ghaffari // Compare reference solution values with current test run for CI 2208fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_SOLVER) { 2213c9e7ad1SJames Wright PetscCall(RegressionTest(user->app_ctx, Q)); 22277841947SLeila Ghaffari } 223ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 22477841947SLeila Ghaffari } 22577841947SLeila Ghaffari 2260de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN = 0xceedf00; // for backwards compatibility 2270de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32; 2280de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64; 2294de8550aSJed Brown 23077841947SLeila Ghaffari // Gather initial Q values in case of continuation of simulation 23177841947SLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) { 23277841947SLeila Ghaffari PetscViewer viewer; 2332b730f8bSJeremy L Thompson 234f17d818dSJames Wright PetscFunctionBeginUser; 2352b730f8bSJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer)); 236530ad8c4SKenneth E. Jansen PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps)); 2372b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 238ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 23977841947SLeila Ghaffari } 24077841947SLeila Ghaffari 241841e4c73SJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes 242841e4c73SJed Brown int FreeContextPetsc(void *data) { 2432b730f8bSJeremy L Thompson if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed"); 244841e4c73SJed Brown return CEED_ERROR_SUCCESS; 245841e4c73SJed Brown } 246ef080ff9SJames Wright 247ef080ff9SJames Wright // Return mass qfunction specification for number of components N 248ef080ff9SJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) { 249ef080ff9SJames Wright PetscFunctionBeginUser; 250ef080ff9SJames Wright switch (N) { 251ef080ff9SJames Wright case 1: 252a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf)); 253ef080ff9SJames Wright break; 254ef080ff9SJames Wright case 5: 255a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf)); 256ef080ff9SJames Wright break; 257052409adSJames Wright case 7: 258a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf)); 259052409adSJames Wright break; 260ef080ff9SJames Wright case 9: 261a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf)); 262ef080ff9SJames Wright break; 263ef080ff9SJames Wright case 22: 264a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf)); 265ef080ff9SJames Wright break; 266ef080ff9SJames Wright default: 26783ae4962SJames Wright SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N); 268ef080ff9SJames Wright } 269ef080ff9SJames Wright 270a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP)); 271a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE)); 272a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP)); 2737f2a9303SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf, N)); 274ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 275ef080ff9SJames Wright } 2760df12fefSJames Wright 277999ff5c7SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) { 278999ff5c7SJames Wright PetscFunctionBeginUser; 279ee4ca9cbSJames Wright if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS); 280999ff5c7SJames Wright 281999ff5c7SJames Wright PetscCall(DMDestroy(&context->dm)); 282999ff5c7SJames Wright PetscCall(KSPDestroy(&context->ksp)); 283999ff5c7SJames Wright 284999ff5c7SJames Wright PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx)); 285999ff5c7SJames Wright 286999ff5c7SJames Wright PetscCall(PetscFree(context)); 287ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 288999ff5c7SJames Wright } 289999ff5c7SJames Wright 2908b892a05SJames Wright /* 2918b892a05SJames Wright * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer 2928b892a05SJames Wright * 2938b892a05SJames Wright * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`. 2948b892a05SJames Wright * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example). 2958b892a05SJames Wright * 2968b892a05SJames 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. 2978b892a05SJames Wright * 2988b892a05SJames Wright * @param[in] comm MPI_Comm for the program 2998b892a05SJames Wright * @param[in] path Path to the file 3008b892a05SJames Wright * @param[in] char_array_len Length of the character array that should contain each line 3018b892a05SJames Wright * @param[out] dims Dimensions of the file, taken from the first line of the file 3028b892a05SJames Wright * @param[out] fp File pointer to the opened file 3038b892a05SJames Wright */ 304cbef7084SJames Wright PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2], 3058b892a05SJames Wright FILE **fp) { 306457e73b2SJames Wright int ndims; 3078b892a05SJames Wright char line[char_array_len]; 3088b892a05SJames Wright char **array; 3098b892a05SJames Wright 3108b892a05SJames Wright PetscFunctionBeginUser; 3118b892a05SJames Wright PetscCall(PetscFOpen(comm, path, "r", fp)); 3128b892a05SJames Wright PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line)); 3138b892a05SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 314457e73b2SJames Wright PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path); 3158b892a05SJames Wright 3168b892a05SJames Wright for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]); 3178b892a05SJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array)); 318ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3198b892a05SJames Wright } 3208b892a05SJames Wright 3218b892a05SJames Wright /* 3228b892a05SJames Wright * @brief Get the number of rows for the PHASTA file at path. 3238b892a05SJames Wright * 3248b892a05SJames 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. 3258b892a05SJames Wright * 3268b892a05SJames Wright * @param[in] comm MPI_Comm for the program 3278b892a05SJames Wright * @param[in] path Path to the file 3288b892a05SJames Wright * @param[out] nrows Number of rows 3298b892a05SJames Wright */ 330cbef7084SJames Wright PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) { 3318b892a05SJames Wright const PetscInt char_array_len = 512; 3328b892a05SJames Wright PetscInt dims[2]; 3338b892a05SJames Wright FILE *fp; 3348b892a05SJames Wright 3358b892a05SJames Wright PetscFunctionBeginUser; 336cbef7084SJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 3378b892a05SJames Wright *nrows = dims[0]; 3388b892a05SJames Wright PetscCall(PetscFClose(comm, fp)); 339ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3408b892a05SJames Wright } 3414cc9442bSJames Wright 342cbef7084SJames Wright PetscErrorCode PhastaDatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) { 343457e73b2SJames Wright PetscInt dims[2]; 3444cc9442bSJames Wright FILE *fp; 3454cc9442bSJames Wright const PetscInt char_array_len = 512; 3464cc9442bSJames Wright char line[char_array_len]; 3474cc9442bSJames Wright 348f17d818dSJames Wright PetscFunctionBeginUser; 349cbef7084SJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 3504cc9442bSJames Wright 3514cc9442bSJames Wright for (PetscInt i = 0; i < dims[0]; i++) { 352*38690fecSJames Wright int ndims; 353*38690fecSJames Wright char **row_array; 354*38690fecSJames Wright 3554cc9442bSJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 3564cc9442bSJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array)); 3570e654f56SJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 358457e73b2SJames Wright "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 3594cc9442bSJames Wright 360*38690fecSJames Wright for (PetscInt j = 0; j < dims[1]; j++) array[i * dims[1] + j] = (PetscReal)atof(row_array[j]); 361*38690fecSJames Wright PetscCall(PetscStrToArrayDestroy(ndims, row_array)); 3624cc9442bSJames Wright } 3634cc9442bSJames Wright 3644cc9442bSJames Wright PetscCall(PetscFClose(comm, fp)); 365ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3664cc9442bSJames Wright } 36775d1798cSJames Wright 3682e31f396SJames Wright // Print information about the given simulation run 3697bc7b61fSJames Wright PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData problem, TS ts) { 370a424bcd0SJames Wright Ceed ceed = user->ceed; 3717bc7b61fSJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)ts); 3727bc7b61fSJames Wright 3732e31f396SJames Wright PetscFunctionBeginUser; 3742e31f396SJames Wright // Header and rank 3752e31f396SJames Wright char host_name[PETSC_MAX_PATH_LEN]; 3762e31f396SJames Wright PetscMPIInt rank, comm_size; 3772e31f396SJames Wright PetscCall(PetscGetHostName(host_name, sizeof host_name)); 3782e31f396SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank)); 3792e31f396SJames Wright PetscCallMPI(MPI_Comm_size(comm, &comm_size)); 3802e31f396SJames Wright PetscCall(PetscPrintf(comm, 3812e31f396SJames Wright "\n-- Navier-Stokes solver - libCEED + PETSc --\n" 3822e31f396SJames Wright " MPI:\n" 3832e31f396SJames Wright " Host Name : %s\n" 3842e31f396SJames Wright " Total ranks : %d\n", 3852e31f396SJames Wright host_name, comm_size)); 3862e31f396SJames Wright 3872e31f396SJames Wright // Problem specific info 388367c849eSJames Wright PetscCall(problem->print_info(user, problem, user->app_ctx)); 3892e31f396SJames Wright 3902e31f396SJames Wright // libCEED 3912e31f396SJames Wright const char *used_resource; 3922e31f396SJames Wright CeedMemType mem_type_backend; 393a424bcd0SJames Wright PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource)); 394a424bcd0SJames Wright PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend)); 3952e31f396SJames Wright PetscCall(PetscPrintf(comm, 3962e31f396SJames Wright " libCEED:\n" 3972e31f396SJames Wright " libCEED Backend : %s\n" 3982e31f396SJames Wright " libCEED Backend MemType : %s\n", 3992e31f396SJames Wright used_resource, CeedMemTypes[mem_type_backend])); 4002e31f396SJames Wright // PETSc 4017bc7b61fSJames Wright VecType vec_type; 4022e31f396SJames Wright char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3"; 4032e31f396SJames Wright if (problem->dim == 2) box_faces_str[3] = '\0'; 4042e31f396SJames Wright PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL)); 4052e31f396SJames Wright PetscCall(DMGetVecType(user->dm, &vec_type)); 4062e31f396SJames Wright PetscCall(PetscPrintf(comm, 4072e31f396SJames Wright " PETSc:\n" 4082e31f396SJames Wright " Box Faces : %s\n" 4092e31f396SJames Wright " DM VecType : %s\n" 4102e31f396SJames Wright " Time Stepping Scheme : %s\n", 4117bc7b61fSJames Wright box_faces_str, vec_type, phys_ctx->implicit ? "implicit" : "explicit")); 4127bc7b61fSJames Wright { 4137bc7b61fSJames Wright char pmat_type_str[PETSC_MAX_PATH_LEN]; 4147bc7b61fSJames Wright MatType amat_type, pmat_type; 4157bc7b61fSJames Wright Mat Amat, Pmat; 4167bc7b61fSJames Wright TSIJacobianFn *ijacob_function; 4177bc7b61fSJames Wright 4187bc7b61fSJames Wright PetscCall(TSGetIJacobian(ts, &Amat, &Pmat, &ijacob_function, NULL)); 4197bc7b61fSJames Wright PetscCall(MatGetType(Amat, &amat_type)); 4207bc7b61fSJames Wright PetscCall(MatGetType(Pmat, &pmat_type)); 4217bc7b61fSJames Wright 4227bc7b61fSJames Wright PetscCall(PetscStrncpy(pmat_type_str, pmat_type, sizeof(pmat_type_str))); 4237bc7b61fSJames Wright if (!strcmp(pmat_type, MATCEED)) { 4247bc7b61fSJames Wright MatType pmat_coo_type; 4257bc7b61fSJames Wright char pmat_coo_type_str[PETSC_MAX_PATH_LEN]; 4267bc7b61fSJames Wright 4277bc7b61fSJames Wright PetscCall(MatCeedGetCOOMatType(Pmat, &pmat_coo_type)); 4287bc7b61fSJames Wright PetscCall(PetscSNPrintf(pmat_coo_type_str, sizeof(pmat_coo_type_str), " (COO MatType: %s)", pmat_coo_type)); 4297bc7b61fSJames Wright PetscCall(PetscStrlcat(pmat_type_str, pmat_coo_type_str, sizeof(pmat_type_str))); 4307bc7b61fSJames Wright } 4317bc7b61fSJames Wright if (ijacob_function) { 4327bc7b61fSJames Wright PetscCall(PetscPrintf(comm, 4337bc7b61fSJames Wright " IJacobian A MatType : %s\n" 4347bc7b61fSJames Wright " IJacobian P MatType : %s\n", 4357bc7b61fSJames Wright amat_type, pmat_type_str)); 4367bc7b61fSJames Wright } 4377bc7b61fSJames Wright } 4382e31f396SJames Wright if (user->app_ctx->cont_steps) { 4392e31f396SJames Wright PetscCall(PetscPrintf(comm, 4402e31f396SJames Wright " Continue:\n" 4412e31f396SJames Wright " Filename: : %s\n" 4422e31f396SJames Wright " Step: : %" PetscInt_FMT "\n" 4432e31f396SJames Wright " Time: : %g\n", 4442e31f396SJames Wright user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time)); 4452e31f396SJames Wright } 4462e31f396SJames Wright // Mesh 4472e31f396SJames Wright const PetscInt num_comp_q = 5; 4482e31f396SJames Wright PetscInt glob_dofs, owned_dofs, local_dofs; 4492e31f396SJames Wright const CeedInt num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra; 4502e31f396SJames Wright PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL)); 4512e31f396SJames Wright PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL)); 4522e31f396SJames Wright PetscCall(PetscPrintf(comm, 4532e31f396SJames Wright " Mesh:\n" 4542e31f396SJames Wright " Number of 1D Basis Nodes (P) : %" CeedInt_FMT "\n" 4552e31f396SJames Wright " Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n" 4562e31f396SJames Wright " Global DoFs : %" PetscInt_FMT "\n" 4572e31f396SJames Wright " DoFs per node : %" PetscInt_FMT "\n" 458ce11f295SJames Wright " Global %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT "\n", 459ce11f295SJames Wright num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q)); 4602e31f396SJames Wright // -- Get Partition Statistics 4612e31f396SJames Wright PetscCall(PetscPrintf(comm, " Partition: (min,max,median,max/median)\n")); 4622e31f396SJames Wright { 4632e31f396SJames Wright PetscInt *gather_buffer = NULL; 464ce11f295SJames Wright PetscInt part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3]; 4652e31f396SJames Wright PetscInt median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1; 4662e31f396SJames Wright if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer)); 4672e31f396SJames Wright 468ce11f295SJames Wright PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 4692e31f396SJames Wright if (!rank) { 4702e31f396SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 4712e31f396SJames Wright part_owned_dofs[0] = gather_buffer[0]; // min 4722e31f396SJames Wright part_owned_dofs[1] = gather_buffer[comm_size - 1]; // max 4732e31f396SJames Wright part_owned_dofs[2] = gather_buffer[median_index]; // median 4742e31f396SJames Wright PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2]; 475ce11f295SJames Wright PetscCall(PetscPrintf( 476ce11f295SJames Wright comm, " Global Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 4772e31f396SJames Wright part_owned_dofs[0] / num_comp_q, part_owned_dofs[1] / num_comp_q, part_owned_dofs[2] / num_comp_q, part_owned_dof_ratio)); 4782e31f396SJames Wright } 4792e31f396SJames Wright 480ce11f295SJames Wright PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 4812e31f396SJames Wright if (!rank) { 4822e31f396SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 4832e31f396SJames Wright part_local_dofs[0] = gather_buffer[0]; // min 4842e31f396SJames Wright part_local_dofs[1] = gather_buffer[comm_size - 1]; // max 4852e31f396SJames Wright part_local_dofs[2] = gather_buffer[median_index]; // median 4862e31f396SJames Wright PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2]; 487ce11f295SJames Wright PetscCall(PetscPrintf( 488ce11f295SJames Wright comm, " Local Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 4892e31f396SJames Wright part_local_dofs[0] / num_comp_q, part_local_dofs[1] / num_comp_q, part_local_dofs[2] / num_comp_q, part_local_dof_ratio)); 4902e31f396SJames Wright } 4912e31f396SJames Wright 49265ba01baSJames Wright if (comm_size != 1) { 493ce11f295SJames Wright PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0; 4942e31f396SJames Wright { 4952e31f396SJames Wright PetscSF sf; 496ce11f295SJames Wright PetscInt nrranks, niranks; 497ce11f295SJames Wright const PetscInt *roffset, *rmine, *rremote, *ioffset, *irootloc; 498ce11f295SJames Wright const PetscMPIInt *rranks, *iranks; 4992e31f396SJames Wright PetscCall(DMGetSectionSF(user->dm, &sf)); 5002e31f396SJames Wright PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote)); 501ce11f295SJames Wright PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc)); 5022e31f396SJames Wright for (PetscInt i = 0; i < nrranks; i++) { 5032e31f396SJames Wright if (rranks[i] == rank) continue; // Ignore same-part global->local transfers 5042e31f396SJames Wright num_remote_roots_total += roffset[i + 1] - roffset[i]; 505ce11f295SJames Wright num_ghost_interface_ranks++; 506ce11f295SJames Wright } 507ce11f295SJames Wright for (PetscInt i = 0; i < niranks; i++) { 508ce11f295SJames Wright if (iranks[i] == rank) continue; 509ce11f295SJames Wright num_remote_leaves_total += ioffset[i + 1] - ioffset[i]; 510ce11f295SJames Wright num_owned_interface_ranks++; 5112e31f396SJames Wright } 5122e31f396SJames Wright } 513ce11f295SJames Wright PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 5142e31f396SJames Wright if (!rank) { 5152e31f396SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 516ce11f295SJames Wright part_boundary_dofs[0] = gather_buffer[0]; // min 517ce11f295SJames Wright part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 518ce11f295SJames Wright part_boundary_dofs[2] = gather_buffer[median_index]; // median 519ce11f295SJames Wright PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 520ce11f295SJames Wright PetscCall(PetscPrintf( 52165ba01baSJames Wright comm, " Ghost Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 52265ba01baSJames Wright num_comp_q, part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q, 52365ba01baSJames Wright part_shared_dof_ratio)); 524ce11f295SJames Wright } 525ce11f295SJames Wright 526ce11f295SJames Wright PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 527ce11f295SJames Wright if (!rank) { 528ce11f295SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 529ce11f295SJames Wright part_neighbors[0] = gather_buffer[0]; // min 530ce11f295SJames Wright part_neighbors[1] = gather_buffer[comm_size - 1]; // max 531ce11f295SJames Wright part_neighbors[2] = gather_buffer[median_index]; // median 532ce11f295SJames Wright PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 533ce11f295SJames Wright PetscCall(PetscPrintf(comm, " Ghost Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 534ce11f295SJames Wright part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 535ce11f295SJames Wright } 536ce11f295SJames Wright 537ce11f295SJames Wright PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 538ce11f295SJames Wright if (!rank) { 539ce11f295SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 540ce11f295SJames Wright part_boundary_dofs[0] = gather_buffer[0]; // min 541ce11f295SJames Wright part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 542ce11f295SJames Wright part_boundary_dofs[2] = gather_buffer[median_index]; // median 543ce11f295SJames Wright PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 544ce11f295SJames Wright PetscCall(PetscPrintf( 54565ba01baSJames Wright comm, " Owned Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 54665ba01baSJames Wright num_comp_q, part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q, 54765ba01baSJames Wright part_shared_dof_ratio)); 548ce11f295SJames Wright } 549ce11f295SJames Wright 550ce11f295SJames Wright PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 551ce11f295SJames Wright if (!rank) { 552ce11f295SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 553ce11f295SJames Wright part_neighbors[0] = gather_buffer[0]; // min 554ce11f295SJames Wright part_neighbors[1] = gather_buffer[comm_size - 1]; // max 555ce11f295SJames Wright part_neighbors[2] = gather_buffer[median_index]; // median 556ce11f295SJames Wright PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 557ce11f295SJames Wright PetscCall(PetscPrintf(comm, " Owned Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 558ce11f295SJames Wright part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 5592e31f396SJames Wright } 56065ba01baSJames Wright } 5612e31f396SJames Wright 5622e31f396SJames Wright if (!rank) PetscCall(PetscFree(gather_buffer)); 5632e31f396SJames Wright } 5642e31f396SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 5652e31f396SJames Wright } 566