Line data Source code
1 : // Copyright Amazon.com, Inc. or its affiliates. All Rights Reserved.
2 : // SPDX-License-Identifier: Apache-2.0
3 :
4 : #include "fem/integrator.hpp"
5 :
6 : #include "fem/libceed/coefficient.hpp"
7 : #include "fem/libceed/integrator.hpp"
8 : #include "utils/diagnostic.hpp"
9 :
10 : PalacePragmaDiagnosticPush
11 : PalacePragmaDiagnosticDisableUnused
12 :
13 : #include "fem/qfunctions/hdiv_qf.h"
14 : #include "fem/qfunctions/l2_qf.h"
15 :
16 : PalacePragmaDiagnosticPop
17 :
18 : namespace palace
19 : {
20 :
21 : using namespace ceed;
22 :
23 995 : void CurlCurlIntegrator::Assemble(Ceed ceed, CeedElemRestriction trial_restr,
24 : CeedElemRestriction test_restr, CeedBasis trial_basis,
25 : CeedBasis test_basis, CeedVector geom_data,
26 : CeedElemRestriction geom_data_restr,
27 : CeedOperator *op) const
28 : {
29 : CeedQFunctionInfo info;
30 995 : info.assemble_q_data = assemble_q_data;
31 :
32 : // Set up QFunctions.
33 : CeedInt dim, space_dim, trial_num_comp, test_num_comp;
34 995 : PalaceCeedCall(ceed, CeedBasisGetDimension(trial_basis, &dim));
35 995 : PalaceCeedCall(ceed, CeedGeometryDataGetSpaceDimension(geom_data_restr, dim, &space_dim));
36 995 : PalaceCeedCall(ceed, CeedBasisGetNumComponents(trial_basis, &trial_num_comp));
37 995 : PalaceCeedCall(ceed, CeedBasisGetNumComponents(test_basis, &test_num_comp));
38 995 : MFEM_VERIFY(trial_num_comp == test_num_comp && trial_num_comp == 1,
39 : "CurlCurlIntegrator requires test and trial spaces with a single component!");
40 995 : switch (10 * space_dim + dim)
41 : {
42 288 : case 22:
43 : // Curl in 2D has a single component.
44 288 : info.apply_qf = assemble_q_data ? f_build_l2_1 : f_apply_l2_1;
45 288 : info.apply_qf_path = PalaceQFunctionRelativePath(assemble_q_data ? f_build_l2_1_loc
46 : : f_apply_l2_1_loc);
47 : break;
48 420 : case 33:
49 420 : info.apply_qf = assemble_q_data ? f_build_hdiv_33 : f_apply_hdiv_33;
50 420 : info.apply_qf_path = PalaceQFunctionRelativePath(
51 : assemble_q_data ? f_build_hdiv_33_loc : f_apply_hdiv_33_loc);
52 : break;
53 287 : case 32:
54 : // Curl in 2D has a single component.
55 287 : info.apply_qf = assemble_q_data ? f_build_l2_1 : f_apply_l2_1;
56 287 : info.apply_qf_path = PalaceQFunctionRelativePath(assemble_q_data ? f_build_l2_1_loc
57 : : f_apply_l2_1_loc);
58 : break;
59 0 : default:
60 0 : MFEM_ABORT("Invalid value of (dim, space_dim) = (" << dim << ", " << space_dim
61 : << ") for CurlCurlIntegrator!");
62 : }
63 995 : info.trial_ops = EvalMode::Curl;
64 995 : info.test_ops = EvalMode::Curl;
65 995 : if (dim < 3)
66 : {
67 575 : info.trial_ops |= EvalMode::Weight;
68 : }
69 :
70 : // Set up the coefficient and assemble.
71 1415 : auto ctx = PopulateCoefficientContext((dim < 3) ? 1 : dim, Q, transpose);
72 995 : AssembleCeedOperator(info, (void *)ctx.data(), ctx.size() * sizeof(CeedIntScalar), ceed,
73 : trial_restr, test_restr, trial_basis, test_basis, geom_data,
74 : geom_data_restr, op);
75 995 : }
76 :
77 : } // namespace palace
|