Skip to content

Commit 1088199

Browse files
authored
Merge pull request #1632 from CEED/jrwrigh/fix_basis_project
fix(basis): Use basis_from dim for projection matrices
2 parents 0de0756 + 706bf52 commit 1088199

File tree

3 files changed

+109
-1
lines changed

3 files changed

+109
-1
lines changed

interface/ceed-basis.c

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -234,7 +234,7 @@ static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, CeedBasis bas
234234
CeedScalar *interp_to_inv, *interp_from;
235235
const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL, *grad_from_source = NULL;
236236

237-
CeedCall(CeedBasisGetDimension(basis_to, &dim));
237+
CeedCall(CeedBasisGetDimension(basis_from, &dim));
238238
if (are_both_tensor) {
239239
CeedCall(CeedBasisGetInterp1D(basis_to, &interp_to_source));
240240
CeedCall(CeedBasisGetInterp1D(basis_from, &interp_from_source));

tests/t319-basis.c

Lines changed: 36 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,7 @@
11
/// @file
22
/// Test projection interp and grad in multiple dimensions
33
/// \test Test projection interp and grad in multiple dimensions
4+
#include "t319-basis.h"
45
#include <ceed.h>
56
#include <math.h>
67
#include <stdio.h>
@@ -201,6 +202,41 @@ int main(int argc, char **argv) {
201202
CeedBasisDestroy(&basis_to_nontensor);
202203
CeedBasisDestroy(&basis_project);
203204
}
205+
206+
// Test projection between basis of different topological dimension
207+
{
208+
CeedInt face_dim = 2, P_1D = 2;
209+
CeedBasis basis_face, basis_cell_to_face, basis_proj;
210+
211+
CeedScalar *q_ref = NULL, *q_weights = NULL;
212+
const CeedScalar *grad, *interp;
213+
CeedInt P, Q;
214+
GetCellToFaceTabulation(CEED_GAUSS, &P, &Q, &interp, &grad);
215+
216+
CeedBasisCreateTensorH1Lagrange(ceed, face_dim, 1, 2, P_1D, CEED_GAUSS, &basis_face);
217+
CeedBasisCreateH1(ceed, CEED_TOPOLOGY_HEX, 1, P, Q, (CeedScalar *)interp, (CeedScalar *)grad, q_ref, q_weights, &basis_cell_to_face);
218+
CeedBasisCreateProjection(basis_cell_to_face, basis_face, &basis_proj);
219+
const CeedScalar *interp_proj, *grad_proj, *interp_proj_ref, *grad_proj_ref;
220+
221+
GetCellToFaceTabulation(CEED_GAUSS_LOBATTO, NULL, NULL, &interp_proj_ref, &grad_proj_ref);
222+
CeedBasisGetInterp(basis_proj, &interp_proj);
223+
CeedBasisGetGrad(basis_proj, &grad_proj);
224+
CeedScalar tol = 100 * CEED_EPSILON;
225+
226+
for (CeedInt i = 0; i < 4 * 8; i++) {
227+
if (fabs(interp_proj[i] - ((CeedScalar *)interp_proj_ref)[i]) > tol)
228+
printf("Mixed Topology Projection: interp[%" CeedInt_FMT "] expected %f, got %f\n", i, interp_proj[i], ((CeedScalar *)interp_proj_ref)[i]);
229+
}
230+
231+
for (CeedInt i = 0; i < 3 * 4 * 8; i++) {
232+
if (fabs(grad_proj[i] - ((CeedScalar *)grad_proj_ref)[i]) > tol)
233+
printf("Mixed Topology Projection: grad[%" CeedInt_FMT "] expected %f, got %f\n", i, grad_proj[i], ((CeedScalar *)grad_proj_ref)[i]);
234+
}
235+
236+
CeedBasisDestroy(&basis_face);
237+
CeedBasisDestroy(&basis_cell_to_face);
238+
CeedBasisDestroy(&basis_proj);
239+
}
204240
CeedDestroy(&ceed);
205241
return 0;
206242
}

tests/t319-basis.h

Lines changed: 72 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,72 @@
1+
// Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2+
// All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3+
//
4+
// SPDX-License-Identifier: BSD-2-Clause
5+
//
6+
// This file is part of CEED: http://github.com/ceed
7+
8+
#include <ceed.h>
9+
10+
// Interpolation matrices for cell-to-face of Q1 hexahedral element onto it's "5" face (in PETSc)
11+
// Nodes are at Gauss-Lobatto points and quadrature points are Gauss, all over [-1,1] domain range
12+
const CeedScalar Q1_interp_gauss[4][8] = {
13+
{0.62200846792814612, 0, 0.16666666666666669, 0, 0.16666666666666669, 0, 0.044658198738520463, 0},
14+
{0.16666666666666669, 0, 0.62200846792814612, 0, 0.044658198738520463, 0, 0.16666666666666669, 0},
15+
{0.16666666666666669, 0, 0.044658198738520463, 0, 0.62200846792814612, 0, 0.16666666666666669, 0},
16+
{0.044658198738520463, 0, 0.16666666666666669, 0, 0.16666666666666669, 0, 0.62200846792814612, 0}
17+
};
18+
const CeedScalar Q1_grad_gauss[3][4][8] = {
19+
{{-0.31100423396407312, 0.31100423396407312, -0.083333333333333343, 0.083333333333333343, -0.083333333333333343, 0.083333333333333343,
20+
-0.022329099369260232, 0.022329099369260232},
21+
{-0.083333333333333343, 0.083333333333333343, -0.31100423396407312, 0.31100423396407312, -0.022329099369260232, 0.022329099369260232,
22+
-0.083333333333333343, 0.083333333333333343},
23+
{-0.083333333333333343, 0.083333333333333343, -0.022329099369260232, 0.022329099369260232, -0.31100423396407312, 0.31100423396407312,
24+
-0.083333333333333343, 0.083333333333333343},
25+
{-0.022329099369260232, 0.022329099369260232, -0.083333333333333343, 0.083333333333333343, -0.083333333333333343, 0.083333333333333343,
26+
-0.31100423396407312, 0.31100423396407312} },
27+
{{-0.39433756729740643, 0, 0.39433756729740643, 0, -0.10566243270259357, 0, 0.10566243270259357, 0},
28+
{-0.39433756729740643, 0, 0.39433756729740643, 0, -0.10566243270259357, 0, 0.10566243270259357, 0},
29+
{-0.10566243270259357, 0, 0.10566243270259357, 0, -0.39433756729740643, 0, 0.39433756729740643, 0},
30+
{-0.10566243270259357, 0, 0.10566243270259357, 0, -0.39433756729740643, 0, 0.39433756729740643, 0}},
31+
{{-0.39433756729740643, 0, -0.10566243270259357, 0, 0.39433756729740643, 0, 0.10566243270259357, 0},
32+
{-0.10566243270259357, 0, -0.39433756729740643, 0, 0.10566243270259357, 0, 0.39433756729740643, 0},
33+
{-0.39433756729740643, 0, -0.10566243270259357, 0, 0.39433756729740643, 0, 0.10566243270259357, 0},
34+
{-0.10566243270259357, 0, -0.39433756729740643, 0, 0.10566243270259357, 0, 0.39433756729740643, 0}}
35+
};
36+
37+
const CeedScalar Q1_interp_gauss_lobatto[4][8] = {
38+
{1, 0, 0, 0, 0, 0, 0, 0},
39+
{0, 0, 1, 0, 0, 0, 0, 0},
40+
{0, 0, 0, 0, 1, 0, 0, 0},
41+
{0, 0, 0, 0, 0, 0, 1, 0}
42+
};
43+
/* clang-format off */
44+
const CeedScalar Q1_grad_gauss_lobatto[3][4][8] = {
45+
{{-0.5, 0.5, 0, 0, 0, 0, 0, 0},
46+
{0, 0, -0.5, 0.5, 0, 0, 0, 0},
47+
{0, 0, 0, 0, -0.5, 0.5, 0, 0},
48+
{0, 0, 0, 0, 0, 0, -0.5, 0.5}},
49+
{{-0.5, 0, 0.5, 0, 0, 0, 0, 0},
50+
{-0.5, 0, 0.5, 0, 0, 0, 0, 0},
51+
{0, 0, 0, 0, -0.5, 0, 0.5, 0},
52+
{0, 0, 0, 0, -0.5, 0, 0.5, 0}},
53+
{{-0.5, 0, 0, 0, 0.5, 0, 0, 0},
54+
{0, 0, -0.5, 0, 0, 0, 0.5, 0},
55+
{-0.5, 0, 0, 0, 0.5, 0, 0, 0},
56+
{0, 0, -0.5, 0, 0, 0, 0.5, 0}}
57+
};
58+
/* clang-format on */
59+
60+
static void GetCellToFaceTabulation(CeedQuadMode quad_mode, CeedInt *P, CeedInt *Q, const CeedScalar **interp, const CeedScalar **grad) {
61+
if (P) *P = 8;
62+
if (Q) *Q = 4;
63+
64+
if (quad_mode == CEED_GAUSS) {
65+
*interp = (const CeedScalar *)Q1_interp_gauss;
66+
*grad = (const CeedScalar *)Q1_grad_gauss;
67+
}
68+
if (quad_mode == CEED_GAUSS_LOBATTO) {
69+
*interp = (const CeedScalar *)Q1_interp_gauss_lobatto;
70+
*grad = (const CeedScalar *)Q1_grad_gauss_lobatto;
71+
}
72+
}

0 commit comments

Comments
 (0)