1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 /// @file 9 /// Miscellaneous utility functions 10 11 #include "../navierstokes.h" 12 #include "../qfunctions/mass.h" 13 14 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) { 15 PetscFunctionBeginUser; 16 17 // --------------------------------------------------------------------------- 18 // Update time for evaluation 19 // --------------------------------------------------------------------------- 20 if (user->phys->ics_time_label) CeedOperatorSetContextDouble(ceed_data->op_ics, user->phys->ics_time_label, &time); 21 22 // --------------------------------------------------------------------------- 23 // ICs 24 // --------------------------------------------------------------------------- 25 // -- CEED Restriction 26 CeedVector q0_ceed; 27 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &q0_ceed, NULL); 28 29 // -- Place PETSc vector in CEED vector 30 PetscMemType q0_mem_type; 31 PetscCall(VecP2C(Q_loc, &q0_mem_type, q0_ceed)); 32 33 // -- Apply CEED Operator 34 CeedOperatorApply(ceed_data->op_ics, ceed_data->x_coord, q0_ceed, CEED_REQUEST_IMMEDIATE); 35 36 // -- Restore vectors 37 PetscCall(VecC2P(q0_ceed, q0_mem_type, Q_loc)); 38 39 // -- Local-to-Global 40 PetscCall(VecZeroEntries(Q)); 41 PetscCall(DMLocalToGlobal(dm, Q_loc, ADD_VALUES, Q)); 42 43 // --------------------------------------------------------------------------- 44 // Fix multiplicity for output of ICs 45 // --------------------------------------------------------------------------- 46 // -- CEED Restriction 47 CeedVector mult_vec; 48 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL); 49 50 // -- Place PETSc vector in CEED vector 51 PetscMemType m_mem_type; 52 Vec multiplicity_loc; 53 PetscCall(DMGetLocalVector(dm, &multiplicity_loc)); 54 PetscCall(VecP2C(multiplicity_loc, &m_mem_type, mult_vec)); 55 56 // -- Get multiplicity 57 CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec); 58 59 // -- Restore vectors 60 PetscCall(VecC2P(mult_vec, m_mem_type, multiplicity_loc)); 61 62 // -- Local-to-Global 63 Vec multiplicity; 64 PetscCall(DMGetGlobalVector(dm, &multiplicity)); 65 PetscCall(VecZeroEntries(multiplicity)); 66 PetscCall(DMLocalToGlobal(dm, multiplicity_loc, ADD_VALUES, multiplicity)); 67 68 // -- Fix multiplicity 69 PetscCall(VecPointwiseDivide(Q, Q, multiplicity)); 70 PetscCall(VecPointwiseDivide(Q_loc, Q_loc, multiplicity_loc)); 71 72 // -- Restore vectors 73 PetscCall(DMRestoreLocalVector(dm, &multiplicity_loc)); 74 PetscCall(DMRestoreGlobalVector(dm, &multiplicity)); 75 76 // Cleanup 77 CeedVectorDestroy(&mult_vec); 78 CeedVectorDestroy(&q0_ceed); 79 80 PetscFunctionReturn(0); 81 } 82 83 PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 84 Vec grad_FVM) { 85 Vec Qbc, boundary_mask; 86 PetscFunctionBegin; 87 88 // Mask (zero) Dirichlet entries 89 PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 90 PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask)); 91 PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 92 93 PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 94 PetscCall(VecAXPY(Q_loc, 1., Qbc)); 95 PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 96 97 PetscFunctionReturn(0); 98 } 99 100 // @brief Load vector from binary file, possibly with embedded solution time and step number 101 PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) { 102 PetscInt token, file_step_number; 103 PetscReal file_time; 104 PetscFunctionBeginUser; 105 106 // Attempt 107 PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT)); 108 if (token == FLUIDS_FILE_TOKEN) { // New style format; we're reading a file with step number and time in the header 109 PetscCall(PetscViewerBinaryRead(viewer, &file_step_number, 1, NULL, PETSC_INT)); 110 PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL)); 111 if (time) *time = file_time; 112 if (step_number) *step_number = file_step_number; 113 } else if (token == VEC_FILE_CLASSID) { // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ] 114 PetscInt length, N; 115 PetscCall(PetscViewerBinaryRead(viewer, &length, 1, NULL, PETSC_INT)); 116 PetscCall(VecGetSize(Q, &N)); 117 PetscCheck(length == N, comm, PETSC_ERR_ARG_INCOMP, "File Vec has length %" PetscInt_FMT " but DM has global Vec size %" PetscInt_FMT, length, N); 118 PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE)); 119 } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file"); 120 121 // Load Q from existent solution 122 PetscCall(VecLoad(Q, viewer)); 123 124 PetscFunctionReturn(0); 125 } 126 127 // Compare reference solution values with current test run for CI 128 PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q) { 129 Vec Qref; 130 PetscViewer viewer; 131 PetscReal error, Qrefnorm; 132 MPI_Comm comm = PetscObjectComm((PetscObject)Q); 133 PetscFunctionBegin; 134 135 // Read reference file 136 PetscCall(VecDuplicate(Q, &Qref)); 137 PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer)); 138 PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL)); 139 140 // Compute error with respect to reference solution 141 PetscCall(VecAXPY(Q, -1.0, Qref)); 142 PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm)); 143 PetscCall(VecScale(Q, 1. / Qrefnorm)); 144 PetscCall(VecNorm(Q, NORM_MAX, &error)); 145 146 // Check error 147 if (error > app_ctx->test_tol) { 148 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error)); 149 } 150 151 // Cleanup 152 PetscCall(PetscViewerDestroy(&viewer)); 153 PetscCall(VecDestroy(&Qref)); 154 155 PetscFunctionReturn(0); 156 } 157 158 // Get error for problems with exact solutions 159 PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) { 160 PetscInt loc_nodes; 161 Vec Q_exact, Q_exact_loc; 162 PetscReal rel_error, norm_error, norm_exact; 163 PetscFunctionBegin; 164 165 // Get exact solution at final time 166 PetscCall(DMCreateGlobalVector(dm, &Q_exact)); 167 PetscCall(DMGetLocalVector(dm, &Q_exact_loc)); 168 PetscCall(VecGetSize(Q_exact_loc, &loc_nodes)); 169 PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time)); 170 171 // Get |exact solution - obtained solution| 172 PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact)); 173 PetscCall(VecAXPY(Q, -1.0, Q_exact)); 174 PetscCall(VecNorm(Q, NORM_1, &norm_error)); 175 176 // Compute relative error 177 rel_error = norm_error / norm_exact; 178 179 // Output relative error 180 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error)); 181 // Cleanup 182 PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc)); 183 PetscCall(VecDestroy(&Q_exact)); 184 185 PetscFunctionReturn(0); 186 } 187 188 // Post-processing 189 PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) { 190 PetscInt steps; 191 TSConvergedReason reason; 192 PetscFunctionBegin; 193 194 // Print relative error 195 if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) { 196 PetscCall(GetError_NS(ceed_data, dm, user, Q, final_time)); 197 } 198 199 // Print final time and number of steps 200 PetscCall(TSGetStepNumber(ts, &steps)); 201 PetscCall(TSGetConvergedReason(ts, &reason)); 202 if (user->app_ctx->test_type == TESTTYPE_NONE) { 203 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason], 204 steps, (double)final_time)); 205 } 206 207 // Output numerical values from command line 208 PetscCall(VecViewFromOptions(Q, NULL, "-vec_view")); 209 210 // Compare reference solution values with current test run for CI 211 if (user->app_ctx->test_type == TESTTYPE_SOLVER) { 212 PetscCall(RegressionTests_NS(user->app_ctx, Q)); 213 } 214 215 PetscFunctionReturn(0); 216 } 217 218 const PetscInt FLUIDS_FILE_TOKEN = 0xceedf00; 219 220 // Gather initial Q values in case of continuation of simulation 221 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) { 222 PetscViewer viewer; 223 224 PetscFunctionBegin; 225 226 PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer)); 227 PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps)); 228 PetscCall(PetscViewerDestroy(&viewer)); 229 230 PetscFunctionReturn(0); 231 } 232 233 // Record boundary values from initial condition 234 PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc) { 235 Vec Qbc, boundary_mask; 236 PetscFunctionBegin; 237 238 PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 239 PetscCall(VecCopy(Q_loc, Qbc)); 240 PetscCall(VecZeroEntries(Q_loc)); 241 PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc)); 242 PetscCall(VecAXPY(Qbc, -1., Q_loc)); 243 PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 244 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS)); 245 246 PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 247 PetscCall(DMGetGlobalVector(dm, &Q)); 248 PetscCall(VecZeroEntries(boundary_mask)); 249 PetscCall(VecSet(Q, 1.0)); 250 PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask)); 251 PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 252 253 PetscFunctionReturn(0); 254 } 255 256 // Free a plain data context that was allocated using PETSc; returning libCEED error codes 257 int FreeContextPetsc(void *data) { 258 if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed"); 259 return CEED_ERROR_SUCCESS; 260 } 261 262 // Return mass qfunction specification for number of components N 263 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) { 264 CeedQFunctionUser qfunction_ptr; 265 const char *qfunction_loc; 266 PetscFunctionBeginUser; 267 268 switch (N) { 269 case 1: 270 qfunction_ptr = Mass_1; 271 qfunction_loc = Mass_1_loc; 272 break; 273 case 5: 274 qfunction_ptr = Mass_5; 275 qfunction_loc = Mass_5_loc; 276 break; 277 case 9: 278 qfunction_ptr = Mass_9; 279 qfunction_loc = Mass_9_loc; 280 break; 281 case 22: 282 qfunction_ptr = Mass_22; 283 qfunction_loc = Mass_22_loc; 284 break; 285 default: 286 SETERRQ(PETSC_COMM_WORLD, -1, "Could not find mass qfunction of size %d", N); 287 } 288 CeedQFunctionCreateInterior(ceed, 1, qfunction_ptr, qfunction_loc, qf); 289 290 CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP); 291 CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE); 292 CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP); 293 PetscFunctionReturn(0); 294 } 295 296 /* @brief L^2 Projection of a source FEM function to a target FEM space 297 * 298 * To solve system using a lumped mass matrix, pass a KSP object with ksp_type=preonly, pc_type=jacobi, pc_jacobi_type=rowsum. 299 * 300 * @param[in] source_vec Global Vec of the source FEM function. NULL indicates using rhs_matop_ctx->X_loc 301 * @param[out] target_vec Global Vec of the target (result) FEM function. NULL indicates using rhs_matop_ctx->Y_loc 302 * @param[in] rhs_matop_ctx MatopApplyContext for performing the RHS evaluation 303 * @param[in] ksp KSP for solving the consistent projection problem 304 */ 305 PetscErrorCode ComputeL2Projection(Vec source_vec, Vec target_vec, MatopApplyContext rhs_matop_ctx, KSP ksp) { 306 PetscFunctionBeginUser; 307 308 PetscCall(ApplyLocal_Ceed(source_vec, target_vec, rhs_matop_ctx)); 309 PetscCall(KSPSolve(ksp, target_vec, target_vec)); 310 311 PetscFunctionReturn(0); 312 } 313