1 /// @file 2 /// Test assembly of mass matrix operator point block diagonal 3 /// \test Test assembly of mass matrix operator point block diagonal 4 #include "t537-operator.h" 5 6 #include <ceed.h> 7 #include <math.h> 8 #include <stdio.h> 9 #include <stdlib.h> 10 11 int main(int argc, char **argv) { 12 Ceed ceed; 13 CeedElemRestriction elem_restriction_x, elem_restriction_u, elem_restriction_q_data; 14 CeedBasis basis_x, basis_u; 15 CeedQFunction qf_setup, qf_mass; 16 CeedOperator op_setup, op_mass; 17 CeedVector q_data, x, assembled, u, v; 18 CeedInt num_elem = 6, p = 3, q = 4, dim = 2, num_comp = 2; 19 CeedInt n_x = 3, n_y = 2; 20 CeedInt num_dofs = (n_x * 2 + 1) * (n_y * 2 + 1), num_qpts = num_elem * q * q; 21 CeedInt ind_x[num_elem * p * p]; 22 CeedInt *rows, *cols; 23 CeedSize num_entries; 24 CeedScalar assembled_true[num_comp * num_comp * num_dofs]; 25 CeedScalar assembled_full_true[num_comp * num_dofs][num_comp * num_dofs]; 26 27 CeedInit(argv[1], &ceed); 28 29 // Vectors 30 CeedVectorCreate(ceed, dim * num_dofs, &x); 31 { 32 CeedScalar x_array[dim * num_dofs]; 33 34 for (CeedInt i = 0; i < n_x * (p - 1) + 1; i++) { 35 for (CeedInt j = 0; j < n_y * (p - 1) + 1; j++) { 36 x_array[i + j * (n_x * (p - 1) + 1) + 0 * num_dofs] = (CeedScalar)i / ((p - 1) * n_x); 37 x_array[i + j * (n_x * (p - 1) + 1) + 1 * num_dofs] = (CeedScalar)j / ((p - 1) * n_y); 38 } 39 } 40 CeedVectorSetArray(x, CEED_MEM_HOST, CEED_COPY_VALUES, x_array); 41 } 42 CeedVectorCreate(ceed, num_qpts, &q_data); 43 44 // Restrictions 45 for (CeedInt i = 0; i < num_elem; i++) { 46 CeedInt col, row, offset; 47 col = i % n_x; 48 row = i / n_x; 49 offset = col * (p - 1) + row * (n_x * 2 + 1) * (p - 1); 50 for (CeedInt j = 0; j < p; j++) { 51 for (CeedInt k = 0; k < p; k++) ind_x[p * (p * i + k) + j] = offset + k * (n_x * 2 + 1) + j; 52 } 53 } 54 CeedElemRestrictionCreate(ceed, num_elem, p * p, dim, num_dofs, dim * num_dofs, CEED_MEM_HOST, CEED_USE_POINTER, ind_x, &elem_restriction_x); 55 CeedElemRestrictionCreate(ceed, num_elem, p * p, num_comp, num_dofs, num_comp * num_dofs, CEED_MEM_HOST, CEED_USE_POINTER, ind_x, 56 &elem_restriction_u); 57 58 CeedInt strides_q_data[3] = {1, q * q, q * q}; 59 CeedElemRestrictionCreateStrided(ceed, num_elem, q * q, 1, num_qpts, strides_q_data, &elem_restriction_q_data); 60 61 // Bases 62 CeedBasisCreateTensorH1Lagrange(ceed, dim, dim, p, q, CEED_GAUSS, &basis_x); 63 CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, p, q, CEED_GAUSS, &basis_u); 64 65 // QFunctions 66 CeedQFunctionCreateInterior(ceed, 1, setup, setup_loc, &qf_setup); 67 CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT); 68 CeedQFunctionAddInput(qf_setup, "dx", dim * dim, CEED_EVAL_GRAD); 69 CeedQFunctionAddOutput(qf_setup, "rho", 1, CEED_EVAL_NONE); 70 71 CeedQFunctionCreateInterior(ceed, 1, mass, mass_loc, &qf_mass); 72 CeedQFunctionAddInput(qf_mass, "rho", 1, CEED_EVAL_NONE); 73 CeedQFunctionAddInput(qf_mass, "u", num_comp, CEED_EVAL_INTERP); 74 CeedQFunctionAddOutput(qf_mass, "v", num_comp, CEED_EVAL_INTERP); 75 76 // Operators 77 CeedOperatorCreate(ceed, qf_setup, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup); 78 CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE); 79 CeedOperatorSetField(op_setup, "dx", elem_restriction_x, basis_x, CEED_VECTOR_ACTIVE); 80 CeedOperatorSetField(op_setup, "rho", elem_restriction_q_data, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE); 81 82 CeedOperatorCreate(ceed, qf_mass, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_mass); 83 CeedOperatorSetField(op_mass, "rho", elem_restriction_q_data, CEED_BASIS_NONE, q_data); 84 CeedOperatorSetField(op_mass, "u", elem_restriction_u, basis_u, CEED_VECTOR_ACTIVE); 85 CeedOperatorSetField(op_mass, "v", elem_restriction_u, basis_u, CEED_VECTOR_ACTIVE); 86 87 // Apply Setup Operator 88 CeedOperatorApply(op_setup, x, q_data, CEED_REQUEST_IMMEDIATE); 89 90 // Assemble diagonal 91 CeedVectorCreate(ceed, num_comp * num_comp * num_dofs, &assembled); 92 CeedOperatorLinearAssemblePointBlockDiagonal(op_mass, assembled, CEED_REQUEST_IMMEDIATE); 93 CeedOperatorLinearAssemblePointBlockDiagonalSymbolic(op_mass, &num_entries, &rows, &cols); 94 95 // Manually assemble point-block diagonal 96 CeedVectorCreate(ceed, num_comp * num_dofs, &u); 97 CeedVectorSetValue(u, 0.0); 98 CeedVectorCreate(ceed, num_comp * num_dofs, &v); 99 for (CeedInt i = 0; i < num_comp * num_dofs; i++) { 100 for (CeedInt j = 0; j < num_comp * num_dofs; j++) { 101 assembled_full_true[i][j] = 0.; 102 } 103 } 104 for (int i = 0; i < num_comp * num_comp * num_dofs; i++) assembled_true[i] = 0.0; 105 CeedInt ind_old = -1; 106 for (int i = 0; i < num_dofs; i++) { 107 for (int j = 0; j < num_comp; j++) { 108 CeedScalar *u_array; 109 const CeedScalar *v_array; 110 111 // Set input 112 CeedVectorGetArray(u, CEED_MEM_HOST, &u_array); 113 CeedInt ind = i + j * num_dofs; 114 u_array[ind] = 1.0; 115 if (ind > 0) u_array[ind_old] = 0.0; 116 ind_old = ind; 117 CeedVectorRestoreArray(u, &u_array); 118 119 // Compute effect of DoF i, comp j 120 CeedOperatorApply(op_mass, u, v, CEED_REQUEST_IMMEDIATE); 121 122 // Retrieve entry 123 CeedVectorGetArrayRead(v, CEED_MEM_HOST, &v_array); 124 for (int k = 0; k < num_comp; k++) assembled_true[i * num_comp * num_comp + k * num_comp + j] += v_array[i + k * num_dofs]; 125 for (CeedInt k = 0; k < num_comp * num_dofs; k++) assembled_full_true[k][i + num_dofs * j] += v_array[k]; 126 for (CeedInt k = 0; k < num_comp; k++) assembled_full_true[i * num_comp + k][i * num_comp + j] = v_array[i + k * num_dofs]; 127 CeedVectorRestoreArrayRead(v, &v_array); 128 } 129 } 130 131 // Check output 132 { 133 const CeedScalar *assembled_array; 134 135 CeedVectorGetArrayRead(assembled, CEED_MEM_HOST, &assembled_array); 136 for (int i = 0; i < num_comp * num_comp * num_dofs; i++) { 137 if (fabs(assembled_array[i] - assembled_true[i]) > 100. * CEED_EPSILON) { 138 // LCOV_EXCL_START 139 printf("[%" CeedInt_FMT "] Error in assembly: %f != %f\n", i, assembled_array[i], assembled_true[i]); 140 // LCOV_EXCL_STOP 141 } 142 } 143 144 for (CeedInt i = 0; i < num_entries; i++) { 145 if (fabs(assembled_array[i] - assembled_full_true[rows[i]][cols[i]]) > 100. * CEED_EPSILON) { 146 // LCOV_EXCL_START 147 printf("[%" CeedInt_FMT "] Error in symbolic assembly: %f != %f for row,col indices (%" CeedInt_FMT ", %" CeedInt_FMT ")\n", i, 148 assembled_array[i], assembled_full_true[rows[i]][cols[i]], rows[i], cols[i]); 149 // LCOV_EXCL_STOP 150 } 151 } 152 CeedVectorRestoreArrayRead(assembled, &assembled_array); 153 } 154 155 // Cleanup 156 CeedVectorDestroy(&x); 157 CeedVectorDestroy(&assembled); 158 CeedVectorDestroy(&q_data); 159 CeedVectorDestroy(&u); 160 CeedVectorDestroy(&v); 161 CeedElemRestrictionDestroy(&elem_restriction_u); 162 CeedElemRestrictionDestroy(&elem_restriction_x); 163 CeedElemRestrictionDestroy(&elem_restriction_q_data); 164 CeedBasisDestroy(&basis_u); 165 CeedBasisDestroy(&basis_x); 166 CeedQFunctionDestroy(&qf_setup); 167 CeedQFunctionDestroy(&qf_mass); 168 CeedOperatorDestroy(&op_setup); 169 CeedOperatorDestroy(&op_mass); 170 CeedDestroy(&ceed); 171 free(rows); 172 free(cols); 173 return 0; 174 } 175