xref: /honee/src/honee-mass.c (revision 64dd23fed3bf7ffd926f9c8ee5c8cd1891ca06e3)
1*64dd23feSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2*64dd23feSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3*64dd23feSJames Wright 
4*64dd23feSJames Wright #include <ceed.h>
5*64dd23feSJames Wright #include <petsc-ceed.h>
6*64dd23feSJames Wright #include <petsc.h>
7*64dd23feSJames Wright #include "../qfunctions/mass.h"
8*64dd23feSJames Wright 
9*64dd23feSJames Wright /**
10*64dd23feSJames Wright   @brief Create mass CeedQFunction for number of components `N`
11*64dd23feSJames Wright 
12*64dd23feSJames Wright   Output QFunction has two inputs ("u", and "qdata") and one output ("v").
13*64dd23feSJames Wright   "qdata" is assumed to have `wdetJ` as it's first component; all other components are ignored.
14*64dd23feSJames Wright 
15*64dd23feSJames Wright   @param[in]  ceed        Ceed object
16*64dd23feSJames Wright   @param[in]  N           Number of components
17*64dd23feSJames Wright   @param[in]  q_data_size Number of components of `CeedVector` holding wdetJ
18*64dd23feSJames Wright   @param[out] qf          CeedQFunction of mass operator
19*64dd23feSJames Wright **/
20*64dd23feSJames Wright PetscErrorCode HoneeMassQFunctionCreate(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
21*64dd23feSJames Wright   PetscFunctionBeginUser;
22*64dd23feSJames Wright   switch (N) {
23*64dd23feSJames Wright     case 1:
24*64dd23feSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf));
25*64dd23feSJames Wright       break;
26*64dd23feSJames Wright     case 2:
27*64dd23feSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 2, Mass_2, Mass_2_loc, qf));
28*64dd23feSJames Wright       break;
29*64dd23feSJames Wright     case 3:
30*64dd23feSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 3, Mass_3, Mass_3_loc, qf));
31*64dd23feSJames Wright       break;
32*64dd23feSJames Wright     case 4:
33*64dd23feSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 4, Mass_4, Mass_4_loc, qf));
34*64dd23feSJames Wright       break;
35*64dd23feSJames Wright     case 5:
36*64dd23feSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf));
37*64dd23feSJames Wright       break;
38*64dd23feSJames Wright     case 7:
39*64dd23feSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf));
40*64dd23feSJames Wright       break;
41*64dd23feSJames Wright     case 9:
42*64dd23feSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf));
43*64dd23feSJames Wright       break;
44*64dd23feSJames Wright     case 12:
45*64dd23feSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_12, Mass_12_loc, qf));
46*64dd23feSJames Wright       break;
47*64dd23feSJames Wright     case 22:
48*64dd23feSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf));
49*64dd23feSJames Wright       break;
50*64dd23feSJames Wright     default:
51*64dd23feSJames Wright       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP,
52*64dd23feSJames Wright               "Mass qfunction of size %" CeedInt_FMT " does not exist.\nA new mass qfunction can be easily added; see source code for pattern.", N);
53*64dd23feSJames Wright   }
54*64dd23feSJames Wright 
55*64dd23feSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP));
56*64dd23feSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE));
57*64dd23feSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP));
58*64dd23feSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf, N));
59*64dd23feSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
60*64dd23feSJames Wright }
61