xref: /libCEED/examples/fluids/src/misc.c (revision b7d6643930f63dbc81dde690883cc613b302935e)
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"
12ef080ff9SJames Wright #include "../qfunctions/mass.h"
1377841947SLeila Ghaffari 
142b730f8bSJeremy L Thompson PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) {
1577841947SLeila Ghaffari   PetscFunctionBeginUser;
1677841947SLeila Ghaffari 
1777841947SLeila Ghaffari   // ---------------------------------------------------------------------------
18a0add3c9SJed Brown   // Update time for evaluation
1977841947SLeila Ghaffari   // ---------------------------------------------------------------------------
202b730f8bSJeremy L Thompson   if (user->phys->ics_time_label) CeedOperatorContextSetDouble(ceed_data->op_ics, user->phys->ics_time_label, &time);
2177841947SLeila Ghaffari 
2277841947SLeila Ghaffari   // ---------------------------------------------------------------------------
2377841947SLeila Ghaffari   // ICs
2477841947SLeila Ghaffari   // ---------------------------------------------------------------------------
2577841947SLeila Ghaffari   // -- CEED Restriction
2677841947SLeila Ghaffari   CeedVector q0_ceed;
2777841947SLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &q0_ceed, NULL);
2877841947SLeila Ghaffari 
2977841947SLeila Ghaffari   // -- Place PETSc vector in CEED vector
3077841947SLeila Ghaffari   CeedScalar  *q0;
3177841947SLeila Ghaffari   PetscMemType q0_mem_type;
322b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(Q_loc, (PetscScalar **)&q0, &q0_mem_type));
3377841947SLeila Ghaffari   CeedVectorSetArray(q0_ceed, MemTypeP2C(q0_mem_type), CEED_USE_POINTER, q0);
3477841947SLeila Ghaffari 
3577841947SLeila Ghaffari   // -- Apply CEED Operator
362b730f8bSJeremy L Thompson   CeedOperatorApply(ceed_data->op_ics, ceed_data->x_coord, q0_ceed, CEED_REQUEST_IMMEDIATE);
3777841947SLeila Ghaffari 
3877841947SLeila Ghaffari   // -- Restore vectors
3977841947SLeila Ghaffari   CeedVectorTakeArray(q0_ceed, MemTypeP2C(q0_mem_type), NULL);
402b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, (const PetscScalar **)&q0));
4177841947SLeila Ghaffari 
4277841947SLeila Ghaffari   // -- Local-to-Global
432b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(Q));
442b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, Q_loc, ADD_VALUES, Q));
4577841947SLeila Ghaffari 
4677841947SLeila Ghaffari   // ---------------------------------------------------------------------------
4777841947SLeila Ghaffari   // Fix multiplicity for output of ICs
4877841947SLeila Ghaffari   // ---------------------------------------------------------------------------
4977841947SLeila Ghaffari   // -- CEED Restriction
5077841947SLeila Ghaffari   CeedVector mult_vec;
5177841947SLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL);
5277841947SLeila Ghaffari 
5377841947SLeila Ghaffari   // -- Place PETSc vector in CEED vector
5477841947SLeila Ghaffari   CeedScalar  *mult;
5577841947SLeila Ghaffari   PetscMemType m_mem_type;
5677841947SLeila Ghaffari   Vec          multiplicity_loc;
572b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &multiplicity_loc));
582b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(multiplicity_loc, (PetscScalar **)&mult, &m_mem_type));
5977841947SLeila Ghaffari   CeedVectorSetArray(mult_vec, MemTypeP2C(m_mem_type), CEED_USE_POINTER, mult);
6077841947SLeila Ghaffari 
6177841947SLeila Ghaffari   // -- Get multiplicity
6277841947SLeila Ghaffari   CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec);
6377841947SLeila Ghaffari 
6477841947SLeila Ghaffari   // -- Restore vectors
6577841947SLeila Ghaffari   CeedVectorTakeArray(mult_vec, MemTypeP2C(m_mem_type), NULL);
662b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(multiplicity_loc, (const PetscScalar **)&mult));
6777841947SLeila Ghaffari 
6877841947SLeila Ghaffari   // -- Local-to-Global
6977841947SLeila Ghaffari   Vec multiplicity;
702b730f8bSJeremy L Thompson   PetscCall(DMGetGlobalVector(dm, &multiplicity));
712b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(multiplicity));
722b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, multiplicity_loc, ADD_VALUES, multiplicity));
7377841947SLeila Ghaffari 
7477841947SLeila Ghaffari   // -- Fix multiplicity
752b730f8bSJeremy L Thompson   PetscCall(VecPointwiseDivide(Q, Q, multiplicity));
762b730f8bSJeremy L Thompson   PetscCall(VecPointwiseDivide(Q_loc, Q_loc, multiplicity_loc));
7777841947SLeila Ghaffari 
7877841947SLeila Ghaffari   // -- Restore vectors
792b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &multiplicity_loc));
802b730f8bSJeremy L Thompson   PetscCall(DMRestoreGlobalVector(dm, &multiplicity));
8177841947SLeila Ghaffari 
8277841947SLeila Ghaffari   // Cleanup
8377841947SLeila Ghaffari   CeedVectorDestroy(&mult_vec);
8477841947SLeila Ghaffari   CeedVectorDestroy(&q0_ceed);
8577841947SLeila Ghaffari 
8677841947SLeila Ghaffari   PetscFunctionReturn(0);
8777841947SLeila Ghaffari }
8877841947SLeila Ghaffari 
892b730f8bSJeremy L Thompson PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
902b730f8bSJeremy L Thompson                                              Vec grad_FVM) {
915571c6fdSJames Wright   Vec Qbc, boundary_mask;
9277841947SLeila Ghaffari   PetscFunctionBegin;
9377841947SLeila Ghaffari 
945571c6fdSJames Wright   // Mask (zero) Dirichlet entries
955571c6fdSJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
965571c6fdSJames Wright   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
975571c6fdSJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
985571c6fdSJames Wright 
992b730f8bSJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
1002b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Q_loc, 1., Qbc));
1012b730f8bSJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
10277841947SLeila Ghaffari 
10377841947SLeila Ghaffari   PetscFunctionReturn(0);
10477841947SLeila Ghaffari }
10577841947SLeila Ghaffari 
10677841947SLeila Ghaffari // Compare reference solution values with current test run for CI
10777841947SLeila Ghaffari PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q) {
10877841947SLeila Ghaffari   Vec         Qref;
10977841947SLeila Ghaffari   PetscViewer viewer;
11077841947SLeila Ghaffari   PetscReal   error, Qrefnorm;
11177841947SLeila Ghaffari   PetscFunctionBegin;
11277841947SLeila Ghaffari 
11377841947SLeila Ghaffari   // Read reference file
1142b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(Q, &Qref));
115*b7d66439SJames Wright   PetscCall(PetscViewerBinaryOpen(PetscObjectComm((PetscObject)Q), app_ctx->test_file_path, FILE_MODE_READ, &viewer));
1162b730f8bSJeremy L Thompson   PetscCall(VecLoad(Qref, viewer));
11777841947SLeila Ghaffari 
11877841947SLeila Ghaffari   // Compute error with respect to reference solution
1192b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Qref));
1202b730f8bSJeremy L Thompson   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
1212b730f8bSJeremy L Thompson   PetscCall(VecScale(Q, 1. / Qrefnorm));
1222b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q, NORM_MAX, &error));
12377841947SLeila Ghaffari 
12477841947SLeila Ghaffari   // Check error
12577841947SLeila Ghaffari   if (error > app_ctx->test_tol) {
1262b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
12777841947SLeila Ghaffari   }
12877841947SLeila Ghaffari 
12977841947SLeila Ghaffari   // Cleanup
1302b730f8bSJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
1312b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&Qref));
13277841947SLeila Ghaffari 
13377841947SLeila Ghaffari   PetscFunctionReturn(0);
13477841947SLeila Ghaffari }
13577841947SLeila Ghaffari 
13677841947SLeila Ghaffari // Get error for problems with exact solutions
1372b730f8bSJeremy L Thompson PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
13877841947SLeila Ghaffari   PetscInt  loc_nodes;
13977841947SLeila Ghaffari   Vec       Q_exact, Q_exact_loc;
14077841947SLeila Ghaffari   PetscReal rel_error, norm_error, norm_exact;
14177841947SLeila Ghaffari   PetscFunctionBegin;
14277841947SLeila Ghaffari 
14377841947SLeila Ghaffari   // Get exact solution at final time
1442b730f8bSJeremy L Thompson   PetscCall(DMCreateGlobalVector(dm, &Q_exact));
1452b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1462b730f8bSJeremy L Thompson   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
1472b730f8bSJeremy L Thompson   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
14877841947SLeila Ghaffari 
14977841947SLeila Ghaffari   // Get |exact solution - obtained solution|
1502b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1512b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_exact));
1522b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q, NORM_1, &norm_error));
15377841947SLeila Ghaffari 
15477841947SLeila Ghaffari   // Compute relative error
15577841947SLeila Ghaffari   rel_error = norm_error / norm_exact;
15677841947SLeila Ghaffari 
15777841947SLeila Ghaffari   // Output relative error
1582b730f8bSJeremy L Thompson   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
15977841947SLeila Ghaffari   // Cleanup
1602b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
1612b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&Q_exact));
16277841947SLeila Ghaffari 
16377841947SLeila Ghaffari   PetscFunctionReturn(0);
16477841947SLeila Ghaffari }
16577841947SLeila Ghaffari 
16677841947SLeila Ghaffari // Post-processing
1672b730f8bSJeremy L Thompson PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) {
16877841947SLeila Ghaffari   PetscInt          steps;
169cf7a0454SJed Brown   TSConvergedReason reason;
17077841947SLeila Ghaffari   PetscFunctionBegin;
17177841947SLeila Ghaffari 
17277841947SLeila Ghaffari   // Print relative error
173841e4c73SJed Brown   if (problem->non_zero_time && !user->app_ctx->test_mode) {
1742b730f8bSJeremy L Thompson     PetscCall(GetError_NS(ceed_data, dm, user, Q, final_time));
17577841947SLeila Ghaffari   }
17677841947SLeila Ghaffari 
17777841947SLeila Ghaffari   // Print final time and number of steps
1782b730f8bSJeremy L Thompson   PetscCall(TSGetStepNumber(ts, &steps));
179cf7a0454SJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
180841e4c73SJed Brown   if (!user->app_ctx->test_mode) {
181cf7a0454SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
182cf7a0454SJed Brown                           steps, (double)final_time));
18377841947SLeila Ghaffari   }
18477841947SLeila Ghaffari 
18577841947SLeila Ghaffari   // Output numerical values from command line
1862b730f8bSJeremy L Thompson   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
18777841947SLeila Ghaffari 
18877841947SLeila Ghaffari   // Compare reference solution values with current test run for CI
189841e4c73SJed Brown   if (user->app_ctx->test_mode) {
1902b730f8bSJeremy L Thompson     PetscCall(RegressionTests_NS(user->app_ctx, Q));
19177841947SLeila Ghaffari   }
19277841947SLeila Ghaffari 
19377841947SLeila Ghaffari   PetscFunctionReturn(0);
19477841947SLeila Ghaffari }
19577841947SLeila Ghaffari 
1964de8550aSJed Brown const PetscInt FLUIDS_FILE_TOKEN = 0xceedf00;
1974de8550aSJed Brown 
19877841947SLeila Ghaffari // Gather initial Q values in case of continuation of simulation
19977841947SLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
20077841947SLeila Ghaffari   PetscViewer viewer;
2014de8550aSJed Brown   PetscInt    token, step_number;
2024de8550aSJed Brown   PetscReal   time;
2032b730f8bSJeremy L Thompson 
20477841947SLeila Ghaffari   PetscFunctionBegin;
20577841947SLeila Ghaffari 
20677841947SLeila Ghaffari   // Read input
2072b730f8bSJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
20877841947SLeila Ghaffari 
2094de8550aSJed Brown   // Attempt
2104de8550aSJed Brown   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT));
2114de8550aSJed Brown   if (token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
2124de8550aSJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &step_number, 1, NULL, PETSC_INT));
2134de8550aSJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &time, 1, NULL, PETSC_REAL));
2144de8550aSJed Brown     app_ctx->cont_steps = step_number;
2154de8550aSJed Brown     app_ctx->cont_time  = time;
2164de8550aSJed Brown   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
2174de8550aSJed Brown     PetscInt length, N;
2184de8550aSJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &length, 1, NULL, PETSC_INT));
2194de8550aSJed Brown     PetscCall(VecGetSize(Q, &N));
2204de8550aSJed 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);
2214de8550aSJed Brown     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
2224de8550aSJed Brown   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
2234de8550aSJed Brown 
22477841947SLeila Ghaffari   // Load Q from existent solution
2252b730f8bSJeremy L Thompson   PetscCall(VecLoad(Q, viewer));
22677841947SLeila Ghaffari 
22777841947SLeila Ghaffari   // Cleanup
2282b730f8bSJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
22977841947SLeila Ghaffari 
23077841947SLeila Ghaffari   PetscFunctionReturn(0);
23177841947SLeila Ghaffari }
23277841947SLeila Ghaffari 
23377841947SLeila Ghaffari // Record boundary values from initial condition
23477841947SLeila Ghaffari PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc) {
2355571c6fdSJames Wright   Vec Qbc, boundary_mask;
23677841947SLeila Ghaffari   PetscFunctionBegin;
23777841947SLeila Ghaffari 
2382b730f8bSJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
2392b730f8bSJeremy L Thompson   PetscCall(VecCopy(Q_loc, Qbc));
2402b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(Q_loc));
2412b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
2422b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Qbc, -1., Q_loc));
2432b730f8bSJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
2442b730f8bSJeremy L Thompson   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS));
24577841947SLeila Ghaffari 
2465571c6fdSJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
2475571c6fdSJames Wright   PetscCall(DMGetGlobalVector(dm, &Q));
2485571c6fdSJames Wright   PetscCall(VecZeroEntries(boundary_mask));
2495571c6fdSJames Wright   PetscCall(VecSet(Q, 1.0));
2505571c6fdSJames Wright   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask));
2515571c6fdSJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
2525571c6fdSJames Wright 
25377841947SLeila Ghaffari   PetscFunctionReturn(0);
25477841947SLeila Ghaffari }
255841e4c73SJed Brown 
256841e4c73SJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
257841e4c73SJed Brown int FreeContextPetsc(void *data) {
2582b730f8bSJeremy L Thompson   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
259841e4c73SJed Brown   return CEED_ERROR_SUCCESS;
260841e4c73SJed Brown }
261ef080ff9SJames Wright 
262ef080ff9SJames Wright // Return mass qfunction specification for number of components N
263ef080ff9SJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
264ef080ff9SJames Wright   CeedQFunctionUser qfunction_ptr;
265ef080ff9SJames Wright   const char       *qfunction_loc;
266ef080ff9SJames Wright   PetscFunctionBeginUser;
267ef080ff9SJames Wright 
268ef080ff9SJames Wright   switch (N) {
269ef080ff9SJames Wright     case 1:
270ef080ff9SJames Wright       qfunction_ptr = Mass_1;
271ef080ff9SJames Wright       qfunction_loc = Mass_1_loc;
272ef080ff9SJames Wright       break;
273ef080ff9SJames Wright     case 5:
274ef080ff9SJames Wright       qfunction_ptr = Mass_5;
275ef080ff9SJames Wright       qfunction_loc = Mass_5_loc;
276ef080ff9SJames Wright       break;
277ef080ff9SJames Wright     case 9:
278ef080ff9SJames Wright       qfunction_ptr = Mass_9;
279ef080ff9SJames Wright       qfunction_loc = Mass_9_loc;
280ef080ff9SJames Wright       break;
281ef080ff9SJames Wright     case 22:
282ef080ff9SJames Wright       qfunction_ptr = Mass_22;
283ef080ff9SJames Wright       qfunction_loc = Mass_22_loc;
284ef080ff9SJames Wright       break;
285ef080ff9SJames Wright     default:
286ef080ff9SJames Wright       SETERRQ(PETSC_COMM_WORLD, -1, "Could not find mass qfunction of size %d", N);
287ef080ff9SJames Wright   }
288ef080ff9SJames Wright   CeedQFunctionCreateInterior(ceed, 1, qfunction_ptr, qfunction_loc, qf);
289ef080ff9SJames Wright 
290ef080ff9SJames Wright   CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP);
291ef080ff9SJames Wright   CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE);
292ef080ff9SJames Wright   CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP);
293ef080ff9SJames Wright   PetscFunctionReturn(0);
294ef080ff9SJames Wright }
295