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 1177841947SLeila Ghaffari #include "../navierstokes.h" 1277841947SLeila Ghaffari 13841e4c73SJed Brown PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, 14841e4c73SJed Brown Vec Q_loc, Vec Q, 1577841947SLeila Ghaffari CeedScalar time) { 1677841947SLeila Ghaffari PetscErrorCode ierr; 1777841947SLeila Ghaffari PetscFunctionBeginUser; 1877841947SLeila Ghaffari 1977841947SLeila Ghaffari // --------------------------------------------------------------------------- 20a0add3c9SJed Brown // Update time for evaluation 2177841947SLeila Ghaffari // --------------------------------------------------------------------------- 22841e4c73SJed Brown if (user->phys->ics_time_label) 23841e4c73SJed Brown CeedOperatorContextSetDouble(ceed_data->op_ics, user->phys->ics_time_label, 24841e4c73SJed Brown &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 3477841947SLeila Ghaffari CeedScalar *q0; 3577841947SLeila Ghaffari PetscMemType q0_mem_type; 3677841947SLeila Ghaffari ierr = VecGetArrayAndMemType(Q_loc, (PetscScalar **)&q0, &q0_mem_type); 3777841947SLeila Ghaffari CHKERRQ(ierr); 3877841947SLeila Ghaffari CeedVectorSetArray(q0_ceed, MemTypeP2C(q0_mem_type), CEED_USE_POINTER, q0); 3977841947SLeila Ghaffari 4077841947SLeila Ghaffari // -- Apply CEED Operator 4177841947SLeila Ghaffari CeedOperatorApply(ceed_data->op_ics, ceed_data->x_coord, q0_ceed, 4277841947SLeila Ghaffari CEED_REQUEST_IMMEDIATE); 4377841947SLeila Ghaffari 4477841947SLeila Ghaffari // -- Restore vectors 4577841947SLeila Ghaffari CeedVectorTakeArray(q0_ceed, MemTypeP2C(q0_mem_type), NULL); 4677841947SLeila Ghaffari ierr = VecRestoreArrayReadAndMemType(Q_loc, (const PetscScalar **)&q0); 4777841947SLeila Ghaffari CHKERRQ(ierr); 4877841947SLeila Ghaffari 4977841947SLeila Ghaffari // -- Local-to-Global 5077841947SLeila Ghaffari ierr = VecZeroEntries(Q); CHKERRQ(ierr); 5177841947SLeila Ghaffari ierr = DMLocalToGlobal(dm, Q_loc, ADD_VALUES, Q); CHKERRQ(ierr); 5277841947SLeila Ghaffari 5377841947SLeila Ghaffari // --------------------------------------------------------------------------- 5477841947SLeila Ghaffari // Fix multiplicity for output of ICs 5577841947SLeila Ghaffari // --------------------------------------------------------------------------- 5677841947SLeila Ghaffari // -- CEED Restriction 5777841947SLeila Ghaffari CeedVector mult_vec; 5877841947SLeila Ghaffari CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL); 5977841947SLeila Ghaffari 6077841947SLeila Ghaffari // -- Place PETSc vector in CEED vector 6177841947SLeila Ghaffari CeedScalar *mult; 6277841947SLeila Ghaffari PetscMemType m_mem_type; 6377841947SLeila Ghaffari Vec multiplicity_loc; 6477841947SLeila Ghaffari ierr = DMGetLocalVector(dm, &multiplicity_loc); CHKERRQ(ierr); 6577841947SLeila Ghaffari ierr = VecGetArrayAndMemType(multiplicity_loc, (PetscScalar **)&mult, 6677841947SLeila Ghaffari &m_mem_type); 6777841947SLeila Ghaffari CHKERRQ(ierr); 6877841947SLeila Ghaffari CeedVectorSetArray(mult_vec, MemTypeP2C(m_mem_type), CEED_USE_POINTER, mult); 6977841947SLeila Ghaffari CHKERRQ(ierr); 7077841947SLeila Ghaffari 7177841947SLeila Ghaffari // -- Get multiplicity 7277841947SLeila Ghaffari CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec); 7377841947SLeila Ghaffari 7477841947SLeila Ghaffari // -- Restore vectors 7577841947SLeila Ghaffari CeedVectorTakeArray(mult_vec, MemTypeP2C(m_mem_type), NULL); 7677841947SLeila Ghaffari ierr = VecRestoreArrayReadAndMemType(multiplicity_loc, 7777841947SLeila Ghaffari (const PetscScalar **)&mult); CHKERRQ(ierr); 7877841947SLeila Ghaffari 7977841947SLeila Ghaffari // -- Local-to-Global 8077841947SLeila Ghaffari Vec multiplicity; 8177841947SLeila Ghaffari ierr = DMGetGlobalVector(dm, &multiplicity); CHKERRQ(ierr); 8277841947SLeila Ghaffari ierr = VecZeroEntries(multiplicity); CHKERRQ(ierr); 8377841947SLeila Ghaffari ierr = DMLocalToGlobal(dm, multiplicity_loc, ADD_VALUES, multiplicity); 8477841947SLeila Ghaffari CHKERRQ(ierr); 8577841947SLeila Ghaffari 8677841947SLeila Ghaffari // -- Fix multiplicity 8777841947SLeila Ghaffari ierr = VecPointwiseDivide(Q, Q, multiplicity); CHKERRQ(ierr); 8877841947SLeila Ghaffari ierr = VecPointwiseDivide(Q_loc, Q_loc, multiplicity_loc); CHKERRQ(ierr); 8977841947SLeila Ghaffari 9077841947SLeila Ghaffari // -- Restore vectors 9177841947SLeila Ghaffari ierr = DMRestoreLocalVector(dm, &multiplicity_loc); CHKERRQ(ierr); 9277841947SLeila Ghaffari ierr = DMRestoreGlobalVector(dm, &multiplicity); CHKERRQ(ierr); 9377841947SLeila Ghaffari 9477841947SLeila Ghaffari // Cleanup 9577841947SLeila Ghaffari CeedVectorDestroy(&mult_vec); 9677841947SLeila Ghaffari CeedVectorDestroy(&q0_ceed); 9777841947SLeila Ghaffari 9877841947SLeila Ghaffari PetscFunctionReturn(0); 9977841947SLeila Ghaffari } 10077841947SLeila Ghaffari 10177841947SLeila Ghaffari PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, 10277841947SLeila Ghaffari PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, 10377841947SLeila Ghaffari Vec cell_geom_FVM, Vec grad_FVM) { 10477841947SLeila Ghaffari 105*5571c6fdSJames Wright Vec Qbc, boundary_mask; 10677841947SLeila Ghaffari PetscErrorCode ierr; 10777841947SLeila Ghaffari PetscFunctionBegin; 10877841947SLeila Ghaffari 109*5571c6fdSJames Wright // Mask (zero) Dirichlet entries 110*5571c6fdSJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 111*5571c6fdSJames Wright PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask)); 112*5571c6fdSJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 113*5571c6fdSJames Wright 11477841947SLeila Ghaffari ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 11577841947SLeila Ghaffari ierr = VecAXPY(Q_loc, 1., Qbc); CHKERRQ(ierr); 11677841947SLeila Ghaffari ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 11777841947SLeila Ghaffari 11877841947SLeila Ghaffari PetscFunctionReturn(0); 11977841947SLeila Ghaffari } 12077841947SLeila Ghaffari 12177841947SLeila Ghaffari // Compare reference solution values with current test run for CI 12277841947SLeila Ghaffari PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q) { 12377841947SLeila Ghaffari 12477841947SLeila Ghaffari Vec Qref; 12577841947SLeila Ghaffari PetscViewer viewer; 12677841947SLeila Ghaffari PetscReal error, Qrefnorm; 12777841947SLeila Ghaffari PetscErrorCode ierr; 12877841947SLeila Ghaffari PetscFunctionBegin; 12977841947SLeila Ghaffari 13077841947SLeila Ghaffari // Read reference file 13177841947SLeila Ghaffari ierr = VecDuplicate(Q, &Qref); CHKERRQ(ierr); 13277841947SLeila Ghaffari ierr = PetscViewerBinaryOpen(PetscObjectComm((PetscObject)Q), 13377841947SLeila Ghaffari app_ctx->file_path, FILE_MODE_READ, 13477841947SLeila Ghaffari &viewer); CHKERRQ(ierr); 13577841947SLeila Ghaffari ierr = VecLoad(Qref, viewer); CHKERRQ(ierr); 13677841947SLeila Ghaffari 13777841947SLeila Ghaffari // Compute error with respect to reference solution 13877841947SLeila Ghaffari ierr = VecAXPY(Q, -1.0, Qref); CHKERRQ(ierr); 13977841947SLeila Ghaffari ierr = VecNorm(Qref, NORM_MAX, &Qrefnorm); CHKERRQ(ierr); 14077841947SLeila Ghaffari ierr = VecScale(Q, 1./Qrefnorm); CHKERRQ(ierr); 14177841947SLeila Ghaffari ierr = VecNorm(Q, NORM_MAX, &error); CHKERRQ(ierr); 14277841947SLeila Ghaffari 14377841947SLeila Ghaffari // Check error 14477841947SLeila Ghaffari if (error > app_ctx->test_tol) { 14577841947SLeila Ghaffari ierr = PetscPrintf(PETSC_COMM_WORLD, 14677841947SLeila Ghaffari "Test failed with error norm %g\n", 14777841947SLeila Ghaffari (double)error); CHKERRQ(ierr); 14877841947SLeila Ghaffari } 14977841947SLeila Ghaffari 15077841947SLeila Ghaffari // Cleanup 15177841947SLeila Ghaffari ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 15277841947SLeila Ghaffari ierr = VecDestroy(&Qref); CHKERRQ(ierr); 15377841947SLeila Ghaffari 15477841947SLeila Ghaffari PetscFunctionReturn(0); 15577841947SLeila Ghaffari } 15677841947SLeila Ghaffari 15777841947SLeila Ghaffari // Get error for problems with exact solutions 158841e4c73SJed Brown PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, 15977841947SLeila Ghaffari PetscScalar final_time) { 16077841947SLeila Ghaffari PetscInt loc_nodes; 16177841947SLeila Ghaffari Vec Q_exact, Q_exact_loc; 16277841947SLeila Ghaffari PetscReal rel_error, norm_error, norm_exact; 16377841947SLeila Ghaffari PetscErrorCode ierr; 16477841947SLeila Ghaffari PetscFunctionBegin; 16577841947SLeila Ghaffari 16677841947SLeila Ghaffari // Get exact solution at final time 16777841947SLeila Ghaffari ierr = DMCreateGlobalVector(dm, &Q_exact); CHKERRQ(ierr); 16877841947SLeila Ghaffari ierr = DMGetLocalVector(dm, &Q_exact_loc); CHKERRQ(ierr); 16977841947SLeila Ghaffari ierr = VecGetSize(Q_exact_loc, &loc_nodes); CHKERRQ(ierr); 170841e4c73SJed Brown ierr = ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, 171841e4c73SJed Brown final_time); 17277841947SLeila Ghaffari CHKERRQ(ierr); 17377841947SLeila Ghaffari 17477841947SLeila Ghaffari // Get |exact solution - obtained solution| 17577841947SLeila Ghaffari ierr = VecNorm(Q_exact, NORM_1, &norm_exact); CHKERRQ(ierr); 17677841947SLeila Ghaffari ierr = VecAXPY(Q, -1.0, Q_exact); CHKERRQ(ierr); 17777841947SLeila Ghaffari ierr = VecNorm(Q, NORM_1, &norm_error); CHKERRQ(ierr); 17877841947SLeila Ghaffari 17977841947SLeila Ghaffari // Compute relative error 18077841947SLeila Ghaffari rel_error = norm_error / norm_exact; 18177841947SLeila Ghaffari 18277841947SLeila Ghaffari // Output relative error 18377841947SLeila Ghaffari ierr = PetscPrintf(PETSC_COMM_WORLD, 18477841947SLeila Ghaffari "Relative Error: %g\n", 18577841947SLeila Ghaffari (double)rel_error); CHKERRQ(ierr); 18677841947SLeila Ghaffari // Cleanup 18777841947SLeila Ghaffari ierr = DMRestoreLocalVector(dm, &Q_exact_loc); CHKERRQ(ierr); 18877841947SLeila Ghaffari ierr = VecDestroy(&Q_exact); CHKERRQ(ierr); 18977841947SLeila Ghaffari 19077841947SLeila Ghaffari PetscFunctionReturn(0); 19177841947SLeila Ghaffari } 19277841947SLeila Ghaffari 19377841947SLeila Ghaffari // Post-processing 19477841947SLeila Ghaffari PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, 195841e4c73SJed Brown ProblemData *problem, User user, 19677841947SLeila Ghaffari Vec Q, PetscScalar final_time) { 19777841947SLeila Ghaffari PetscInt steps; 19877841947SLeila Ghaffari PetscErrorCode ierr; 19977841947SLeila Ghaffari PetscFunctionBegin; 20077841947SLeila Ghaffari 20177841947SLeila Ghaffari // Print relative error 202841e4c73SJed Brown if (problem->non_zero_time && !user->app_ctx->test_mode) { 203841e4c73SJed Brown ierr = GetError_NS(ceed_data, dm, user, Q, final_time); CHKERRQ(ierr); 20477841947SLeila Ghaffari } 20577841947SLeila Ghaffari 20677841947SLeila Ghaffari // Print final time and number of steps 20777841947SLeila Ghaffari ierr = TSGetStepNumber(ts, &steps); CHKERRQ(ierr); 208841e4c73SJed Brown if (!user->app_ctx->test_mode) { 20977841947SLeila Ghaffari ierr = PetscPrintf(PETSC_COMM_WORLD, 21008140895SJed Brown "Time integrator took %" PetscInt_FMT " time steps to reach final time %g\n", 21177841947SLeila Ghaffari steps, (double)final_time); CHKERRQ(ierr); 21277841947SLeila Ghaffari } 21377841947SLeila Ghaffari 21477841947SLeila Ghaffari // Output numerical values from command line 21577841947SLeila Ghaffari ierr = VecViewFromOptions(Q, NULL, "-vec_view"); CHKERRQ(ierr); 21677841947SLeila Ghaffari 21777841947SLeila Ghaffari // Compare reference solution values with current test run for CI 218841e4c73SJed Brown if (user->app_ctx->test_mode) { 219841e4c73SJed Brown ierr = RegressionTests_NS(user->app_ctx, Q); CHKERRQ(ierr); 22077841947SLeila Ghaffari } 22177841947SLeila Ghaffari 22277841947SLeila Ghaffari PetscFunctionReturn(0); 22377841947SLeila Ghaffari } 22477841947SLeila Ghaffari 22577841947SLeila Ghaffari // Gather initial Q values in case of continuation of simulation 22677841947SLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) { 22777841947SLeila Ghaffari 22877841947SLeila Ghaffari PetscViewer viewer; 22977841947SLeila Ghaffari char file_path[PETSC_MAX_PATH_LEN]; 23077841947SLeila Ghaffari PetscErrorCode ierr; 23177841947SLeila Ghaffari PetscFunctionBegin; 23277841947SLeila Ghaffari 23377841947SLeila Ghaffari // Read input 23477841947SLeila Ghaffari ierr = PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", 23577841947SLeila Ghaffari app_ctx->output_dir); CHKERRQ(ierr); 23677841947SLeila Ghaffari ierr = PetscViewerBinaryOpen(comm, file_path, FILE_MODE_READ, &viewer); 23777841947SLeila Ghaffari CHKERRQ(ierr); 23877841947SLeila Ghaffari 23977841947SLeila Ghaffari // Load Q from existent solution 24077841947SLeila Ghaffari ierr = VecLoad(Q, viewer); CHKERRQ(ierr); 24177841947SLeila Ghaffari 24277841947SLeila Ghaffari // Cleanup 24377841947SLeila Ghaffari ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 24477841947SLeila Ghaffari 24577841947SLeila Ghaffari PetscFunctionReturn(0); 24677841947SLeila Ghaffari } 24777841947SLeila Ghaffari 24877841947SLeila Ghaffari // Record boundary values from initial condition 24977841947SLeila Ghaffari PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc) { 25077841947SLeila Ghaffari 251*5571c6fdSJames Wright Vec Qbc, boundary_mask; 25277841947SLeila Ghaffari PetscErrorCode ierr; 25377841947SLeila Ghaffari PetscFunctionBegin; 25477841947SLeila Ghaffari 25577841947SLeila Ghaffari ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 25677841947SLeila Ghaffari ierr = VecCopy(Q_loc, Qbc); CHKERRQ(ierr); 25777841947SLeila Ghaffari ierr = VecZeroEntries(Q_loc); CHKERRQ(ierr); 25877841947SLeila Ghaffari ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc); CHKERRQ(ierr); 25977841947SLeila Ghaffari ierr = VecAXPY(Qbc, -1., Q_loc); CHKERRQ(ierr); 26077841947SLeila Ghaffari ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 26177841947SLeila Ghaffari ierr = PetscObjectComposeFunction((PetscObject)dm, 26277841947SLeila Ghaffari "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS); 26377841947SLeila Ghaffari CHKERRQ(ierr); 26477841947SLeila Ghaffari 265*5571c6fdSJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 266*5571c6fdSJames Wright PetscCall(DMGetGlobalVector(dm, &Q)); 267*5571c6fdSJames Wright PetscCall(VecZeroEntries(boundary_mask)); 268*5571c6fdSJames Wright PetscCall(VecSet(Q, 1.0)); 269*5571c6fdSJames Wright PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask)); 270*5571c6fdSJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 271*5571c6fdSJames Wright 27277841947SLeila Ghaffari PetscFunctionReturn(0); 27377841947SLeila Ghaffari } 274841e4c73SJed Brown 275841e4c73SJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes 276841e4c73SJed Brown int FreeContextPetsc(void *data) { 277841e4c73SJed Brown if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, 278841e4c73SJed Brown "PetscFree failed"); 279841e4c73SJed Brown return CEED_ERROR_SUCCESS; 280841e4c73SJed Brown } 281