2 * Jalview - A Sequence Alignment Editor and Viewer ($$Version-Rel$$)
3 * Copyright (C) $$Year-Rel$$ The Jalview Authors
5 * This file is part of Jalview.
7 * Jalview is free software: you can redistribute it and/or
8 * modify it under the terms of the GNU General Public License
9 * as published by the Free Software Foundation, either version 3
10 * of the License, or (at your option) any later version.
12 * Jalview is distributed in the hope that it will be useful, but
13 * WITHOUT ANY WARRANTY; without even the implied warranty
14 * of MERCHANTABILITY or FITNESS FOR A PARTICULAR
15 * PURPOSE. See the GNU General Public License for more details.
17 * You should have received a copy of the GNU General Public License
18 * along with Jalview. If not, see <http://www.gnu.org/licenses/>.
19 * The Jalview Authors are detailed in the 'AUTHORS' file.
21 package jalview.analysis;
23 import jalview.analysis.scoremodels.PIDModel;
24 import jalview.analysis.scoremodels.SimilarityParams;
25 import jalview.datamodel.AlignmentAnnotation;
26 import jalview.datamodel.AlignmentI;
27 import jalview.datamodel.AlignmentOrder;
28 import jalview.datamodel.SequenceFeature;
29 import jalview.datamodel.SequenceGroup;
30 import jalview.datamodel.SequenceI;
31 import jalview.datamodel.SequenceNode;
32 import jalview.util.QuickSort;
34 import java.util.ArrayList;
35 import java.util.Collections;
36 import java.util.Iterator;
37 import java.util.List;
40 * Routines for manipulating the order of a multiple sequence alignment TODO:
41 * this class retains some global states concerning sort-order which should be
42 * made attributes for the caller's alignment visualization. TODO: refactor to
43 * allow a subset of selected sequences to be sorted within the context of a
44 * whole alignment. Sort method template is: SequenceI[] tobesorted, [ input
45 * data mapping to each tobesorted element to use ], Alignment context of
46 * tobesorted that are to be re-ordered, boolean sortinplace, [special data - ie
47 * seuqence to be sorted w.r.t.]) sortinplace implies that the sorted vector
48 * resulting from applying the operation to tobesorted should be mapped back to
49 * the original positions in alignment. Otherwise, normal behaviour is to re
50 * order alignment so that tobesorted is sorted and grouped together starting
51 * from the first tobesorted position in the alignment. e.g. (a,tb2,b,tb1,c,tb3
52 * becomes a,tb1,tb2,tb3,b,c)
54 public class AlignmentSorter
57 * todo: refactor searches to follow a basic pattern: (search property, last
58 * search state, current sort direction)
60 static boolean sortIdAscending = true;
62 static int lastGroupHash = 0;
64 static boolean sortGroupAscending = true;
66 static AlignmentOrder lastOrder = null;
68 static boolean sortOrderAscending = true;
70 static TreeModel lastTree = null;
72 static boolean sortTreeAscending = true;
75 * last Annotation Label used for sort by Annotation score
77 private static String lastSortByAnnotation;
80 * string hash of last arguments to sortByFeature
81 * (sort order toggles if this is unchanged between sorts)
83 private static String sortByFeatureCriteria;
85 private static boolean sortByFeatureAscending = true;
87 private static boolean sortLengthAscending;
90 * Sorts sequences in the alignment by Percentage Identity with the given
91 * reference sequence, sorting the highest identity to the top
99 public static void sortByPID(AlignmentI align, SequenceI s)
101 int nSeq = align.getHeight();
103 float[] scores = new float[nSeq];
104 SequenceI[] seqs = new SequenceI[nSeq];
105 String refSeq = s.getSequenceAsString();
107 SimilarityParams pidParams = new SimilarityParams(true, true, true,
109 for (int i = 0; i < nSeq; i++)
111 scores[i] = (float) PIDModel.computePID(
112 align.getSequenceAt(i).getSequenceAsString(), refSeq,
114 seqs[i] = align.getSequenceAt(i);
117 QuickSort.sort(scores, seqs);
119 setReverseOrder(align, seqs);
123 * Reverse the order of the sort
130 private static void setReverseOrder(AlignmentI align, SequenceI[] seqs)
132 int nSeq = seqs.length;
142 len = (nSeq + 1) / 2;
145 // NOTE: DO NOT USE align.setSequenceAt() here - it will NOT work
146 List<SequenceI> asq = align.getSequences();
149 for (int i = 0; i < len; i++)
151 // SequenceI tmp = seqs[i];
152 asq.set(i, seqs[nSeq - i - 1]);
153 asq.set(nSeq - i - 1, seqs[i]);
159 * Sets the Alignment object with the given sequences
162 * Alignment object to be updated
164 * sequences as a vector
166 private static void setOrder(AlignmentI align, List<SequenceI> tmp)
168 setOrder(align, vectorSubsetToArray(tmp, align.getSequences()));
172 * Sets the Alignment object with the given sequences
177 * sequences as an array
179 public static void setOrder(AlignmentI align, SequenceI[] seqs)
181 // NOTE: DO NOT USE align.setSequenceAt() here - it will NOT work
182 List<SequenceI> algn = align.getSequences();
185 List<SequenceI> tmp = new ArrayList<>();
187 for (int i = 0; i < seqs.length; i++)
189 if (algn.contains(seqs[i]))
196 // User may have hidden seqs, then clicked undo or redo
197 for (int i = 0; i < tmp.size(); i++)
199 algn.add(tmp.get(i));
205 * Sorts by ID. Numbers are sorted before letters.
208 * The alignment object to sort
210 public static void sortByID(AlignmentI align)
212 int nSeq = align.getHeight();
214 String[] ids = new String[nSeq];
215 SequenceI[] seqs = new SequenceI[nSeq];
217 for (int i = 0; i < nSeq; i++)
219 ids[i] = align.getSequenceAt(i).getName();
220 seqs[i] = align.getSequenceAt(i);
223 QuickSort.sort(ids, seqs);
227 setReverseOrder(align, seqs);
231 setOrder(align, seqs);
234 sortIdAscending = !sortIdAscending;
238 * Sorts by sequence length
241 * The alignment object to sort
243 public static void sortByLength(AlignmentI align)
245 int nSeq = align.getHeight();
247 float[] length = new float[nSeq];
248 SequenceI[] seqs = new SequenceI[nSeq];
250 for (int i = 0; i < nSeq; i++)
252 seqs[i] = align.getSequenceAt(i);
253 length[i] = (seqs[i].getEnd() - seqs[i].getStart());
256 QuickSort.sort(length, seqs);
258 if (sortLengthAscending)
260 setReverseOrder(align, seqs);
264 setOrder(align, seqs);
267 sortLengthAscending = !sortLengthAscending;
271 * Sorts the alignment by size of group. <br>
272 * Maintains the order of sequences in each group by order in given alignment
276 * sorts the given alignment object by group
278 public static void sortByGroup(AlignmentI align)
280 // MAINTAINS ORIGNAL SEQUENCE ORDER,
281 // ORDERS BY GROUP SIZE
282 List<SequenceGroup> groups = new ArrayList<>();
284 if (groups.hashCode() != lastGroupHash)
286 sortGroupAscending = true;
287 lastGroupHash = groups.hashCode();
291 sortGroupAscending = !sortGroupAscending;
294 // SORTS GROUPS BY SIZE
295 // ////////////////////
296 for (SequenceGroup sg : align.getGroups())
298 for (int j = 0; j < groups.size(); j++)
300 SequenceGroup sg2 = groups.get(j);
302 if (sg.getSize() > sg2.getSize())
310 if (!groups.contains(sg))
316 // NOW ADD SEQUENCES MAINTAINING ALIGNMENT ORDER
317 // /////////////////////////////////////////////
318 List<SequenceI> seqs = new ArrayList<>();
320 for (int i = 0; i < groups.size(); i++)
322 SequenceGroup sg = groups.get(i);
323 SequenceI[] orderedseqs = sg.getSequencesInOrder(align);
325 for (int j = 0; j < orderedseqs.length; j++)
327 seqs.add(orderedseqs[j]);
331 if (sortGroupAscending)
333 setOrder(align, seqs);
337 setReverseOrder(align,
338 vectorSubsetToArray(seqs, align.getSequences()));
343 * Select sequences in order from tmp that is present in mask, and any
344 * remaining sequences in mask not in tmp
347 * thread safe collection of sequences
349 * thread safe collection of sequences
351 * @return intersect(tmp,mask)+intersect(complement(tmp),mask)
353 private static SequenceI[] vectorSubsetToArray(List<SequenceI> tmp,
354 List<SequenceI> mask)
357 // tmp2 = tmp.retainAll(mask);
358 // return tmp2.addAll(mask.removeAll(tmp2))
360 ArrayList<SequenceI> seqs = new ArrayList<>();
362 boolean[] tmask = new boolean[mask.size()];
364 for (i = 0; i < mask.size(); i++)
369 for (i = 0; i < tmp.size(); i++)
371 SequenceI sq = tmp.get(i);
372 idx = mask.indexOf(sq);
373 if (idx > -1 && tmask[idx])
380 for (i = 0; i < tmask.length; i++)
384 seqs.add(mask.get(i));
388 return seqs.toArray(new SequenceI[seqs.size()]);
392 * Sorts by a given AlignmentOrder object
397 * specified order for alignment
399 public static void sortBy(AlignmentI align, AlignmentOrder order)
401 // Get an ordered vector of sequences which may also be present in align
402 List<SequenceI> tmp = order.getOrder();
404 if (lastOrder == order)
406 sortOrderAscending = !sortOrderAscending;
410 sortOrderAscending = true;
413 if (sortOrderAscending)
415 setOrder(align, tmp);
419 setReverseOrder(align,
420 vectorSubsetToArray(tmp, align.getSequences()));
432 * @return DOCUMENT ME!
434 private static List<SequenceI> getOrderByTree(AlignmentI align,
437 int nSeq = align.getHeight();
439 List<SequenceI> tmp = new ArrayList<>();
441 tmp = _sortByTree(tree.getTopNode(), tmp, align.getSequences());
443 if (tmp.size() != nSeq)
445 // TODO: JBPNote - decide if this is always an error
446 // (eg. not when a tree is associated to another alignment which has more
448 if (tmp.size() != nSeq)
450 addStrays(align, tmp);
453 if (tmp.size() != nSeq)
455 System.err.println("WARNING: tmp.size()=" + tmp.size() + " != nseq="
457 + " in getOrderByTree - tree contains sequences not in alignment");
465 * Sorts the alignment by a given tree
472 public static void sortByTree(AlignmentI align, TreeModel tree)
474 List<SequenceI> tmp = getOrderByTree(align, tree);
476 // tmp should properly permute align with tree.
477 if (lastTree != tree)
479 sortTreeAscending = true;
484 sortTreeAscending = !sortTreeAscending;
487 if (sortTreeAscending)
489 setOrder(align, tmp);
493 setReverseOrder(align,
494 vectorSubsetToArray(tmp, align.getSequences()));
506 private static void addStrays(AlignmentI align, List<SequenceI> tmp)
508 int nSeq = align.getHeight();
510 for (int i = 0; i < nSeq; i++)
512 if (!tmp.contains(align.getSequenceAt(i)))
514 tmp.add(align.getSequenceAt(i));
518 if (nSeq != tmp.size())
521 .println("ERROR: Size still not right even after addStrays");
535 * @return DOCUMENT ME!
537 private static List<SequenceI> _sortByTree(SequenceNode node,
538 List<SequenceI> tmp, List<SequenceI> seqset)
545 SequenceNode left = (SequenceNode) node.left();
546 SequenceNode right = (SequenceNode) node.right();
548 if ((left == null) && (right == null))
550 if (!node.isPlaceholder() && (node.element() != null))
552 if (node.element() instanceof SequenceI)
554 if (!tmp.contains(node.element())) // && (seqset==null ||
555 // seqset.size()==0 ||
556 // seqset.contains(tmp)))
558 tmp.add((SequenceI) node.element());
567 _sortByTree(left, tmp, seqset);
568 _sortByTree(right, tmp, seqset);
575 // Alignment.sortBy(OrderObj) - sequence of sequence pointer refs in
580 * recover the order of sequences given by the safe numbering scheme introducd
581 * SeqsetUtils.uniquify.
583 public static void recoverOrder(SequenceI[] alignment)
585 float[] ids = new float[alignment.length];
587 for (int i = 0; i < alignment.length; i++)
589 ids[i] = (Float.valueOf(alignment[i].getName().substring(8)))
593 jalview.util.QuickSort.sort(ids, alignment);
597 * Sort sequence in order of increasing score attribute for annotation with a
598 * particular scoreLabel. Or reverse if same label was used previously
601 * exact label for sequence associated AlignmentAnnotation scores to
604 * sequences to be sorted
606 public static void sortByAnnotationScore(String scoreLabel,
607 AlignmentI alignment)
609 SequenceI[] seqs = alignment.getSequencesArray();
610 boolean[] hasScore = new boolean[seqs.length]; // per sequence score
612 int hasScores = 0; // number of scores present on set
613 double[] scores = new double[seqs.length];
614 double min = 0, max = 0;
615 for (int i = 0; i < seqs.length; i++)
617 AlignmentAnnotation[] scoreAnn = seqs[i].getAnnotation(scoreLabel);
618 if (scoreAnn != null)
622 scores[i] = scoreAnn[0].getScore(); // take the first instance of this
626 max = min = scores[i];
647 return; // do nothing - no scores present to sort by.
649 if (hasScores < seqs.length)
651 for (int i = 0; i < seqs.length; i++)
655 scores[i] = (max + i + 1.0);
660 jalview.util.QuickSort.sort(scores, seqs);
661 if (lastSortByAnnotation != scoreLabel)
663 lastSortByAnnotation = scoreLabel;
664 setOrder(alignment, seqs);
668 setReverseOrder(alignment, seqs);
673 * types of feature ordering: Sort by score : average score - or total score -
674 * over all features in region Sort by feature label text: (or if null -
675 * feature type text) - numerical or alphabetical Sort by feature density:
676 * based on counts - ignoring individual text or scores for each feature
678 public static String FEATURE_SCORE = "average_score";
680 public static String FEATURE_LABEL = "text";
682 public static String FEATURE_DENSITY = "density";
685 * Sort sequences by feature score or density, optionally restricted by
686 * feature types, feature groups, or alignment start/end positions.
688 * If the sort is repeated for the same combination of types and groups, sort
691 * @param featureTypes
692 * a list of feature types to include (or null for all)
694 * a list of feature groups to include (or null for all)
696 * start column position to include (base zero)
698 * end column position to include (base zero)
700 * the alignment to be sorted
702 * either "average_score" or "density" ("text" not yet implemented)
704 public static void sortByFeature(List<String> featureTypes,
705 List<String> groups, final int startCol, final int endCol,
706 AlignmentI alignment, String method)
708 if (method != FEATURE_SCORE && method != FEATURE_LABEL
709 && method != FEATURE_DENSITY)
712 .format("Implementation Error - sortByFeature method must be either '%s' or '%s'",
713 FEATURE_SCORE, FEATURE_DENSITY);
714 System.err.println(msg);
718 flipFeatureSortIfUnchanged(method, featureTypes, groups, startCol, endCol);
720 SequenceI[] seqs = alignment.getSequencesArray();
722 boolean[] hasScore = new boolean[seqs.length]; // per sequence score
724 int hasScores = 0; // number of scores present on set
725 double[] scores = new double[seqs.length];
726 int[] seqScores = new int[seqs.length];
727 Object[][] feats = new Object[seqs.length][];
731 for (int i = 0; i < seqs.length; i++)
734 * get sequence residues overlapping column region
735 * and features for residue positions and specified types
737 String[] types = featureTypes == null ? null : featureTypes
738 .toArray(new String[featureTypes.size()]);
739 List<SequenceFeature> sfs = seqs[i].findFeatures(startCol + 1,
745 Iterator<SequenceFeature> it = sfs.listIterator();
748 SequenceFeature sf = it.next();
751 * accept all features with null or empty group, otherwise
752 * check group is one of the currently visible groups
754 String featureGroup = sf.getFeatureGroup();
755 if (groups != null && featureGroup != null
756 && !"".equals(featureGroup)
757 && !groups.contains(featureGroup))
763 float score = sf.getScore();
764 if (FEATURE_SCORE.equals(method) && !Float.isNaN(score))
766 if (seqScores[i] == 0)
773 // take the first instance of this score // ??
778 feats[i] = sfs.toArray(new SequenceFeature[sfs.size()]);
781 if (method == FEATURE_LABEL)
783 // order the labels by alphabet (not yet implemented)
784 String[] labs = new String[sfs.size()];
785 for (int l = 0; l < sfs.size(); l++)
787 SequenceFeature sf = sfs.get(l);
788 String description = sf.getDescription();
789 labs[l] = (description != null ? description : sf.getType());
791 QuickSort.sort(labs, feats[i]);
796 // compute average score
797 scores[i] /= seqScores[i];
798 // update the score bounds.
806 max = Math.max(max, scores[i]);
807 min = Math.min(min, scores[i]);
812 if (FEATURE_SCORE.equals(method))
816 return; // do nothing - no scores present to sort by.
819 if (hasScores < seqs.length)
821 for (int i = 0; i < seqs.length; i++)
825 scores[i] = (max + 1 + i);
829 // int nf = (feats[i] == null) ? 0
830 // : ((SequenceFeature[]) feats[i]).length;
831 // // System.err.println("Sorting on Score: seq " +
833 // + " Feats: " + nf + " Score : " + scores[i]);
837 QuickSort.sortByDouble(scores, seqs, sortByFeatureAscending);
839 else if (FEATURE_DENSITY.equals(method))
841 for (int i = 0; i < seqs.length; i++)
843 int featureCount = feats[i] == null ? 0
844 : ((SequenceFeature[]) feats[i]).length;
845 scores[i] = featureCount;
846 // System.err.println("Sorting on Density: seq "+seqs[i].getName()+
847 // " Feats: "+featureCount+" Score : "+scores[i]);
849 QuickSort.sortByDouble(scores, seqs, sortByFeatureAscending);
852 setOrder(alignment, seqs);
856 * Builds a string hash of criteria for sorting, and if unchanged from last
857 * time, reverse the sort order
860 * @param featureTypes
865 protected static void flipFeatureSortIfUnchanged(String method,
866 List<String> featureTypes, List<String> groups,
867 final int startCol, final int endCol)
869 StringBuilder sb = new StringBuilder(64);
870 sb.append(startCol).append(method).append(endCol);
871 if (featureTypes != null)
873 Collections.sort(featureTypes);
874 sb.append(featureTypes.toString());
878 Collections.sort(groups);
879 sb.append(groups.toString());
881 String scoreCriteria = sb.toString();
884 * if resorting on the same criteria, toggle sort order
886 if (sortByFeatureCriteria == null
887 || !scoreCriteria.equals(sortByFeatureCriteria))
889 sortByFeatureAscending = true;
893 sortByFeatureAscending = !sortByFeatureAscending;
895 sortByFeatureCriteria = scoreCriteria;