xref: /honee/src/misc.c (revision b7f03f12333883f854e9fffaf41e1cf2f63d5bb5)
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 
1315a3537eSJed Brown PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user,
1415a3537eSJed Brown                                    Vec Q_loc, Vec Q,
15a515125bSLeila Ghaffari                                    CeedScalar time) {
16a515125bSLeila Ghaffari   PetscErrorCode ierr;
17a515125bSLeila Ghaffari   PetscFunctionBeginUser;
18a515125bSLeila Ghaffari 
19a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
20*b7f03f12SJed Brown   // Update time for evaluation
21a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
2215a3537eSJed Brown   if (user->phys->ics_time_label)
2315a3537eSJed Brown     CeedOperatorContextSetDouble(ceed_data->op_ics, user->phys->ics_time_label,
2415a3537eSJed Brown                                  &time);
25a515125bSLeila Ghaffari 
26a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
27a515125bSLeila Ghaffari   // ICs
28a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
29a515125bSLeila Ghaffari   // -- CEED Restriction
30a515125bSLeila Ghaffari   CeedVector q0_ceed;
31a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &q0_ceed, NULL);
32a515125bSLeila Ghaffari 
33a515125bSLeila Ghaffari   // -- Place PETSc vector in CEED vector
34a515125bSLeila Ghaffari   CeedScalar *q0;
35a515125bSLeila Ghaffari   PetscMemType q0_mem_type;
36a515125bSLeila Ghaffari   ierr = VecGetArrayAndMemType(Q_loc, (PetscScalar **)&q0, &q0_mem_type);
37a515125bSLeila Ghaffari   CHKERRQ(ierr);
38a515125bSLeila Ghaffari   CeedVectorSetArray(q0_ceed, MemTypeP2C(q0_mem_type), CEED_USE_POINTER, q0);
39a515125bSLeila Ghaffari 
40a515125bSLeila Ghaffari   // -- Apply CEED Operator
41a515125bSLeila Ghaffari   CeedOperatorApply(ceed_data->op_ics, ceed_data->x_coord, q0_ceed,
42a515125bSLeila Ghaffari                     CEED_REQUEST_IMMEDIATE);
43a515125bSLeila Ghaffari 
44a515125bSLeila Ghaffari   // -- Restore vectors
45a515125bSLeila Ghaffari   CeedVectorTakeArray(q0_ceed, MemTypeP2C(q0_mem_type), NULL);
46a515125bSLeila Ghaffari   ierr = VecRestoreArrayReadAndMemType(Q_loc, (const PetscScalar **)&q0);
47a515125bSLeila Ghaffari   CHKERRQ(ierr);
48a515125bSLeila Ghaffari 
49a515125bSLeila Ghaffari   // -- Local-to-Global
50a515125bSLeila Ghaffari   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
51a515125bSLeila Ghaffari   ierr = DMLocalToGlobal(dm, Q_loc, ADD_VALUES, Q); CHKERRQ(ierr);
52a515125bSLeila Ghaffari 
53a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
54a515125bSLeila Ghaffari   // Fix multiplicity for output of ICs
55a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
56a515125bSLeila Ghaffari   // -- CEED Restriction
57a515125bSLeila Ghaffari   CeedVector mult_vec;
58a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL);
59a515125bSLeila Ghaffari 
60a515125bSLeila Ghaffari   // -- Place PETSc vector in CEED vector
61a515125bSLeila Ghaffari   CeedScalar *mult;
62a515125bSLeila Ghaffari   PetscMemType m_mem_type;
63a515125bSLeila Ghaffari   Vec multiplicity_loc;
64a515125bSLeila Ghaffari   ierr = DMGetLocalVector(dm, &multiplicity_loc); CHKERRQ(ierr);
65a515125bSLeila Ghaffari   ierr = VecGetArrayAndMemType(multiplicity_loc, (PetscScalar **)&mult,
66a515125bSLeila Ghaffari                                &m_mem_type);
67a515125bSLeila Ghaffari   CHKERRQ(ierr);
68a515125bSLeila Ghaffari   CeedVectorSetArray(mult_vec, MemTypeP2C(m_mem_type), CEED_USE_POINTER, mult);
69a515125bSLeila Ghaffari   CHKERRQ(ierr);
70a515125bSLeila Ghaffari 
71a515125bSLeila Ghaffari   // -- Get multiplicity
72a515125bSLeila Ghaffari   CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec);
73a515125bSLeila Ghaffari 
74a515125bSLeila Ghaffari   // -- Restore vectors
75a515125bSLeila Ghaffari   CeedVectorTakeArray(mult_vec, MemTypeP2C(m_mem_type), NULL);
76a515125bSLeila Ghaffari   ierr = VecRestoreArrayReadAndMemType(multiplicity_loc,
77a515125bSLeila Ghaffari                                        (const PetscScalar **)&mult); CHKERRQ(ierr);
78a515125bSLeila Ghaffari 
79a515125bSLeila Ghaffari   // -- Local-to-Global
80a515125bSLeila Ghaffari   Vec multiplicity;
81a515125bSLeila Ghaffari   ierr = DMGetGlobalVector(dm, &multiplicity); CHKERRQ(ierr);
82a515125bSLeila Ghaffari   ierr = VecZeroEntries(multiplicity); CHKERRQ(ierr);
83a515125bSLeila Ghaffari   ierr = DMLocalToGlobal(dm, multiplicity_loc, ADD_VALUES, multiplicity);
84a515125bSLeila Ghaffari   CHKERRQ(ierr);
85a515125bSLeila Ghaffari 
86a515125bSLeila Ghaffari   // -- Fix multiplicity
87a515125bSLeila Ghaffari   ierr = VecPointwiseDivide(Q, Q, multiplicity); CHKERRQ(ierr);
88a515125bSLeila Ghaffari   ierr = VecPointwiseDivide(Q_loc, Q_loc, multiplicity_loc); CHKERRQ(ierr);
89a515125bSLeila Ghaffari 
90a515125bSLeila Ghaffari   // -- Restore vectors
91a515125bSLeila Ghaffari   ierr = DMRestoreLocalVector(dm, &multiplicity_loc); CHKERRQ(ierr);
92a515125bSLeila Ghaffari   ierr = DMRestoreGlobalVector(dm, &multiplicity); CHKERRQ(ierr);
93a515125bSLeila Ghaffari 
94a515125bSLeila Ghaffari   // Cleanup
95a515125bSLeila Ghaffari   CeedVectorDestroy(&mult_vec);
96a515125bSLeila Ghaffari   CeedVectorDestroy(&q0_ceed);
97a515125bSLeila Ghaffari 
98a515125bSLeila Ghaffari   PetscFunctionReturn(0);
99a515125bSLeila Ghaffari }
100a515125bSLeila Ghaffari 
101a515125bSLeila Ghaffari PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm,
102a515125bSLeila Ghaffari     PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM,
103a515125bSLeila Ghaffari     Vec cell_geom_FVM, Vec grad_FVM) {
104a515125bSLeila Ghaffari 
105a515125bSLeila Ghaffari   Vec            Qbc;
106a515125bSLeila Ghaffari   PetscErrorCode ierr;
107a515125bSLeila Ghaffari   PetscFunctionBegin;
108a515125bSLeila Ghaffari 
109a515125bSLeila Ghaffari   ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
110a515125bSLeila Ghaffari   ierr = VecAXPY(Q_loc, 1., Qbc); CHKERRQ(ierr);
111a515125bSLeila Ghaffari   ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
112a515125bSLeila Ghaffari 
113a515125bSLeila Ghaffari   PetscFunctionReturn(0);
114a515125bSLeila Ghaffari }
115a515125bSLeila Ghaffari 
116a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI
117a515125bSLeila Ghaffari PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q) {
118a515125bSLeila Ghaffari 
119a515125bSLeila Ghaffari   Vec            Qref;
120a515125bSLeila Ghaffari   PetscViewer    viewer;
121a515125bSLeila Ghaffari   PetscReal      error, Qrefnorm;
122a515125bSLeila Ghaffari   PetscErrorCode ierr;
123a515125bSLeila Ghaffari   PetscFunctionBegin;
124a515125bSLeila Ghaffari 
125a515125bSLeila Ghaffari   // Read reference file
126a515125bSLeila Ghaffari   ierr = VecDuplicate(Q, &Qref); CHKERRQ(ierr);
127a515125bSLeila Ghaffari   ierr = PetscViewerBinaryOpen(PetscObjectComm((PetscObject)Q),
128a515125bSLeila Ghaffari                                app_ctx->file_path, FILE_MODE_READ,
129a515125bSLeila Ghaffari                                &viewer); CHKERRQ(ierr);
130a515125bSLeila Ghaffari   ierr = VecLoad(Qref, viewer); CHKERRQ(ierr);
131a515125bSLeila Ghaffari 
132a515125bSLeila Ghaffari   // Compute error with respect to reference solution
133a515125bSLeila Ghaffari   ierr = VecAXPY(Q, -1.0, Qref); CHKERRQ(ierr);
134a515125bSLeila Ghaffari   ierr = VecNorm(Qref, NORM_MAX, &Qrefnorm); CHKERRQ(ierr);
135a515125bSLeila Ghaffari   ierr = VecScale(Q, 1./Qrefnorm); CHKERRQ(ierr);
136a515125bSLeila Ghaffari   ierr = VecNorm(Q, NORM_MAX, &error); CHKERRQ(ierr);
137a515125bSLeila Ghaffari 
138a515125bSLeila Ghaffari   // Check error
139a515125bSLeila Ghaffari   if (error > app_ctx->test_tol) {
140a515125bSLeila Ghaffari     ierr = PetscPrintf(PETSC_COMM_WORLD,
141a515125bSLeila Ghaffari                        "Test failed with error norm %g\n",
142a515125bSLeila Ghaffari                        (double)error); CHKERRQ(ierr);
143a515125bSLeila Ghaffari   }
144a515125bSLeila Ghaffari 
145a515125bSLeila Ghaffari   // Cleanup
146a515125bSLeila Ghaffari   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
147a515125bSLeila Ghaffari   ierr = VecDestroy(&Qref); CHKERRQ(ierr);
148a515125bSLeila Ghaffari 
149a515125bSLeila Ghaffari   PetscFunctionReturn(0);
150a515125bSLeila Ghaffari }
151a515125bSLeila Ghaffari 
152a515125bSLeila Ghaffari // Get error for problems with exact solutions
15315a3537eSJed Brown PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q,
154a515125bSLeila Ghaffari                            PetscScalar final_time) {
155a515125bSLeila Ghaffari   PetscInt       loc_nodes;
156a515125bSLeila Ghaffari   Vec            Q_exact, Q_exact_loc;
157a515125bSLeila Ghaffari   PetscReal      rel_error, norm_error, norm_exact;
158a515125bSLeila Ghaffari   PetscErrorCode ierr;
159a515125bSLeila Ghaffari   PetscFunctionBegin;
160a515125bSLeila Ghaffari 
161a515125bSLeila Ghaffari   // Get exact solution at final time
162a515125bSLeila Ghaffari   ierr = DMCreateGlobalVector(dm, &Q_exact); CHKERRQ(ierr);
163a515125bSLeila Ghaffari   ierr = DMGetLocalVector(dm, &Q_exact_loc); CHKERRQ(ierr);
164a515125bSLeila Ghaffari   ierr = VecGetSize(Q_exact_loc, &loc_nodes); CHKERRQ(ierr);
16515a3537eSJed Brown   ierr = ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact,
16615a3537eSJed Brown                              final_time);
167a515125bSLeila Ghaffari   CHKERRQ(ierr);
168a515125bSLeila Ghaffari 
169a515125bSLeila Ghaffari   // Get |exact solution - obtained solution|
170a515125bSLeila Ghaffari   ierr = VecNorm(Q_exact, NORM_1, &norm_exact); CHKERRQ(ierr);
171a515125bSLeila Ghaffari   ierr = VecAXPY(Q, -1.0, Q_exact);  CHKERRQ(ierr);
172a515125bSLeila Ghaffari   ierr = VecNorm(Q, NORM_1, &norm_error); CHKERRQ(ierr);
173a515125bSLeila Ghaffari 
174a515125bSLeila Ghaffari   // Compute relative error
175a515125bSLeila Ghaffari   rel_error = norm_error / norm_exact;
176a515125bSLeila Ghaffari 
177a515125bSLeila Ghaffari   // Output relative error
178a515125bSLeila Ghaffari   ierr = PetscPrintf(PETSC_COMM_WORLD,
179a515125bSLeila Ghaffari                      "Relative Error: %g\n",
180a515125bSLeila Ghaffari                      (double)rel_error); CHKERRQ(ierr);
181a515125bSLeila Ghaffari   // Cleanup
182a515125bSLeila Ghaffari   ierr = DMRestoreLocalVector(dm, &Q_exact_loc); CHKERRQ(ierr);
183a515125bSLeila Ghaffari   ierr = VecDestroy(&Q_exact); CHKERRQ(ierr);
184a515125bSLeila Ghaffari 
185a515125bSLeila Ghaffari   PetscFunctionReturn(0);
186a515125bSLeila Ghaffari }
187a515125bSLeila Ghaffari 
188a515125bSLeila Ghaffari // Post-processing
189a515125bSLeila Ghaffari PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm,
19015a3537eSJed Brown                               ProblemData *problem, User user,
191a515125bSLeila Ghaffari                               Vec Q, PetscScalar final_time) {
192a515125bSLeila Ghaffari   PetscInt       steps;
193a515125bSLeila Ghaffari   PetscErrorCode ierr;
194a515125bSLeila Ghaffari   PetscFunctionBegin;
195a515125bSLeila Ghaffari 
196a515125bSLeila Ghaffari   // Print relative error
19715a3537eSJed Brown   if (problem->non_zero_time && !user->app_ctx->test_mode) {
19815a3537eSJed Brown     ierr = GetError_NS(ceed_data, dm, user, Q, final_time); CHKERRQ(ierr);
199a515125bSLeila Ghaffari   }
200a515125bSLeila Ghaffari 
201a515125bSLeila Ghaffari   // Print final time and number of steps
202a515125bSLeila Ghaffari   ierr = TSGetStepNumber(ts, &steps); CHKERRQ(ierr);
20315a3537eSJed Brown   if (!user->app_ctx->test_mode) {
204a515125bSLeila Ghaffari     ierr = PetscPrintf(PETSC_COMM_WORLD,
205d940ca4eSJed Brown                        "Time integrator took %" PetscInt_FMT " time steps to reach final time %g\n",
206a515125bSLeila Ghaffari                        steps, (double)final_time); CHKERRQ(ierr);
207a515125bSLeila Ghaffari   }
208a515125bSLeila Ghaffari 
209a515125bSLeila Ghaffari   // Output numerical values from command line
210a515125bSLeila Ghaffari   ierr = VecViewFromOptions(Q, NULL, "-vec_view"); CHKERRQ(ierr);
211a515125bSLeila Ghaffari 
212a515125bSLeila Ghaffari   // Compare reference solution values with current test run for CI
21315a3537eSJed Brown   if (user->app_ctx->test_mode) {
21415a3537eSJed Brown     ierr = RegressionTests_NS(user->app_ctx, Q); CHKERRQ(ierr);
215a515125bSLeila Ghaffari   }
216a515125bSLeila Ghaffari 
217a515125bSLeila Ghaffari   PetscFunctionReturn(0);
218a515125bSLeila Ghaffari }
219a515125bSLeila Ghaffari 
220a515125bSLeila Ghaffari // Gather initial Q values in case of continuation of simulation
221a515125bSLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
222a515125bSLeila Ghaffari 
223a515125bSLeila Ghaffari   PetscViewer    viewer;
224a515125bSLeila Ghaffari   char           file_path[PETSC_MAX_PATH_LEN];
225a515125bSLeila Ghaffari   PetscErrorCode ierr;
226a515125bSLeila Ghaffari   PetscFunctionBegin;
227a515125bSLeila Ghaffari 
228a515125bSLeila Ghaffari   // Read input
229a515125bSLeila Ghaffari   ierr = PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin",
230a515125bSLeila Ghaffari                        app_ctx->output_dir); CHKERRQ(ierr);
231a515125bSLeila Ghaffari   ierr = PetscViewerBinaryOpen(comm, file_path, FILE_MODE_READ, &viewer);
232a515125bSLeila Ghaffari   CHKERRQ(ierr);
233a515125bSLeila Ghaffari 
234a515125bSLeila Ghaffari   // Load Q from existent solution
235a515125bSLeila Ghaffari   ierr = VecLoad(Q, viewer); CHKERRQ(ierr);
236a515125bSLeila Ghaffari 
237a515125bSLeila Ghaffari   // Cleanup
238a515125bSLeila Ghaffari   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
239a515125bSLeila Ghaffari 
240a515125bSLeila Ghaffari   PetscFunctionReturn(0);
241a515125bSLeila Ghaffari }
242a515125bSLeila Ghaffari 
243a515125bSLeila Ghaffari // Record boundary values from initial condition
244a515125bSLeila Ghaffari PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc) {
245a515125bSLeila Ghaffari 
246a515125bSLeila Ghaffari   Vec            Qbc;
247a515125bSLeila Ghaffari   PetscErrorCode ierr;
248a515125bSLeila Ghaffari   PetscFunctionBegin;
249a515125bSLeila Ghaffari 
250a515125bSLeila Ghaffari   ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
251a515125bSLeila Ghaffari   ierr = VecCopy(Q_loc, Qbc); CHKERRQ(ierr);
252a515125bSLeila Ghaffari   ierr = VecZeroEntries(Q_loc); CHKERRQ(ierr);
253a515125bSLeila Ghaffari   ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc); CHKERRQ(ierr);
254a515125bSLeila Ghaffari   ierr = VecAXPY(Qbc, -1., Q_loc); CHKERRQ(ierr);
255a515125bSLeila Ghaffari   ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
256a515125bSLeila Ghaffari   ierr = PetscObjectComposeFunction((PetscObject)dm,
257a515125bSLeila Ghaffari                                     "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS);
258a515125bSLeila Ghaffari   CHKERRQ(ierr);
259a515125bSLeila Ghaffari 
260a515125bSLeila Ghaffari   PetscFunctionReturn(0);
261a515125bSLeila Ghaffari }
26215a3537eSJed Brown 
26315a3537eSJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
26415a3537eSJed Brown int FreeContextPetsc(void *data) {
26515a3537eSJed Brown   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS,
26615a3537eSJed Brown                                           "PetscFree failed");
26715a3537eSJed Brown   return CEED_ERROR_SUCCESS;
26815a3537eSJed Brown }
269