xref: /libCEED/examples/solids/src/setup-libceed.c (revision 2b730f8b5a9c809740a0b3b302db43a719c636b1)
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.
33d8e8822SJeremy L Thompson //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
53d8e8822SJeremy L Thompson //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
73d8e8822SJeremy L Thompson 
85754ecacSJeremy L Thompson /// @file
95754ecacSJeremy L Thompson /// libCEED setup for solid mechanics example using PETSc
105754ecacSJeremy L Thompson 
115754ecacSJeremy L Thompson #include "../include/setup-libceed.h"
12*2b730f8bSJeremy L Thompson 
135754ecacSJeremy L Thompson #include "../include/structs.h"
145754ecacSJeremy L Thompson #include "../include/utils.h"
155754ecacSJeremy L Thompson #include "../qfunctions/constant-force.h"      // Constant forcing function
165754ecacSJeremy L Thompson #include "../qfunctions/manufactured-force.h"  // Manufactured solution forcing
17*2b730f8bSJeremy L Thompson #include "../qfunctions/traction-boundary.h"   // Traction boundaries
185754ecacSJeremy L Thompson 
195754ecacSJeremy L Thompson #if PETSC_VERSION_LT(3, 14, 0)
205754ecacSJeremy L Thompson #define DMPlexGetClosureIndices(a, b, c, d, e, f, g, h, i) DMPlexGetClosureIndices(a, b, c, d, f, g, i)
215754ecacSJeremy L Thompson #define DMPlexRestoreClosureIndices(a, b, c, d, e, f, g, h, i) DMPlexRestoreClosureIndices(a, b, c, d, f, g, i)
225754ecacSJeremy L Thompson #endif
235754ecacSJeremy L Thompson 
245754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
255754ecacSJeremy L Thompson // Problem options
265754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
275754ecacSJeremy L Thompson // Forcing function data
285754ecacSJeremy L Thompson forcingData forcing_options[3] = {
29*2b730f8bSJeremy L Thompson     [FORCE_NONE]  = {.setup_forcing = NULL,               .setup_forcing_loc = NULL                  },
30*2b730f8bSJeremy L Thompson     [FORCE_CONST] = {.setup_forcing = SetupConstantForce, .setup_forcing_loc = SetupConstantForce_loc},
31*2b730f8bSJeremy L Thompson     [FORCE_MMS]   = {.setup_forcing = SetupMMSForce,      .setup_forcing_loc = SetupMMSForce_loc     }
325754ecacSJeremy L Thompson };
335754ecacSJeremy L Thompson 
345754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
355754ecacSJeremy L Thompson // libCEED Functions
365754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
375754ecacSJeremy L Thompson // Destroy libCEED objects
385754ecacSJeremy L Thompson PetscErrorCode CeedDataDestroy(CeedInt level, CeedData data) {
395754ecacSJeremy L Thompson   PetscFunctionBegin;
405754ecacSJeremy L Thompson 
415754ecacSJeremy L Thompson   // Vectors
425754ecacSJeremy L Thompson   CeedVectorDestroy(&data->x_ceed);
435754ecacSJeremy L Thompson   CeedVectorDestroy(&data->y_ceed);
445754ecacSJeremy L Thompson   CeedVectorDestroy(&data->geo_data);
45*2b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < SOLIDS_MAX_NUMBER_FIELDS; i++) CeedVectorDestroy(&data->stored_fields[i]);
465754ecacSJeremy L Thompson   CeedVectorDestroy(&data->geo_data_diagnostic);
475754ecacSJeremy L Thompson   CeedVectorDestroy(&data->true_soln);
485754ecacSJeremy L Thompson   // Restrictions
495754ecacSJeremy L Thompson   CeedElemRestrictionDestroy(&data->elem_restr_x);
505754ecacSJeremy L Thompson   CeedElemRestrictionDestroy(&data->elem_restr_u);
515754ecacSJeremy L Thompson   CeedElemRestrictionDestroy(&data->elem_restr_geo_data_i);
52*2b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < SOLIDS_MAX_NUMBER_FIELDS; i++) CeedElemRestrictionDestroy(&data->elem_restr_stored_fields_i[i]);
535754ecacSJeremy L Thompson   CeedElemRestrictionDestroy(&data->elem_restr_energy);
545754ecacSJeremy L Thompson   CeedElemRestrictionDestroy(&data->elem_restr_diagnostic);
555754ecacSJeremy L Thompson   CeedElemRestrictionDestroy(&data->elem_restr_geo_data_diagnostic_i);
565754ecacSJeremy L Thompson   // Bases
575754ecacSJeremy L Thompson   CeedBasisDestroy(&data->basis_x);
585754ecacSJeremy L Thompson   CeedBasisDestroy(&data->basis_u);
595754ecacSJeremy L Thompson   CeedBasisDestroy(&data->basis_energy);
605754ecacSJeremy L Thompson   CeedBasisDestroy(&data->basis_diagnostic);
615754ecacSJeremy L Thompson   // QFunctions
625754ecacSJeremy L Thompson   CeedQFunctionDestroy(&data->qf_residual);
635754ecacSJeremy L Thompson   CeedQFunctionDestroy(&data->qf_jacobian);
645754ecacSJeremy L Thompson   CeedQFunctionDestroy(&data->qf_energy);
655754ecacSJeremy L Thompson   CeedQFunctionDestroy(&data->qf_diagnostic);
665754ecacSJeremy L Thompson   // Operators
675754ecacSJeremy L Thompson   CeedOperatorDestroy(&data->op_residual);
685754ecacSJeremy L Thompson   CeedOperatorDestroy(&data->op_jacobian);
695754ecacSJeremy L Thompson   CeedOperatorDestroy(&data->op_energy);
705754ecacSJeremy L Thompson   CeedOperatorDestroy(&data->op_diagnostic);
715754ecacSJeremy L Thompson   // Restriction and Prolongation data
725754ecacSJeremy L Thompson   CeedBasisDestroy(&data->basis_c_to_f);
735754ecacSJeremy L Thompson   CeedOperatorDestroy(&data->op_prolong);
745754ecacSJeremy L Thompson   CeedOperatorDestroy(&data->op_restrict);
755754ecacSJeremy L Thompson 
76*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(data));
775754ecacSJeremy L Thompson 
785754ecacSJeremy L Thompson   PetscFunctionReturn(0);
795754ecacSJeremy L Thompson };
805754ecacSJeremy L Thompson 
815754ecacSJeremy L Thompson // Utility function to create local CEED restriction from DMPlex
82*2b730f8bSJeremy L Thompson PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr) {
837ed3e4cdSJeremy L Thompson   PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets;
845754ecacSJeremy L Thompson 
855754ecacSJeremy L Thompson   PetscFunctionBeginUser;
86*2b730f8bSJeremy L Thompson   PetscCall(DMPlexGetLocalOffsets(dm, domain_label, value, height, 0, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets));
875754ecacSJeremy L Thompson 
88*2b730f8bSJeremy L Thompson   CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, elem_restr_offsets, elem_restr);
89*2b730f8bSJeremy L Thompson   PetscCall(PetscFree(elem_restr_offsets));
905754ecacSJeremy L Thompson 
915754ecacSJeremy L Thompson   PetscFunctionReturn(0);
925754ecacSJeremy L Thompson };
935754ecacSJeremy L Thompson 
945754ecacSJeremy L Thompson // Utility function to get Ceed Restriction for each domain
95*2b730f8bSJeremy L Thompson PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, PetscInt value, CeedInt Q, CeedInt q_data_size,
96*2b730f8bSJeremy L Thompson                                        CeedElemRestriction *elem_restr_q, CeedElemRestriction *elem_restr_x, CeedElemRestriction *elem_restr_qd_i) {
975754ecacSJeremy L Thompson   DM      dm_coord;
985754ecacSJeremy L Thompson   CeedInt dim, num_local_elem;
995754ecacSJeremy L Thompson   CeedInt Q_dim;
1005754ecacSJeremy L Thompson 
1015754ecacSJeremy L Thompson   PetscFunctionBeginUser;
1025754ecacSJeremy L Thompson 
103*2b730f8bSJeremy L Thompson   PetscCall(DMGetDimension(dm, &dim));
1045754ecacSJeremy L Thompson   dim -= height;
1055754ecacSJeremy L Thompson   Q_dim = CeedIntPow(Q, dim);
106*2b730f8bSJeremy L Thompson   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
107*2b730f8bSJeremy L Thompson   PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL));
1085754ecacSJeremy L Thompson   if (elem_restr_q) {
109*2b730f8bSJeremy L Thompson     PetscCall(CreateRestrictionFromPlex(ceed, dm, height, domain_label, value, elem_restr_q));
1105754ecacSJeremy L Thompson   }
1115754ecacSJeremy L Thompson   if (elem_restr_x) {
112*2b730f8bSJeremy L Thompson     PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, height, domain_label, value, elem_restr_x));
1135754ecacSJeremy L Thompson   }
1145754ecacSJeremy L Thompson   if (elem_restr_qd_i) {
1155754ecacSJeremy L Thompson     CeedElemRestrictionGetNumElements(*elem_restr_q, &num_local_elem);
116*2b730f8bSJeremy L Thompson     CeedElemRestrictionCreateStrided(ceed, num_local_elem, Q_dim, q_data_size, q_data_size * num_local_elem * Q_dim, CEED_STRIDES_BACKEND,
117*2b730f8bSJeremy L Thompson                                      elem_restr_qd_i);
1185754ecacSJeremy L Thompson   }
1195754ecacSJeremy L Thompson 
1205754ecacSJeremy L Thompson   PetscFunctionReturn(0);
1215754ecacSJeremy L Thompson };
1225754ecacSJeremy L Thompson 
1235754ecacSJeremy L Thompson // Set up libCEED on the fine grid for a given degree
124*2b730f8bSJeremy L Thompson PetscErrorCode SetupLibceedFineLevel(DM dm, DM dm_energy, DM dm_diagnostic, Ceed ceed, AppCtx app_ctx, CeedQFunctionContext phys_ctx,
125*2b730f8bSJeremy L Thompson                                      ProblemData problem_data, PetscInt fine_level, PetscInt num_comp_u, PetscInt U_g_size, PetscInt U_loc_size,
126*2b730f8bSJeremy L Thompson                                      CeedVector force_ceed, CeedVector neumann_ceed, CeedData *data) {
1275754ecacSJeremy L Thompson   CeedInt            P = app_ctx->level_degrees[fine_level] + 1;
1285754ecacSJeremy L Thompson   CeedInt            Q = app_ctx->level_degrees[fine_level] + 1 + app_ctx->q_extra;
1295754ecacSJeremy L Thompson   CeedInt            dim, num_comp_x, num_comp_e = 1, num_comp_d = 5;
1305754ecacSJeremy L Thompson   CeedInt            num_qpts;
131a61c78d6SJeremy L Thompson   CeedInt            q_data_size    = problem_data.q_data_size;
1325754ecacSJeremy L Thompson   forcingType        forcing_choice = app_ctx->forcing_choice;
1335754ecacSJeremy L Thompson   DM                 dm_coord;
1345754ecacSJeremy L Thompson   Vec                coords;
1355754ecacSJeremy L Thompson   PetscInt           c_start, c_end, num_elem;
1365754ecacSJeremy L Thompson   const PetscScalar *coordArray;
1375754ecacSJeremy L Thompson   CeedVector         x_coord;
1385754ecacSJeremy L Thompson   CeedQFunction      qf_setup_geo, qf_residual, qf_jacobian, qf_energy, qf_diagnostic;
1395754ecacSJeremy L Thompson   CeedOperator       op_setup_geo, op_residual, op_jacobian, op_energy, op_diagnostic;
1405754ecacSJeremy L Thompson 
1415754ecacSJeremy L Thompson   PetscFunctionBeginUser;
1425754ecacSJeremy L Thompson 
1435754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
1445754ecacSJeremy L Thompson   // libCEED bases
1455754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
146*2b730f8bSJeremy L Thompson   PetscCall(DMGetDimension(dm, &dim));
1475754ecacSJeremy L Thompson   num_comp_x = dim;
1485754ecacSJeremy L Thompson   // -- Coordinate basis
149*2b730f8bSJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, problem_data.quadrature_mode, &data[fine_level]->basis_x);
1505754ecacSJeremy L Thompson   // -- Solution basis
151*2b730f8bSJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_u, P, Q, problem_data.quadrature_mode, &data[fine_level]->basis_u);
1525754ecacSJeremy L Thompson   // -- Energy basis
153*2b730f8bSJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_e, P, Q, problem_data.quadrature_mode, &data[fine_level]->basis_energy);
1545754ecacSJeremy L Thompson   // -- Diagnostic output basis
155*2b730f8bSJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_u, P, P, CEED_GAUSS_LOBATTO, &data[fine_level]->basis_diagnostic);
1565754ecacSJeremy L Thompson 
1575754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
1585754ecacSJeremy L Thompson   // libCEED restrictions
1595754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
160*2b730f8bSJeremy L Thompson   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
161*2b730f8bSJeremy L Thompson   PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL));
1625754ecacSJeremy L Thompson 
1635754ecacSJeremy L Thompson   // -- Coordinate restriction
164*2b730f8bSJeremy L Thompson   PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, 0, 0, 0, &(data[fine_level]->elem_restr_x)));
1655754ecacSJeremy L Thompson   // -- Solution restriction
166*2b730f8bSJeremy L Thompson   PetscCall(CreateRestrictionFromPlex(ceed, dm, 0, 0, 0, &data[fine_level]->elem_restr_u));
1675754ecacSJeremy L Thompson   // -- Energy restriction
168*2b730f8bSJeremy L Thompson   PetscCall(CreateRestrictionFromPlex(ceed, dm_energy, 0, 0, 0, &data[fine_level]->elem_restr_energy));
1695754ecacSJeremy L Thompson   // -- Diagnostic data restriction
170*2b730f8bSJeremy L Thompson   PetscCall(CreateRestrictionFromPlex(ceed, dm_diagnostic, 0, 0, 0, &data[fine_level]->elem_restr_diagnostic));
1715754ecacSJeremy L Thompson 
1725754ecacSJeremy L Thompson   // -- Stored data at quadrature points
173*2b730f8bSJeremy L Thompson   PetscCall(DMPlexGetHeightStratum(dm, 0, &c_start, &c_end));
1745754ecacSJeremy L Thompson   num_elem = c_end - c_start;
1755754ecacSJeremy L Thompson   CeedBasisGetNumQuadraturePoints(data[fine_level]->basis_u, &num_qpts);
1765754ecacSJeremy L Thompson   // ---- Geometric data restriction, residual and Jacobian operators
177*2b730f8bSJeremy L Thompson   CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, q_data_size, num_elem * num_qpts * q_data_size, CEED_STRIDES_BACKEND,
1785754ecacSJeremy L Thompson                                    &data[fine_level]->elem_restr_geo_data_i);
1795754ecacSJeremy L Thompson   // ---- Stored field restrictions
1805754ecacSJeremy L Thompson   for (CeedInt i = 0; i < problem_data.number_fields_stored; i++) {
181*2b730f8bSJeremy L Thompson     CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, problem_data.field_sizes[i], num_elem * num_qpts * problem_data.field_sizes[i],
182*2b730f8bSJeremy L Thompson                                      CEED_STRIDES_BACKEND, &data[fine_level]->elem_restr_stored_fields_i[i]);
1835754ecacSJeremy L Thompson   }
1845754ecacSJeremy L Thompson   // ---- Geometric data restriction, diagnostic operator
185*2b730f8bSJeremy L Thompson   CeedElemRestrictionCreateStrided(ceed, num_elem, P * P * P, q_data_size, num_elem * P * P * P * q_data_size, CEED_STRIDES_BACKEND,
1865754ecacSJeremy L Thompson                                    &data[fine_level]->elem_restr_geo_data_diagnostic_i);
1875754ecacSJeremy L Thompson 
1885754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
1895754ecacSJeremy L Thompson   // Element coordinates
1905754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
191*2b730f8bSJeremy L Thompson   PetscCall(DMGetCoordinatesLocal(dm, &coords));
192*2b730f8bSJeremy L Thompson   PetscCall(VecGetArrayRead(coords, &coordArray));
1935754ecacSJeremy L Thompson 
1945754ecacSJeremy L Thompson   CeedElemRestrictionCreateVector(data[fine_level]->elem_restr_x, &x_coord, NULL);
195*2b730f8bSJeremy L Thompson   CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, (PetscScalar *)coordArray);
196*2b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayRead(coords, &coordArray));
1975754ecacSJeremy L Thompson 
1985754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
1995754ecacSJeremy L Thompson   // Persistent libCEED vectors
2005754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
2015754ecacSJeremy L Thompson   // -- Operator action variables
2025754ecacSJeremy L Thompson   CeedVectorCreate(ceed, U_loc_size, &data[fine_level]->x_ceed);
2035754ecacSJeremy L Thompson   CeedVectorCreate(ceed, U_loc_size, &data[fine_level]->y_ceed);
2045754ecacSJeremy L Thompson   // -- Geometric data vector
205*2b730f8bSJeremy L Thompson   CeedVectorCreate(ceed, num_elem * num_qpts * q_data_size, &data[fine_level]->geo_data);
2065754ecacSJeremy L Thompson   // -- Stored field vectors
2075754ecacSJeremy L Thompson   for (CeedInt i = 0; i < problem_data.number_fields_stored; i++) {
208*2b730f8bSJeremy L Thompson     CeedVectorCreate(ceed, num_elem * num_qpts * problem_data.field_sizes[i], &data[fine_level]->stored_fields[i]);
2095754ecacSJeremy L Thompson   }
2105754ecacSJeremy L Thompson   // -- Collocated geometric data vector
211*2b730f8bSJeremy L Thompson   CeedVectorCreate(ceed, num_elem * P * P * P * q_data_size, &data[fine_level]->geo_data_diagnostic);
2125754ecacSJeremy L Thompson 
2135754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
2145754ecacSJeremy L Thompson   // Geometric factor computation
2155754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
2165754ecacSJeremy L Thompson   // Create the QFunction and Operator that computes the quadrature data
2175754ecacSJeremy L Thompson   //   geo_data returns dXdx_i,j and w * det.
2185754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
2195754ecacSJeremy L Thompson   // -- QFunction
220*2b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem_data.setup_geo, problem_data.setup_geo_loc, &qf_setup_geo);
2215754ecacSJeremy L Thompson   CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x * dim, CEED_EVAL_GRAD);
2225754ecacSJeremy L Thompson   CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT);
223a61c78d6SJeremy L Thompson   CeedQFunctionAddOutput(qf_setup_geo, "qdata", q_data_size, CEED_EVAL_NONE);
2245754ecacSJeremy L Thompson   // -- Operator
225*2b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_setup_geo, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup_geo);
226*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_geo, "dx", data[fine_level]->elem_restr_x, data[fine_level]->basis_x, CEED_VECTOR_ACTIVE);
227*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, data[fine_level]->basis_x, CEED_VECTOR_NONE);
228*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_geo, "qdata", data[fine_level]->elem_restr_geo_data_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
2295754ecacSJeremy L Thompson   // -- Compute the quadrature data
230*2b730f8bSJeremy L Thompson   CeedOperatorApply(op_setup_geo, x_coord, data[fine_level]->geo_data, CEED_REQUEST_IMMEDIATE);
2315754ecacSJeremy L Thompson   // -- Cleanup
2325754ecacSJeremy L Thompson   CeedQFunctionDestroy(&qf_setup_geo);
2335754ecacSJeremy L Thompson   CeedOperatorDestroy(&op_setup_geo);
2345754ecacSJeremy L Thompson 
2355754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
2365754ecacSJeremy L Thompson   // Local residual evaluator
2375754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
2385754ecacSJeremy L Thompson   // Create the QFunction and Operator that computes the residual of the
2395754ecacSJeremy L Thompson   //   non-linear PDE.
2405754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
2415754ecacSJeremy L Thompson   // -- QFunction
242*2b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem_data.residual, problem_data.residual_loc, &qf_residual);
2435754ecacSJeremy L Thompson   CeedQFunctionAddInput(qf_residual, "du", num_comp_u * dim, CEED_EVAL_GRAD);
244a61c78d6SJeremy L Thompson   CeedQFunctionAddInput(qf_residual, "qdata", q_data_size, CEED_EVAL_NONE);
2455754ecacSJeremy L Thompson   CeedQFunctionAddOutput(qf_residual, "dv", num_comp_u * dim, CEED_EVAL_GRAD);
2465754ecacSJeremy L Thompson   for (CeedInt i = 0; i < problem_data.number_fields_stored; i++) {
247*2b730f8bSJeremy L Thompson     CeedQFunctionAddOutput(qf_residual, problem_data.field_names[i], problem_data.field_sizes[i], CEED_EVAL_NONE);
2485754ecacSJeremy L Thompson   }
2495754ecacSJeremy L Thompson   CeedQFunctionSetContext(qf_residual, phys_ctx);
2505754ecacSJeremy L Thompson   // -- Operator
251*2b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_residual, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_residual);
252*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_residual, "du", data[fine_level]->elem_restr_u, data[fine_level]->basis_u, CEED_VECTOR_ACTIVE);
253*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_residual, "qdata", data[fine_level]->elem_restr_geo_data_i, CEED_BASIS_COLLOCATED, data[fine_level]->geo_data);
254*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_residual, "dv", data[fine_level]->elem_restr_u, data[fine_level]->basis_u, CEED_VECTOR_ACTIVE);
2555754ecacSJeremy L Thompson   for (CeedInt i = 0; i < problem_data.number_fields_stored; i++) {
256*2b730f8bSJeremy L Thompson     CeedOperatorSetField(op_residual, problem_data.field_names[i], data[fine_level]->elem_restr_stored_fields_i[i], CEED_BASIS_COLLOCATED,
2575754ecacSJeremy L Thompson                          data[fine_level]->stored_fields[i]);
2585754ecacSJeremy L Thompson   }
2595754ecacSJeremy L Thompson   // -- Save libCEED data
2605754ecacSJeremy L Thompson   data[fine_level]->qf_residual = qf_residual;
2615754ecacSJeremy L Thompson   data[fine_level]->op_residual = op_residual;
2625754ecacSJeremy L Thompson 
2635754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
2645754ecacSJeremy L Thompson   // Jacobian evaluator
2655754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
2665754ecacSJeremy L Thompson   // Create the QFunction and Operator that computes the action of the
2675754ecacSJeremy L Thompson   //   Jacobian for each linear solve.
2685754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
2695754ecacSJeremy L Thompson   // -- QFunction
270*2b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem_data.jacobian, problem_data.jacobian_loc, &qf_jacobian);
2715754ecacSJeremy L Thompson   CeedQFunctionAddInput(qf_jacobian, "delta du", num_comp_u * dim, CEED_EVAL_GRAD);
272a61c78d6SJeremy L Thompson   CeedQFunctionAddInput(qf_jacobian, "qdata", q_data_size, CEED_EVAL_NONE);
2735754ecacSJeremy L Thompson   for (CeedInt i = 0; i < problem_data.number_fields_stored; i++) {
274*2b730f8bSJeremy L Thompson     CeedQFunctionAddInput(qf_jacobian, problem_data.field_names[i], problem_data.field_sizes[i], CEED_EVAL_NONE);
2755754ecacSJeremy L Thompson   }
2765754ecacSJeremy L Thompson   CeedQFunctionAddOutput(qf_jacobian, "delta dv", num_comp_u * dim, CEED_EVAL_GRAD);
2775754ecacSJeremy L Thompson   CeedQFunctionSetContext(qf_jacobian, phys_ctx);
2785754ecacSJeremy L Thompson   // -- Operator
279*2b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_jacobian, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_jacobian);
280*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_jacobian, "delta du", data[fine_level]->elem_restr_u, data[fine_level]->basis_u, CEED_VECTOR_ACTIVE);
281*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_jacobian, "qdata", data[fine_level]->elem_restr_geo_data_i, CEED_BASIS_COLLOCATED, data[fine_level]->geo_data);
282*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_jacobian, "delta dv", data[fine_level]->elem_restr_u, data[fine_level]->basis_u, CEED_VECTOR_ACTIVE);
2835754ecacSJeremy L Thompson   for (CeedInt i = 0; i < problem_data.number_fields_stored; i++) {
284*2b730f8bSJeremy L Thompson     CeedOperatorSetField(op_jacobian, problem_data.field_names[i], data[fine_level]->elem_restr_stored_fields_i[i], CEED_BASIS_COLLOCATED,
2855754ecacSJeremy L Thompson                          data[fine_level]->stored_fields[i]);
2865754ecacSJeremy L Thompson   }
2875754ecacSJeremy L Thompson   // -- Save libCEED data
2885754ecacSJeremy L Thompson   data[fine_level]->qf_jacobian = qf_jacobian;
2895754ecacSJeremy L Thompson   data[fine_level]->op_jacobian = op_jacobian;
2905754ecacSJeremy L Thompson 
2915754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
2925754ecacSJeremy L Thompson   // Traction boundary conditions, if needed
2935754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
2945754ecacSJeremy L Thompson   if (app_ctx->bc_traction_count > 0) {
2955754ecacSJeremy L Thompson     // -- Setup
2965754ecacSJeremy L Thompson     DMLabel domain_label;
297*2b730f8bSJeremy L Thompson     PetscCall(DMGetLabel(dm, "Face Sets", &domain_label));
298*2b730f8bSJeremy L Thompson     PetscCall(DMGetDimension(dm, &dim));
2995754ecacSJeremy L Thompson 
3005754ecacSJeremy L Thompson     // -- Basis
3015754ecacSJeremy L Thompson     CeedInt   height = 1;
3025754ecacSJeremy L Thompson     CeedBasis basis_x_face, basis_u_face;
303*2b730f8bSJeremy L Thompson     CeedBasisCreateTensorH1Lagrange(ceed, dim - height, num_comp_x, 2, Q, problem_data.quadrature_mode, &basis_x_face);
304*2b730f8bSJeremy L Thompson     CeedBasisCreateTensorH1Lagrange(ceed, dim - height, num_comp_u, P, Q, problem_data.quadrature_mode, &basis_u_face);
3055754ecacSJeremy L Thompson     // -- QFunction
3065754ecacSJeremy L Thompson     CeedQFunction        qf_traction;
3075754ecacSJeremy L Thompson     CeedQFunctionContext traction_ctx;
308*2b730f8bSJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, SetupTractionBCs, SetupTractionBCs_loc, &qf_traction);
3095754ecacSJeremy L Thompson     CeedQFunctionContextCreate(ceed, &traction_ctx);
3105754ecacSJeremy L Thompson     CeedQFunctionSetContext(qf_traction, traction_ctx);
311*2b730f8bSJeremy L Thompson     CeedQFunctionAddInput(qf_traction, "dx", num_comp_x * (num_comp_x - height), CEED_EVAL_GRAD);
3125754ecacSJeremy L Thompson     CeedQFunctionAddInput(qf_traction, "weight", 1, CEED_EVAL_WEIGHT);
3135754ecacSJeremy L Thompson     CeedQFunctionAddOutput(qf_traction, "v", num_comp_u, CEED_EVAL_INTERP);
3145754ecacSJeremy L Thompson 
3155754ecacSJeremy L Thompson     // -- Compute contribution on each boundary face
3165754ecacSJeremy L Thompson     for (CeedInt i = 0; i < app_ctx->bc_traction_count; i++) {
3175754ecacSJeremy L Thompson       CeedElemRestriction elem_restr_x_face, elem_restr_u_face;
3185754ecacSJeremy L Thompson       CeedOperator        op_traction;
319*2b730f8bSJeremy L Thompson       CeedQFunctionContextSetData(traction_ctx, CEED_MEM_HOST, CEED_USE_POINTER, 3 * sizeof(CeedScalar), app_ctx->bc_traction_vector[i]);
3205754ecacSJeremy L Thompson       // Setup restriction
321*2b730f8bSJeremy L Thompson       PetscCall(
322*2b730f8bSJeremy L Thompson           GetRestrictionForDomain(ceed, dm, 1, domain_label, app_ctx->bc_traction_faces[i], Q, 0, &elem_restr_u_face, &elem_restr_x_face, NULL));
3235754ecacSJeremy L Thompson       // ---- Create boundary Operator
3245754ecacSJeremy L Thompson       CeedOperatorCreate(ceed, qf_traction, NULL, NULL, &op_traction);
325*2b730f8bSJeremy L Thompson       CeedOperatorSetField(op_traction, "dx", elem_restr_x_face, basis_x_face, CEED_VECTOR_ACTIVE);
326*2b730f8bSJeremy L Thompson       CeedOperatorSetField(op_traction, "weight", CEED_ELEMRESTRICTION_NONE, basis_x_face, CEED_VECTOR_NONE);
327*2b730f8bSJeremy L Thompson       CeedOperatorSetField(op_traction, "v", elem_restr_u_face, basis_u_face, CEED_VECTOR_ACTIVE);
3285754ecacSJeremy L Thompson       // ---- Compute traction on face
329*2b730f8bSJeremy L Thompson       CeedOperatorApplyAdd(op_traction, x_coord, neumann_ceed, CEED_REQUEST_IMMEDIATE);
3305754ecacSJeremy L Thompson       // ---- Cleanup
3315754ecacSJeremy L Thompson       CeedElemRestrictionDestroy(&elem_restr_x_face);
3325754ecacSJeremy L Thompson       CeedElemRestrictionDestroy(&elem_restr_u_face);
3335754ecacSJeremy L Thompson       CeedOperatorDestroy(&op_traction);
3345754ecacSJeremy L Thompson     }
3355754ecacSJeremy L Thompson     // -- Cleanup
3365754ecacSJeremy L Thompson     CeedBasisDestroy(&basis_x_face);
3375754ecacSJeremy L Thompson     CeedBasisDestroy(&basis_u_face);
3385754ecacSJeremy L Thompson     CeedQFunctionDestroy(&qf_traction);
3395754ecacSJeremy L Thompson     CeedQFunctionContextDestroy(&traction_ctx);
3405754ecacSJeremy L Thompson   }
3415754ecacSJeremy L Thompson 
3425754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
3435754ecacSJeremy L Thompson   // Forcing term, if needed
3445754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
3455754ecacSJeremy L Thompson   // Create the QFunction and Operator that computes the forcing term (RHS)
3465754ecacSJeremy L Thompson   //   for the non-linear PDE.
3475754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
3485754ecacSJeremy L Thompson   if (forcing_choice != FORCE_NONE) {
3495754ecacSJeremy L Thompson     CeedQFunction qf_setup_force;
3505754ecacSJeremy L Thompson     CeedOperator  op_setup_force;
3515754ecacSJeremy L Thompson 
3525754ecacSJeremy L Thompson     // -- QFunction
353*2b730f8bSJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, forcing_options[forcing_choice].setup_forcing, forcing_options[forcing_choice].setup_forcing_loc,
3545754ecacSJeremy L Thompson                                 &qf_setup_force);
3555754ecacSJeremy L Thompson     CeedQFunctionAddInput(qf_setup_force, "x", num_comp_x, CEED_EVAL_INTERP);
356*2b730f8bSJeremy L Thompson     CeedQFunctionAddInput(qf_setup_force, "qdata", q_data_size, CEED_EVAL_NONE);
3575754ecacSJeremy L Thompson     CeedQFunctionAddOutput(qf_setup_force, "force", num_comp_u, CEED_EVAL_INTERP);
3585754ecacSJeremy L Thompson     if (forcing_choice == FORCE_MMS) {
3595754ecacSJeremy L Thompson       CeedQFunctionSetContext(qf_setup_force, phys_ctx);
3605754ecacSJeremy L Thompson     } else {
3615754ecacSJeremy L Thompson       CeedQFunctionContext ctxForcing;
3625754ecacSJeremy L Thompson       CeedQFunctionContextCreate(ceed, &ctxForcing);
363*2b730f8bSJeremy L Thompson       CeedQFunctionContextSetData(ctxForcing, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*app_ctx->forcing_vector), app_ctx->forcing_vector);
3645754ecacSJeremy L Thompson       CeedQFunctionSetContext(qf_setup_force, ctxForcing);
3655754ecacSJeremy L Thompson       CeedQFunctionContextDestroy(&ctxForcing);
3665754ecacSJeremy L Thompson     }
3675754ecacSJeremy L Thompson     // -- Operator
368*2b730f8bSJeremy L Thompson     CeedOperatorCreate(ceed, qf_setup_force, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup_force);
369*2b730f8bSJeremy L Thompson     CeedOperatorSetField(op_setup_force, "x", data[fine_level]->elem_restr_x, data[fine_level]->basis_x, CEED_VECTOR_ACTIVE);
370*2b730f8bSJeremy L Thompson     CeedOperatorSetField(op_setup_force, "qdata", data[fine_level]->elem_restr_geo_data_i, CEED_BASIS_COLLOCATED, data[fine_level]->geo_data);
371*2b730f8bSJeremy L Thompson     CeedOperatorSetField(op_setup_force, "force", data[fine_level]->elem_restr_u, data[fine_level]->basis_u, CEED_VECTOR_ACTIVE);
3725754ecacSJeremy L Thompson     // -- Compute forcing term
3735754ecacSJeremy L Thompson     CeedOperatorApply(op_setup_force, x_coord, force_ceed, CEED_REQUEST_IMMEDIATE);
3745754ecacSJeremy L Thompson     // -- Cleanup
3755754ecacSJeremy L Thompson     CeedQFunctionDestroy(&qf_setup_force);
3765754ecacSJeremy L Thompson     CeedOperatorDestroy(&op_setup_force);
3775754ecacSJeremy L Thompson   }
3785754ecacSJeremy L Thompson 
3795754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
3805754ecacSJeremy L Thompson   // True solution, for MMS
3815754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
3825754ecacSJeremy L Thompson   // Create the QFunction and Operator that computes the true solution at
3835754ecacSJeremy L Thompson   //   the mesh nodes for validation with the manufactured solution.
3845754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
3855754ecacSJeremy L Thompson   if (problem_data.true_soln) {
3865754ecacSJeremy L Thompson     CeedScalar       *true_array;
3875754ecacSJeremy L Thompson     const CeedScalar *mult_array;
3885754ecacSJeremy L Thompson     CeedVector        mult_vec;
3895754ecacSJeremy L Thompson     CeedBasis         basis_x_true;
3905754ecacSJeremy L Thompson     CeedQFunction     qf_true;
3915754ecacSJeremy L Thompson     CeedOperator      op_true;
3925754ecacSJeremy L Thompson 
3935754ecacSJeremy L Thompson     // -- Solution vector
3945754ecacSJeremy L Thompson     CeedVectorCreate(ceed, U_loc_size, &(data[fine_level]->true_soln));
3955754ecacSJeremy L Thompson     // -- Basis
396*2b730f8bSJeremy L Thompson     CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, P, CEED_GAUSS_LOBATTO, &basis_x_true);
3975754ecacSJeremy L Thompson     // QFunction
398*2b730f8bSJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, problem_data.true_soln, problem_data.true_soln_loc, &qf_true);
3995754ecacSJeremy L Thompson     CeedQFunctionAddInput(qf_true, "x", num_comp_x, CEED_EVAL_INTERP);
400a61c78d6SJeremy L Thompson     CeedQFunctionAddOutput(qf_true, "true solution", num_comp_u, CEED_EVAL_NONE);
4015754ecacSJeremy L Thompson     // Operator
402*2b730f8bSJeremy L Thompson     CeedOperatorCreate(ceed, qf_true, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_true);
403*2b730f8bSJeremy L Thompson     CeedOperatorSetField(op_true, "x", data[fine_level]->elem_restr_x, basis_x_true, CEED_VECTOR_ACTIVE);
404*2b730f8bSJeremy L Thompson     CeedOperatorSetField(op_true, "true solution", data[fine_level]->elem_restr_u, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
4055754ecacSJeremy L Thompson     // -- Compute true solution
406*2b730f8bSJeremy L Thompson     CeedOperatorApply(op_true, x_coord, data[fine_level]->true_soln, CEED_REQUEST_IMMEDIATE);
4075754ecacSJeremy L Thompson     // -- Multiplicity calculation
408*2b730f8bSJeremy L Thompson     CeedElemRestrictionCreateVector(data[fine_level]->elem_restr_u, &mult_vec, NULL);
4095754ecacSJeremy L Thompson     CeedVectorSetValue(mult_vec, 0.);
4105754ecacSJeremy L Thompson     CeedElemRestrictionGetMultiplicity(data[fine_level]->elem_restr_u, mult_vec);
4115754ecacSJeremy L Thompson     // -- Multiplicity correction
4125754ecacSJeremy L Thompson     CeedVectorGetArray(data[fine_level]->true_soln, CEED_MEM_HOST, &true_array);
4135754ecacSJeremy L Thompson     CeedVectorGetArrayRead(mult_vec, CEED_MEM_HOST, &mult_array);
414*2b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < U_loc_size; i++) true_array[i] /= mult_array[i];
4155754ecacSJeremy L Thompson     CeedVectorRestoreArray(data[fine_level]->true_soln, &true_array);
4165754ecacSJeremy L Thompson     CeedVectorRestoreArrayRead(mult_vec, &mult_array);
4175754ecacSJeremy L Thompson     // -- Cleanup
4185754ecacSJeremy L Thompson     CeedVectorDestroy(&mult_vec);
4195754ecacSJeremy L Thompson     CeedBasisDestroy(&basis_x_true);
4205754ecacSJeremy L Thompson     CeedQFunctionDestroy(&qf_true);
4215754ecacSJeremy L Thompson     CeedOperatorDestroy(&op_true);
4225754ecacSJeremy L Thompson   }
4235754ecacSJeremy L Thompson 
4245754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
4255754ecacSJeremy L Thompson   // Local energy computation
4265754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
4275754ecacSJeremy L Thompson   // Create the QFunction and Operator that computes the strain energy
4285754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
4295754ecacSJeremy L Thompson   // -- QFunction
430*2b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem_data.energy, problem_data.energy_loc, &qf_energy);
4315754ecacSJeremy L Thompson   CeedQFunctionAddInput(qf_energy, "du", num_comp_u * dim, CEED_EVAL_GRAD);
432a61c78d6SJeremy L Thompson   CeedQFunctionAddInput(qf_energy, "qdata", q_data_size, CEED_EVAL_NONE);
4335754ecacSJeremy L Thompson   CeedQFunctionAddOutput(qf_energy, "energy", num_comp_e, CEED_EVAL_INTERP);
4345754ecacSJeremy L Thompson   CeedQFunctionSetContext(qf_energy, phys_ctx);
4355754ecacSJeremy L Thompson   // -- Operator
436*2b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_energy, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_energy);
437*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_energy, "du", data[fine_level]->elem_restr_u, data[fine_level]->basis_u, CEED_VECTOR_ACTIVE);
438*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_energy, "qdata", data[fine_level]->elem_restr_geo_data_i, CEED_BASIS_COLLOCATED, data[fine_level]->geo_data);
439*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_energy, "energy", data[fine_level]->elem_restr_energy, data[fine_level]->basis_energy, CEED_VECTOR_ACTIVE);
4405754ecacSJeremy L Thompson   // -- Save libCEED data
4415754ecacSJeremy L Thompson   data[fine_level]->qf_energy = qf_energy;
4425754ecacSJeremy L Thompson   data[fine_level]->op_energy = op_energy;
4435754ecacSJeremy L Thompson 
4445754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
4455754ecacSJeremy L Thompson   // Diagnostic value computation
4465754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
4475754ecacSJeremy L Thompson   // Create the QFunction and Operator that computes nodal diagnostic quantities
4485754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
4495754ecacSJeremy L Thompson   // Geometric factors
4505754ecacSJeremy L Thompson   // -- Coordinate basis
4515754ecacSJeremy L Thompson   CeedBasis basis_x;
452*2b730f8bSJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS_LOBATTO, &basis_x);
4535754ecacSJeremy L Thompson   // -- QFunction
454*2b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem_data.setup_geo, problem_data.setup_geo_loc, &qf_setup_geo);
4555754ecacSJeremy L Thompson   CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x * dim, CEED_EVAL_GRAD);
4565754ecacSJeremy L Thompson   CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT);
457a61c78d6SJeremy L Thompson   CeedQFunctionAddOutput(qf_setup_geo, "qdata", q_data_size, CEED_EVAL_NONE);
4585754ecacSJeremy L Thompson   // -- Operator
459*2b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_setup_geo, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup_geo);
460*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_geo, "dx", data[fine_level]->elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
461*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE);
462*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_geo, "qdata", data[fine_level]->elem_restr_geo_data_diagnostic_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
4635754ecacSJeremy L Thompson   // -- Compute the quadrature data
464*2b730f8bSJeremy L Thompson   CeedOperatorApply(op_setup_geo, x_coord, data[fine_level]->geo_data_diagnostic, CEED_REQUEST_IMMEDIATE);
4655754ecacSJeremy L Thompson   // -- Cleanup
4665754ecacSJeremy L Thompson   CeedBasisDestroy(&basis_x);
4675754ecacSJeremy L Thompson   CeedQFunctionDestroy(&qf_setup_geo);
4685754ecacSJeremy L Thompson   CeedOperatorDestroy(&op_setup_geo);
4695754ecacSJeremy L Thompson 
4705754ecacSJeremy L Thompson   // Diagnostic quantities
4715754ecacSJeremy L Thompson   // -- QFunction
472*2b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem_data.diagnostic, problem_data.diagnostic_loc, &qf_diagnostic);
4735754ecacSJeremy L Thompson   CeedQFunctionAddInput(qf_diagnostic, "u", num_comp_u, CEED_EVAL_INTERP);
4745754ecacSJeremy L Thompson   CeedQFunctionAddInput(qf_diagnostic, "du", num_comp_u * dim, CEED_EVAL_GRAD);
475a61c78d6SJeremy L Thompson   CeedQFunctionAddInput(qf_diagnostic, "qdata", q_data_size, CEED_EVAL_NONE);
476*2b730f8bSJeremy L Thompson   CeedQFunctionAddOutput(qf_diagnostic, "diagnostic values", num_comp_u + num_comp_d, CEED_EVAL_NONE);
4775754ecacSJeremy L Thompson   CeedQFunctionSetContext(qf_diagnostic, phys_ctx);
4785754ecacSJeremy L Thompson   // -- Operator
479*2b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, qf_diagnostic, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_diagnostic);
480*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_diagnostic, "u", data[fine_level]->elem_restr_u, data[fine_level]->basis_diagnostic, CEED_VECTOR_ACTIVE);
481*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_diagnostic, "du", data[fine_level]->elem_restr_u, data[fine_level]->basis_diagnostic, CEED_VECTOR_ACTIVE);
482*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_diagnostic, "qdata", data[fine_level]->elem_restr_geo_data_diagnostic_i, CEED_BASIS_COLLOCATED,
483*2b730f8bSJeremy L Thompson                        data[fine_level]->geo_data_diagnostic);
484*2b730f8bSJeremy L Thompson   CeedOperatorSetField(op_diagnostic, "diagnostic values", data[fine_level]->elem_restr_diagnostic, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
4855754ecacSJeremy L Thompson   // -- Save libCEED data
4865754ecacSJeremy L Thompson   data[fine_level]->qf_diagnostic = qf_diagnostic;
4875754ecacSJeremy L Thompson   data[fine_level]->op_diagnostic = op_diagnostic;
4885754ecacSJeremy L Thompson 
4895754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
4905754ecacSJeremy L Thompson   // Cleanup
4915754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
4925754ecacSJeremy L Thompson   CeedVectorDestroy(&x_coord);
4935754ecacSJeremy L Thompson 
4945754ecacSJeremy L Thompson   PetscFunctionReturn(0);
4955754ecacSJeremy L Thompson };
4965754ecacSJeremy L Thompson 
4975754ecacSJeremy L Thompson // Set up libCEED multigrid level for a given degree
4985754ecacSJeremy L Thompson //   Prolongation and Restriction are between level and level+1
499*2b730f8bSJeremy L Thompson PetscErrorCode SetupLibceedLevel(DM dm, Ceed ceed, AppCtx app_ctx, ProblemData problem_data, PetscInt level, PetscInt num_comp_u, PetscInt U_g_size,
500*2b730f8bSJeremy L Thompson                                  PetscInt U_loc_size, CeedVector fine_mult, CeedData *data) {
5015754ecacSJeremy L Thompson   CeedInt      fine_level = app_ctx->num_levels - 1;
5025754ecacSJeremy L Thompson   CeedInt      P          = app_ctx->level_degrees[level] + 1;
5035754ecacSJeremy L Thompson   CeedInt      Q          = app_ctx->level_degrees[fine_level] + 1 + app_ctx->q_extra;
5045754ecacSJeremy L Thompson   CeedInt      dim;
5055754ecacSJeremy L Thompson   CeedOperator op_jacobian, op_prolong, op_restrict;
5065754ecacSJeremy L Thompson 
5075754ecacSJeremy L Thompson   PetscFunctionBeginUser;
5085754ecacSJeremy L Thompson 
509*2b730f8bSJeremy L Thompson   PetscCall(DMGetDimension(dm, &dim));
5105754ecacSJeremy L Thompson 
5115754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
5125754ecacSJeremy L Thompson   // libCEED restrictions
5135754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
5145754ecacSJeremy L Thompson   // -- Solution restriction
515*2b730f8bSJeremy L Thompson   PetscCall(CreateRestrictionFromPlex(ceed, dm, 0, 0, 0, &data[level]->elem_restr_u));
5165754ecacSJeremy L Thompson 
5175754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
5185754ecacSJeremy L Thompson   // libCEED bases
5195754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
5205754ecacSJeremy L Thompson   // -- Solution basis
521*2b730f8bSJeremy L Thompson   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_u, P, Q, problem_data.quadrature_mode, &data[level]->basis_u);
5225754ecacSJeremy L Thompson 
5235754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
5245754ecacSJeremy L Thompson   // Persistent libCEED vectors
5255754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
5265754ecacSJeremy L Thompson   CeedVectorCreate(ceed, U_loc_size, &data[level]->x_ceed);
5275754ecacSJeremy L Thompson   CeedVectorCreate(ceed, U_loc_size, &data[level]->y_ceed);
5285754ecacSJeremy L Thompson 
5295754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
5305754ecacSJeremy L Thompson   // Coarse Grid, Prolongation, and Restriction Operators
5315754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
5325754ecacSJeremy L Thompson   // Create the Operators that compute the prolongation and
5335754ecacSJeremy L Thompson   //   restriction between the p-multigrid levels and the coarse grid eval.
5345754ecacSJeremy L Thompson   // ---------------------------------------------------------------------------
535*2b730f8bSJeremy L Thompson   CeedOperatorMultigridLevelCreate(data[level + 1]->op_jacobian, fine_mult, data[level]->elem_restr_u, data[level]->basis_u, &op_jacobian,
536*2b730f8bSJeremy L Thompson                                    &op_prolong, &op_restrict);
5375754ecacSJeremy L Thompson 
5385754ecacSJeremy L Thompson   // -- Save libCEED data
5395754ecacSJeremy L Thompson   data[level]->op_jacobian     = op_jacobian;
5405754ecacSJeremy L Thompson   data[level + 1]->op_prolong  = op_prolong;
5415754ecacSJeremy L Thompson   data[level + 1]->op_restrict = op_restrict;
5425754ecacSJeremy L Thompson 
5435754ecacSJeremy L Thompson   PetscFunctionReturn(0);
5445754ecacSJeremy L Thompson };
545