xref: /libCEED/tests/t596-operator.c (revision debdbcfd6d1bf75f737a96117abce91cc9b450a2)
1*1b95d8c6SJeremy L Thompson /// @file
2*1b95d8c6SJeremy L Thompson /// Test full assembly of mass matrix operator
3*1b95d8c6SJeremy L Thompson /// \test Test full assembly of mass matrix operator AtPoints
4*1b95d8c6SJeremy L Thompson #include <ceed.h>
5*1b95d8c6SJeremy L Thompson #include <math.h>
6*1b95d8c6SJeremy L Thompson #include <stdio.h>
7*1b95d8c6SJeremy L Thompson #include <stdlib.h>
8*1b95d8c6SJeremy L Thompson 
9*1b95d8c6SJeremy L Thompson #include "t596-operator.h"
10*1b95d8c6SJeremy L Thompson 
main(int argc,char ** argv)11*1b95d8c6SJeremy L Thompson int main(int argc, char **argv) {
12*1b95d8c6SJeremy L Thompson   Ceed ceed;
13*1b95d8c6SJeremy L Thompson 
14*1b95d8c6SJeremy L Thompson   CeedInit(argv[1], &ceed);
15*1b95d8c6SJeremy L Thompson 
16*1b95d8c6SJeremy L Thompson   for (CeedInt num_comp = 1; num_comp <= 3; num_comp++) {
17*1b95d8c6SJeremy L Thompson     CeedElemRestriction elem_restriction_x, elem_restriction_x_points, elem_restriction_u, elem_restriction_q_data;
18*1b95d8c6SJeremy L Thompson     CeedBasis           basis_x, basis_u;
19*1b95d8c6SJeremy L Thompson     CeedQFunction       qf_setup, qf_mass;
20*1b95d8c6SJeremy L Thompson     CeedOperator        op_setup, op_mass;
21*1b95d8c6SJeremy L Thompson     CeedVector          q_data, x, x_points, u, v;
22*1b95d8c6SJeremy L Thompson     CeedInt             p = 3, q = 4, dim = 2;
23*1b95d8c6SJeremy L Thompson     CeedInt             n_x = 3, n_y = 2;
24*1b95d8c6SJeremy L Thompson     CeedInt             num_elem = n_x * n_y;
25*1b95d8c6SJeremy L Thompson     CeedInt             num_dofs = (n_x * 2 + 1) * (n_y * 2 + 1), num_points_per_elem = 4, num_points = num_elem * num_points_per_elem;
26*1b95d8c6SJeremy L Thompson     CeedInt             ind_x[num_elem * p * p];
27*1b95d8c6SJeremy L Thompson     CeedScalar          assembled_values[num_comp * num_comp * num_dofs * num_dofs];
28*1b95d8c6SJeremy L Thompson     CeedScalar          assembled_true[num_comp * num_comp * num_dofs * num_dofs];
29*1b95d8c6SJeremy L Thompson 
30*1b95d8c6SJeremy L Thompson     // Points
31*1b95d8c6SJeremy L Thompson     CeedVectorCreate(ceed, dim * num_points, &x_points);
32*1b95d8c6SJeremy L Thompson     {
33*1b95d8c6SJeremy L Thompson       CeedScalar x_array[dim * num_points];
34*1b95d8c6SJeremy L Thompson 
35*1b95d8c6SJeremy L Thompson       for (CeedInt e = 0; e < num_elem; e++) {
36*1b95d8c6SJeremy L Thompson         for (CeedInt d = 0; d < dim; d++) {
37*1b95d8c6SJeremy L Thompson           x_array[num_points_per_elem * (e * dim + d) + 0] = 0.25;
38*1b95d8c6SJeremy L Thompson           x_array[num_points_per_elem * (e * dim + d) + 1] = d == 0 ? -0.25 : 0.25;
39*1b95d8c6SJeremy L Thompson           x_array[num_points_per_elem * (e * dim + d) + 2] = d == 0 ? 0.25 : -0.25;
40*1b95d8c6SJeremy L Thompson           x_array[num_points_per_elem * (e * dim + d) + 3] = 0.25;
41*1b95d8c6SJeremy L Thompson         }
42*1b95d8c6SJeremy L Thompson       }
43*1b95d8c6SJeremy L Thompson       CeedVectorSetArray(x_points, CEED_MEM_HOST, CEED_COPY_VALUES, x_array);
44*1b95d8c6SJeremy L Thompson     }
45*1b95d8c6SJeremy L Thompson     {
46*1b95d8c6SJeremy L Thompson       CeedInt ind_x[num_elem + 1 + num_points];
47*1b95d8c6SJeremy L Thompson 
48*1b95d8c6SJeremy L Thompson       for (CeedInt i = 0; i <= num_elem; i++) ind_x[i] = num_elem + 1 + i * num_points_per_elem;
49*1b95d8c6SJeremy L Thompson       for (CeedInt i = 0; i < num_points; i++) ind_x[num_elem + 1 + i] = i;
50*1b95d8c6SJeremy L Thompson       CeedElemRestrictionCreateAtPoints(ceed, num_elem, num_points, dim, num_points * dim, CEED_MEM_HOST, CEED_COPY_VALUES, ind_x,
51*1b95d8c6SJeremy L Thompson                                         &elem_restriction_x_points);
52*1b95d8c6SJeremy L Thompson       CeedElemRestrictionCreateAtPoints(ceed, num_elem, num_points, 1, num_points, CEED_MEM_HOST, CEED_COPY_VALUES, ind_x, &elem_restriction_q_data);
53*1b95d8c6SJeremy L Thompson     }
54*1b95d8c6SJeremy L Thompson 
55*1b95d8c6SJeremy L Thompson     // Vectors
56*1b95d8c6SJeremy L Thompson     CeedVectorCreate(ceed, dim * num_dofs, &x);
57*1b95d8c6SJeremy L Thompson     {
58*1b95d8c6SJeremy L Thompson       CeedScalar x_array[dim * num_dofs];
59*1b95d8c6SJeremy L Thompson 
60*1b95d8c6SJeremy L Thompson       for (CeedInt i = 0; i < n_x * 2 + 1; i++) {
61*1b95d8c6SJeremy L Thompson         for (CeedInt j = 0; j < n_y * 2 + 1; j++) {
62*1b95d8c6SJeremy L Thompson           x_array[i + j * (n_x * 2 + 1) + 0 * num_dofs] = (CeedScalar)i / (2 * n_x);
63*1b95d8c6SJeremy L Thompson           x_array[i + j * (n_x * 2 + 1) + 1 * num_dofs] = (CeedScalar)j / (2 * n_y);
64*1b95d8c6SJeremy L Thompson         }
65*1b95d8c6SJeremy L Thompson       }
66*1b95d8c6SJeremy L Thompson       CeedVectorSetArray(x, CEED_MEM_HOST, CEED_COPY_VALUES, x_array);
67*1b95d8c6SJeremy L Thompson     }
68*1b95d8c6SJeremy L Thompson     CeedVectorCreate(ceed, num_comp * num_dofs, &u);
69*1b95d8c6SJeremy L Thompson     CeedVectorCreate(ceed, num_comp * num_dofs, &v);
70*1b95d8c6SJeremy L Thompson     CeedVectorCreate(ceed, num_points, &q_data);
71*1b95d8c6SJeremy L Thompson 
72*1b95d8c6SJeremy L Thompson     // Restrictions
73*1b95d8c6SJeremy L Thompson     for (CeedInt i = 0; i < num_elem; i++) {
74*1b95d8c6SJeremy L Thompson       CeedInt col, row, offset;
75*1b95d8c6SJeremy L Thompson       col    = i % n_x;
76*1b95d8c6SJeremy L Thompson       row    = i / n_x;
77*1b95d8c6SJeremy L Thompson       offset = col * (p - 1) + row * (n_x * 2 + 1) * (p - 1);
78*1b95d8c6SJeremy L Thompson       for (CeedInt j = 0; j < p; j++) {
79*1b95d8c6SJeremy L Thompson         for (CeedInt k = 0; k < p; k++) ind_x[p * (p * i + k) + j] = offset + k * (n_x * 2 + 1) + j;
80*1b95d8c6SJeremy L Thompson       }
81*1b95d8c6SJeremy L Thompson     }
82*1b95d8c6SJeremy L Thompson     CeedElemRestrictionCreate(ceed, num_elem, p * p, dim, num_dofs, dim * num_dofs, CEED_MEM_HOST, CEED_USE_POINTER, ind_x, &elem_restriction_x);
83*1b95d8c6SJeremy L Thompson     CeedElemRestrictionCreate(ceed, num_elem, p * p, num_comp, num_dofs, num_comp * num_dofs, CEED_MEM_HOST, CEED_USE_POINTER, ind_x,
84*1b95d8c6SJeremy L Thompson                               &elem_restriction_u);
85*1b95d8c6SJeremy L Thompson 
86*1b95d8c6SJeremy L Thompson     // Bases
87*1b95d8c6SJeremy L Thompson     CeedBasisCreateTensorH1Lagrange(ceed, dim, dim, p, q, CEED_GAUSS, &basis_x);
88*1b95d8c6SJeremy L Thompson     CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, p, q, CEED_GAUSS, &basis_u);
89*1b95d8c6SJeremy L Thompson 
90*1b95d8c6SJeremy L Thompson     // QFunctions
91*1b95d8c6SJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, setup, setup_loc, &qf_setup);
92*1b95d8c6SJeremy L Thompson     CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT);
93*1b95d8c6SJeremy L Thompson     CeedQFunctionAddInput(qf_setup, "dx", dim * dim, CEED_EVAL_GRAD);
94*1b95d8c6SJeremy L Thompson     CeedQFunctionAddOutput(qf_setup, "rho", 1, CEED_EVAL_NONE);
95*1b95d8c6SJeremy L Thompson 
96*1b95d8c6SJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, mass, mass_loc, &qf_mass);
97*1b95d8c6SJeremy L Thompson     CeedQFunctionAddInput(qf_mass, "rho", 1, CEED_EVAL_NONE);
98*1b95d8c6SJeremy L Thompson     CeedQFunctionAddInput(qf_mass, "u", num_comp, CEED_EVAL_INTERP);
99*1b95d8c6SJeremy L Thompson     CeedQFunctionAddOutput(qf_mass, "v", num_comp, CEED_EVAL_INTERP);
100*1b95d8c6SJeremy L Thompson     {
101*1b95d8c6SJeremy L Thompson       CeedQFunctionContext qf_context;
102*1b95d8c6SJeremy L Thompson 
103*1b95d8c6SJeremy L Thompson       CeedQFunctionContextCreate(ceed, &qf_context);
104*1b95d8c6SJeremy L Thompson       CeedQFunctionContextSetData(qf_context, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof(CeedInt), &num_comp);
105*1b95d8c6SJeremy L Thompson       CeedQFunctionSetContext(qf_mass, qf_context);
106*1b95d8c6SJeremy L Thompson       CeedQFunctionContextDestroy(&qf_context);
107*1b95d8c6SJeremy L Thompson     }
108*1b95d8c6SJeremy L Thompson 
109*1b95d8c6SJeremy L Thompson     // Operators
110*1b95d8c6SJeremy L Thompson     CeedOperatorCreateAtPoints(ceed, qf_setup, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup);
111*1b95d8c6SJeremy L Thompson     CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE);
112*1b95d8c6SJeremy L Thompson     CeedOperatorSetField(op_setup, "dx", elem_restriction_x, basis_x, CEED_VECTOR_ACTIVE);
113*1b95d8c6SJeremy L Thompson     CeedOperatorSetField(op_setup, "rho", elem_restriction_q_data, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
114*1b95d8c6SJeremy L Thompson     CeedOperatorAtPointsSetPoints(op_setup, elem_restriction_x_points, x_points);
115*1b95d8c6SJeremy L Thompson 
116*1b95d8c6SJeremy L Thompson     CeedOperatorCreateAtPoints(ceed, qf_mass, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_mass);
117*1b95d8c6SJeremy L Thompson     CeedOperatorSetField(op_mass, "rho", elem_restriction_q_data, CEED_BASIS_NONE, q_data);
118*1b95d8c6SJeremy L Thompson     CeedOperatorSetField(op_mass, "u", elem_restriction_u, basis_u, CEED_VECTOR_ACTIVE);
119*1b95d8c6SJeremy L Thompson     CeedOperatorSetField(op_mass, "v", elem_restriction_u, basis_u, CEED_VECTOR_ACTIVE);
120*1b95d8c6SJeremy L Thompson     CeedOperatorAtPointsSetPoints(op_mass, elem_restriction_x_points, x_points);
121*1b95d8c6SJeremy L Thompson 
122*1b95d8c6SJeremy L Thompson     // Apply Setup Operator
123*1b95d8c6SJeremy L Thompson     CeedOperatorApply(op_setup, x, q_data, CEED_REQUEST_IMMEDIATE);
124*1b95d8c6SJeremy L Thompson 
125*1b95d8c6SJeremy L Thompson     // Fully assemble operator
126*1b95d8c6SJeremy L Thompson     CeedSize   num_entries;
127*1b95d8c6SJeremy L Thompson     CeedInt   *rows;
128*1b95d8c6SJeremy L Thompson     CeedInt   *cols;
129*1b95d8c6SJeremy L Thompson     CeedVector assembled;
130*1b95d8c6SJeremy L Thompson 
131*1b95d8c6SJeremy L Thompson     for (CeedInt k = 0; k < num_comp * num_comp * num_dofs * num_dofs; ++k) {
132*1b95d8c6SJeremy L Thompson       assembled_values[k] = 0.0;
133*1b95d8c6SJeremy L Thompson       assembled_true[k]   = 0.0;
134*1b95d8c6SJeremy L Thompson     }
135*1b95d8c6SJeremy L Thompson     CeedOperatorLinearAssembleSymbolic(op_mass, &num_entries, &rows, &cols);
136*1b95d8c6SJeremy L Thompson     CeedVectorCreate(ceed, num_entries, &assembled);
137*1b95d8c6SJeremy L Thompson     CeedOperatorLinearAssemble(op_mass, assembled);
138*1b95d8c6SJeremy L Thompson     {
139*1b95d8c6SJeremy L Thompson       const CeedScalar *assembled_array;
140*1b95d8c6SJeremy L Thompson 
141*1b95d8c6SJeremy L Thompson       CeedVectorGetArrayRead(assembled, CEED_MEM_HOST, &assembled_array);
142*1b95d8c6SJeremy L Thompson       for (CeedInt k = 0; k < num_entries; k++) {
143*1b95d8c6SJeremy L Thompson         assembled_values[rows[k] * num_comp * num_dofs + cols[k]] += assembled_array[k];
144*1b95d8c6SJeremy L Thompson       }
145*1b95d8c6SJeremy L Thompson       CeedVectorRestoreArrayRead(assembled, &assembled_array);
146*1b95d8c6SJeremy L Thompson     }
147*1b95d8c6SJeremy L Thompson 
148*1b95d8c6SJeremy L Thompson     // Manually assemble operator
149*1b95d8c6SJeremy L Thompson     CeedVectorSetValue(u, 0.0);
150*1b95d8c6SJeremy L Thompson     for (CeedInt j = 0; j < num_comp * num_dofs; j++) {
151*1b95d8c6SJeremy L Thompson       CeedScalar       *u_array;
152*1b95d8c6SJeremy L Thompson       const CeedScalar *v_array;
153*1b95d8c6SJeremy L Thompson 
154*1b95d8c6SJeremy L Thompson       // Set input
155*1b95d8c6SJeremy L Thompson       CeedVectorGetArray(u, CEED_MEM_HOST, &u_array);
156*1b95d8c6SJeremy L Thompson       u_array[j] = 1.0;
157*1b95d8c6SJeremy L Thompson       if (j) u_array[j - 1] = 0.0;
158*1b95d8c6SJeremy L Thompson       CeedVectorRestoreArray(u, &u_array);
159*1b95d8c6SJeremy L Thompson 
160*1b95d8c6SJeremy L Thompson       // Compute entries for column j
161*1b95d8c6SJeremy L Thompson       CeedOperatorApply(op_mass, u, v, CEED_REQUEST_IMMEDIATE);
162*1b95d8c6SJeremy L Thompson 
163*1b95d8c6SJeremy L Thompson       CeedVectorGetArrayRead(v, CEED_MEM_HOST, &v_array);
164*1b95d8c6SJeremy L Thompson       for (CeedInt i = 0; i < num_comp * num_dofs; i++) assembled_true[i * num_comp * num_dofs + j] = v_array[i];
165*1b95d8c6SJeremy L Thompson       CeedVectorRestoreArrayRead(v, &v_array);
166*1b95d8c6SJeremy L Thompson     }
167*1b95d8c6SJeremy L Thompson 
168*1b95d8c6SJeremy L Thompson     // Check output
169*1b95d8c6SJeremy L Thompson     for (CeedInt i = 0; i < num_comp * num_dofs; i++) {
170*1b95d8c6SJeremy L Thompson       for (CeedInt j = 0; j < num_comp * num_dofs; j++) {
171*1b95d8c6SJeremy L Thompson         if (fabs(assembled_values[i * num_dofs * num_comp + j] - assembled_true[i * num_dofs * num_comp + j]) > 100. * CEED_EPSILON) {
172*1b95d8c6SJeremy L Thompson           // LCOV_EXCL_START
173*1b95d8c6SJeremy L Thompson           printf("[%" CeedInt_FMT ", %" CeedInt_FMT "] Error in assembly: %f != %f\n", i, j, assembled_values[i * num_dofs * num_comp + j],
174*1b95d8c6SJeremy L Thompson                  assembled_true[i * num_dofs * num_comp + j]);
175*1b95d8c6SJeremy L Thompson           // LCOV_EXCL_STOP
176*1b95d8c6SJeremy L Thompson         }
177*1b95d8c6SJeremy L Thompson       }
178*1b95d8c6SJeremy L Thompson     }
179*1b95d8c6SJeremy L Thompson 
180*1b95d8c6SJeremy L Thompson     // Cleanup
181*1b95d8c6SJeremy L Thompson     free(rows);
182*1b95d8c6SJeremy L Thompson     free(cols);
183*1b95d8c6SJeremy L Thompson     CeedVectorDestroy(&x);
184*1b95d8c6SJeremy L Thompson     CeedVectorDestroy(&x_points);
185*1b95d8c6SJeremy L Thompson     CeedVectorDestroy(&q_data);
186*1b95d8c6SJeremy L Thompson     CeedVectorDestroy(&u);
187*1b95d8c6SJeremy L Thompson     CeedVectorDestroy(&v);
188*1b95d8c6SJeremy L Thompson     CeedVectorDestroy(&assembled);
189*1b95d8c6SJeremy L Thompson     CeedElemRestrictionDestroy(&elem_restriction_u);
190*1b95d8c6SJeremy L Thompson     CeedElemRestrictionDestroy(&elem_restriction_x);
191*1b95d8c6SJeremy L Thompson     CeedElemRestrictionDestroy(&elem_restriction_x_points);
192*1b95d8c6SJeremy L Thompson     CeedElemRestrictionDestroy(&elem_restriction_q_data);
193*1b95d8c6SJeremy L Thompson     CeedBasisDestroy(&basis_u);
194*1b95d8c6SJeremy L Thompson     CeedBasisDestroy(&basis_x);
195*1b95d8c6SJeremy L Thompson     CeedQFunctionDestroy(&qf_setup);
196*1b95d8c6SJeremy L Thompson     CeedQFunctionDestroy(&qf_mass);
197*1b95d8c6SJeremy L Thompson     CeedOperatorDestroy(&op_setup);
198*1b95d8c6SJeremy L Thompson     CeedOperatorDestroy(&op_mass);
199*1b95d8c6SJeremy L Thompson   }
200*1b95d8c6SJeremy L Thompson   CeedDestroy(&ceed);
201*1b95d8c6SJeremy L Thompson   return 0;
202*1b95d8c6SJeremy L Thompson }
203