public
Created

My first Open-MPI program: reads the reads in one or more BAM files

  • Download Gist
Makefile
Makefile
1 2 3 4 5 6 7 8 9 10
SAMDIR=/commun/data/packages/samtools-0.1.18
OPENMPI.dir=/commun/data/packages/openmpi
CC=$(OPENMPI.dir)/bin/ortecc
.PHONY: all
all: ./ompibam
 
./ompibam: ompibam.c
$(CC) -o $@ -O3 -I$(SAMDIR) \
-Wl,-rpath -Wl,$(OPENMPI.dir)/lib \
-L$(SAMDIR) -Wall $< -lmpi -lbam -lz
ompibam.c
C
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236
/**
* Author : Pierre Lindenbaum PhD
* Date : 2013
* WWW: http://plindenbaum.blogspot.com
* Motivation:
* my first Open-MPI source. Takes a list of BAM files
* and count the number of reads.
*/
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <mpi.h>
#include "sam.h"
 
/** message send to the node: "count the read in a bam" */
#define TAG_COUNT_BAM 1
/** message send to the node: "shutdown" */
#define TAG_RELEASE_SLAVE 2
 
 
/**
* fixed-size structure send to the nodes. contains the BAMs and the counts of reads
*/
typedef struct Params_t
{
char filename[FILENAME_MAX];
int is_error;
long count_total;
long count_mapped;
long count_properly_paired;
long count_dup;
}Params;
 
#define IS_FLAG_SET(b,f) ((b->core.flag & (f))!=0)
 
 
/**
* count the reads in a BAM file using the SAMTOOLS API
* fills the structure Params*
*/
static void count_reads(Params* p)
{
int ret;
 
bam1_t* b;
bamFile fp=bam_open(p->filename, "r") ;
bam_header_t *header;
if(fp==0)
{
fprintf(stderr,"Cannot open %s\n",p->filename);
p->is_error=1;
return;
}
header= bam_header_read(fp);
b=bam_init1();
while ((ret = bam_read1(fp, b)) >= 0)
{
p->count_total++;
if(IS_FLAG_SET(b,BAM_FUNMAP)) continue;
p->count_mapped++;
if(IS_FLAG_SET(b,BAM_FPROPER_PAIR))
{
p->count_properly_paired++;
}
if(IS_FLAG_SET(b,BAM_FDUP))
{
p->count_dup++;
}
}
bam_destroy1(b);
bam_header_destroy(header);
bam_close(fp);
}
 
/**
* This is the slave method, it receive a message 'Params' from the master method
* if the status is TAG_RELEASE_SLAVE the slave is released
* else the bam in param.filename is analyszed and the message is sent back to the master
*/
static void slave()
{
for(;;)
{
MPI_Status status;
Params param;
/* Receive a message from the master */
MPI_Recv(
&param,
sizeof(Params),
MPI_CHAR,
0,
MPI_ANY_TAG,
MPI_COMM_WORLD,
&status);
if(status.MPI_TAG==TAG_RELEASE_SLAVE)
{
break;
}
count_reads(&param);
/* Send the result back */
MPI_Send(
&param,
sizeof(Params),
MPI_CHAR,
0,
0,
MPI_COMM_WORLD
);
}
}
 
/**
* this is the master method, while there is
* a bam in the queue, send the analysis to the slaves.
* At the end, close the slaves
*/
static void master(const int n_bams,char** argv)
{
int rank=0;
int optind=0;
int proc_count=0;
/* get number of available processor proc_count */
if(MPI_Comm_size(MPI_COMM_WORLD, &proc_count)!=MPI_SUCCESS )
{
return;
}
printf("#filename\ttotal\tmapped\tproperly_paired\tdup\n");
/* loop over the bam files */
while(optind< n_bams)
{
int i=0;
/* number of message sent */
int count_sent=0;
/* loop over the available procs */
for (rank = 1;
rank <proc_count && optind+count_sent < n_bams;
++rank)
{
/* create a new param */
Params param;
/** fill to '0' */
memset((void*)&param,0,sizeof(Params));
/* set the BAM filename */
strcpy(param.filename,argv[count_sent+optind]);
/* Send it to each rank */
MPI_Send(
&param, /* message buffer */
sizeof(Params), /* one data item */
MPI_CHAR, /* data item is an array of bytes */
rank, /* destination process rank */
TAG_COUNT_BAM, /* user chosen message tag */
MPI_COMM_WORLD /* default communicator */
);
++count_sent;
}
/* now we get the results */
for(i=0;i< count_sent;++i)
{
/* get a message from the slave */
Params param;
MPI_Status status;
MPI_Recv(&param, /* message buffer */
sizeof(Params), /* one data item */
MPI_CHAR, /* data item is an array of bytes */
MPI_ANY_SOURCE, /* receive from any sender */
MPI_ANY_TAG, /* any type of message */
MPI_COMM_WORLD, /* default communicator */
&status); /* info about the received message */
if(param.is_error!=0)
{
fprintf(stderr,"Error for %s\n",param.filename);
}
else /* print the result */
{
printf("%s\t%ld\t%ld\t%ld\t%ld\n",
param.filename,
param.count_total,
param.count_mapped,
param.count_properly_paired,
param.count_dup
);
}
}
if(count_sent==0)
{
break;
}
optind+=count_sent;
}
/* shutdown the slave */
for (rank = 1; rank < proc_count; ++rank)
{
Params empty;
MPI_Send(
&empty, /* message buffer */
sizeof(Params), /* one data item */
MPI_CHAR, /* data item is an array of bytes */
rank, /* destination process rank */
TAG_RELEASE_SLAVE, /* SHutdown */
MPI_COMM_WORLD
); /* default communicator */
}
}
 
int main(int argc, char *argv[])
{
int rank=0;
/* Initialisation */
MPI_Init(&argc, &argv);
/* Find out my identity in the default communicator */
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
/* rand is zero = I'm the master */
if (rank == 0)
{
master(argc-1,&argv[1]);
}
else /* I'm a slave */
{
slave(argc-1,&argv[1]);
}
/* Shut down MPI */
MPI_Finalize();
return EXIT_SUCCESS;
}

Please sign in to comment on this gist.

Something went wrong with that request. Please try again.