Skip to content
Snippets Groups Projects
Commit d1fa9aeb authored by Moritz Hanke's avatar Moritz Hanke
Browse files

improves error messages in interp_method_avg

parent 2d1f481e
No related branches found
No related tags found
No related merge requests found
......@@ -522,7 +522,8 @@ static int compute_weights_dist_no(
// respect to the given triangle
static void compute_barycentric_coords(
double * barycentric_coords, size_t triangle_indices[3],
yac_const_coordinate_pointer grid_coords) {
yac_const_coordinate_pointer grid_coords,
char const * caller) {
double A[3][3];
lapack_int n = 3, nrhs = 1, lda = n, ldx = n, ipiv[3];
......@@ -535,11 +536,25 @@ static void compute_barycentric_coords(
// where: A is the matrix consisting of the vertex coordinates of the three
// corners of the triangle
// we compute b by solving this linear system using LAPACK
YAC_ASSERT(
YAC_ASSERT_F(
!LAPACKE_dgesv(
LAPACK_COL_MAJOR, n, nrhs, &A[0][0], lda,
ipiv, barycentric_coords, ldx),
"ERROR: internal error (could not solve linear 3x3 system)")
LAPACK_COL_MAJOR, n, nrhs, &A[0][0], lda, ipiv, barycentric_coords, ldx),
"ERROR(%s::compute_sb_coords): "
"internal error (could not solve linear 3x3 system)\n"
"(vector: (% .6e;% .6e;% .6e)\n"
" triangle: ((% .6e;% .6e;% .6e),\n"
" (% .6e;% .6e;% .6e),\n"
" (% .6e;% .6e;% .6e))",
caller, barycentric_coords[0], barycentric_coords[1], barycentric_coords[2],
grid_coords[triangle_indices[0]][0],
grid_coords[triangle_indices[0]][1],
grid_coords[triangle_indices[0]][2],
grid_coords[triangle_indices[1]][0],
grid_coords[triangle_indices[1]][1],
grid_coords[triangle_indices[1]][2],
grid_coords[triangle_indices[2]][0],
grid_coords[triangle_indices[2]][1],
grid_coords[triangle_indices[2]][2])
}
// returns 0 if edge_direction_vec points south; 1 otherwise
......@@ -729,7 +744,8 @@ static int compute_weights_bary_yes(
// compute barycentric coordinates for the current triangle
compute_barycentric_coords(
temp_barycentric_coords, triangle_indices[i], src_coords);
temp_barycentric_coords, triangle_indices[i], src_coords,
"compute_weights_bary_yes");
double curr_min_barycentric_coord =
MIN(temp_barycentric_coords[0],
......@@ -832,7 +848,8 @@ static int compute_weights_bary_no(
// compute barycentric coordinates for the current triangle
compute_barycentric_coords(
temp_barycentric_coords, triangle_indices[i], src_coords);
temp_barycentric_coords, triangle_indices[i], src_coords,
"compute_weights_bary_no");
double curr_min_barycentric_coord =
MIN(temp_barycentric_coords[0],
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment