Skip to content

Commit c11aec6

Browse files
Fix memory leak
1 parent f8f9dae commit c11aec6

File tree

1 file changed

+24
-25
lines changed

1 file changed

+24
-25
lines changed

interface/ceed-basis.c

+24-25
Original file line numberDiff line numberDiff line change
@@ -845,8 +845,7 @@ int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, Ceed
845845
@param[in] num_nodes Total number of nodes (dofs per element)
846846
@param[in] num_qpts Total number of quadrature points
847847
@param[in] interp Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points
848-
@param[in] curl Row-major (cdim * num_qpts * num_nodes, cdim = 1 if dim < 3 else dim) matrix expressing curl of basis functions at quadrature
849-
points
848+
@param[in] curl Row-major (cdim * num_qpts * num_nodes, cdim = 1 if dim < 3 else dim) matrix expressing curl of basis functions at quadrature points
850849
@param[in] q_ref Array of length num_qpts * dim holding the locations of quadrature points on the reference element
851850
@param[in] q_weight Array of length num_qpts holding the quadrature weights on the reference element
852851
@param[out] basis Address of the variable where the newly created CeedBasis will be stored.
@@ -1095,7 +1094,7 @@ int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode,
10951094
break;
10961095
case CEED_EVAL_NONE:
10971096
case CEED_EVAL_INTERP:
1098-
qdim = (fe_space == CEED_FE_SPACE_H1) ? 1 : dim;
1097+
qdim = (fe_space == CEED_FE_SPACE_H1) ? 1 : dim;
10991098
bad_dims = ((t_mode == CEED_TRANSPOSE && (u_length < num_elem * num_comp * num_qpts * qdim || v_length < num_elem * num_comp * num_nodes)) ||
11001099
(t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp * qdim || u_length < num_elem * num_comp * num_nodes)));
11011100
break;
@@ -1108,7 +1107,7 @@ int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode,
11081107
(t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp || u_length < num_elem * num_comp * num_nodes)));
11091108
break;
11101109
case CEED_EVAL_CURL:
1111-
cdim = (dim < 3) ? 1 : dim;
1110+
cdim = (dim < 3) ? 1 : dim;
11121111
bad_dims = ((t_mode == CEED_TRANSPOSE && (u_length < num_elem * num_comp * num_qpts * cdim || v_length < num_elem * num_comp * num_nodes)) ||
11131112
(t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp * cdim || u_length < num_elem * num_comp * num_nodes)));
11141113
break;
@@ -1451,13 +1450,14 @@ int CeedBasisDestroy(CeedBasis *basis) {
14511450
if (!*basis || --(*basis)->ref_count > 0) return CEED_ERROR_SUCCESS;
14521451
if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis));
14531452
if ((*basis)->contract) CeedCall(CeedTensorContractDestroy(&(*basis)->contract));
1453+
CeedCall(CeedFree(&(*basis)->q_ref_1d));
1454+
CeedCall(CeedFree(&(*basis)->q_weight_1d));
14541455
CeedCall(CeedFree(&(*basis)->interp));
14551456
CeedCall(CeedFree(&(*basis)->interp_1d));
14561457
CeedCall(CeedFree(&(*basis)->grad));
1457-
CeedCall(CeedFree(&(*basis)->div));
14581458
CeedCall(CeedFree(&(*basis)->grad_1d));
1459-
CeedCall(CeedFree(&(*basis)->q_ref_1d));
1460-
CeedCall(CeedFree(&(*basis)->q_weight_1d));
1459+
CeedCall(CeedFree(&(*basis)->div));
1460+
CeedCall(CeedFree(&(*basis)->curl));
14611461
CeedCall(CeedDestroy(&(*basis)->ceed));
14621462
CeedCall(CeedFree(basis));
14631463
return CEED_ERROR_SUCCESS;
@@ -1872,27 +1872,26 @@ CeedPragmaOptimizeOn
18721872
}
18731873
CeedPragmaOptimizeOn
18741874

1875-
/**
1876-
@brief Apply Householder Q matrix
1875+
/**
1876+
@brief Apply Householder Q matrix
18771877
1878-
Compute A = Q A, where Q is mxm and A is mxn.
1878+
Compute A = Q A, where Q is mxm and A is mxn.
18791879
1880-
@param[in,out] A Matrix to apply Householder Q to, in place
1881-
@param[in] Q Householder Q matrix
1882-
@param[in] tau Householder scaling factors
1883-
@param[in] t_mode Transpose mode for application
1884-
@param[in] m Number of rows in A
1885-
@param[in] n Number of columns in A
1886-
@param[in] k Number of elementary reflectors in Q, k<m
1887-
@param[in] row Row stride in A
1888-
@param[in] col Col stride in A
1880+
@param[in,out] A Matrix to apply Householder Q to, in place
1881+
@param[in] Q Householder Q matrix
1882+
@param[in] tau Householder scaling factors
1883+
@param[in] t_mode Transpose mode for application
1884+
@param[in] m Number of rows in A
1885+
@param[in] n Number of columns in A
1886+
@param[in] k Number of elementary reflectors in Q, k<m
1887+
@param[in] row Row stride in A
1888+
@param[in] col Col stride in A
18891889
1890-
@return An error code: 0 - success, otherwise - failure
1890+
@return An error code: 0 - success, otherwise - failure
18911891
1892-
@ref Developer
1893-
**/
1894-
int
1895-
CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n, CeedInt k,
1892+
@ref Developer
1893+
**/
1894+
int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n, CeedInt k,
18961895
CeedInt row, CeedInt col) {
18971896
CeedScalar *v;
18981897
CeedCall(CeedMalloc(m, &v));
@@ -1933,4 +1932,4 @@ int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScala
19331932
return CEED_ERROR_SUCCESS;
19341933
}
19351934

1936-
/// @}
1935+
/// @}

0 commit comments

Comments
 (0)