I agree it could be a bit weird, but the merge method will catch any
inconsistency and throw an error, so I think it's fairly safe. And as you
say, the common case should be that the files all have the same seqinfo
anyway.
On Apr 25, 2014 10:32 PM, "Martin Morgan" <mtmor...@fhcrc.org> wrote:

> Hi Ryan --
>
> I implemented your suggestion in Rsamtools 1.17.8; it's a little weird to
> come up with such a synthetic Seqinfo, but I think the usual use case will
> be that the bam files are actually identical with respect to their seqinfo.
> And it's better than the previous behavior.
>
> Thanks for the suggestion,
>
> Martin
>
> On 04/25/2014 03:58 PM, Ryan C. Thompson wrote:
>
>> Hi all,
>>
>> I noticed that the seqinfo works on BamFile objects, but not on
>> BamFileList
>> objects. For BamFileList, it does not throw an error, but rather uses the
>> inherited method for "List", which does not return a useful result for
>> BamFileList. I suggest the following implementation of a useful seqinfo
>> function
>> for BamFileList, along with some code demonstrating the problem:
>>
>> fl <- BamFile(system.file("extdata", "ex1.bam", package="Rsamtools",
>> mustWork=TRUE))
>> # This works.
>> seqinfo(fl)
>>
>> fll <- BamFileList(fl, fl)
>> ## This works, but it uses the generic for "List" which gives a useless
>> result.
>> seqinfo(fll)
>>
>> ## Now add a method
>> setMethod("seqinfo", signature=list(x="BamFileList"), function (x)
>> {
>>      Reduce(merge, lapply(x, seqinfo))
>> }
>> )
>>
>> ## Now this returns a good result
>> seqinfo(fll)
>> ## So does this
>> seqlengths(fll)
>>
>>
>>
>> -Ryan Thompson
>>
>> _______________________________________________
>> Bioc-devel@r-project.org mailing list
>> https://stat.ethz.ch/mailman/listinfo/bioc-devel
>>
>
>
> --
> Computational Biology / Fred Hutchinson Cancer Research Center
> 1100 Fairview Ave. N.
> PO Box 19024 Seattle, WA 98109
>
> Location: Arnold Building M1 B861
> Phone: (206) 667-2793
>

        [[alternative HTML version deleted]]

_______________________________________________
Bioc-devel@r-project.org mailing list
https://stat.ethz.ch/mailman/listinfo/bioc-devel

Reply via email to