xref: /honee/src/misc.c (revision 9f59f36e7a2e19e53ad001515b81028b86a73a2c)
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"
12*9f59f36eSJames Wright #include "../qfunctions/mass.h"
13a515125bSLeila Ghaffari 
142b916ea7SJeremy L Thompson PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) {
15a515125bSLeila Ghaffari   PetscFunctionBeginUser;
16a515125bSLeila Ghaffari 
17a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
18b7f03f12SJed Brown   // Update time for evaluation
19a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
202b916ea7SJeremy L Thompson   if (user->phys->ics_time_label) CeedOperatorContextSetDouble(ceed_data->op_ics, user->phys->ics_time_label, &time);
21a515125bSLeila Ghaffari 
22a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
23a515125bSLeila Ghaffari   // ICs
24a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
25a515125bSLeila Ghaffari   // -- CEED Restriction
26a515125bSLeila Ghaffari   CeedVector q0_ceed;
27a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &q0_ceed, NULL);
28a515125bSLeila Ghaffari 
29a515125bSLeila Ghaffari   // -- Place PETSc vector in CEED vector
30a515125bSLeila Ghaffari   CeedScalar  *q0;
31a515125bSLeila Ghaffari   PetscMemType q0_mem_type;
322b916ea7SJeremy L Thompson   PetscCall(VecGetArrayAndMemType(Q_loc, (PetscScalar **)&q0, &q0_mem_type));
33a515125bSLeila Ghaffari   CeedVectorSetArray(q0_ceed, MemTypeP2C(q0_mem_type), CEED_USE_POINTER, q0);
34a515125bSLeila Ghaffari 
35a515125bSLeila Ghaffari   // -- Apply CEED Operator
362b916ea7SJeremy L Thompson   CeedOperatorApply(ceed_data->op_ics, ceed_data->x_coord, q0_ceed, CEED_REQUEST_IMMEDIATE);
37a515125bSLeila Ghaffari 
38a515125bSLeila Ghaffari   // -- Restore vectors
39a515125bSLeila Ghaffari   CeedVectorTakeArray(q0_ceed, MemTypeP2C(q0_mem_type), NULL);
402b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, (const PetscScalar **)&q0));
41a515125bSLeila Ghaffari 
42a515125bSLeila Ghaffari   // -- Local-to-Global
432b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(Q));
442b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, Q_loc, ADD_VALUES, Q));
45a515125bSLeila Ghaffari 
46a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
47a515125bSLeila Ghaffari   // Fix multiplicity for output of ICs
48a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
49a515125bSLeila Ghaffari   // -- CEED Restriction
50a515125bSLeila Ghaffari   CeedVector mult_vec;
51a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL);
52a515125bSLeila Ghaffari 
53a515125bSLeila Ghaffari   // -- Place PETSc vector in CEED vector
54a515125bSLeila Ghaffari   CeedScalar  *mult;
55a515125bSLeila Ghaffari   PetscMemType m_mem_type;
56a515125bSLeila Ghaffari   Vec          multiplicity_loc;
572b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &multiplicity_loc));
582b916ea7SJeremy L Thompson   PetscCall(VecGetArrayAndMemType(multiplicity_loc, (PetscScalar **)&mult, &m_mem_type));
59a515125bSLeila Ghaffari   CeedVectorSetArray(mult_vec, MemTypeP2C(m_mem_type), CEED_USE_POINTER, mult);
60a515125bSLeila Ghaffari 
61a515125bSLeila Ghaffari   // -- Get multiplicity
62a515125bSLeila Ghaffari   CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec);
63a515125bSLeila Ghaffari 
64a515125bSLeila Ghaffari   // -- Restore vectors
65a515125bSLeila Ghaffari   CeedVectorTakeArray(mult_vec, MemTypeP2C(m_mem_type), NULL);
662b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(multiplicity_loc, (const PetscScalar **)&mult));
67a515125bSLeila Ghaffari 
68a515125bSLeila Ghaffari   // -- Local-to-Global
69a515125bSLeila Ghaffari   Vec multiplicity;
702b916ea7SJeremy L Thompson   PetscCall(DMGetGlobalVector(dm, &multiplicity));
712b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(multiplicity));
722b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, multiplicity_loc, ADD_VALUES, multiplicity));
73a515125bSLeila Ghaffari 
74a515125bSLeila Ghaffari   // -- Fix multiplicity
752b916ea7SJeremy L Thompson   PetscCall(VecPointwiseDivide(Q, Q, multiplicity));
762b916ea7SJeremy L Thompson   PetscCall(VecPointwiseDivide(Q_loc, Q_loc, multiplicity_loc));
77a515125bSLeila Ghaffari 
78a515125bSLeila Ghaffari   // -- Restore vectors
792b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &multiplicity_loc));
802b916ea7SJeremy L Thompson   PetscCall(DMRestoreGlobalVector(dm, &multiplicity));
81a515125bSLeila Ghaffari 
82a515125bSLeila Ghaffari   // Cleanup
83a515125bSLeila Ghaffari   CeedVectorDestroy(&mult_vec);
84a515125bSLeila Ghaffari   CeedVectorDestroy(&q0_ceed);
85a515125bSLeila Ghaffari 
86a515125bSLeila Ghaffari   PetscFunctionReturn(0);
87a515125bSLeila Ghaffari }
88a515125bSLeila Ghaffari 
892b916ea7SJeremy L Thompson PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
902b916ea7SJeremy L Thompson                                              Vec grad_FVM) {
919d437337SJames Wright   Vec Qbc, boundary_mask;
92a515125bSLeila Ghaffari   PetscFunctionBegin;
93a515125bSLeila Ghaffari 
949d437337SJames Wright   // Mask (zero) Dirichlet entries
959d437337SJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
969d437337SJames Wright   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
979d437337SJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
989d437337SJames Wright 
992b916ea7SJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
1002b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q_loc, 1., Qbc));
1012b916ea7SJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
102a515125bSLeila Ghaffari 
103a515125bSLeila Ghaffari   PetscFunctionReturn(0);
104a515125bSLeila Ghaffari }
105a515125bSLeila Ghaffari 
106a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI
107a515125bSLeila Ghaffari PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q) {
108a515125bSLeila Ghaffari   Vec         Qref;
109a515125bSLeila Ghaffari   PetscViewer viewer;
110a515125bSLeila Ghaffari   PetscReal   error, Qrefnorm;
111a515125bSLeila Ghaffari   PetscFunctionBegin;
112a515125bSLeila Ghaffari 
113a515125bSLeila Ghaffari   // Read reference file
1142b916ea7SJeremy L Thompson   PetscCall(VecDuplicate(Q, &Qref));
1152b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(PetscObjectComm((PetscObject)Q), app_ctx->file_path, FILE_MODE_READ, &viewer));
1162b916ea7SJeremy L Thompson   PetscCall(VecLoad(Qref, viewer));
117a515125bSLeila Ghaffari 
118a515125bSLeila Ghaffari   // Compute error with respect to reference solution
1192b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Qref));
1202b916ea7SJeremy L Thompson   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
1212b916ea7SJeremy L Thompson   PetscCall(VecScale(Q, 1. / Qrefnorm));
1222b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_MAX, &error));
123a515125bSLeila Ghaffari 
124a515125bSLeila Ghaffari   // Check error
125a515125bSLeila Ghaffari   if (error > app_ctx->test_tol) {
1262b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
127a515125bSLeila Ghaffari   }
128a515125bSLeila Ghaffari 
129a515125bSLeila Ghaffari   // Cleanup
1302b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
1312b916ea7SJeremy L Thompson   PetscCall(VecDestroy(&Qref));
132a515125bSLeila Ghaffari 
133a515125bSLeila Ghaffari   PetscFunctionReturn(0);
134a515125bSLeila Ghaffari }
135a515125bSLeila Ghaffari 
136a515125bSLeila Ghaffari // Get error for problems with exact solutions
1372b916ea7SJeremy L Thompson PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
138a515125bSLeila Ghaffari   PetscInt  loc_nodes;
139a515125bSLeila Ghaffari   Vec       Q_exact, Q_exact_loc;
140a515125bSLeila Ghaffari   PetscReal rel_error, norm_error, norm_exact;
141a515125bSLeila Ghaffari   PetscFunctionBegin;
142a515125bSLeila Ghaffari 
143a515125bSLeila Ghaffari   // Get exact solution at final time
1442b916ea7SJeremy L Thompson   PetscCall(DMCreateGlobalVector(dm, &Q_exact));
1452b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1462b916ea7SJeremy L Thompson   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
1472b916ea7SJeremy L Thompson   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
148a515125bSLeila Ghaffari 
149a515125bSLeila Ghaffari   // Get |exact solution - obtained solution|
1502b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1512b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_exact));
1522b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_1, &norm_error));
153a515125bSLeila Ghaffari 
154a515125bSLeila Ghaffari   // Compute relative error
155a515125bSLeila Ghaffari   rel_error = norm_error / norm_exact;
156a515125bSLeila Ghaffari 
157a515125bSLeila Ghaffari   // Output relative error
1582b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
159a515125bSLeila Ghaffari   // Cleanup
1602b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
1612b916ea7SJeremy L Thompson   PetscCall(VecDestroy(&Q_exact));
162a515125bSLeila Ghaffari 
163a515125bSLeila Ghaffari   PetscFunctionReturn(0);
164a515125bSLeila Ghaffari }
165a515125bSLeila Ghaffari 
166a515125bSLeila Ghaffari // Post-processing
1672b916ea7SJeremy L Thompson PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) {
168a515125bSLeila Ghaffari   PetscInt          steps;
169f0784ed3SJed Brown   TSConvergedReason reason;
170a515125bSLeila Ghaffari   PetscFunctionBegin;
171a515125bSLeila Ghaffari 
172a515125bSLeila Ghaffari   // Print relative error
17315a3537eSJed Brown   if (problem->non_zero_time && !user->app_ctx->test_mode) {
1742b916ea7SJeremy L Thompson     PetscCall(GetError_NS(ceed_data, dm, user, Q, final_time));
175a515125bSLeila Ghaffari   }
176a515125bSLeila Ghaffari 
177a515125bSLeila Ghaffari   // Print final time and number of steps
1782b916ea7SJeremy L Thompson   PetscCall(TSGetStepNumber(ts, &steps));
179f0784ed3SJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
18015a3537eSJed Brown   if (!user->app_ctx->test_mode) {
181f0784ed3SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
182f0784ed3SJed Brown                           steps, (double)final_time));
183a515125bSLeila Ghaffari   }
184a515125bSLeila Ghaffari 
185a515125bSLeila Ghaffari   // Output numerical values from command line
1862b916ea7SJeremy L Thompson   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
187a515125bSLeila Ghaffari 
188a515125bSLeila Ghaffari   // Compare reference solution values with current test run for CI
18915a3537eSJed Brown   if (user->app_ctx->test_mode) {
1902b916ea7SJeremy L Thompson     PetscCall(RegressionTests_NS(user->app_ctx, Q));
191a515125bSLeila Ghaffari   }
192a515125bSLeila Ghaffari 
193a515125bSLeila Ghaffari   PetscFunctionReturn(0);
194a515125bSLeila Ghaffari }
195a515125bSLeila Ghaffari 
1969293eaa1SJed Brown const PetscInt FLUIDS_FILE_TOKEN = 0xceedf00;
1979293eaa1SJed Brown 
198a515125bSLeila Ghaffari // Gather initial Q values in case of continuation of simulation
199a515125bSLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
200a515125bSLeila Ghaffari   PetscViewer viewer;
2019293eaa1SJed Brown   PetscInt    token, step_number;
2029293eaa1SJed Brown   PetscReal   time;
2032b916ea7SJeremy L Thompson 
204a515125bSLeila Ghaffari   PetscFunctionBegin;
205a515125bSLeila Ghaffari 
206a515125bSLeila Ghaffari   // Read input
2072b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
208a515125bSLeila Ghaffari 
2099293eaa1SJed Brown   // Attempt
2109293eaa1SJed Brown   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT));
2119293eaa1SJed Brown   if (token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
2129293eaa1SJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &step_number, 1, NULL, PETSC_INT));
2139293eaa1SJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &time, 1, NULL, PETSC_REAL));
2149293eaa1SJed Brown     app_ctx->cont_steps = step_number;
2159293eaa1SJed Brown     app_ctx->cont_time  = time;
2169293eaa1SJed Brown   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
2179293eaa1SJed Brown     PetscInt length, N;
2189293eaa1SJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &length, 1, NULL, PETSC_INT));
2199293eaa1SJed Brown     PetscCall(VecGetSize(Q, &N));
2209293eaa1SJed 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);
2219293eaa1SJed Brown     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
2229293eaa1SJed Brown   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
2239293eaa1SJed Brown 
224a515125bSLeila Ghaffari   // Load Q from existent solution
2252b916ea7SJeremy L Thompson   PetscCall(VecLoad(Q, viewer));
226a515125bSLeila Ghaffari 
227a515125bSLeila Ghaffari   // Cleanup
2282b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
229a515125bSLeila Ghaffari 
230a515125bSLeila Ghaffari   PetscFunctionReturn(0);
231a515125bSLeila Ghaffari }
232a515125bSLeila Ghaffari 
233a515125bSLeila Ghaffari // Record boundary values from initial condition
234a515125bSLeila Ghaffari PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc) {
2359d437337SJames Wright   Vec Qbc, boundary_mask;
236a515125bSLeila Ghaffari   PetscFunctionBegin;
237a515125bSLeila Ghaffari 
2382b916ea7SJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
2392b916ea7SJeremy L Thompson   PetscCall(VecCopy(Q_loc, Qbc));
2402b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(Q_loc));
2412b916ea7SJeremy L Thompson   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
2422b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Qbc, -1., Q_loc));
2432b916ea7SJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
2442b916ea7SJeremy L Thompson   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS));
245a515125bSLeila Ghaffari 
2469d437337SJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
2479d437337SJames Wright   PetscCall(DMGetGlobalVector(dm, &Q));
2489d437337SJames Wright   PetscCall(VecZeroEntries(boundary_mask));
2499d437337SJames Wright   PetscCall(VecSet(Q, 1.0));
2509d437337SJames Wright   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask));
2519d437337SJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
2529d437337SJames Wright 
253a515125bSLeila Ghaffari   PetscFunctionReturn(0);
254a515125bSLeila Ghaffari }
25515a3537eSJed Brown 
25615a3537eSJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
25715a3537eSJed Brown int FreeContextPetsc(void *data) {
2582b916ea7SJeremy L Thompson   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
25915a3537eSJed Brown   return CEED_ERROR_SUCCESS;
26015a3537eSJed Brown }
261*9f59f36eSJames Wright 
262*9f59f36eSJames Wright // Return mass qfunction specification for number of components N
263*9f59f36eSJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
264*9f59f36eSJames Wright   CeedQFunctionUser qfunction_ptr;
265*9f59f36eSJames Wright   const char       *qfunction_loc;
266*9f59f36eSJames Wright   PetscFunctionBeginUser;
267*9f59f36eSJames Wright 
268*9f59f36eSJames Wright   switch (N) {
269*9f59f36eSJames Wright     case 1:
270*9f59f36eSJames Wright       qfunction_ptr = Mass_1;
271*9f59f36eSJames Wright       qfunction_loc = Mass_1_loc;
272*9f59f36eSJames Wright       break;
273*9f59f36eSJames Wright     case 5:
274*9f59f36eSJames Wright       qfunction_ptr = Mass_5;
275*9f59f36eSJames Wright       qfunction_loc = Mass_5_loc;
276*9f59f36eSJames Wright       break;
277*9f59f36eSJames Wright     case 9:
278*9f59f36eSJames Wright       qfunction_ptr = Mass_9;
279*9f59f36eSJames Wright       qfunction_loc = Mass_9_loc;
280*9f59f36eSJames Wright       break;
281*9f59f36eSJames Wright     case 22:
282*9f59f36eSJames Wright       qfunction_ptr = Mass_22;
283*9f59f36eSJames Wright       qfunction_loc = Mass_22_loc;
284*9f59f36eSJames Wright       break;
285*9f59f36eSJames Wright     default:
286*9f59f36eSJames Wright       SETERRQ(PETSC_COMM_WORLD, -1, "Could not find mass qfunction of size %d", N);
287*9f59f36eSJames Wright   }
288*9f59f36eSJames Wright   CeedQFunctionCreateInterior(ceed, 1, qfunction_ptr, qfunction_loc, qf);
289*9f59f36eSJames Wright 
290*9f59f36eSJames Wright   CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP);
291*9f59f36eSJames Wright   CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE);
292*9f59f36eSJames Wright   CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP);
293*9f59f36eSJames Wright   PetscFunctionReturn(0);
294*9f59f36eSJames Wright }
295