xref: /honee/src/misc.c (revision f0784ed3459207699c8e6c42e22de54252c413f5)
1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Miscellaneous utility functions
10a515125bSLeila Ghaffari 
11a515125bSLeila Ghaffari #include "../navierstokes.h"
12a515125bSLeila Ghaffari 
132b916ea7SJeremy L Thompson PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) {
14a515125bSLeila Ghaffari   PetscFunctionBeginUser;
15a515125bSLeila Ghaffari 
16a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
17b7f03f12SJed Brown   // Update time for evaluation
18a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
192b916ea7SJeremy L Thompson   if (user->phys->ics_time_label) CeedOperatorContextSetDouble(ceed_data->op_ics, user->phys->ics_time_label, &time);
20a515125bSLeila Ghaffari 
21a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
22a515125bSLeila Ghaffari   // ICs
23a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
24a515125bSLeila Ghaffari   // -- CEED Restriction
25a515125bSLeila Ghaffari   CeedVector q0_ceed;
26a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &q0_ceed, NULL);
27a515125bSLeila Ghaffari 
28a515125bSLeila Ghaffari   // -- Place PETSc vector in CEED vector
29a515125bSLeila Ghaffari   CeedScalar  *q0;
30a515125bSLeila Ghaffari   PetscMemType q0_mem_type;
312b916ea7SJeremy L Thompson   PetscCall(VecGetArrayAndMemType(Q_loc, (PetscScalar **)&q0, &q0_mem_type));
32a515125bSLeila Ghaffari   CeedVectorSetArray(q0_ceed, MemTypeP2C(q0_mem_type), CEED_USE_POINTER, q0);
33a515125bSLeila Ghaffari 
34a515125bSLeila Ghaffari   // -- Apply CEED Operator
352b916ea7SJeremy L Thompson   CeedOperatorApply(ceed_data->op_ics, ceed_data->x_coord, q0_ceed, CEED_REQUEST_IMMEDIATE);
36a515125bSLeila Ghaffari 
37a515125bSLeila Ghaffari   // -- Restore vectors
38a515125bSLeila Ghaffari   CeedVectorTakeArray(q0_ceed, MemTypeP2C(q0_mem_type), NULL);
392b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, (const PetscScalar **)&q0));
40a515125bSLeila Ghaffari 
41a515125bSLeila Ghaffari   // -- Local-to-Global
422b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(Q));
432b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, Q_loc, ADD_VALUES, Q));
44a515125bSLeila Ghaffari 
45a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
46a515125bSLeila Ghaffari   // Fix multiplicity for output of ICs
47a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
48a515125bSLeila Ghaffari   // -- CEED Restriction
49a515125bSLeila Ghaffari   CeedVector mult_vec;
50a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL);
51a515125bSLeila Ghaffari 
52a515125bSLeila Ghaffari   // -- Place PETSc vector in CEED vector
53a515125bSLeila Ghaffari   CeedScalar  *mult;
54a515125bSLeila Ghaffari   PetscMemType m_mem_type;
55a515125bSLeila Ghaffari   Vec          multiplicity_loc;
562b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &multiplicity_loc));
572b916ea7SJeremy L Thompson   PetscCall(VecGetArrayAndMemType(multiplicity_loc, (PetscScalar **)&mult, &m_mem_type));
58a515125bSLeila Ghaffari   CeedVectorSetArray(mult_vec, MemTypeP2C(m_mem_type), CEED_USE_POINTER, mult);
59a515125bSLeila Ghaffari 
60a515125bSLeila Ghaffari   // -- Get multiplicity
61a515125bSLeila Ghaffari   CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec);
62a515125bSLeila Ghaffari 
63a515125bSLeila Ghaffari   // -- Restore vectors
64a515125bSLeila Ghaffari   CeedVectorTakeArray(mult_vec, MemTypeP2C(m_mem_type), NULL);
652b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(multiplicity_loc, (const PetscScalar **)&mult));
66a515125bSLeila Ghaffari 
67a515125bSLeila Ghaffari   // -- Local-to-Global
68a515125bSLeila Ghaffari   Vec multiplicity;
692b916ea7SJeremy L Thompson   PetscCall(DMGetGlobalVector(dm, &multiplicity));
702b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(multiplicity));
712b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, multiplicity_loc, ADD_VALUES, multiplicity));
72a515125bSLeila Ghaffari 
73a515125bSLeila Ghaffari   // -- Fix multiplicity
742b916ea7SJeremy L Thompson   PetscCall(VecPointwiseDivide(Q, Q, multiplicity));
752b916ea7SJeremy L Thompson   PetscCall(VecPointwiseDivide(Q_loc, Q_loc, multiplicity_loc));
76a515125bSLeila Ghaffari 
77a515125bSLeila Ghaffari   // -- Restore vectors
782b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &multiplicity_loc));
792b916ea7SJeremy L Thompson   PetscCall(DMRestoreGlobalVector(dm, &multiplicity));
80a515125bSLeila Ghaffari 
81a515125bSLeila Ghaffari   // Cleanup
82a515125bSLeila Ghaffari   CeedVectorDestroy(&mult_vec);
83a515125bSLeila Ghaffari   CeedVectorDestroy(&q0_ceed);
84a515125bSLeila Ghaffari 
85a515125bSLeila Ghaffari   PetscFunctionReturn(0);
86a515125bSLeila Ghaffari }
87a515125bSLeila Ghaffari 
882b916ea7SJeremy L Thompson PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
892b916ea7SJeremy L Thompson                                              Vec grad_FVM) {
909d437337SJames Wright   Vec Qbc, boundary_mask;
91a515125bSLeila Ghaffari   PetscFunctionBegin;
92a515125bSLeila Ghaffari 
939d437337SJames Wright   // Mask (zero) Dirichlet entries
949d437337SJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
959d437337SJames Wright   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
969d437337SJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
979d437337SJames Wright 
982b916ea7SJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
992b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q_loc, 1., Qbc));
1002b916ea7SJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
101a515125bSLeila Ghaffari 
102a515125bSLeila Ghaffari   PetscFunctionReturn(0);
103a515125bSLeila Ghaffari }
104a515125bSLeila Ghaffari 
105a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI
106a515125bSLeila Ghaffari PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q) {
107a515125bSLeila Ghaffari   Vec         Qref;
108a515125bSLeila Ghaffari   PetscViewer viewer;
109a515125bSLeila Ghaffari   PetscReal   error, Qrefnorm;
110a515125bSLeila Ghaffari   PetscFunctionBegin;
111a515125bSLeila Ghaffari 
112a515125bSLeila Ghaffari   // Read reference file
1132b916ea7SJeremy L Thompson   PetscCall(VecDuplicate(Q, &Qref));
1142b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(PetscObjectComm((PetscObject)Q), app_ctx->file_path, FILE_MODE_READ, &viewer));
1152b916ea7SJeremy L Thompson   PetscCall(VecLoad(Qref, viewer));
116a515125bSLeila Ghaffari 
117a515125bSLeila Ghaffari   // Compute error with respect to reference solution
1182b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Qref));
1192b916ea7SJeremy L Thompson   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
1202b916ea7SJeremy L Thompson   PetscCall(VecScale(Q, 1. / Qrefnorm));
1212b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_MAX, &error));
122a515125bSLeila Ghaffari 
123a515125bSLeila Ghaffari   // Check error
124a515125bSLeila Ghaffari   if (error > app_ctx->test_tol) {
1252b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
126a515125bSLeila Ghaffari   }
127a515125bSLeila Ghaffari 
128a515125bSLeila Ghaffari   // Cleanup
1292b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
1302b916ea7SJeremy L Thompson   PetscCall(VecDestroy(&Qref));
131a515125bSLeila Ghaffari 
132a515125bSLeila Ghaffari   PetscFunctionReturn(0);
133a515125bSLeila Ghaffari }
134a515125bSLeila Ghaffari 
135a515125bSLeila Ghaffari // Get error for problems with exact solutions
1362b916ea7SJeremy L Thompson PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
137a515125bSLeila Ghaffari   PetscInt  loc_nodes;
138a515125bSLeila Ghaffari   Vec       Q_exact, Q_exact_loc;
139a515125bSLeila Ghaffari   PetscReal rel_error, norm_error, norm_exact;
140a515125bSLeila Ghaffari   PetscFunctionBegin;
141a515125bSLeila Ghaffari 
142a515125bSLeila Ghaffari   // Get exact solution at final time
1432b916ea7SJeremy L Thompson   PetscCall(DMCreateGlobalVector(dm, &Q_exact));
1442b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1452b916ea7SJeremy L Thompson   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
1462b916ea7SJeremy L Thompson   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
147a515125bSLeila Ghaffari 
148a515125bSLeila Ghaffari   // Get |exact solution - obtained solution|
1492b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1502b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_exact));
1512b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_1, &norm_error));
152a515125bSLeila Ghaffari 
153a515125bSLeila Ghaffari   // Compute relative error
154a515125bSLeila Ghaffari   rel_error = norm_error / norm_exact;
155a515125bSLeila Ghaffari 
156a515125bSLeila Ghaffari   // Output relative error
1572b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
158a515125bSLeila Ghaffari   // Cleanup
1592b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
1602b916ea7SJeremy L Thompson   PetscCall(VecDestroy(&Q_exact));
161a515125bSLeila Ghaffari 
162a515125bSLeila Ghaffari   PetscFunctionReturn(0);
163a515125bSLeila Ghaffari }
164a515125bSLeila Ghaffari 
165a515125bSLeila Ghaffari // Post-processing
1662b916ea7SJeremy L Thompson PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) {
167a515125bSLeila Ghaffari   PetscInt          steps;
168*f0784ed3SJed Brown   TSConvergedReason reason;
169a515125bSLeila Ghaffari   PetscFunctionBegin;
170a515125bSLeila Ghaffari 
171a515125bSLeila Ghaffari   // Print relative error
17215a3537eSJed Brown   if (problem->non_zero_time && !user->app_ctx->test_mode) {
1732b916ea7SJeremy L Thompson     PetscCall(GetError_NS(ceed_data, dm, user, Q, final_time));
174a515125bSLeila Ghaffari   }
175a515125bSLeila Ghaffari 
176a515125bSLeila Ghaffari   // Print final time and number of steps
1772b916ea7SJeremy L Thompson   PetscCall(TSGetStepNumber(ts, &steps));
178*f0784ed3SJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
17915a3537eSJed Brown   if (!user->app_ctx->test_mode) {
180*f0784ed3SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
181*f0784ed3SJed Brown                           steps, (double)final_time));
182a515125bSLeila Ghaffari   }
183a515125bSLeila Ghaffari 
184a515125bSLeila Ghaffari   // Output numerical values from command line
1852b916ea7SJeremy L Thompson   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
186a515125bSLeila Ghaffari 
187a515125bSLeila Ghaffari   // Compare reference solution values with current test run for CI
18815a3537eSJed Brown   if (user->app_ctx->test_mode) {
1892b916ea7SJeremy L Thompson     PetscCall(RegressionTests_NS(user->app_ctx, Q));
190a515125bSLeila Ghaffari   }
191a515125bSLeila Ghaffari 
192a515125bSLeila Ghaffari   PetscFunctionReturn(0);
193a515125bSLeila Ghaffari }
194a515125bSLeila Ghaffari 
1959293eaa1SJed Brown const PetscInt FLUIDS_FILE_TOKEN = 0xceedf00;
1969293eaa1SJed Brown 
197a515125bSLeila Ghaffari // Gather initial Q values in case of continuation of simulation
198a515125bSLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
199a515125bSLeila Ghaffari   PetscViewer viewer;
2009293eaa1SJed Brown   PetscInt    token, step_number;
2019293eaa1SJed Brown   PetscReal   time;
2022b916ea7SJeremy L Thompson 
203a515125bSLeila Ghaffari   PetscFunctionBegin;
204a515125bSLeila Ghaffari 
205a515125bSLeila Ghaffari   // Read input
2062b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
207a515125bSLeila Ghaffari 
2089293eaa1SJed Brown   // Attempt
2099293eaa1SJed Brown   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT));
2109293eaa1SJed Brown   if (token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
2119293eaa1SJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &step_number, 1, NULL, PETSC_INT));
2129293eaa1SJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &time, 1, NULL, PETSC_REAL));
2139293eaa1SJed Brown     app_ctx->cont_steps = step_number;
2149293eaa1SJed Brown     app_ctx->cont_time  = time;
2159293eaa1SJed Brown   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
2169293eaa1SJed Brown     PetscInt length, N;
2179293eaa1SJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &length, 1, NULL, PETSC_INT));
2189293eaa1SJed Brown     PetscCall(VecGetSize(Q, &N));
2199293eaa1SJed Brown     PetscCheck(length == N, comm, PETSC_ERR_ARG_INCOMP, "File Vec has length %" PetscInt_FMT " but DM has global Vec size %" PetscInt_FMT, length, N);
2209293eaa1SJed Brown     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
2219293eaa1SJed Brown   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
2229293eaa1SJed Brown 
223a515125bSLeila Ghaffari   // Load Q from existent solution
2242b916ea7SJeremy L Thompson   PetscCall(VecLoad(Q, viewer));
225a515125bSLeila Ghaffari 
226a515125bSLeila Ghaffari   // Cleanup
2272b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
228a515125bSLeila Ghaffari 
229a515125bSLeila Ghaffari   PetscFunctionReturn(0);
230a515125bSLeila Ghaffari }
231a515125bSLeila Ghaffari 
232a515125bSLeila Ghaffari // Record boundary values from initial condition
233a515125bSLeila Ghaffari PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc) {
2349d437337SJames Wright   Vec Qbc, boundary_mask;
235a515125bSLeila Ghaffari   PetscFunctionBegin;
236a515125bSLeila Ghaffari 
2372b916ea7SJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
2382b916ea7SJeremy L Thompson   PetscCall(VecCopy(Q_loc, Qbc));
2392b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(Q_loc));
2402b916ea7SJeremy L Thompson   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
2412b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Qbc, -1., Q_loc));
2422b916ea7SJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
2432b916ea7SJeremy L Thompson   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS));
244a515125bSLeila Ghaffari 
2459d437337SJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
2469d437337SJames Wright   PetscCall(DMGetGlobalVector(dm, &Q));
2479d437337SJames Wright   PetscCall(VecZeroEntries(boundary_mask));
2489d437337SJames Wright   PetscCall(VecSet(Q, 1.0));
2499d437337SJames Wright   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask));
2509d437337SJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
2519d437337SJames Wright 
252a515125bSLeila Ghaffari   PetscFunctionReturn(0);
253a515125bSLeila Ghaffari }
25415a3537eSJed Brown 
25515a3537eSJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
25615a3537eSJed Brown int FreeContextPetsc(void *data) {
2572b916ea7SJeremy L Thompson   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
25815a3537eSJed Brown   return CEED_ERROR_SUCCESS;
25915a3537eSJed Brown }
260