xref: /libCEED/examples/solids/problems/linear.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 #include "../qfunctions/linear.h"
9*2b730f8bSJeremy L Thompson 
10*2b730f8bSJeremy L Thompson #include <ceed.h>
11*2b730f8bSJeremy L Thompson 
12*2b730f8bSJeremy L Thompson #include "../include/setup-libceed.h"
13*2b730f8bSJeremy L Thompson #include "../include/structs.h"
14*2b730f8bSJeremy L Thompson #include "../problems/neo-hookean.h"
15*2b730f8bSJeremy L Thompson #include "../problems/problems.h"
16*2b730f8bSJeremy L Thompson #include "../qfunctions/common.h"
175754ecacSJeremy L Thompson #include "../qfunctions/manufactured-true.h"
185754ecacSJeremy L Thompson 
195754ecacSJeremy L Thompson ProblemData linear_elasticity = {
205754ecacSJeremy L Thompson     .setup_geo            = SetupGeo,
215754ecacSJeremy L Thompson     .setup_geo_loc        = SetupGeo_loc,
22a61c78d6SJeremy L Thompson     .q_data_size          = 10,
235754ecacSJeremy L Thompson     .quadrature_mode      = CEED_GAUSS,
245754ecacSJeremy L Thompson     .residual             = ElasLinearF,
255754ecacSJeremy L Thompson     .residual_loc         = ElasLinearF_loc,
265754ecacSJeremy L Thompson     .number_fields_stored = 0,
275754ecacSJeremy L Thompson     .jacobian             = ElasLineardF,
285754ecacSJeremy L Thompson     .jacobian_loc         = ElasLineardF_loc,
295754ecacSJeremy L Thompson     .energy               = ElasLinearEnergy,
305754ecacSJeremy L Thompson     .energy_loc           = ElasLinearEnergy_loc,
315754ecacSJeremy L Thompson     .diagnostic           = ElasLinearDiagnostic,
325754ecacSJeremy L Thompson     .diagnostic_loc       = ElasLinearDiagnostic_loc,
335754ecacSJeremy L Thompson     .true_soln            = MMSTrueSoln,
345754ecacSJeremy L Thompson     .true_soln_loc        = MMSTrueSoln_loc,
355754ecacSJeremy L Thompson };
365754ecacSJeremy L Thompson 
37*2b730f8bSJeremy L Thompson PetscErrorCode SetupLibceedFineLevel_ElasLinear(DM dm, DM dm_energy, DM dm_diagnostic, Ceed ceed, AppCtx app_ctx, CeedQFunctionContext phys_ctx,
38*2b730f8bSJeremy L Thompson                                                 PetscInt fine_level, PetscInt num_comp_u, PetscInt U_g_size, PetscInt U_loc_size,
39*2b730f8bSJeremy L Thompson                                                 CeedVector force_ceed, CeedVector neumann_ceed, CeedData *data) {
405754ecacSJeremy L Thompson   PetscFunctionBegin;
415754ecacSJeremy L Thompson 
42*2b730f8bSJeremy L Thompson   PetscCall(SetupLibceedFineLevel(dm, dm_energy, dm_diagnostic, ceed, app_ctx, phys_ctx, linear_elasticity, fine_level, num_comp_u, U_g_size,
43*2b730f8bSJeremy L Thompson                                   U_loc_size, force_ceed, neumann_ceed, data));
445754ecacSJeremy L Thompson 
455754ecacSJeremy L Thompson   PetscFunctionReturn(0);
465754ecacSJeremy L Thompson };
475754ecacSJeremy L Thompson 
48*2b730f8bSJeremy L Thompson PetscErrorCode SetupLibceedLevel_ElasLinear(DM dm, Ceed ceed, AppCtx app_ctx, PetscInt level, PetscInt num_comp_u, PetscInt U_g_size,
49*2b730f8bSJeremy L Thompson                                             PetscInt U_loc_size, CeedVector fine_mult, CeedData *data) {
505754ecacSJeremy L Thompson   PetscFunctionBegin;
515754ecacSJeremy L Thompson 
52*2b730f8bSJeremy L Thompson   PetscCall(SetupLibceedLevel(dm, ceed, app_ctx, linear_elasticity, level, num_comp_u, U_g_size, U_loc_size, fine_mult, data));
535754ecacSJeremy L Thompson 
545754ecacSJeremy L Thompson   PetscFunctionReturn(0);
555754ecacSJeremy L Thompson };
56