960615aac8e418fca4bc1a55ba90a49f6d51795c
[jalview.git] / src / jalview / io / vcf / VCFLoader.java
1 package jalview.io.vcf;
2
3 import htsjdk.samtools.util.CloseableIterator;
4 import htsjdk.variant.variantcontext.Allele;
5 import htsjdk.variant.variantcontext.VariantContext;
6 import htsjdk.variant.vcf.VCFHeader;
7 import htsjdk.variant.vcf.VCFHeaderLine;
8 import htsjdk.variant.vcf.VCFHeaderLineCount;
9 import htsjdk.variant.vcf.VCFHeaderLineType;
10 import htsjdk.variant.vcf.VCFInfoHeaderLine;
11
12 import jalview.analysis.AlignmentUtils;
13 import jalview.analysis.Dna;
14 import jalview.api.AlignViewControllerGuiI;
15 import jalview.bin.Cache;
16 import jalview.datamodel.AlignmentI;
17 import jalview.datamodel.DBRefEntry;
18 import jalview.datamodel.GeneLociI;
19 import jalview.datamodel.Mapping;
20 import jalview.datamodel.SequenceFeature;
21 import jalview.datamodel.SequenceI;
22 import jalview.datamodel.features.FeatureAttributeType;
23 import jalview.datamodel.features.FeatureSource;
24 import jalview.datamodel.features.FeatureSources;
25 import jalview.ext.ensembl.EnsemblMap;
26 import jalview.ext.htsjdk.VCFReader;
27 import jalview.io.gff.Gff3Helper;
28 import jalview.io.gff.SequenceOntologyI;
29 import jalview.util.MapList;
30 import jalview.util.MappingUtils;
31 import jalview.util.MessageManager;
32
33 import java.io.IOException;
34 import java.util.ArrayList;
35 import java.util.HashMap;
36 import java.util.List;
37 import java.util.Map;
38 import java.util.Map.Entry;
39 import java.util.regex.Pattern;
40
41 /**
42  * A class to read VCF data (using the htsjdk) and add variants as sequence
43  * features on dna and any related protein product sequences
44  * 
45  * @author gmcarstairs
46  */
47 public class VCFLoader
48 {
49   /*
50    * Lookup keys, and default values, for Preference entries that describe
51    * patterns for VCF and VEP fields to capture 
52    */
53   private static final String VEP_FIELDS_PREF = "VEP_FIELDS";
54
55   private static final String VCF_FIELDS_PREF = "VCF_FIELDS";
56
57   private static final String DEFAULT_VCF_FIELDS = "AF,AC*";
58
59   private static final String DEFAULT_VEP_FIELDS = "Allele,Consequence,IMPACT,SWISSPROT,SIFT,PolyPhen,CLIN_SIG";
60
61   /*
62    * keys to fields of VEP CSQ consequence data
63    * see https://www.ensembl.org/info/docs/tools/vep/vep_formats.html
64    */
65   private static final String ALLELE_KEY = "Allele";
66
67   private static final String ALLELE_NUM_KEY = "ALLELE_NUM"; // 0 (ref), 1...
68   private static final String FEATURE_KEY = "Feature"; // Ensembl stable id
69
70   /*
71    * default VCF INFO key for VEP consequence data
72    * NB this can be overridden running VEP with --vcf_info_field
73    * - we don't handle this case (require identifier to be CSQ)
74    */
75   private static final String CSQ_FIELD = "CSQ";
76
77   /*
78    * separator for fields in consequence data is '|'
79    */
80   private static final String PIPE_REGEX = "\\|";
81
82   /*
83    * key for Allele Frequency output by VEP
84    * see http://www.ensembl.org/info/docs/tools/vep/vep_formats.html
85    */
86   private static final String ALLELE_FREQUENCY_KEY = "AF";
87
88   /*
89    * delimiter that separates multiple consequence data blocks
90    */
91   private static final String COMMA = ",";
92
93   /*
94    * the feature group assigned to a VCF variant in Jalview
95    */
96   private static final String FEATURE_GROUP_VCF = "VCF";
97
98   /*
99    * internal delimiter used to build keys for assemblyMappings
100    * 
101    */
102   private static final String EXCL = "!";
103
104   /*
105    * the alignment we are associating VCF data with
106    */
107   private AlignmentI al;
108
109   /*
110    * mappings between VCF and sequence reference assembly regions, as 
111    * key = "species!chromosome!fromAssembly!toAssembly
112    * value = Map{fromRange, toRange}
113    */
114   private Map<String, Map<int[], int[]>> assemblyMappings;
115
116   /*
117    * holds details of the VCF header lines (metadata)
118    */
119   private VCFHeader header;
120
121   /*
122    * the position (0...) of field in each block of
123    * CSQ (consequence) data (if declared in the VCF INFO header for CSQ)
124    * see http://www.ensembl.org/info/docs/tools/vep/vep_formats.html
125    */
126   private int csqAlleleFieldIndex = -1;
127   private int csqAlleleNumberFieldIndex = -1;
128   private int csqFeatureFieldIndex = -1;
129
130   /*
131    * a unique identifier under which to save metadata about feature
132    * attributes (selected INFO field data)
133    */
134   private String sourceId;
135
136   /*
137    * The INFO IDs of data that is both present in the VCF file, and
138    * also matched by any filters for data of interest
139    */
140   List<String> vcfFieldsOfInterest;
141
142   /*
143    * The field offsets and identifiers for VEP (CSQ) data that is both present
144    * in the VCF file, and also matched by any filters for data of interest
145    * for example 0 -> Allele, 1 -> Consequence, ..., 36 -> SIFT, ...
146    */
147   Map<Integer, String> vepFieldsOfInterest;
148
149   /**
150    * Constructor given an alignment context
151    * 
152    * @param alignment
153    */
154   public VCFLoader(AlignmentI alignment)
155   {
156     al = alignment;
157
158     // map of species!chromosome!fromAssembly!toAssembly to {fromRange, toRange}
159     assemblyMappings = new HashMap<String, Map<int[], int[]>>();
160   }
161
162   /**
163    * Starts a new thread to query and load VCF variant data on to the alignment
164    * <p>
165    * This method is not thread safe - concurrent threads should use separate
166    * instances of this class.
167    * 
168    * @param filePath
169    * @param gui
170    */
171   public void loadVCF(final String filePath,
172           final AlignViewControllerGuiI gui)
173   {
174     if (gui != null)
175     {
176       gui.setStatus(MessageManager.getString("label.searching_vcf"));
177     }
178
179     new Thread()
180     {
181
182       @Override
183       public void run()
184       {
185         VCFLoader.this.doLoad(filePath, gui);
186       }
187
188     }.start();
189   }
190
191   /**
192    * Loads VCF on to an alignment - provided it can be related to one or more
193    * sequence's chromosomal coordinates
194    * 
195    * @param filePath
196    * @param gui
197    *          optional callback handler for messages
198    */
199   protected void doLoad(String filePath, AlignViewControllerGuiI gui)
200   {
201     VCFReader reader = null;
202     try
203     {
204       // long start = System.currentTimeMillis();
205       reader = new VCFReader(filePath);
206
207       header = reader.getFileHeader();
208
209       sourceId = filePath;
210
211       saveMetadata(sourceId);
212
213       /*
214        * get offset of CSQ ALLELE_NUM and Feature if declared
215        */
216       parseCsqHeader();
217
218       VCFHeaderLine ref = header
219               .getOtherHeaderLine(VCFHeader.REFERENCE_KEY);
220       String vcfAssembly = ref.getValue();
221
222       int varCount = 0;
223       int seqCount = 0;
224
225       /*
226        * query for VCF overlapping each sequence in turn
227        */
228       for (SequenceI seq : al.getSequences())
229       {
230         int added = loadSequenceVCF(seq, reader, vcfAssembly);
231         if (added > 0)
232         {
233           seqCount++;
234           varCount += added;
235           transferAddedFeatures(seq);
236         }
237       }
238       if (gui != null)
239       {
240         // long elapsed = System.currentTimeMillis() - start;
241         String msg = MessageManager.formatMessage("label.added_vcf",
242                 varCount, seqCount);
243         gui.setStatus(msg);
244         if (gui.getFeatureSettingsUI() != null)
245         {
246           gui.getFeatureSettingsUI().discoverAllFeatureData();
247         }
248       }
249     } catch (Throwable e)
250     {
251       System.err.println("Error processing VCF: " + e.getMessage());
252       e.printStackTrace();
253       if (gui != null)
254       {
255         gui.setStatus("Error occurred - see console for details");
256       }
257     } finally
258     {
259       if (reader != null)
260       {
261         try
262         {
263           reader.close();
264         } catch (IOException e)
265         {
266           // ignore
267         }
268       }
269     }
270   }
271
272   /**
273    * Reads metadata (such as INFO field descriptions and datatypes) and saves
274    * them for future reference
275    * 
276    * @param theSourceId
277    */
278   void saveMetadata(String theSourceId)
279   {
280     List<Pattern> vcfFieldPatterns = getFieldMatchers(VCF_FIELDS_PREF,
281             DEFAULT_VCF_FIELDS);
282     vcfFieldsOfInterest = new ArrayList<>();
283
284     FeatureSource metadata = new FeatureSource(theSourceId);
285
286     for (VCFInfoHeaderLine info : header.getInfoHeaderLines())
287     {
288       String attributeId = info.getID();
289       String desc = info.getDescription();
290       VCFHeaderLineType type = info.getType();
291       FeatureAttributeType attType = null;
292       switch (type)
293       {
294       case Character:
295         attType = FeatureAttributeType.Character;
296         break;
297       case Flag:
298         attType = FeatureAttributeType.Flag;
299         break;
300       case Float:
301         attType = FeatureAttributeType.Float;
302         break;
303       case Integer:
304         attType = FeatureAttributeType.Integer;
305         break;
306       case String:
307         attType = FeatureAttributeType.String;
308         break;
309       }
310       metadata.setAttributeName(attributeId, desc);
311       metadata.setAttributeType(attributeId, attType);
312
313       if (isFieldWanted(attributeId, vcfFieldPatterns))
314       {
315         vcfFieldsOfInterest.add(attributeId);
316       }
317     }
318
319     FeatureSources.getInstance().addSource(theSourceId, metadata);
320   }
321
322   /**
323    * Answers true if the field id is matched by any of the filter patterns, else
324    * false. Matching is against regular expression patterns, and is not
325    * case-sensitive.
326    * 
327    * @param id
328    * @param filters
329    * @return
330    */
331   private boolean isFieldWanted(String id, List<Pattern> filters)
332   {
333     for (Pattern p : filters)
334     {
335       if (p.matcher(id.toUpperCase()).matches())
336       {
337         return true;
338       }
339     }
340     return false;
341   }
342
343   /**
344    * Records 'wanted' fields defined in the CSQ INFO header (if there is one).
345    * Also records the position of selected fields (Allele, ALLELE_NUM, Feature)
346    * required for processing.
347    * <p>
348    * CSQ fields are declared in the CSQ INFO Description e.g.
349    * <p>
350    * Description="Consequence ...from ... VEP. Format: Allele|Consequence|...
351    */
352   protected void parseCsqHeader()
353   {
354     List<Pattern> vepFieldFilters = getFieldMatchers(VEP_FIELDS_PREF,
355             DEFAULT_VEP_FIELDS);
356     vepFieldsOfInterest = new HashMap<>();
357
358     VCFInfoHeaderLine csqInfo = header.getInfoHeaderLine(CSQ_FIELD);
359     if (csqInfo == null)
360     {
361       return;
362     }
363
364     /*
365      * parse out the pipe-separated list of CSQ fields; we assume here that
366      * these form the last part of the description, and contain no spaces
367      */
368     String desc = csqInfo.getDescription();
369     int spacePos = desc.lastIndexOf(" ");
370     desc = desc.substring(spacePos + 1);
371
372     if (desc != null)
373     {
374       String[] format = desc.split(PIPE_REGEX);
375       int index = 0;
376       for (String field : format)
377       {
378         if (ALLELE_NUM_KEY.equals(field))
379         {
380           csqAlleleNumberFieldIndex = index;
381         }
382         if (ALLELE_KEY.equals(field))
383         {
384           csqAlleleFieldIndex = index;
385         }
386         if (FEATURE_KEY.equals(field))
387         {
388           csqFeatureFieldIndex = index;
389         }
390
391         if (isFieldWanted(field, vepFieldFilters))
392         {
393           vepFieldsOfInterest.put(index, field);
394         }
395
396         index++;
397       }
398     }
399   }
400
401   /**
402    * Reads the Preference value for the given key, with default specified if no
403    * preference set. The value is interpreted as a comma-separated list of
404    * regular expressions, and converted into a list of compiled patterns ready
405    * for matching. Patterns are forced to upper-case for non-case-sensitive
406    * matching.
407    * <p>
408    * This supports user-defined filters for fields of interest to capture while
409    * processing data. For example, VCF_FIELDS = AF,AC* would mean that VCF INFO
410    * fields with an ID of AF, or starting with AC, would be matched.
411    * 
412    * @param key
413    * @param def
414    * @return
415    */
416   private List<Pattern> getFieldMatchers(String key, String def)
417   {
418     String pref = Cache.getDefault(key, def);
419     List<Pattern> patterns = new ArrayList<>();
420     String[] tokens = pref.split(",");
421     for (String token : tokens)
422     {
423       patterns.add(Pattern.compile(token.toUpperCase()));
424     }
425     return patterns;
426   }
427
428   /**
429    * Transfers VCF features to sequences to which this sequence has a mapping.
430    * If the mapping is 3:1, computes peptide variants from nucleotide variants.
431    * 
432    * @param seq
433    */
434   protected void transferAddedFeatures(SequenceI seq)
435   {
436     DBRefEntry[] dbrefs = seq.getDBRefs();
437     if (dbrefs == null)
438     {
439       return;
440     }
441     for (DBRefEntry dbref : dbrefs)
442     {
443       Mapping mapping = dbref.getMap();
444       if (mapping == null || mapping.getTo() == null)
445       {
446         continue;
447       }
448
449       SequenceI mapTo = mapping.getTo();
450       MapList map = mapping.getMap();
451       if (map.getFromRatio() == 3)
452       {
453         /*
454          * dna-to-peptide product mapping
455          */
456         AlignmentUtils.computeProteinFeatures(seq, mapTo, map);
457       }
458       else
459       {
460         /*
461          * nucleotide-to-nucleotide mapping e.g. transcript to CDS
462          */
463         List<SequenceFeature> features = seq.getFeatures()
464                 .getPositionalFeatures(SequenceOntologyI.SEQUENCE_VARIANT);
465         for (SequenceFeature sf : features)
466         {
467           if (FEATURE_GROUP_VCF.equals(sf.getFeatureGroup()))
468           {
469             transferFeature(sf, mapTo, map);
470           }
471         }
472       }
473     }
474   }
475
476   /**
477    * Tries to add overlapping variants read from a VCF file to the given
478    * sequence, and returns the number of variant features added. Note that this
479    * requires the sequence to hold information as to its species, chromosomal
480    * positions and reference assembly, in order to be able to map the VCF
481    * variants to the sequence (or not)
482    * 
483    * @param seq
484    * @param reader
485    * @param vcfAssembly
486    * @return
487    */
488   protected int loadSequenceVCF(SequenceI seq, VCFReader reader,
489           String vcfAssembly)
490   {
491     int count = 0;
492     GeneLociI seqCoords = seq.getGeneLoci();
493     if (seqCoords == null)
494     {
495       System.out.println(String.format(
496               "Can't query VCF for %s as chromosome coordinates not known",
497               seq.getName()));
498       return 0;
499     }
500
501     if (!vcfSpeciesMatchesSequence(vcfAssembly, seqCoords.getSpeciesId()))
502     {
503       return 0;
504     }
505
506     List<int[]> seqChromosomalContigs = seqCoords.getMap().getToRanges();
507     for (int[] range : seqChromosomalContigs)
508     {
509       count += addVcfVariants(seq, reader, range, vcfAssembly);
510     }
511
512     return count;
513   }
514
515   /**
516    * Answers true if the species inferred from the VCF reference identifier
517    * matches that for the sequence
518    * 
519    * @param vcfAssembly
520    * @param speciesId
521    * @return
522    */
523   boolean vcfSpeciesMatchesSequence(String vcfAssembly, String speciesId)
524   {
525     // PROBLEM 1
526     // there are many aliases for species - how to equate one with another?
527     // PROBLEM 2
528     // VCF ##reference header is an unstructured URI - how to extract species?
529     // perhaps check if ref includes any (Ensembl) alias of speciesId??
530     // TODO ask the user to confirm this??
531
532     if (vcfAssembly.contains("Homo_sapiens") // gnomAD exome data example
533             && "HOMO_SAPIENS".equals(speciesId)) // Ensembl species id
534     {
535       return true;
536     }
537
538     if (vcfAssembly.contains("c_elegans") // VEP VCF response example
539             && "CAENORHABDITIS_ELEGANS".equals(speciesId)) // Ensembl
540     {
541       return true;
542     }
543
544     // this is not a sustainable solution...
545
546     return false;
547   }
548
549   /**
550    * Queries the VCF reader for any variants that overlap the given chromosome
551    * region of the sequence, and adds as variant features. Returns the number of
552    * overlapping variants found.
553    * 
554    * @param seq
555    * @param reader
556    * @param range
557    *          start-end range of a sequence region in its chromosomal
558    *          coordinates
559    * @param vcfAssembly
560    *          the '##reference' identifier for the VCF reference assembly
561    * @return
562    */
563   protected int addVcfVariants(SequenceI seq, VCFReader reader,
564           int[] range, String vcfAssembly)
565   {
566     GeneLociI seqCoords = seq.getGeneLoci();
567
568     String chromosome = seqCoords.getChromosomeId();
569     String seqRef = seqCoords.getAssemblyId();
570     String species = seqCoords.getSpeciesId();
571
572     /*
573      * map chromosomal coordinates from sequence to VCF if the VCF
574      * data has a different reference assembly to the sequence
575      */
576     // TODO generalise for non-human species
577     // - or get the user to choose in a dialog
578
579     int offset = 0;
580     if ("GRCh38".equalsIgnoreCase(seqRef) // Ensembl
581             && vcfAssembly.contains("Homo_sapiens_assembly19")) // gnomAD
582     {
583       String toRef = "GRCh37";
584       int[] newRange = mapReferenceRange(range, chromosome, "human",
585               seqRef, toRef);
586       if (newRange == null)
587       {
588         System.err.println(String.format(
589                 "Failed to map %s:%s:%s:%d:%d to %s", species, chromosome,
590                 seqRef, range[0], range[1], toRef));
591         return 0;
592       }
593       offset = newRange[0] - range[0];
594       range = newRange;
595     }
596
597     boolean forwardStrand = range[0] <= range[1];
598
599     /*
600      * query the VCF for overlaps
601      * (convert a reverse strand range to forwards)
602      */
603     int count = 0;
604     MapList mapping = seqCoords.getMap();
605
606     int fromLocus = Math.min(range[0], range[1]);
607     int toLocus = Math.max(range[0], range[1]);
608     CloseableIterator<VariantContext> variants = reader.query(chromosome,
609             fromLocus, toLocus);
610     while (variants.hasNext())
611     {
612       /*
613        * get variant location in sequence chromosomal coordinates
614        */
615       VariantContext variant = variants.next();
616
617       int start = variant.getStart() - offset;
618       int end = variant.getEnd() - offset;
619
620       /*
621        * convert chromosomal location to sequence coordinates
622        * - may be reverse strand (convert to forward for sequence feature)
623        * - null if a partially overlapping feature
624        */
625       int[] seqLocation = mapping.locateInFrom(start, end);
626       if (seqLocation != null)
627       {
628         int featureStart = Math.min(seqLocation[0], seqLocation[1]);
629         int featureEnd = Math.max(seqLocation[0], seqLocation[1]);
630         count += addAlleleFeatures(seq, variant, featureStart, featureEnd,
631                 forwardStrand);
632       }
633     }
634
635     variants.close();
636
637     return count;
638   }
639
640   /**
641    * A convenience method to get the AF value for the given alternate allele
642    * index
643    * 
644    * @param variant
645    * @param alleleIndex
646    * @return
647    */
648   protected float getAlleleFrequency(VariantContext variant, int alleleIndex)
649   {
650     float score = 0f;
651     String attributeValue = getAttributeValue(variant,
652             ALLELE_FREQUENCY_KEY, alleleIndex);
653     if (attributeValue != null)
654     {
655       try
656       {
657         score = Float.parseFloat(attributeValue);
658       } catch (NumberFormatException e)
659       {
660         // leave as 0
661       }
662     }
663
664     return score;
665   }
666
667   /**
668    * A convenience method to get an attribute value for an alternate allele
669    * 
670    * @param variant
671    * @param attributeName
672    * @param alleleIndex
673    * @return
674    */
675   protected String getAttributeValue(VariantContext variant,
676           String attributeName, int alleleIndex)
677   {
678     Object att = variant.getAttribute(attributeName);
679
680     if (att instanceof String)
681     {
682       return (String) att;
683     }
684     else if (att instanceof ArrayList)
685     {
686       return ((List<String>) att).get(alleleIndex);
687     }
688
689     return null;
690   }
691
692   /**
693    * Adds one variant feature for each allele in the VCF variant record, and
694    * returns the number of features added.
695    * 
696    * @param seq
697    * @param variant
698    * @param featureStart
699    * @param featureEnd
700    * @param forwardStrand
701    * @return
702    */
703   protected int addAlleleFeatures(SequenceI seq, VariantContext variant,
704           int featureStart, int featureEnd, boolean forwardStrand)
705   {
706     int added = 0;
707
708     /*
709      * Javadoc says getAlternateAlleles() imposes no order on the list returned
710      * so we proceed defensively to get them in strict order
711      */
712     int altAlleleCount = variant.getAlternateAlleles().size();
713     for (int i = 0; i < altAlleleCount; i++)
714     {
715       added += addAlleleFeature(seq, variant, i, featureStart, featureEnd,
716               forwardStrand);
717     }
718     return added;
719   }
720
721   /**
722    * Inspects one allele and attempts to add a variant feature for it to the
723    * sequence. We extract as much as possible of the additional data associated
724    * with this allele to store in the feature's key-value map. Answers the
725    * number of features added (0 or 1).
726    * 
727    * @param seq
728    * @param variant
729    * @param altAlleleIndex
730    *          (0, 1..)
731    * @param featureStart
732    * @param featureEnd
733    * @param forwardStrand
734    * @return
735    */
736   protected int addAlleleFeature(SequenceI seq, VariantContext variant,
737           int altAlleleIndex, int featureStart, int featureEnd,
738           boolean forwardStrand)
739   {
740     String reference = variant.getReference().getBaseString();
741     Allele alt = variant.getAlternateAllele(altAlleleIndex);
742     String allele = alt.getBaseString();
743
744     /*
745      * build the ref,alt allele description e.g. "G,A", using the base
746      * complement if the sequence is on the reverse strand
747      */
748     // TODO check how structural variants are shown on reverse strand
749     StringBuilder sb = new StringBuilder();
750     sb.append(forwardStrand ? reference : Dna.reverseComplement(reference));
751     sb.append(COMMA);
752     sb.append(forwardStrand ? allele : Dna.reverseComplement(allele));
753     String alleles = sb.toString(); // e.g. G,A
754
755     String type = SequenceOntologyI.SEQUENCE_VARIANT;
756     float score = getAlleleFrequency(variant, altAlleleIndex);
757
758     SequenceFeature sf = new SequenceFeature(type, alleles, featureStart,
759             featureEnd, score, FEATURE_GROUP_VCF);
760     sf.setSource(sourceId);
761
762     sf.setValue(Gff3Helper.ALLELES, alleles);
763
764     addAlleleProperties(variant, seq, sf, altAlleleIndex);
765
766     seq.addSequenceFeature(sf);
767
768     return 1;
769   }
770
771   /**
772    * Add any allele-specific VCF key-value data to the sequence feature
773    * 
774    * @param variant
775    * @param seq
776    * @param sf
777    * @param altAlelleIndex
778    *          (0, 1..)
779    */
780   protected void addAlleleProperties(VariantContext variant, SequenceI seq,
781           SequenceFeature sf, final int altAlelleIndex)
782   {
783     Map<String, Object> atts = variant.getAttributes();
784
785     for (Entry<String, Object> att : atts.entrySet())
786     {
787       String key = att.getKey();
788
789       /*
790        * extract Consequence data (if present) that we are able to
791        * associated with the allele for this variant feature
792        */
793       if (CSQ_FIELD.equals(key))
794       {
795         addConsequences(variant, seq, sf, altAlelleIndex);
796         continue;
797       }
798
799       /*
800        * filter out fields we don't want to capture
801        */
802       if (!vcfFieldsOfInterest.contains(key))
803       {
804         continue;
805       }
806
807       /*
808        * we extract values for other data which are allele-specific; 
809        * these may be per alternate allele (INFO[key].Number = 'A') 
810        * or per allele including reference (INFO[key].Number = 'R') 
811        */
812       VCFInfoHeaderLine infoHeader = header.getInfoHeaderLine(key);
813       if (infoHeader == null)
814       {
815         /*
816          * can't be sure what data belongs to this allele, so
817          * play safe and don't take any
818          */
819         continue;
820       }
821
822       VCFHeaderLineCount number = infoHeader.getCountType();
823       int index = altAlelleIndex;
824       if (number == VCFHeaderLineCount.R)
825       {
826         /*
827          * one value per allele including reference, so bump index
828          * e.g. the 3rd value is for the  2nd alternate allele
829          */
830         index++;
831       }
832       else if (number != VCFHeaderLineCount.A)
833       {
834         /*
835          * don't save other values as not allele-related
836          */
837         continue;
838       }
839
840       /*
841        * take the index'th value
842        */
843       String value = getAttributeValue(variant, key, index);
844       if (value != null)
845       {
846         sf.setValue(key, value);
847       }
848     }
849   }
850
851   /**
852    * Inspects CSQ data blocks (consequences) and adds attributes on the sequence
853    * feature for the current allele (and transcript if applicable)
854    * <p>
855    * Allele matching: if field ALLELE_NUM is present, it must match
856    * altAlleleIndex. If not present, then field Allele value must match the VCF
857    * Allele.
858    * <p>
859    * Transcript matching: if sequence name can be identified to at least one of
860    * the consequences' Feature values, then select only consequences that match
861    * the value (i.e. consequences for the current transcript sequence). If not,
862    * take all consequences (this is the case when adding features to the gene
863    * sequence).
864    * 
865    * @param variant
866    * @param seq
867    * @param sf
868    * @param altAlelleIndex
869    *          (0, 1..)
870    */
871   protected void addConsequences(VariantContext variant, SequenceI seq,
872           SequenceFeature sf, int altAlelleIndex)
873   {
874     Object value = variant.getAttribute(CSQ_FIELD);
875
876     if (value == null || !(value instanceof ArrayList<?>))
877     {
878       return;
879     }
880
881     List<String> consequences = (List<String>) value;
882
883     /*
884      * if CSQ data includes 'Feature', and any value matches the sequence name,
885      * then restrict consequence data to only the matching value (transcript)
886      * i.e. just pick out consequences for the transcript the variant feature is on
887      */
888     String seqName = seq.getName()== null ? "" : seq.getName().toLowerCase();
889     String matchFeature = null;
890     if (csqFeatureFieldIndex > -1)
891     {
892       for (String consequence : consequences)
893       {
894         String[] csqFields = consequence.split(PIPE_REGEX);
895         if (csqFields.length > csqFeatureFieldIndex)
896         {
897           String featureIdentifier = csqFields[csqFeatureFieldIndex];
898           if (featureIdentifier.length() > 4
899                   && seqName.indexOf(featureIdentifier.toLowerCase()) > -1)
900           {
901             matchFeature = featureIdentifier;
902           }
903         }
904       }
905     }
906
907     /*
908      * inspect CSQ consequences; where possible restrict to the consequence
909      * associated with the current transcript (Feature)
910      */
911     Map<String, String> csqValues = new HashMap<>();
912
913     for (String consequence : consequences)
914     {
915       String[] csqFields = consequence.split(PIPE_REGEX);
916
917       if (includeConsequence(csqFields, matchFeature, variant,
918               altAlelleIndex))
919       {
920         /*
921          * inspect individual fields of this consequence, copying non-null
922          * values which are 'fields of interest'
923          */
924         int i = 0;
925         for (String field : csqFields)
926         {
927           if (field != null && field.length() > 0)
928           {
929             String id = vepFieldsOfInterest.get(i);
930             if (id != null)
931             {
932               csqValues.put(id, field);
933             }
934           }
935           i++;
936         }
937       }
938     }
939
940     if (!csqValues.isEmpty())
941     {
942       sf.setValue(CSQ_FIELD, csqValues);
943     }
944   }
945
946   /**
947    * Answers true if we want to associate this block of consequence data with
948    * the specified alternate allele of the VCF variant.
949    * <p>
950    * If consequence data includes the ALLELE_NUM field, then this has to match
951    * altAlleleIndex. Otherwise the Allele field of the consequence data has to
952    * match the allele value.
953    * <p>
954    * Optionally (if matchFeature is not null), restrict to only include
955    * consequences whose Feature value matches. This allows us to attach
956    * consequences to their respective transcripts.
957    * 
958    * @param csqFields
959    * @param matchFeature
960    * @param variant
961    * @param altAlelleIndex
962    *          (0, 1..)
963    * @return
964    */
965   protected boolean includeConsequence(String[] csqFields,
966           String matchFeature, VariantContext variant, int altAlelleIndex)
967   {
968     /*
969      * check consequence is for the current transcript
970      */
971     if (matchFeature != null)
972     {
973       if (csqFields.length <= csqFeatureFieldIndex)
974       {
975         return false;
976       }
977       String featureIdentifier = csqFields[csqFeatureFieldIndex];
978       if (!featureIdentifier.equals(matchFeature))
979       {
980         return false; // consequence is for a different transcript
981       }
982     }
983
984     /*
985      * if ALLELE_NUM is present, it must match altAlleleIndex
986      * NB first alternate allele is 1 for ALLELE_NUM, 0 for altAlleleIndex
987      */
988     if (csqAlleleNumberFieldIndex > -1)
989     {
990       if (csqFields.length <= csqAlleleNumberFieldIndex)
991       {
992         return false;
993       }
994       String alleleNum = csqFields[csqAlleleNumberFieldIndex];
995       return String.valueOf(altAlelleIndex + 1).equals(alleleNum);
996     }
997
998     /*
999      * else consequence allele must match variant allele
1000      */
1001     if (csqAlleleFieldIndex > -1 && csqFields.length > csqAlleleFieldIndex)
1002     {
1003       String csqAllele = csqFields[csqAlleleFieldIndex];
1004       String vcfAllele = variant.getAlternateAllele(altAlelleIndex)
1005               .getBaseString();
1006       return csqAllele.equals(vcfAllele);
1007     }
1008
1009     return false;
1010   }
1011
1012   /**
1013    * A convenience method to complement a dna base and return the string value
1014    * of its complement
1015    * 
1016    * @param reference
1017    * @return
1018    */
1019   protected String complement(byte[] reference)
1020   {
1021     return String.valueOf(Dna.getComplement((char) reference[0]));
1022   }
1023
1024   /**
1025    * Determines the location of the query range (chromosome positions) in a
1026    * different reference assembly.
1027    * <p>
1028    * If the range is just a subregion of one for which we already have a mapping
1029    * (for example, an exon sub-region of a gene), then the mapping is just
1030    * computed arithmetically.
1031    * <p>
1032    * Otherwise, calls the Ensembl REST service that maps from one assembly
1033    * reference's coordinates to another's
1034    * 
1035    * @param queryRange
1036    *          start-end chromosomal range in 'fromRef' coordinates
1037    * @param chromosome
1038    * @param species
1039    * @param fromRef
1040    *          assembly reference for the query coordinates
1041    * @param toRef
1042    *          assembly reference we wish to translate to
1043    * @return the start-end range in 'toRef' coordinates
1044    */
1045   protected int[] mapReferenceRange(int[] queryRange, String chromosome,
1046           String species, String fromRef, String toRef)
1047   {
1048     /*
1049      * first try shorcut of computing the mapping as a subregion of one
1050      * we already have (e.g. for an exon, if we have the gene mapping)
1051      */
1052     int[] mappedRange = findSubsumedRangeMapping(queryRange, chromosome,
1053             species, fromRef, toRef);
1054     if (mappedRange != null)
1055     {
1056       return mappedRange;
1057     }
1058
1059     /*
1060      * call (e.g.) http://rest.ensembl.org/map/human/GRCh38/17:45051610..45109016:1/GRCh37
1061      */
1062     EnsemblMap mapper = new EnsemblMap();
1063     int[] mapping = mapper.getAssemblyMapping(species, chromosome, fromRef,
1064             toRef, queryRange);
1065
1066     if (mapping == null)
1067     {
1068       // mapping service failure
1069       return null;
1070     }
1071
1072     /*
1073      * save mapping for possible future re-use
1074      */
1075     String key = makeRangesKey(chromosome, species, fromRef, toRef);
1076     if (!assemblyMappings.containsKey(key))
1077     {
1078       assemblyMappings.put(key, new HashMap<int[], int[]>());
1079     }
1080
1081     assemblyMappings.get(key).put(queryRange, mapping);
1082
1083     return mapping;
1084   }
1085
1086   /**
1087    * If we already have a 1:1 contiguous mapping which subsumes the given query
1088    * range, this method just calculates and returns the subset of that mapping,
1089    * else it returns null. In practical terms, if a gene has a contiguous
1090    * mapping between (for example) GRCh37 and GRCh38, then we assume that its
1091    * subsidiary exons occupy unchanged relative positions, and just compute
1092    * these as offsets, rather than do another lookup of the mapping.
1093    * <p>
1094    * If in future these assumptions prove invalid (e.g. for bacterial dna?!),
1095    * simply remove this method or let it always return null.
1096    * <p>
1097    * Warning: many rapid calls to the /map service map result in a 429 overload
1098    * error response
1099    * 
1100    * @param queryRange
1101    * @param chromosome
1102    * @param species
1103    * @param fromRef
1104    * @param toRef
1105    * @return
1106    */
1107   protected int[] findSubsumedRangeMapping(int[] queryRange, String chromosome,
1108           String species, String fromRef, String toRef)
1109   {
1110     String key = makeRangesKey(chromosome, species, fromRef, toRef);
1111     if (assemblyMappings.containsKey(key))
1112     {
1113       Map<int[], int[]> mappedRanges = assemblyMappings.get(key);
1114       for (Entry<int[], int[]> mappedRange : mappedRanges.entrySet())
1115       {
1116         int[] fromRange = mappedRange.getKey();
1117         int[] toRange = mappedRange.getValue();
1118         if (fromRange[1] - fromRange[0] == toRange[1] - toRange[0])
1119         {
1120           /*
1121            * mapping is 1:1 in length, so we trust it to have no discontinuities
1122            */
1123           if (MappingUtils.rangeContains(fromRange, queryRange))
1124           {
1125             /*
1126              * fromRange subsumes our query range
1127              */
1128             int offset = queryRange[0] - fromRange[0];
1129             int mappedRangeFrom = toRange[0] + offset;
1130             int mappedRangeTo = mappedRangeFrom + (queryRange[1] - queryRange[0]);
1131             return new int[] { mappedRangeFrom, mappedRangeTo };
1132           }
1133         }
1134       }
1135     }
1136     return null;
1137   }
1138
1139   /**
1140    * Transfers the sequence feature to the target sequence, locating its start
1141    * and end range based on the mapping. Features which do not overlap the
1142    * target sequence are ignored.
1143    * 
1144    * @param sf
1145    * @param targetSequence
1146    * @param mapping
1147    *          mapping from the feature's coordinates to the target sequence
1148    */
1149   protected void transferFeature(SequenceFeature sf,
1150           SequenceI targetSequence, MapList mapping)
1151   {
1152     int[] mappedRange = mapping.locateInTo(sf.getBegin(), sf.getEnd());
1153   
1154     if (mappedRange != null)
1155     {
1156       String group = sf.getFeatureGroup();
1157       int newBegin = Math.min(mappedRange[0], mappedRange[1]);
1158       int newEnd = Math.max(mappedRange[0], mappedRange[1]);
1159       SequenceFeature copy = new SequenceFeature(sf, newBegin, newEnd,
1160               group, sf.getScore());
1161       targetSequence.addSequenceFeature(copy);
1162     }
1163   }
1164
1165   /**
1166    * Formats a ranges map lookup key
1167    * 
1168    * @param chromosome
1169    * @param species
1170    * @param fromRef
1171    * @param toRef
1172    * @return
1173    */
1174   protected static String makeRangesKey(String chromosome, String species,
1175           String fromRef, String toRef)
1176   {
1177     return species + EXCL + chromosome + EXCL + fromRef + EXCL
1178             + toRef;
1179   }
1180 }