xref: /libCEED/examples/fluids/src/misc.c (revision 0df12fefc12866395131d5961861faf43d192251)
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   PetscMemType q0_mem_type;
31c798d55aSJames Wright   PetscCall(VecP2C(Q_loc, &q0_mem_type, q0_ceed));
3277841947SLeila Ghaffari 
3377841947SLeila Ghaffari   // -- Apply CEED Operator
342b730f8bSJeremy L Thompson   CeedOperatorApply(ceed_data->op_ics, ceed_data->x_coord, q0_ceed, CEED_REQUEST_IMMEDIATE);
3577841947SLeila Ghaffari 
3677841947SLeila Ghaffari   // -- Restore vectors
37c798d55aSJames Wright   PetscCall(VecC2P(q0_ceed, q0_mem_type, Q_loc));
3877841947SLeila Ghaffari 
3977841947SLeila Ghaffari   // -- Local-to-Global
402b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(Q));
412b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, Q_loc, ADD_VALUES, Q));
4277841947SLeila Ghaffari 
4377841947SLeila Ghaffari   // ---------------------------------------------------------------------------
4477841947SLeila Ghaffari   // Fix multiplicity for output of ICs
4577841947SLeila Ghaffari   // ---------------------------------------------------------------------------
4677841947SLeila Ghaffari   // -- CEED Restriction
4777841947SLeila Ghaffari   CeedVector mult_vec;
4877841947SLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL);
4977841947SLeila Ghaffari 
5077841947SLeila Ghaffari   // -- Place PETSc vector in CEED vector
5177841947SLeila Ghaffari   PetscMemType m_mem_type;
5277841947SLeila Ghaffari   Vec          multiplicity_loc;
532b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &multiplicity_loc));
54c798d55aSJames Wright   PetscCall(VecP2C(multiplicity_loc, &m_mem_type, mult_vec));
5577841947SLeila Ghaffari 
5677841947SLeila Ghaffari   // -- Get multiplicity
5777841947SLeila Ghaffari   CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec);
5877841947SLeila Ghaffari 
5977841947SLeila Ghaffari   // -- Restore vectors
60c798d55aSJames Wright   PetscCall(VecC2P(mult_vec, m_mem_type, multiplicity_loc));
6177841947SLeila Ghaffari 
6277841947SLeila Ghaffari   // -- Local-to-Global
6377841947SLeila Ghaffari   Vec multiplicity;
642b730f8bSJeremy L Thompson   PetscCall(DMGetGlobalVector(dm, &multiplicity));
652b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(multiplicity));
662b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, multiplicity_loc, ADD_VALUES, multiplicity));
6777841947SLeila Ghaffari 
6877841947SLeila Ghaffari   // -- Fix multiplicity
692b730f8bSJeremy L Thompson   PetscCall(VecPointwiseDivide(Q, Q, multiplicity));
702b730f8bSJeremy L Thompson   PetscCall(VecPointwiseDivide(Q_loc, Q_loc, multiplicity_loc));
7177841947SLeila Ghaffari 
7277841947SLeila Ghaffari   // -- Restore vectors
732b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &multiplicity_loc));
742b730f8bSJeremy L Thompson   PetscCall(DMRestoreGlobalVector(dm, &multiplicity));
7577841947SLeila Ghaffari 
7677841947SLeila Ghaffari   // Cleanup
7777841947SLeila Ghaffari   CeedVectorDestroy(&mult_vec);
7877841947SLeila Ghaffari   CeedVectorDestroy(&q0_ceed);
7977841947SLeila Ghaffari 
8077841947SLeila Ghaffari   PetscFunctionReturn(0);
8177841947SLeila Ghaffari }
8277841947SLeila Ghaffari 
832b730f8bSJeremy L Thompson PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
842b730f8bSJeremy L Thompson                                              Vec grad_FVM) {
855571c6fdSJames Wright   Vec Qbc, boundary_mask;
8677841947SLeila Ghaffari   PetscFunctionBegin;
8777841947SLeila Ghaffari 
885571c6fdSJames Wright   // Mask (zero) Dirichlet entries
895571c6fdSJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
905571c6fdSJames Wright   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
915571c6fdSJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
925571c6fdSJames Wright 
932b730f8bSJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
942b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Q_loc, 1., Qbc));
952b730f8bSJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
9677841947SLeila Ghaffari 
9777841947SLeila Ghaffari   PetscFunctionReturn(0);
9877841947SLeila Ghaffari }
9977841947SLeila Ghaffari 
10077841947SLeila Ghaffari // Compare reference solution values with current test run for CI
10177841947SLeila Ghaffari PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q) {
10277841947SLeila Ghaffari   Vec         Qref;
10377841947SLeila Ghaffari   PetscViewer viewer;
10477841947SLeila Ghaffari   PetscReal   error, Qrefnorm;
10577841947SLeila Ghaffari   PetscFunctionBegin;
10677841947SLeila Ghaffari 
10777841947SLeila Ghaffari   // Read reference file
1082b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(Q, &Qref));
109b7d66439SJames Wright   PetscCall(PetscViewerBinaryOpen(PetscObjectComm((PetscObject)Q), app_ctx->test_file_path, FILE_MODE_READ, &viewer));
1102b730f8bSJeremy L Thompson   PetscCall(VecLoad(Qref, viewer));
11177841947SLeila Ghaffari 
11277841947SLeila Ghaffari   // Compute error with respect to reference solution
1132b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Qref));
1142b730f8bSJeremy L Thompson   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
1152b730f8bSJeremy L Thompson   PetscCall(VecScale(Q, 1. / Qrefnorm));
1162b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q, NORM_MAX, &error));
11777841947SLeila Ghaffari 
11877841947SLeila Ghaffari   // Check error
11977841947SLeila Ghaffari   if (error > app_ctx->test_tol) {
1202b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
12177841947SLeila Ghaffari   }
12277841947SLeila Ghaffari 
12377841947SLeila Ghaffari   // Cleanup
1242b730f8bSJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
1252b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&Qref));
12677841947SLeila Ghaffari 
12777841947SLeila Ghaffari   PetscFunctionReturn(0);
12877841947SLeila Ghaffari }
12977841947SLeila Ghaffari 
13077841947SLeila Ghaffari // Get error for problems with exact solutions
1312b730f8bSJeremy L Thompson PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
13277841947SLeila Ghaffari   PetscInt  loc_nodes;
13377841947SLeila Ghaffari   Vec       Q_exact, Q_exact_loc;
13477841947SLeila Ghaffari   PetscReal rel_error, norm_error, norm_exact;
13577841947SLeila Ghaffari   PetscFunctionBegin;
13677841947SLeila Ghaffari 
13777841947SLeila Ghaffari   // Get exact solution at final time
1382b730f8bSJeremy L Thompson   PetscCall(DMCreateGlobalVector(dm, &Q_exact));
1392b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1402b730f8bSJeremy L Thompson   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
1412b730f8bSJeremy L Thompson   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
14277841947SLeila Ghaffari 
14377841947SLeila Ghaffari   // Get |exact solution - obtained solution|
1442b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1452b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_exact));
1462b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q, NORM_1, &norm_error));
14777841947SLeila Ghaffari 
14877841947SLeila Ghaffari   // Compute relative error
14977841947SLeila Ghaffari   rel_error = norm_error / norm_exact;
15077841947SLeila Ghaffari 
15177841947SLeila Ghaffari   // Output relative error
1522b730f8bSJeremy L Thompson   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
15377841947SLeila Ghaffari   // Cleanup
1542b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
1552b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&Q_exact));
15677841947SLeila Ghaffari 
15777841947SLeila Ghaffari   PetscFunctionReturn(0);
15877841947SLeila Ghaffari }
15977841947SLeila Ghaffari 
16077841947SLeila Ghaffari // Post-processing
1612b730f8bSJeremy L Thompson PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) {
16277841947SLeila Ghaffari   PetscInt          steps;
163cf7a0454SJed Brown   TSConvergedReason reason;
16477841947SLeila Ghaffari   PetscFunctionBegin;
16577841947SLeila Ghaffari 
16677841947SLeila Ghaffari   // Print relative error
1678fb33541SJames Wright   if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) {
1682b730f8bSJeremy L Thompson     PetscCall(GetError_NS(ceed_data, dm, user, Q, final_time));
16977841947SLeila Ghaffari   }
17077841947SLeila Ghaffari 
17177841947SLeila Ghaffari   // Print final time and number of steps
1722b730f8bSJeremy L Thompson   PetscCall(TSGetStepNumber(ts, &steps));
173cf7a0454SJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
1748fb33541SJames Wright   if (user->app_ctx->test_type == TESTTYPE_NONE) {
175cf7a0454SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
176cf7a0454SJed Brown                           steps, (double)final_time));
17777841947SLeila Ghaffari   }
17877841947SLeila Ghaffari 
17977841947SLeila Ghaffari   // Output numerical values from command line
1802b730f8bSJeremy L Thompson   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
18177841947SLeila Ghaffari 
18277841947SLeila Ghaffari   // Compare reference solution values with current test run for CI
1838fb33541SJames Wright   if (user->app_ctx->test_type == TESTTYPE_SOLVER) {
1842b730f8bSJeremy L Thompson     PetscCall(RegressionTests_NS(user->app_ctx, Q));
18577841947SLeila Ghaffari   }
18677841947SLeila Ghaffari 
18777841947SLeila Ghaffari   PetscFunctionReturn(0);
18877841947SLeila Ghaffari }
18977841947SLeila Ghaffari 
1904de8550aSJed Brown const PetscInt FLUIDS_FILE_TOKEN = 0xceedf00;
1914de8550aSJed Brown 
19277841947SLeila Ghaffari // Gather initial Q values in case of continuation of simulation
19377841947SLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
19477841947SLeila Ghaffari   PetscViewer viewer;
1954de8550aSJed Brown   PetscInt    token, step_number;
1964de8550aSJed Brown   PetscReal   time;
1972b730f8bSJeremy L Thompson 
19877841947SLeila Ghaffari   PetscFunctionBegin;
19977841947SLeila Ghaffari 
20077841947SLeila Ghaffari   // Read input
2012b730f8bSJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
20277841947SLeila Ghaffari 
2034de8550aSJed Brown   // Attempt
2044de8550aSJed Brown   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT));
2054de8550aSJed Brown   if (token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
2064de8550aSJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &step_number, 1, NULL, PETSC_INT));
2074de8550aSJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &time, 1, NULL, PETSC_REAL));
2084de8550aSJed Brown     app_ctx->cont_steps = step_number;
2094de8550aSJed Brown     app_ctx->cont_time  = time;
2104de8550aSJed Brown   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
2114de8550aSJed Brown     PetscInt length, N;
2124de8550aSJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &length, 1, NULL, PETSC_INT));
2134de8550aSJed Brown     PetscCall(VecGetSize(Q, &N));
2144de8550aSJed 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);
2154de8550aSJed Brown     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
2164de8550aSJed Brown   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
2174de8550aSJed Brown 
21877841947SLeila Ghaffari   // Load Q from existent solution
2192b730f8bSJeremy L Thompson   PetscCall(VecLoad(Q, viewer));
22077841947SLeila Ghaffari 
22177841947SLeila Ghaffari   // Cleanup
2222b730f8bSJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
22377841947SLeila Ghaffari 
22477841947SLeila Ghaffari   PetscFunctionReturn(0);
22577841947SLeila Ghaffari }
22677841947SLeila Ghaffari 
22777841947SLeila Ghaffari // Record boundary values from initial condition
22877841947SLeila Ghaffari PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc) {
2295571c6fdSJames Wright   Vec Qbc, boundary_mask;
23077841947SLeila Ghaffari   PetscFunctionBegin;
23177841947SLeila Ghaffari 
2322b730f8bSJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
2332b730f8bSJeremy L Thompson   PetscCall(VecCopy(Q_loc, Qbc));
2342b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(Q_loc));
2352b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
2362b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Qbc, -1., Q_loc));
2372b730f8bSJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
2382b730f8bSJeremy L Thompson   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS));
23977841947SLeila Ghaffari 
2405571c6fdSJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
2415571c6fdSJames Wright   PetscCall(DMGetGlobalVector(dm, &Q));
2425571c6fdSJames Wright   PetscCall(VecZeroEntries(boundary_mask));
2435571c6fdSJames Wright   PetscCall(VecSet(Q, 1.0));
2445571c6fdSJames Wright   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask));
2455571c6fdSJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
2465571c6fdSJames Wright 
24777841947SLeila Ghaffari   PetscFunctionReturn(0);
24877841947SLeila Ghaffari }
249841e4c73SJed Brown 
250841e4c73SJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
251841e4c73SJed Brown int FreeContextPetsc(void *data) {
2522b730f8bSJeremy L Thompson   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
253841e4c73SJed Brown   return CEED_ERROR_SUCCESS;
254841e4c73SJed Brown }
255ef080ff9SJames Wright 
256ef080ff9SJames Wright // Return mass qfunction specification for number of components N
257ef080ff9SJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
258ef080ff9SJames Wright   CeedQFunctionUser qfunction_ptr;
259ef080ff9SJames Wright   const char       *qfunction_loc;
260ef080ff9SJames Wright   PetscFunctionBeginUser;
261ef080ff9SJames Wright 
262ef080ff9SJames Wright   switch (N) {
263ef080ff9SJames Wright     case 1:
264ef080ff9SJames Wright       qfunction_ptr = Mass_1;
265ef080ff9SJames Wright       qfunction_loc = Mass_1_loc;
266ef080ff9SJames Wright       break;
267ef080ff9SJames Wright     case 5:
268ef080ff9SJames Wright       qfunction_ptr = Mass_5;
269ef080ff9SJames Wright       qfunction_loc = Mass_5_loc;
270ef080ff9SJames Wright       break;
271ef080ff9SJames Wright     case 9:
272ef080ff9SJames Wright       qfunction_ptr = Mass_9;
273ef080ff9SJames Wright       qfunction_loc = Mass_9_loc;
274ef080ff9SJames Wright       break;
275ef080ff9SJames Wright     case 22:
276ef080ff9SJames Wright       qfunction_ptr = Mass_22;
277ef080ff9SJames Wright       qfunction_loc = Mass_22_loc;
278ef080ff9SJames Wright       break;
279ef080ff9SJames Wright     default:
280ef080ff9SJames Wright       SETERRQ(PETSC_COMM_WORLD, -1, "Could not find mass qfunction of size %d", N);
281ef080ff9SJames Wright   }
282ef080ff9SJames Wright   CeedQFunctionCreateInterior(ceed, 1, qfunction_ptr, qfunction_loc, qf);
283ef080ff9SJames Wright 
284ef080ff9SJames Wright   CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP);
285ef080ff9SJames Wright   CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE);
286ef080ff9SJames Wright   CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP);
287ef080ff9SJames Wright   PetscFunctionReturn(0);
288ef080ff9SJames Wright }
289*0df12fefSJames Wright 
290*0df12fefSJames Wright /* @brief L^2 Projection of a source FEM function to a target FEM space
291*0df12fefSJames Wright  *
292*0df12fefSJames Wright  * To solve system using a lumped mass matrix, pass a KSP object with ksp_type=preonly, pc_type=jacobi, pc_jacobi_type=rowsum.
293*0df12fefSJames Wright  *
294*0df12fefSJames Wright  * @param[in]  source_vec    Global Vec of the source FEM function. NULL indicates using rhs_matop_ctx->X_loc
295*0df12fefSJames Wright  * @param[out] target_vec    Global Vec of the target (result) FEM function. NULL indicates using rhs_matop_ctx->Y_loc
296*0df12fefSJames Wright  * @param[in]  rhs_matop_ctx MatopApplyContext for performing the RHS evaluation
297*0df12fefSJames Wright  * @param[in]  ksp           KSP for solving the consistent projection problem
298*0df12fefSJames Wright  */
299*0df12fefSJames Wright PetscErrorCode ComputeL2Projection(Vec source_vec, Vec target_vec, MatopApplyContext rhs_matop_ctx, KSP ksp) {
300*0df12fefSJames Wright   PetscFunctionBeginUser;
301*0df12fefSJames Wright 
302*0df12fefSJames Wright   PetscCall(ApplyLocal_Ceed(source_vec, target_vec, rhs_matop_ctx));
303*0df12fefSJames Wright   PetscCall(KSPSolve(ksp, target_vec, target_vec));
304*0df12fefSJames Wright 
305*0df12fefSJames Wright   PetscFunctionReturn(0);
306*0df12fefSJames Wright }
307