164dd23feSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 264dd23feSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 364dd23feSJames Wright 464dd23feSJames Wright #include <ceed.h> 564dd23feSJames Wright #include <petsc-ceed.h> 664dd23feSJames Wright #include <petsc.h> 764dd23feSJames Wright #include "../qfunctions/mass.h" 864dd23feSJames Wright 964dd23feSJames Wright /** 1064dd23feSJames Wright @brief Create mass CeedQFunction for number of components `N` 1164dd23feSJames Wright 1264dd23feSJames Wright Output QFunction has two inputs ("u", and "qdata") and one output ("v"). 1364dd23feSJames Wright "qdata" is assumed to have `wdetJ` as it's first component; all other components are ignored. 1464dd23feSJames Wright 1564dd23feSJames Wright @param[in] ceed Ceed object 1664dd23feSJames Wright @param[in] N Number of components 1764dd23feSJames Wright @param[in] q_data_size Number of components of `CeedVector` holding wdetJ 1864dd23feSJames Wright @param[out] qf CeedQFunction of mass operator 1964dd23feSJames Wright **/ 2064dd23feSJames Wright PetscErrorCode HoneeMassQFunctionCreate(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) { 2164dd23feSJames Wright PetscFunctionBeginUser; 2264dd23feSJames Wright switch (N) { 2364dd23feSJames Wright case 1: 2464dd23feSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf)); 2564dd23feSJames Wright break; 2664dd23feSJames Wright case 2: 27*bb4806caSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_2, Mass_2_loc, qf)); 2864dd23feSJames Wright break; 2964dd23feSJames Wright case 3: 30*bb4806caSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_3, Mass_3_loc, qf)); 3164dd23feSJames Wright break; 3264dd23feSJames Wright case 4: 33*bb4806caSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_4, Mass_4_loc, qf)); 3464dd23feSJames Wright break; 3564dd23feSJames Wright case 5: 3664dd23feSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf)); 3764dd23feSJames Wright break; 3864dd23feSJames Wright case 7: 3964dd23feSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf)); 4064dd23feSJames Wright break; 4164dd23feSJames Wright case 9: 4264dd23feSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf)); 4364dd23feSJames Wright break; 4464dd23feSJames Wright case 12: 4564dd23feSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_12, Mass_12_loc, qf)); 4664dd23feSJames Wright break; 4764dd23feSJames Wright case 22: 4864dd23feSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf)); 4964dd23feSJames Wright break; 5064dd23feSJames Wright default: 5164dd23feSJames Wright SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, 5264dd23feSJames Wright "Mass qfunction of size %" CeedInt_FMT " does not exist.\nA new mass qfunction can be easily added; see source code for pattern.", N); 5364dd23feSJames Wright } 5464dd23feSJames Wright 5564dd23feSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP)); 5664dd23feSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE)); 5764dd23feSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP)); 5864dd23feSJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf, N)); 5964dd23feSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 6064dd23feSJames Wright } 61