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