summaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorMagnus Lundborg <lundborg.magnus@gmail.com>2013-10-11 08:26:12 (GMT)
committerMagnus Lundborg <lundborg.magnus@gmail.com>2013-10-11 08:26:12 (GMT)
commit5a1a62f1e62db1afbac3eb138671a3f1222eb888 (patch)
tree2927355ae8d4ef89c6a2f221e14177a475399149 /src
parent437ebc84d8d5baa2564d9ffe9a10800b15e961e3 (diff)
Fixing many llvm clang warnings.
Diffstat (limited to 'src')
-rw-r--r--src/compression/bwt.c14
-rw-r--r--src/compression/coder.c8
-rw-r--r--src/compression/tng_compress.c22
-rw-r--r--src/compression/xtc3.c32
-rw-r--r--src/lib/tng_io.c215
-rw-r--r--src/tests/tng_io_testing.c38
-rw-r--r--src/tests/tng_parallel_read.c31
7 files changed, 206 insertions, 154 deletions
diff --git a/src/compression/bwt.c b/src/compression/bwt.c
index ea75199..e5cb5ce 100644
--- a/src/compression/bwt.c
+++ b/src/compression/bwt.c
@@ -37,8 +37,8 @@ static int compare_index(int i1,int i2,int nvals,unsigned int *vals,unsigned int
int k1=(int)(nrepeat[i1]&0xFFU);
int repeat2=(int)(nrepeat[i2]>>8);
int k2=(int)(nrepeat[i2]&0xFFU);
-
- if ((repeat1>1) && (repeat1>1) && (k1==k2))
+
+ if ((repeat1>1) && (repeat2>1) && (k1==k2))
{
int maxskip=0;
/* Yes. Compare the repeating patterns. */
@@ -88,7 +88,7 @@ void Ptngc_bwt_merge_sort_inner(int *indices, int nvals,unsigned int *vals,
middle=start+(end-start)/2;
#if 0
printf("For start %d end %d obtained new middle: %d\n",start,end,middle);
-#endif
+#endif
Ptngc_bwt_merge_sort_inner(indices,nvals,vals,
start,middle,
nrepeat,
@@ -126,7 +126,7 @@ void Ptngc_bwt_merge_sort_inner(int *indices, int nvals,unsigned int *vals,
if (compare_index(indices[ileft],indices[iright],nvals,vals,nrepeat)>0)
{
workarray[i]=indices[iright];
- iright++;
+ iright++;
}
else
{
@@ -221,7 +221,7 @@ void Ptngc_comp_to_bwt(unsigned int *vals, int nvals,
pattern. */
#ifdef SHOWIT
printf("Best j and k is now %d and %d\n",good_j,good_k);
-#endif
+#endif
}
}
else
@@ -257,7 +257,7 @@ void Ptngc_comp_to_bwt(unsigned int *vals, int nvals,
nrepeat[i+m]=((unsigned int) (good_k)) | (((unsigned int) (repeat))<<8);
}
/* If no repetition was found for this value signal that here. */
- if (!nrepeat[i])
+ if (!nrepeat[i])
nrepeat[i+m]=257U; /* This is 1<<8 | 1 */
}
}
@@ -334,4 +334,4 @@ void Ptngc_comp_from_bwt(unsigned int *input, int nvals, int index,
free(p);
free(c);
}
-
+
diff --git a/src/compression/coder.c b/src/compression/coder.c
index 270e4d3..ca3de86 100644
--- a/src/compression/coder.c
+++ b/src/compression/coder.c
@@ -171,7 +171,6 @@ static int pack_stopbits_item(struct coder *coder_inst,int item, unsigned char *
else if (item<0)
s=2+(-item-1)*2;
return write_stop_bit_code(coder_inst,s,coding_parameter,output);
- return 0;
}
static int pack_triplet(struct coder *coder_inst,unsigned int *s, unsigned char **output, int coding_parameter,
@@ -194,7 +193,6 @@ static int pack_triplet(struct coder *coder_inst,unsigned int *s, unsigned char
{
if (this_base>max_base)
return 1;
- this_base=max_base;
bits_per_value=maxbits;
jbase=3;
}
@@ -332,7 +330,7 @@ unsigned char *Ptngc_pack_array(struct coder *coder_inst,int *input, int *length
static int unpack_array_stop_bits(struct coder *coder_inst,unsigned char *packed,int *output, int length, int coding_parameter)
{
- (void)coder_inst;
+ (void) coder_inst;
int i,j;
unsigned int extract_mask=0x80;
unsigned char *ptr=packed;
@@ -385,7 +383,7 @@ static int unpack_array_stop_bits(struct coder *coder_inst,unsigned char *packed
static int unpack_array_triplet(struct coder *coder_inst,unsigned char *packed,int *output, int length, int coding_parameter)
{
- (void)coder_inst;
+ (void) coder_inst;
int i,j;
unsigned int extract_mask=0x80;
unsigned char *ptr=packed;
@@ -457,7 +455,7 @@ static int unpack_array_triplet(struct coder *coder_inst,unsigned char *packed,i
static int unpack_array_bwlzh(struct coder *coder_inst,unsigned char *packed,int *output, int length, int natoms)
{
- (void)coder_inst;
+ (void) coder_inst;
int i,j,k,n=length;
unsigned int *pval=warnmalloc(n*sizeof *pval);
int nframes=n/natoms/3;
diff --git a/src/compression/tng_compress.c b/src/compression/tng_compress.c
index ebd495b..2578336 100644
--- a/src/compression/tng_compress.c
+++ b/src/compression/tng_compress.c
@@ -481,13 +481,14 @@ static void compress_quantized_vel(int *quant, int *quant_inter,
if (data)
bufferfix((unsigned char*)data+bufloc,prec_hi,4);
bufloc+=4;
+
+ length=natoms*3;
/* The initial frame */
if ((initial_coding==TNG_COMPRESS_ALGO_VEL_STOPBIT_ONETOONE) ||
(initial_coding==TNG_COMPRESS_ALGO_VEL_TRIPLET_ONETOONE) ||
(initial_coding==TNG_COMPRESS_ALGO_VEL_BWLZH_ONETOONE))
{
struct coder *coder=Ptngc_coder_init();
- length=natoms*3;
datablock=(char*)Ptngc_pack_array(coder,quant,&length,
initial_coding,initial_coding_parameter,natoms,speed);
Ptngc_coder_deinit(coder);
@@ -497,10 +498,12 @@ static void compress_quantized_vel(int *quant, int *quant_inter,
bufferfix((unsigned char*)data+bufloc,(fix_t)length,4);
bufloc+=4;
/* The actual data block. */
- if (data)
- memcpy(data+bufloc,datablock,length);
- free(datablock);
- bufloc+=length;
+ if (data && datablock)
+ {
+ memcpy(data+bufloc,datablock,length);
+ free(datablock);
+ bufloc+=length;
+ }
/* The remaining frames */
if (nframes>1)
{
@@ -683,7 +686,6 @@ static void determine_best_pos_initial_coding(int *quant, int *quant_intra, int
{
best_coding=current_coding;
best_coding_parameter=current_coding_parameter;
- best_code_size=current_code_size;
}
}
*initial_coding=best_coding;
@@ -846,7 +848,6 @@ static void determine_best_pos_coding(int *quant, int *quant_inter, int *quant_i
{
best_coding=current_coding;
best_coding_parameter=current_coding_parameter;
- best_code_size=current_code_size;
}
}
*coding=best_coding;
@@ -950,7 +951,6 @@ static void determine_best_vel_initial_coding(int *quant, int natoms, int speed,
{
best_coding=current_coding;
best_coding_parameter=current_coding_parameter;
- best_code_size=current_code_size;
}
}
*initial_coding=best_coding;
@@ -1090,7 +1090,6 @@ static void determine_best_vel_coding(int *quant, int *quant_inter, int natoms,
if (current_code_size<best_code_size)
{
best_coding=current_coding;
- best_code_size=current_code_size;
best_coding_parameter=current_coding_parameter;
}
}
@@ -1492,7 +1491,6 @@ int DECLSPECDLLEXPORT tng_compress_inquire(char *data,int *vel, int *natoms,
prec_lo=readbufferfix((unsigned char *)data+bufloc,4);
bufloc+=4;
prec_hi=readbufferfix((unsigned char *)data+bufloc,4);
- bufloc+=4;
*precision=PRECISION(prec_hi, prec_lo);
algo[0]=initial_coding;
algo[1]=initial_coding_parameter;
@@ -1583,8 +1581,6 @@ static int tng_compress_uncompress_pos_gen(char *data,double *posd,float *posf,i
/* The remaining frames. */
if (nframes>1)
{
- /* The data block length. */
- length=(int)readbufferfix((unsigned char *)data+bufloc,4);
bufloc+=4;
coder=Ptngc_coder_init();
rval=Ptngc_unpack_array(coder,(unsigned char *)data+bufloc,quant+natoms*3,(nframes-1)*natoms*3,
@@ -1721,8 +1717,6 @@ static int tng_compress_uncompress_vel_gen(char *data,double *veld,float *velf,i
/* The remaining frames. */
if (nframes>1)
{
- /* The data block length. */
- length=(int)readbufferfix((unsigned char *)data+bufloc,4);
bufloc+=4;
coder=Ptngc_coder_init();
rval=Ptngc_unpack_array(coder,(unsigned char *)data+bufloc,quant+natoms*3,(nframes-1)*natoms*3,
diff --git a/src/compression/xtc3.c b/src/compression/xtc3.c
index 197cc12..297e948 100644
--- a/src/compression/xtc3.c
+++ b/src/compression/xtc3.c
@@ -512,7 +512,6 @@ static void buffer_large(struct xtc3_context *xtc3_context, int *input, int inpd
thislen=compute_intlen(interdelta);
if (thislen*TRESHOLD_INTRA_INTER_DIRECT<minlen)
{
- minlen=thislen;
best_type=2; /* Inter delta */
}
}
@@ -730,7 +729,6 @@ static void base_decompress(unsigned char *input, int len, unsigned int *output)
while (nvals_left)
{
int n;
- int numints;
if (basegiven==0)
{
base=(unsigned int)(input[0])|
@@ -755,24 +753,26 @@ static void base_decompress(unsigned char *input, int len, unsigned int *output)
#ifdef SHOWIT
fprintf(stderr,"Reading largeint: ");
#endif
- for (j=0; j<numbytes; j++)
- {
- int ilarge=j/4;
- int ibyte=j%4;
- largeint[ilarge]|=((unsigned int)input[j])<<(ibyte*8);
+ if (numbytes/4 < maxbasevals+1)
+ {
+ for (j=0; j<numbytes; j++)
+ {
+ int ilarge=j/4;
+ int ibyte=j%4;
+ largeint[ilarge]|=((unsigned int)input[j])<<(ibyte*8);
#ifdef SHOWIT
- fprintf(stderr,"%02x",(unsigned int)input[j]);
+ fprintf(stderr,"%02x",(unsigned int)input[j]);
#endif
+ }
}
#ifdef SHOWIT
- fprintf(stderr,"\n");
+ fprintf(stderr,"\n");
#endif
input+=numbytes;
/* Do the long division required to get the output values. */
n=maxbasevals;
if (n>nvals_left)
n=nvals_left;
- numints=(numbytes+3)/4;
for (i=n-1; i>=0; i--)
{
output[outvals+i*3]=Ptngc_largeint_div(base,largeint,largeint_tmp,maxbasevals+1);
@@ -1699,21 +1699,21 @@ static void unpack_one_large(struct xtc3_context *xtc3_context,
int natoms, int current_large_type)
{
int large_ints[3]={0,0,0};
- if (current_large_type==0)
+ if (current_large_type==0 && xtc3_context->large_direct)
{
large_ints[0]=(int)xtc3_context->large_direct[(*ilargedir)]+minint[0];
large_ints[1]=(int)xtc3_context->large_direct[(*ilargedir)+1]+minint[1];
large_ints[2]=(int)xtc3_context->large_direct[(*ilargedir)+2]+minint[2];
(*ilargedir)+=3;
}
- else if (current_large_type==1)
+ else if (current_large_type==1 && xtc3_context->large_intra_delta)
{
large_ints[0]=unpositive_int(xtc3_context->large_intra_delta[(*ilargeintra)])+prevcoord[0];
large_ints[1]=unpositive_int(xtc3_context->large_intra_delta[(*ilargeintra)+1])+prevcoord[1];
large_ints[2]=unpositive_int(xtc3_context->large_intra_delta[(*ilargeintra)+2])+prevcoord[2];
(*ilargeintra)+=3;
}
- else
+ else if (xtc3_context->large_inter_delta)
{
large_ints[0]=unpositive_int(xtc3_context->large_inter_delta[(*ilargeinter)])
+output[outdata-natoms*3+didswap*3];
@@ -1838,7 +1838,7 @@ int Ptngc_unpack_array_xtc3(unsigned char *packed,int *output, int length, int n
prevcoord[1]=minint[1];
prevcoord[2]=minint[2];
- while (ntriplets_left>0)
+ while (ntriplets_left>0 && iinstr<xtc3_context.ninstr)
{
int instr=xtc3_context.instructions[iinstr++];
#ifdef SHOWIT
@@ -1895,7 +1895,7 @@ int Ptngc_unpack_array_xtc3(unsigned char *packed,int *output, int length, int n
outdata+=runlength*3;
}
}
- else if (instr==INSTR_LARGE_RLE)
+ else if (instr==INSTR_LARGE_RLE && irle<xtc3_context.nrle)
{
int large_rle=xtc3_context.rle[irle++];
#ifdef SHOWIT
@@ -1910,7 +1910,7 @@ int Ptngc_unpack_array_xtc3(unsigned char *packed,int *output, int length, int n
outdata+=3;
}
}
- else if (instr==INSTR_SMALL_RUNLENGTH)
+ else if (instr==INSTR_SMALL_RUNLENGTH && irle<xtc3_context.nrle)
{
runlength=xtc3_context.rle[irle++];
#ifdef SHOWIT
diff --git a/src/lib/tng_io.c b/src/lib/tng_io.c
index 9acd857..725e4b6 100644
--- a/src/lib/tng_io.c
+++ b/src/lib/tng_io.c
@@ -965,7 +965,6 @@ static tng_function_status tng_block_header_read
__FILE__, __LINE__);
}
}
- offset += sizeof(block->block_version);
return(TNG_SUCCESS);
}
@@ -992,7 +991,7 @@ static tng_function_status tng_block_header_read
// __FILE__, __LINE__);
// return(TNG_CRITICAL);
// }
-//
+//
// if(!block->block_contents)
// {
// printf("No block data to write. %s: %d\n",
@@ -1133,7 +1132,6 @@ static tng_function_status tng_block_header_write
__FILE__, __LINE__);
}
}
- offset += sizeof(block->block_version);
if(fwrite(block->header_contents, block->header_contents_size,
1, tng_data->output_file) != 1)
@@ -2362,7 +2360,8 @@ static tng_function_status tng_molecules_block_read
tng_chain_data_read(tng_data, block, chain, &offset);
- chain->residues = residue;
+ chain->residues = molecule->residues;
+ residue = chain->residues;
/* Read the residues of the chain */
for(k=chain->n_residues; k--;)
@@ -2436,64 +2435,67 @@ static tng_function_status tng_molecules_block_read
}
offset += sizeof(molecule->n_bonds);
- molecule->bonds = malloc(molecule->n_bonds *
- sizeof(struct tng_bond));
- if(!molecule->bonds)
+ if(molecule->n_bonds > 0)
{
- printf("Cannot allocate memory (%"PRId64" bytes). %s: %d\n",
- molecule->n_bonds * sizeof(struct tng_bond),
- __FILE__, __LINE__);
- if(molecule->chains)
- {
- free(molecule->chains);
- molecule->chains = 0;
- }
- if(molecule->residues)
- {
- free(molecule->residues);
- molecule->residues = 0;
- }
- if(molecule->bonds)
+ molecule->bonds = malloc(molecule->n_bonds *
+ sizeof(struct tng_bond));
+ if(!molecule->bonds)
{
- free(molecule->bonds);
- molecule->bonds = 0;
+ printf("Cannot allocate memory (%"PRId64" bytes). %s: %d\n",
+ molecule->n_bonds * sizeof(struct tng_bond),
+ __FILE__, __LINE__);
+ if(molecule->chains)
+ {
+ free(molecule->chains);
+ molecule->chains = 0;
+ }
+ if(molecule->residues)
+ {
+ free(molecule->residues);
+ molecule->residues = 0;
+ }
+ if(molecule->atoms)
+ {
+ free(molecule->atoms);
+ molecule->atoms = 0;
+ }
+ return(TNG_CRITICAL);
}
- return(TNG_CRITICAL);
- }
- bond = molecule->bonds;
+ bond = molecule->bonds;
- for(j=molecule->n_bonds; j--;)
- {
- memcpy(&bond->from_atom_id, block->block_contents+offset,
- sizeof(bond->from_atom_id));
- if(tng_data->input_endianness_swap_func_64)
+ for(j=molecule->n_bonds; j--;)
{
- if(tng_data->input_endianness_swap_func_64(tng_data,
- &bond->from_atom_id)
- != TNG_SUCCESS)
+ memcpy(&bond->from_atom_id, block->block_contents+offset,
+ sizeof(bond->from_atom_id));
+ if(tng_data->input_endianness_swap_func_64)
{
- printf("Cannot swap byte order. %s: %d\n",
- __FILE__, __LINE__);
+ if(tng_data->input_endianness_swap_func_64(tng_data,
+ &bond->from_atom_id)
+ != TNG_SUCCESS)
+ {
+ printf("Cannot swap byte order. %s: %d\n",
+ __FILE__, __LINE__);
+ }
}
- }
- offset += sizeof(bond->from_atom_id);
+ offset += sizeof(bond->from_atom_id);
- memcpy(&bond->to_atom_id, block->block_contents+offset,
- sizeof(bond->to_atom_id));
- if(tng_data->input_endianness_swap_func_64)
- {
- if(tng_data->input_endianness_swap_func_64(tng_data,
- &bond->to_atom_id)
- != TNG_SUCCESS)
+ memcpy(&bond->to_atom_id, block->block_contents+offset,
+ sizeof(bond->to_atom_id));
+ if(tng_data->input_endianness_swap_func_64)
{
- printf("Cannot swap byte order. %s: %d\n",
- __FILE__, __LINE__);
+ if(tng_data->input_endianness_swap_func_64(tng_data,
+ &bond->to_atom_id)
+ != TNG_SUCCESS)
+ {
+ printf("Cannot swap byte order. %s: %d\n",
+ __FILE__, __LINE__);
+ }
}
- }
- offset += sizeof(bond->to_atom_id);
+ offset += sizeof(bond->to_atom_id);
- bond++;
+ bond++;
+ }
}
}
@@ -2621,7 +2623,6 @@ static tng_function_status tng_molecules_block_write
atom++;
}
- bond = molecule->bonds;
for(j = molecule->n_bonds; j--;)
{
len += sizeof(bond->from_atom_id) + sizeof(bond->to_atom_id);
@@ -3191,7 +3192,6 @@ static tng_function_status tng_frame_set_block_read
__FILE__, __LINE__);
}
}
- offset += sizeof(tng_data->time_per_frame);
}
else
{
@@ -3429,7 +3429,6 @@ static tng_function_status tng_frame_set_block_write
__FILE__, __LINE__);
}
}
- offset += sizeof(tng_data->time_per_frame);
if(tng_block_header_write(tng_data, block, hash_mode) != TNG_SUCCESS)
{
@@ -3760,8 +3759,6 @@ static tng_function_status tng_particle_data_block_create
return(TNG_CRITICAL);
}
frame_set->tr_particle_data = data;
- data = &frame_set->tr_particle_data[frame_set->
- n_particle_data_blocks - 1];
}
else
{
@@ -3779,8 +3776,6 @@ static tng_function_status tng_particle_data_block_create
return(TNG_CRITICAL);
}
tng_data->non_tr_particle_data = data;
- data = &tng_data->non_tr_particle_data[tng_data->
- n_particle_data_blocks - 1];
}
return(TNG_SUCCESS);
@@ -4209,6 +4204,11 @@ static tng_function_status tng_allocate_particle_data_mem
void ***values;
int64_t i, j, k, size, frame_alloc;
+ if(n_particles == 0 || n_values_per_frame == 0)
+ {
+ return(TNG_FAILURE);
+ }
+
if(data->strings && data->datatype == TNG_CHAR_DATA)
{
for(i = data->n_frames; i--;)
@@ -5208,7 +5208,6 @@ static tng_function_status tng_data_block_create
return(TNG_CRITICAL);
}
frame_set->tr_data = data;
- data = &frame_set->tr_data[frame_set->n_data_blocks - 1];
}
else
{
@@ -5224,7 +5223,6 @@ static tng_function_status tng_data_block_create
return(TNG_CRITICAL);
}
tng_data->non_tr_data = data;
- data = &tng_data->non_tr_data[tng_data->n_data_blocks - 1];
}
return(TNG_SUCCESS);
@@ -6888,6 +6886,7 @@ tng_function_status DECLSPECDLLEXPORT tng_molecule_w_id_add
sizeof(int64_t) * (tng_data->n_molecules + 1),
__FILE__, __LINE__);
free(tng_data->molecule_cnt_list);
+ free(new_molecules);
return(TNG_CRITICAL);
}
@@ -7206,7 +7205,7 @@ tng_function_status DECLSPECDLLEXPORT tng_chain_residue_w_id_add
int curr_index;
tng_residue_t new_residues, temp_residue, last_residue;
tng_molecule_t molecule = chain->molecule;
- tng_function_status stat;
+ tng_function_status stat = TNG_SUCCESS;
if(chain->n_residues)
{
@@ -7342,7 +7341,7 @@ tng_function_status DECLSPECDLLEXPORT tng_residue_atom_w_id_add
int64_t i;
tng_atom_t new_atoms;
tng_molecule_t molecule = residue->chain->molecule;
- tng_function_status stat;
+ tng_function_status stat = TNG_SUCCESS;
if(!residue->n_atoms)
{
@@ -9320,8 +9319,6 @@ tng_function_status DECLSPECDLLEXPORT tng_num_frame_sets_get
++cnt;
- file_pos = tng_data->current_trajectory_frame_set_input_file_pos;
-
long_stride_length = tng_data->long_stride_length;
medium_stride_length = tng_data->medium_stride_length;
@@ -9444,6 +9441,11 @@ tng_function_status DECLSPECDLLEXPORT tng_frame_set_nr_find
stat = tng_num_frame_sets_get(tng_data, &n_frame_sets);
+ if(stat != TNG_SUCCESS)
+ {
+ return(stat);
+ }
+
if(nr >= n_frame_sets)
{
return(TNG_FAILURE);
@@ -11972,6 +11974,11 @@ static tng_function_status tng_data_values_alloc
int64_t i;
tng_function_status stat;
+ if(n_frames == 0 || n_values_per_frame == 0)
+ {
+ return(TNG_FAILURE);
+ }
+
if(*values)
{
stat = tng_data_values_free(tng_data, *values, n_frames,
@@ -11984,7 +11991,7 @@ static tng_function_status tng_data_values_alloc
return(stat);
}
}
- *values = malloc(sizeof(union data_values **) * n_frames);
+ *values = malloc(sizeof(union data_values *) * n_frames);
if(!*values)
{
printf("Cannot allocate memory (%"PRId64" bytes). %s: %d\n",
@@ -12060,6 +12067,11 @@ static tng_function_status tng_particle_data_values_alloc
int64_t i, j;
tng_function_status stat;
+ if(n_particles == 0 || n_values_per_frame == 0)
+ {
+ return(TNG_FAILURE);
+ }
+
if(*values)
{
stat = tng_particle_data_values_free(tng_data, *values, n_frames,
@@ -12387,6 +12399,7 @@ tng_function_status tng_data_vector_get(tng_trajectory_t tng_data,
printf("Cannot allocate memory (%"PRId64" bytes). %s: %d\n",
data_size, __FILE__, __LINE__);
free(*values);
+ *values = 0;
return(TNG_CRITICAL);
}
@@ -12694,6 +12707,12 @@ tng_function_status DECLSPECDLLEXPORT tng_data_vector_interval_get
stat = tng_frame_set_read_next(tng_data, hash_mode);
if(stat != TNG_SUCCESS)
{
+ if(current_values)
+ {
+ free(current_values);
+ }
+ free(*values);
+ *values = 0;
return(stat);
}
@@ -12729,7 +12748,10 @@ tng_function_status DECLSPECDLLEXPORT tng_data_vector_interval_get
}
}
- free(current_values);
+ if(current_values)
+ {
+ free(current_values);
+ }
return(TNG_SUCCESS);
}
@@ -13361,7 +13383,7 @@ tng_function_status DECLSPECDLLEXPORT tng_particle_data_vector_interval_get
&n_frames, stride_length, n_particles,
n_values_per_frame, type);
- if(stat != TNG_SUCCESS)
+ if(stat != TNG_SUCCESS || *n_particles == 0)
{
if(current_values)
{
@@ -13443,6 +13465,12 @@ tng_function_status DECLSPECDLLEXPORT tng_particle_data_vector_interval_get
stat = tng_frame_set_read_next(tng_data, hash_mode);
if(stat != TNG_SUCCESS)
{
+ if(current_values)
+ {
+ free(current_values);
+ }
+ free(*values);
+ *values = 0;
return(stat);
}
@@ -13478,7 +13506,10 @@ tng_function_status DECLSPECDLLEXPORT tng_particle_data_vector_interval_get
}
}
- free(current_values);
+ if(current_values)
+ {
+ free(current_values);
+ }
return(TNG_SUCCESS);
}
@@ -13553,14 +13584,14 @@ tng_function_status DECLSPECDLLEXPORT tng_util_trajectory_close
tng_function_status DECLSPECDLLEXPORT tng_util_trajectory_molecules_get
(tng_trajectory_t tng_data,
int64_t *n_mols,
- int64_t *molecule_cnt_list,
+ int64_t **molecule_cnt_list,
tng_molecule_t *mols)
{
/* FIXME: This should return a copy of the molecules instead */
*n_mols = tng_data->n_molecules;
- molecule_cnt_list = tng_data->molecule_cnt_list;
- mols = &tng_data->molecules;
+ *molecule_cnt_list = tng_data->molecule_cnt_list;
+ *mols = tng_data->molecules;
return(TNG_SUCCESS);
}
@@ -13601,10 +13632,10 @@ tng_function_status DECLSPECDLLEXPORT tng_util_molecule_particles_get
*n_particles = mol->n_atoms;
- **names = malloc(sizeof(char *) * *n_particles);
- **types = malloc(sizeof(char *) * *n_particles);
- **res_names = malloc(sizeof(char *) * *n_particles);
- **chain_names = malloc(sizeof(char *) * *n_particles);
+ *names = malloc(sizeof(char *) * *n_particles);
+ *types = malloc(sizeof(char *) * *n_particles);
+ *res_names = malloc(sizeof(char *) * *n_particles);
+ *chain_names = malloc(sizeof(char *) * *n_particles);
*res_ids = malloc(sizeof(int64_t) * *n_particles);
*chain_ids = malloc(sizeof(int64_t) * *n_particles);
@@ -13613,16 +13644,16 @@ tng_function_status DECLSPECDLLEXPORT tng_util_molecule_particles_get
atom = &mol->atoms[i];
res = atom->residue;
chain = res->chain;
- *names[i] = malloc(strlen(atom->name));
+ (*names)[i] = malloc(strlen(atom->name));
strcpy(*names[i], atom->name);
- *types[i] = malloc(strlen(atom->atom_type));
+ (*types)[i] = malloc(strlen(atom->atom_type));
strcpy(*types[i], atom->atom_type);
- *res_names[i] = malloc(strlen(res->name));
+ (*res_names)[i] = malloc(strlen(res->name));
strcpy(*res_names[i], res->name);
- *chain_names[i] = malloc(strlen(chain->name));
+ (*chain_names)[i] = malloc(strlen(chain->name));
strcpy(*chain_names[i], chain->name);
- *res_ids[i] = res->id;
- *chain_ids[i] = chain->id;
+ (*res_ids)[i] = res->id;
+ (*chain_ids)[i] = chain->id;
}
return(TNG_SUCCESS);
@@ -13923,6 +13954,12 @@ tng_function_status DECLSPECDLLEXPORT tng_util_generic_write_frequency_set
stat = tng_allocate_particle_data_mem(tng_data, p_data, n_frames,
f, n_particles,
n_values_per_frame);
+ if(stat != TNG_SUCCESS)
+ {
+ printf("Error allocating particle data memory. %s: %d\n",
+ __FILE__, __LINE__);
+ return(stat);
+ }
}
else
{
@@ -13947,6 +13984,12 @@ tng_function_status DECLSPECDLLEXPORT tng_util_generic_write_frequency_set
n_data_blocks - 1];
stat = tng_allocate_data_mem(tng_data, np_data, n_frames,
f, n_values_per_frame);
+ if(stat != TNG_SUCCESS)
+ {
+ printf("Error allocating particle data memory. %s: %d\n",
+ __FILE__, __LINE__);
+ return(stat);
+ }
}
else
{
@@ -14098,6 +14141,12 @@ tng_function_status DECLSPECDLLEXPORT tng_util_generic_write
stat = tng_allocate_particle_data_mem(tng_data, p_data, n_frames,
stride_length, n_particles,
n_values_per_frame);
+ if(stat != TNG_SUCCESS)
+ {
+ printf("Error allocating particle data memory. %s: %d\n",
+ __FILE__, __LINE__);
+ return(stat);
+ }
}
if(block_type_flag == TNG_TRAJECTORY_BLOCK)
@@ -14142,6 +14191,12 @@ tng_function_status DECLSPECDLLEXPORT tng_util_generic_write
}
stat = tng_allocate_data_mem(tng_data, np_data, n_frames,
stride_length, n_values_per_frame);
+ if(stat != TNG_SUCCESS)
+ {
+ printf("Error allocating particle data memory. %s: %d\n",
+ __FILE__, __LINE__);
+ return(stat);
+ }
}
if(block_type_flag == TNG_TRAJECTORY_BLOCK)
diff --git a/src/tests/tng_io_testing.c b/src/tests/tng_io_testing.c
index 96b53f7..0f67ccc 100644
--- a/src/tests/tng_io_testing.c
+++ b/src/tests/tng_io_testing.c
@@ -136,7 +136,7 @@ static tng_function_status tng_test_read_and_write_file
}
stat = tng_frame_set_write(traj, TNG_USE_HASH);
}
-
+
return(stat);
}
@@ -148,11 +148,11 @@ static tng_function_status tng_test_write_and_read_traj(tng_trajectory_t traj)
int64_t frame_nr;
double box_shape[9];
char atom_type[16], annotation[128];
- tng_function_status stat;
+ tng_function_status stat = TNG_SUCCESS;
tng_medium_stride_length_set(traj, 10);
tng_long_stride_length_set(traj, 100);
-
+
/* Create molecules */
if(tng_test_setup_molecules(traj) == TNG_CRITICAL)
{
@@ -200,7 +200,7 @@ static tng_function_status tng_test_write_and_read_traj(tng_trajectory_t traj)
printf("Failed setting partial charges.\n");
return(TNG_CRITICAL);
}
-
+
stat = tng_particle_data_block_add(traj, TNG_TRAJ_PARTIAL_CHARGES, "PARTIAL CHARGES",
TNG_FLOAT_DATA, TNG_NON_TRAJECTORY_BLOCK,
1, 1, 1, 0, n_particles,
@@ -230,9 +230,9 @@ static tng_function_status tng_test_write_and_read_traj(tng_trajectory_t traj)
printf("Cannot write file headers.\n");
}
-
+
tng_num_frames_per_frame_set_get(traj, &n_frames_per_frame_set);
-
+
data = malloc(sizeof(float) * n_particles *
n_frames_per_frame_set * 3);
if(!data)
@@ -366,9 +366,9 @@ static tng_function_status tng_test_write_and_read_traj(tng_trajectory_t traj)
return(TNG_CRITICAL);
}
}
-
+
/* Write two more frame sets one frame at a time */
- cnt = 0;
+
/* Make a new frame set - if always using the same mapping blocks
* it is not necessary to explicitly add a new frame set - it will
* be added automatically when adding data for a frame */
@@ -383,7 +383,7 @@ static tng_function_status tng_test_write_and_read_traj(tng_trajectory_t traj)
}
frame_nr = i * n_frames_per_frame_set;
-
+
for(k=0; k<300; k++)
{
mapping[k]=k;
@@ -470,7 +470,7 @@ static tng_function_status tng_test_write_and_read_traj(tng_trajectory_t traj)
}
}
}
-
+
free(molpos);
free(data);
@@ -481,7 +481,7 @@ static tng_function_status tng_test_write_and_read_traj(tng_trajectory_t traj)
#endif
stat = tng_file_headers_read(traj, TNG_SKIP_HASH);
-
+
while(stat == TNG_SUCCESS)
{
stat = tng_frame_set_read_next(traj, TNG_SKIP_HASH);
@@ -582,7 +582,7 @@ int main()
tng_trajectory_t traj;
tng_function_status stat;
char time_str[TNG_MAX_DATE_STR_LEN];
-
+
if(tng_trajectory_init(&traj) != TNG_SUCCESS)
{
tng_trajectory_destroy(&traj);
@@ -614,7 +614,7 @@ int main()
{
printf("Test Read and write file:\t\t\tSucceeded.\n");
}
-
+
if(tng_test_get_box_data(traj) != TNG_SUCCESS)
{
printf("Test Get data:\t\t\t\t\tFailed. %s: %d\n",
@@ -635,7 +635,7 @@ int main()
{
printf("Test Destroy and init trajectory:\t\tSucceeded.\n");
}
-
+
#ifdef EXAMPLE_FILES_DIR
tng_output_file_set(traj, EXAMPLE_FILES_DIR "tng_test.tng");
@@ -652,7 +652,7 @@ int main()
{
printf("Test Write and read file:\t\t\tSucceeded.\n");
}
-
+
if(tng_test_get_positions_data(traj) != TNG_SUCCESS)
{
printf("Test Get particle data:\t\t\t\tFailed. %s: %d\n",
@@ -674,7 +674,7 @@ int main()
printf("Test Destroy trajectory:\t\t\tSucceeded.\n");
}
-
+
#ifdef EXAMPLE_FILES_DIR
stat = tng_util_trajectory_open(EXAMPLE_FILES_DIR "tng_test.tng", 'r', &traj);
#else
@@ -691,7 +691,7 @@ int main()
printf("Test Utility function open:\t\t\tSucceeded.\n");
}
- stat = tng_util_trajectory_close(&traj);
+ stat = tng_util_trajectory_close(&traj);
if(stat != TNG_SUCCESS)
{
printf("Test Utility function close:\t\t\tFailed. %s: %d.\n",
@@ -702,8 +702,8 @@ int main()
{
printf("Test Utility function close:\t\t\tSucceeded.\n");
}
-
+
printf("Tests finished\n");
-
+
exit(0);
} \ No newline at end of file
diff --git a/src/tests/tng_parallel_read.c b/src/tests/tng_parallel_read.c
index d0aaa15..58863b0 100644
--- a/src/tests/tng_parallel_read.c
+++ b/src/tests/tng_parallel_read.c
@@ -21,7 +21,7 @@
/* N.B. this code is for testing parallel reading of trajectory frame sets. The
* performance is not improved very much and is to a large extent limited by
* disk i/o. It can however be used as inspiration for writing parallel code
- * using the TNG library. */
+ * using the TNG library. The code is NOT fully tested and may behave strangely. */
int main(int argc, char **argv)
{
@@ -29,7 +29,7 @@ int main(int argc, char **argv)
union data_values ***local_positions = 0; // A 3-dimensional array to be populated
union data_values **particle_pos = 0;
int64_t n_particles, n_values_per_frame, n_frame_sets, n_frames;
- int64_t n_frames_per_frame_set;
+ int64_t n_frames_per_frame_set, tot_n_frames = 0;
char data_type;
int i, j;
int64_t particle = 0, local_first_frame, local_last_frame;
@@ -66,7 +66,7 @@ int main(int argc, char **argv)
tng_num_frame_sets_get(traj, &n_frame_sets);
tng_num_frames_per_frame_set_get(traj, &n_frames_per_frame_set);
- particle_pos = malloc(sizeof(union data_values **) * n_frame_sets *
+ particle_pos = malloc(sizeof(union data_values *) * n_frame_sets *
n_frames_per_frame_set);
for(i = n_frame_sets * n_frames_per_frame_set; i--;)
{
@@ -94,7 +94,7 @@ int main(int argc, char **argv)
private (n_frames, n_particles, n_values_per_frame, \
local_first_frame, local_last_frame, j) \
firstprivate (local_traj, local_positions, frame_set)\
-shared(data_type, traj, n_frame_sets, particle_pos, particle, i)\
+shared(data_type, traj, n_frame_sets, particle_pos, particle, i, tot_n_frames)\
default(none)
{
/* Each tng_trajectory_t keeps its own file pointers and i/o positions.
@@ -106,17 +106,22 @@ default(none)
if(tng_frame_set_nr_find(local_traj, i) != TNG_SUCCESS)
{
printf("FAILED finding frame set %d!\n", i);
+ tot_n_frames = 0;
+ break;
}
if(tng_particle_data_get(local_traj, TNG_TRAJ_POSITIONS, &local_positions,
&n_frames, &n_particles, &n_values_per_frame,
&data_type) != TNG_SUCCESS)
{
printf("FAILED getting particle data\n");
+ tot_n_frames = 0;
+ break;
}
tng_current_frame_set_get(local_traj, &frame_set);
tng_frame_set_frame_range_get(local_traj, frame_set, &local_first_frame, &local_last_frame);
// printf("Frame %"PRId64"-%"PRId64":\n", local_first_frame, local_last_frame);
// printf("%"PRId64" %"PRId64" %"PRId64"\n", n_frames, n_particles, n_values_per_frame);
+ tot_n_frames += n_frames;
for(j = 0; j < n_frames; j++)
{
particle_pos[local_first_frame + j][0] = local_positions[j][particle][0];
@@ -136,24 +141,24 @@ default(none)
switch(data_type)
{
case TNG_INT_DATA:
- for(i = 0; i < n_frame_sets * n_frames_per_frame_set; i++)
+ for(j = 0; j < tot_n_frames; j++)
{
- printf("\t%"PRId64"\t%"PRId64"\t%"PRId64"\n", particle_pos[i][0].i,
- particle_pos[i][1].i, particle_pos[i][2].i);
+ printf("\t%"PRId64"\t%"PRId64"\t%"PRId64"\n", particle_pos[j][0].i,
+ particle_pos[j][1].i, particle_pos[j][2].i);
}
break;
case TNG_FLOAT_DATA:
- for(i = 0; i < n_frame_sets * n_frames_per_frame_set; i++)
+ for(j = 0; j < tot_n_frames; j++)
{
- printf("\t%f\t%f\t%f\n", particle_pos[i][0].f,
- particle_pos[i][1].f, particle_pos[i][2].f);
+ printf("\t%f\t%f\t%f\n", particle_pos[j][0].f,
+ particle_pos[j][1].f, particle_pos[j][2].f);
}
break;
case TNG_DOUBLE_DATA:
- for(i = 0; i < n_frame_sets * n_frames_per_frame_set; i++)
+ for(j = 0; j < tot_n_frames; j++)
{
- printf("\t%f\t%f\t%f\n", particle_pos[i][0].d,
- particle_pos[i][1].d, particle_pos[i][2].d);
+ printf("\t%f\t%f\t%f\n", particle_pos[j][0].d,
+ particle_pos[j][1].d, particle_pos[j][2].d);
}
break;
default:
contact: Jan Huwald // Impressum