RangedSummarizedExperiment.js

  1. import * as generics from "./AllGenerics.js";
  2. import * as gr from "./GRanges.js";
  3. import * as ggr from "./GroupedGRanges.js";
  4. import * as se from "./SummarizedExperiment.js";
  5. import * as utils from "./utils.js";
  6. import * as cutils from "./clone-utils.js";
  7. /**
  8. * A RangedSummarizedExperiment is a {@linkplain SummarizedExperiment} subclass where each row represents a genomic interval.
  9. * As such, it stores an additional {@linkplain GRanges} or {@linkplain GroupedGRanges} of length equal to the number of rows,
  10. * where each element represents the genomic range(s) for the corresponding row of the SummarizedExperiment.
  11. *
  12. * The RangedSummarizedExperiment supports the same set of generics as the {@linkplain SummarizedExperiment}.
  13. * Each method will call the base method, with the following extensions:
  14. *
  15. * - {@linkcode SLICE_2D} will additionally slice the supplied genomic ranges by the desired `rows`.
  16. * - {@linkcode COMBINE_ROWS} will combine genomic ranges across objects.
  17. * If some objects contain a GroupedGRanges and other objects contain GRanges, the latter will be coerced to a GroupedGRanges (where each group contains one range) before combining.
  18. * If any object is a base SummarizedExperiment, a GroupedGRanges containing zero-length groups will be automatically constructed to attempt combining.
  19. * - {@linkcode COMBINE_COLUMNS} will use the genomic ranges from the first object.
  20. *
  21. * @extends SummarizedExperiment
  22. */
  23. export class RangedSummarizedExperiment extends se.SummarizedExperiment {
  24. #check_rowRanges(x) {
  25. if (!(x instanceof gr.GRanges) && !(x instanceof ggr.GroupedGRanges)) {
  26. throw new Error("'rowRanges' should be a 'GRanges' or 'GroupedGRanges' instance");
  27. }
  28. if (generics.LENGTH(x) !== this._rowData.numberOfRows()) {
  29. throw utils.formatLengthError("'rowRanges'", "the number of rows");
  30. }
  31. }
  32. /**
  33. * @param {Object} assays - Object where keys are the assay names and values are multi-dimensional arrays of experimental data.
  34. * All arrays should have the same number of rows and columns.
  35. * @param {?(GRanges|GroupedGRanges)} rowRanges - Genomic ranges corresponding to each row.
  36. *
  37. * Alternatively, each row may correspond to a group of genomic ranges.
  38. *
  39. * If `null`, a {@linkplain GroupedGRanges} is constructed where each row corresponds to one group of ranges of zero length.
  40. * @param {Object} [options={}] - Optional parameters, including those used in the {@linkplain SummarizedExperiment} constructor.
  41. */
  42. constructor(assays, rowRanges, options = {}) {
  43. if (arguments.length == 0) {
  44. super();
  45. return;
  46. }
  47. super(assays, options);
  48. if (rowRanges === null) {
  49. rowRanges = ggr.GroupedGRanges.empty(this.numberOfRows());
  50. } else {
  51. this.#check_rowRanges(rowRanges);
  52. }
  53. this._rowRanges = rowRanges;
  54. return;
  55. }
  56. /**************************************************************************
  57. **************************************************************************
  58. **************************************************************************/
  59. /**
  60. * @return {GRanges} Genomic ranges corresponding to each row.
  61. */
  62. rowRanges() {
  63. return this._rowRanges;
  64. }
  65. /**************************************************************************
  66. **************************************************************************
  67. **************************************************************************/
  68. /**
  69. * @param {GRanges} value - Genomic ranges corresponding to each row.
  70. * This should have length equal to the number of rows in this RangedSummarizedExperiment.
  71. * @param {Object} [options={}] - Optional parameters.
  72. * @param {boolean} [options.inPlace=false] - Whether to mutate this Annotated instance in place.
  73. * If `false`, a new instance is returned.
  74. *
  75. * @return {RangedSummarizedExperiment} The RangedSummarizedExperiment after modifying its `rowRanges`.
  76. * If `inPlace = true`, this is a reference to the current instance, otherwise a new instance is created and returned.
  77. */
  78. setRowRanges(value, { inPlace = false } = {}) {
  79. this.#check_rowRanges(value);
  80. let target = cutils.setterTarget(this, inPlace);
  81. target._rowRanges = value;
  82. return target;
  83. }
  84. /**
  85. * @param {GRanges} value - Genomic ranges corresponding to each row.
  86. * This should have length equal to the number of rows in this RangedSummarizedExperiment.
  87. * @return {RangedSummarizedExperiment} A reference to this RangedSummarizedExperiment after modifying its `rowRanges`.
  88. */
  89. $setRowRanges(value) {
  90. return this.setRowRanges(value, { inPlace: true });
  91. }
  92. /**************************************************************************
  93. **************************************************************************
  94. **************************************************************************/
  95. _bioconductor_SLICE_2D(output, rows, columns, { allowView = false }) {
  96. super._bioconductor_SLICE_2D(output, rows, columns, { allowView });
  97. if (rows !== null) {
  98. output._rowRanges = generics.SLICE(this._rowRanges, rows);
  99. } else {
  100. output._rowRanges = this._rowRanges;
  101. }
  102. }
  103. _bioconductor_COMBINE_ROWS(output, objects) {
  104. super._bioconductor_COMBINE_ROWS(output, objects);
  105. let collected = [];
  106. let has_empty = false;
  107. let has_ggr = false;
  108. for (var i = 0; i < objects.length; i++) {
  109. let x = objects[i];
  110. if (x instanceof RangedSummarizedExperiment) {
  111. let y = x._rowRanges;
  112. if (y instanceof ggr.GroupedGRanges) {
  113. has_ggr = true;
  114. }
  115. collected.push(y);
  116. } else if (x instanceof se.SummarizedExperiment) {
  117. has_empty = true;
  118. collected.push(null);
  119. } else {
  120. throw new Error("objects to be combined must be SummarizedExperiments (failing for object " + String(i) + ")");
  121. }
  122. }
  123. // Promoting nulls and GRanges to GroupedGRanges, if necessary.
  124. if (has_empty || has_ggr) {
  125. for (var i = 0; i < collected.length; i++) {
  126. let current = collected[i];
  127. if (current instanceof gr.GRanges) {
  128. let widths = new Int32Array(generics.LENGTH(current));
  129. widths.fill(1);
  130. let options = {
  131. rangeLengths: widths,
  132. names: current.names(),
  133. elementMetadata: current.elementMetadata(),
  134. metadata: current.metadata()
  135. };
  136. if (options.names !== null) {
  137. current = current.setNames(null);
  138. }
  139. if (options.elementMetadata.metadata().size > 0 || options.elementMetadata.numberOfColumns() > 0) {
  140. current = current.setElementMetadata(null);
  141. }
  142. if (options.metadata.size > 0) {
  143. current = current.setMetadata(new Map);
  144. }
  145. collected[i] = new ggr.GroupedGRanges(current, options);
  146. } else if (current === null){
  147. collected[i] = ggr.GroupedGRanges.empty(objects[i].numberOfRows());
  148. }
  149. }
  150. }
  151. output._rowRanges = generics.COMBINE(collected);
  152. return;
  153. }
  154. _bioconductor_COMBINE_COLUMNS(output, objects) {
  155. super._bioconductor_COMBINE_COLUMNS(output, objects);
  156. output._rowRanges = objects[0]._rowRanges;
  157. return;
  158. }
  159. _bioconductor_CLONE(output, { deepCopy }) {
  160. super._bioconductor_CLONE(output, { deepCopy });
  161. output._rowRanges = cutils.cloneField(this._rowRanges, deepCopy);
  162. return;
  163. }
  164. }