aboutsummaryrefslogtreecommitdiff
path: root/libcrystfel
diff options
context:
space:
mode:
authorThomas White <taw@physics.org>2021-07-21 11:51:56 +0200
committerThomas White <taw@physics.org>2021-07-21 15:34:31 +0200
commit6d94d4115c254d344bbb927596a7141ef39fd298 (patch)
tree652d07c4d6a9521c9359aa13bbda27309f53e6f8 /libcrystfel
parent50cf52b0bc0baa683b9508568131a3f6281bf4ff (diff)
Add missing cleanup on error paths
Diffstat (limited to 'libcrystfel')
-rw-r--r--libcrystfel/src/cell-utils.c5
-rw-r--r--libcrystfel/src/fom.c58
-rw-r--r--libcrystfel/src/image-cbf.c11
-rw-r--r--libcrystfel/src/image-hdf5.c10
-rw-r--r--libcrystfel/src/image.c10
-rw-r--r--libcrystfel/src/index.c6
-rw-r--r--libcrystfel/src/indexers/asdf.c12
-rw-r--r--libcrystfel/src/indexers/felix.c2
-rw-r--r--libcrystfel/src/indexers/fromfile.c5
-rw-r--r--libcrystfel/src/indexers/mosflm.c2
-rw-r--r--libcrystfel/src/indexers/taketwo.c5
-rw-r--r--libcrystfel/src/integration.c17
-rw-r--r--libcrystfel/src/spectrum.c2
-rw-r--r--libcrystfel/src/stream.c5
-rw-r--r--libcrystfel/src/utils.c1
15 files changed, 124 insertions, 27 deletions
diff --git a/libcrystfel/src/cell-utils.c b/libcrystfel/src/cell-utils.c
index 30744bd4..2f3dfc6b 100644
--- a/libcrystfel/src/cell-utils.c
+++ b/libcrystfel/src/cell-utils.c
@@ -1588,7 +1588,10 @@ static Rational *find_candidates(double len, double *a, double *b, double *c,
if ( cands == NULL ) return NULL;
rat = rtnl_list(-5, 5, 1, csl ? 4 : 1, &nrat);
- if ( rat == NULL ) return NULL;
+ if ( rat == NULL ) {
+ free(cands);
+ return NULL;
+ }
for ( ia=0; ia<nrat; ia++ ) {
for ( ib=0; ib<nrat; ib++ ) {
diff --git a/libcrystfel/src/fom.c b/libcrystfel/src/fom.c
index 33887bd8..2f74e945 100644
--- a/libcrystfel/src/fom.c
+++ b/libcrystfel/src/fom.c
@@ -88,10 +88,15 @@ static struct fom_context *init_fom(enum fom_type fom, int nmax, int nshells)
fctx->cts[i] = 0;
}
- fctx->num = NULL;
- fctx->den = NULL;
fctx->num2 = NULL;
fctx->den2 = NULL;
+ fctx->num = NULL;
+ fctx->den = NULL;
+ fctx->n_meas = NULL;
+ fctx->vec1 = NULL;
+ fctx->vec2 = NULL;
+ fctx->n = NULL;
+ fctx->n_within = NULL;
fctx->possible = NULL;
switch ( fctx->fom ) {
@@ -99,7 +104,7 @@ static struct fom_context *init_fom(enum fom_type fom, int nmax, int nshells)
case FOM_RANORSPLIT :
fctx->num2 = malloc(nshells*sizeof(double));
fctx->den2 = malloc(nshells*sizeof(double));
- if ( (fctx->num2 == NULL) || (fctx->den2 == NULL) ) return NULL;
+ if ( (fctx->num2 == NULL) || (fctx->den2 == NULL) ) goto out;
for ( i=0; i<nshells; i++ ) {
fctx->num2[i] = 0.0;
fctx->den2[i] = 0.0;
@@ -116,7 +121,7 @@ static struct fom_context *init_fom(enum fom_type fom, int nmax, int nshells)
case FOM_REDUNDANCY :
fctx->num = malloc(nshells*sizeof(double));
fctx->den = malloc(nshells*sizeof(double));
- if ( (fctx->num == NULL) || (fctx->den == NULL) ) return NULL;
+ if ( (fctx->num == NULL) || (fctx->den == NULL) ) goto out;
for ( i=0; i<nshells; i++ ) {
fctx->num[i] = 0.0;
fctx->den[i] = 0.0;
@@ -124,12 +129,12 @@ static struct fom_context *init_fom(enum fom_type fom, int nmax, int nshells)
break;
case FOM_COMPLETENESS :
- /* Uses 'cts' and 'possible' only */
+ /* Uses 'cts' and 'possible' only - see calculate_possible() */
break;
case FOM_NUM_MEASUREMENTS :
fctx->n_meas = calloc(nshells, sizeof(long int));
- if ( fctx->n_meas == NULL ) return NULL;
+ if ( fctx->n_meas == NULL ) goto out;
break;
case FOM_CC :
@@ -138,15 +143,19 @@ static struct fom_context *init_fom(enum fom_type fom, int nmax, int nshells)
case FOM_CRDANO :
fctx->vec1 = malloc(nshells*sizeof(double *));
fctx->vec2 = malloc(nshells*sizeof(double *));
- if ( (fctx->vec1 == NULL) || (fctx->vec2 == NULL) ) return NULL;
+ if ( (fctx->vec1 == NULL) || (fctx->vec2 == NULL) ) goto out;
+ for ( i=0; i<nshells; i++ ) {
+ fctx->vec1[i] = NULL;
+ fctx->vec2[i] = NULL;
+ }
for ( i=0; i<nshells; i++ ) {
fctx->vec1[i] = malloc(nmax*sizeof(double));
- if ( fctx->vec1[i] == NULL ) return NULL;
+ if ( fctx->vec1[i] == NULL ) goto out;
fctx->vec2[i] = malloc(nmax*sizeof(double));
- if ( fctx->vec2[i] == NULL ) return NULL;
- fctx->n = malloc(nshells*sizeof(int));
- if ( fctx->n == NULL ) return NULL;
+ if ( fctx->vec2[i] == NULL ) goto out;
}
+ fctx->n = malloc(nshells*sizeof(int));
+ if ( fctx->n == NULL ) goto out;
for ( i=0; i<nshells; i++ ) {
fctx->n[i] = 0;
}
@@ -156,7 +165,7 @@ static struct fom_context *init_fom(enum fom_type fom, int nmax, int nshells)
case FOM_D1SIG :
case FOM_D2SIG :
fctx->n_within = malloc(nshells*sizeof(int));
- if ( fctx->n_within == NULL ) return NULL;
+ if ( fctx->n_within == NULL ) goto out;
for ( i=0; i<nshells; i++ ) {
fctx->n_within[i] = 0;
}
@@ -165,6 +174,29 @@ static struct fom_context *init_fom(enum fom_type fom, int nmax, int nshells)
}
return fctx;
+
+out:
+ free(fctx->num2);
+ free(fctx->den2);
+ free(fctx->num);
+ free(fctx->den);
+ free(fctx->n_meas);
+ if ( fctx->vec1 != NULL ) {
+ for ( i=0; i<nshells; i++ ) {
+ free(fctx->vec1[i]);
+ }
+ free(fctx->vec1);
+ }
+ if ( fctx->vec2 != NULL ) {
+ for ( i=0; i<nshells; i++ ) {
+ free(fctx->vec2[i]);
+ }
+ free(fctx->vec2);
+ }
+ free(fctx->n);
+ free(fctx->n_within);
+ free(fctx);
+ return NULL;
}
@@ -588,6 +620,8 @@ struct fom_shells *fom_make_resolution_shells(double rmin, double rmax,
if ( (s->rmins==NULL) || (s->rmaxs==NULL) ) {
ERROR("Couldn't allocate memory for resolution shells.\n");
+ free(s->rmins);
+ free(s->rmaxs);
free(s);
return NULL;
}
diff --git a/libcrystfel/src/image-cbf.c b/libcrystfel/src/image-cbf.c
index 859ef731..1e5659f4 100644
--- a/libcrystfel/src/image-cbf.c
+++ b/libcrystfel/src/image-cbf.c
@@ -314,7 +314,10 @@ static float *read_cbf_data(const char *filename, int gz, int *w, int *h)
} while ( len_read == bufinc );
fh = fmemopen(buf, len, "rb");
- if ( fh == NULL ) return NULL;
+ if ( fh == NULL ) {
+ free(buf);
+ return NULL;
+ }
gzclose(gzfh);
@@ -507,7 +510,11 @@ signed int is_cbf_file(const char *filename)
fh = fopen(filename, "r");
if ( fh == NULL ) return -1;
- if ( fgets(line, 1024, fh) == NULL ) return 0;
+ if ( fgets(line, 1024, fh) == NULL ) {
+ fclose(fh);
+ return 0;
+ }
+
fclose(fh);
if ( strstr(line, "CBF") == NULL ) {
diff --git a/libcrystfel/src/image-hdf5.c b/libcrystfel/src/image-hdf5.c
index 38ebf42b..9cad0a33 100644
--- a/libcrystfel/src/image-hdf5.c
+++ b/libcrystfel/src/image-hdf5.c
@@ -94,7 +94,10 @@ char **read_path_parts(const char *ev_orig, int *pn_plvals)
}
ev = strdup(ev_orig);
- if ( ev == NULL ) return NULL;
+ if ( ev == NULL ) {
+ free(plvals);
+ return NULL;
+ }
start = ev;
do {
@@ -166,7 +169,10 @@ int *read_dim_parts(const char *ev_orig, int *pn_dvals)
}
ev = strdup(ev+2);
- if ( ev == NULL ) return NULL;
+ if ( ev == NULL ) {
+ free(dvals);
+ return NULL;
+ }
start = ev;
done = 0;
diff --git a/libcrystfel/src/image.c b/libcrystfel/src/image.c
index 74e6d8d3..3c5af479 100644
--- a/libcrystfel/src/image.c
+++ b/libcrystfel/src/image.c
@@ -603,7 +603,10 @@ int create_detgeom(struct image *image, const DataTemplate *dtempl)
if ( detgeom == NULL ) return 1;
detgeom->panels = malloc(dtempl->n_panels*sizeof(struct detgeom_panel));
- if ( detgeom->panels == NULL ) return 1;
+ if ( detgeom->panels == NULL ) {
+ free(detgeom);
+ return 1;
+ }
detgeom->n_panels = dtempl->n_panels;
@@ -1493,7 +1496,10 @@ char **image_expand_frames(const DataTemplate *dtempl,
list = malloc(sizeof(char *));
if ( list == NULL ) return NULL;
list[0] = strdup("//");
- if ( list[0] == NULL ) return NULL;
+ if ( list[0] == NULL ) {
+ free(list);
+ return NULL;
+ }
*n_frames = 1;
return list;
}
diff --git a/libcrystfel/src/index.c b/libcrystfel/src/index.c
index ce341197..d8a135fc 100644
--- a/libcrystfel/src/index.c
+++ b/libcrystfel/src/index.c
@@ -410,11 +410,15 @@ IndexingPrivate *setup_indexing(const char *method_list,
ttopts,
fromfile_opts);
- if ( ipriv->engine_private[i] == NULL ) return NULL;
+ if ( ipriv->engine_private[i] == NULL ) {
+ free(ipriv->engine_private);
+ return NULL;
+ }
for ( j=0; j<i; j++ ) {
if ( methods[i] == methods[j] ) {
ERROR("Duplicate indexing method.\n");
+ free(ipriv->engine_private);
return NULL;
}
}
diff --git a/libcrystfel/src/indexers/asdf.c b/libcrystfel/src/indexers/asdf.c
index a8943688..7030a3e7 100644
--- a/libcrystfel/src/indexers/asdf.c
+++ b/libcrystfel/src/indexers/asdf.c
@@ -180,11 +180,18 @@ static struct asdf_cell *asdf_cell_new(int n)
if (c->reflections == NULL) return NULL;
c->indices = malloc(sizeof(double *) * n);
- if (c->indices == NULL) return NULL;
+ if (c->indices == NULL) {
+ free(c->reflections);
+ return NULL;
+ }
for ( i = 0; i < n; i++ ) {
c->indices[i] = malloc(sizeof(double) * 3);
- if (c->indices[i] == NULL) return NULL;
+ if (c->indices[i] == NULL) {
+ free(c->reflections);
+ free(c->indices);
+ return NULL;
+ }
}
c->n = 0;
@@ -973,6 +980,7 @@ static int index_refls(gsl_vector **reflections, int N_reflections,
struct tvector *tvectors = malloc(N_triplets * sizeof(struct tvector));
if (tvectors == NULL) {
ERROR("Failed to allocate tvectors in index_refls!\n");
+ free(fits);
return 0;
}
diff --git a/libcrystfel/src/indexers/felix.c b/libcrystfel/src/indexers/felix.c
index db9cccf8..b7d2a8a3 100644
--- a/libcrystfel/src/indexers/felix.c
+++ b/libcrystfel/src/indexers/felix.c
@@ -143,6 +143,7 @@ static int read_felix(struct felix_private *gp, struct image *image,
/* Read and discard first line */
if ( fgets( line, 1024, fh ) == NULL ) {
ERROR("Failed to read *.felix file.\n");
+ fclose(fh);
return 0;
}
@@ -639,6 +640,7 @@ void *felix_prepare(IndexingMethod *indm, UnitCell *cell,
if ( gp->spacegroup == 0 ) {
ERROR("Couldn't determine representative space group for your cell.\n");
ERROR("Try again with a more conventional cell.\n");
+ free(gp);
return NULL;
}
diff --git a/libcrystfel/src/indexers/fromfile.c b/libcrystfel/src/indexers/fromfile.c
index f0372489..d7dfbd09 100644
--- a/libcrystfel/src/indexers/fromfile.c
+++ b/libcrystfel/src/indexers/fromfile.c
@@ -211,7 +211,10 @@ void *fromfile_prepare(IndexingMethod *indm, struct fromfile_options *opts)
}
dp = malloc(sizeof(struct fromfile_private));
- if ( dp == NULL ) return NULL;
+ if ( dp == NULL ) {
+ fclose(fh);
+ return NULL;
+ }
dp->sol_hash = NULL;
diff --git a/libcrystfel/src/indexers/mosflm.c b/libcrystfel/src/indexers/mosflm.c
index c56f3ae6..129e6859 100644
--- a/libcrystfel/src/indexers/mosflm.c
+++ b/libcrystfel/src/indexers/mosflm.c
@@ -244,6 +244,7 @@ static int read_newmat(struct mosflm_data *mosflm, const char *filename,
n += fscanf(fh, "%f %f %f\n", &asz, &bsz, &csz);
if ( n != 9 ) {
STATUS("Fewer than 9 parameters found in NEWMAT file.\n");
+ fclose(fh);
return 1;
}
@@ -371,6 +372,7 @@ static void write_img(struct image *image, const char *filename)
fh = fopen(filename, "w");
if ( !fh ) {
ERROR("Couldn't open temporary file '%s'\n", filename);
+ free(intimage);
return;
}
diff --git a/libcrystfel/src/indexers/taketwo.c b/libcrystfel/src/indexers/taketwo.c
index 60bb6fed..cd4d87a6 100644
--- a/libcrystfel/src/indexers/taketwo.c
+++ b/libcrystfel/src/indexers/taketwo.c
@@ -1769,6 +1769,10 @@ static int match_obs_to_cell_vecs(struct TheoryVec *cell_vecs, int cell_vec_coun
match_array = &(obs_vecs[i].matches);
match_count = &(obs_vecs[i].match_num);
+ if ( for_sort == NULL ) {
+ return 0;
+ }
+
/* Sort in order to get most agreeable matches first */
qsort(for_sort, count, sizeof(struct sortme), sort_theory_distances);
*match_array = malloc(count*sizeof(struct TheoryVec));
@@ -1831,6 +1835,7 @@ static int gen_observed_vecs(struct rvec *rlps, int rlp_count,
spot_vec.match_num = 0;
spot_vec.her_rlp = &rlps[i];
spot_vec.his_rlp = &rlps[j];
+ spot_vec.in_network = 0;
cell->obs_vecs[count - 1] = spot_vec;
}
diff --git a/libcrystfel/src/integration.c b/libcrystfel/src/integration.c
index b2b6ccfa..326e3183 100644
--- a/libcrystfel/src/integration.c
+++ b/libcrystfel/src/integration.c
@@ -509,6 +509,7 @@ struct intcontext *intcontext_new(struct image *image,
ic->bm = malloc(ic->w * ic->w * sizeof(enum boxmask_val));
if ( ic->bm == NULL ) {
ERROR("Failed to allocate box mask.\n");
+ free(ic);
return NULL;
}
@@ -516,12 +517,21 @@ struct intcontext *intcontext_new(struct image *image,
ic->n_reference_profiles = 1;
ic->reference_profiles = calloc(ic->n_reference_profiles,
sizeof(double *));
- if ( ic->reference_profiles == NULL ) return NULL;
+ if ( ic->reference_profiles == NULL ) {
+ free(ic);
+ return NULL;
+ }
ic->reference_den = calloc(ic->n_reference_profiles, sizeof(double *));
- if ( ic->reference_den == NULL ) return NULL;
+ if ( ic->reference_den == NULL ) {
+ free(ic);
+ return NULL;
+ }
ic->n_profiles_in_reference = calloc(ic->n_reference_profiles,
sizeof(int));
- if ( ic->n_profiles_in_reference == NULL ) return NULL;
+ if ( ic->n_profiles_in_reference == NULL ) {
+ free(ic);
+ return NULL;
+ }
for ( i=0; i<ic->n_reference_profiles; i++ ) {
ic->reference_profiles[i] = malloc(ic->w*ic->w*sizeof(double));
if ( ic->reference_profiles[i] == NULL ) return NULL;
@@ -534,6 +544,7 @@ struct intcontext *intcontext_new(struct image *image,
ic->n_boxes = 0;
ic->max_boxes = 0;
if ( alloc_boxes(ic, 32) ) {
+ free(ic);
return NULL;
}
diff --git a/libcrystfel/src/spectrum.c b/libcrystfel/src/spectrum.c
index 78ad2e6e..bf38a73a 100644
--- a/libcrystfel/src/spectrum.c
+++ b/libcrystfel/src/spectrum.c
@@ -435,6 +435,7 @@ Spectrum *spectrum_load(const char *filename)
if ( fgets(line, 1024, fh) != line ) {
ERROR("Failed to read '%s'\n", filename);
spectrum_free(s);
+ fclose(fh);
return NULL;
}
@@ -444,6 +445,7 @@ Spectrum *spectrum_load(const char *filename)
ERROR("Failed to read ESRF spectrum from %s\n",
filename);
spectrum_free(s);
+ fclose(fh);
return NULL;
}
} else {
diff --git a/libcrystfel/src/stream.c b/libcrystfel/src/stream.c
index cd64a485..8ffaed1a 100644
--- a/libcrystfel/src/stream.c
+++ b/libcrystfel/src/stream.c
@@ -1568,7 +1568,10 @@ static void add_index_record(StreamIndex *index,
new_ptrs = realloc(index->ptrs,
new_max_keys*sizeof(long int));
- if ( new_ptrs == NULL ) return;
+ if ( new_ptrs == NULL ) {
+ free(new_keys);
+ return;
+ }
index->keys = new_keys;
index->ptrs = new_ptrs;
diff --git a/libcrystfel/src/utils.c b/libcrystfel/src/utils.c
index 864b1a30..4aecfc37 100644
--- a/libcrystfel/src/utils.c
+++ b/libcrystfel/src/utils.c
@@ -820,6 +820,7 @@ char *load_entire_file(const char *filename)
if ( fread(contents, 1, statbuf.st_size, fh) != statbuf.st_size ) {
ERROR("Failed to read file '%s'\n", filename);
+ fclose(fh);
free(contents);
return NULL;
}