@@ -678,6 +678,141 @@ def get_arrays_from_intent(self, intent):
678678 it = intent_codes .code [intent ]
679679 return [x for x in self .darrays if x .intent == it ]
680680
681+ def agg_data (self , intent_code = None ):
682+ """
683+ Aggregate GIFTI data arrays into an ndarray or tuple of ndarray
684+
685+ In the general case, the numpy data array is extracted from each ``GiftiDataArray``
686+ object and returned in a ``tuple``, in the order they are found in the GIFTI image.
687+
688+ If all ``GiftiDataArray`` s have ``intent`` of 2001 (``NIFTI_INTENT_TIME_SERIES``),
689+ then the data arrays are concatenated as columns, producing a vertex-by-time array.
690+ If an ``intent_code`` is passed, data arrays are filtered by the selected intents,
691+ before being aggregated.
692+ This may be useful for images containing several intents, or ensuring an expected
693+ data type in an image of uncertain provenance.
694+ If ``intent_code`` is a ``tuple``, then a ``tuple`` will be returned with the result of
695+ ``agg_data`` for each element, in order.
696+ This may be useful for ensuring that expected data arrives in a consistent order.
697+
698+ Parameters
699+ ----------
700+ intent_code : None, string, integer or tuple of strings or integers, optional
701+ code(s) specifying nifti intent
702+
703+ Returns
704+ -------
705+ tuple of ndarrays or ndarray
706+ If the input is a tuple, the returned tuple will match the order.
707+
708+ Examples
709+ --------
710+
711+ Consider a surface GIFTI file:
712+
713+ >>> import nibabel as nib
714+ >>> from nibabel.testing import test_data
715+ >>> surf_img = nib.load(test_data('gifti', 'ascii.gii'))
716+
717+ The coordinate data, which is indicated by the ``NIFTI_INTENT_POINTSET``
718+ intent code, may be retrieved using any of the following equivalent
719+ calls:
720+
721+ >>> coords = surf_img.agg_data('NIFTI_INTENT_POINTSET')
722+ >>> coords_2 = surf_img.agg_data('pointset')
723+ >>> coords_3 = surf_img.agg_data(1008) # Numeric code for pointset
724+ >>> print(np.array2string(coords, precision=3))
725+ [[-16.072 -66.188 21.267]
726+ [-16.706 -66.054 21.233]
727+ [-17.614 -65.402 21.071]]
728+ >>> np.array_equal(coords, coords_2)
729+ True
730+ >>> np.array_equal(coords, coords_3)
731+ True
732+
733+ Similarly, the triangle mesh can be retrieved using various intent
734+ specifiers:
735+
736+ >>> triangles = surf_img.agg_data('NIFTI_INTENT_TRIANGLE')
737+ >>> triangles_2 = surf_img.agg_data('triangle')
738+ >>> triangles_3 = surf_img.agg_data(1009) # Numeric code for pointset
739+ >>> print(np.array2string(triangles))
740+ [0 1 2]
741+ >>> np.array_equal(triangles, triangles_2)
742+ True
743+ >>> np.array_equal(triangles, triangles_3)
744+ True
745+
746+ All arrays can be retrieved as a ``tuple`` by omitting the intent
747+ code:
748+
749+ >>> coords_4, triangles_4 = surf_img.agg_data()
750+ >>> np.array_equal(coords, coords_4)
751+ True
752+ >>> np.array_equal(triangles, triangles_4)
753+ True
754+
755+ Finally, a tuple of intent codes may be passed in order to select
756+ the arrays in a specific order:
757+
758+ >>> triangles_5, coords_5 = surf_img.agg_data(('triangle', 'pointset'))
759+ >>> np.array_equal(triangles, triangles_5)
760+ True
761+ >>> np.array_equal(coords, coords_5)
762+ True
763+
764+ The following image is a GIFTI file with ten (10) data arrays of the same
765+ size, and with intent code 2001 (``NIFTI_INTENT_TIME_SERIES``):
766+
767+ >>> func_img = nib.load(test_data('gifti', 'task.func.gii'))
768+
769+ When aggregating time series data, these arrays are concatenated into
770+ a single, vertex-by-timestep array:
771+
772+ >>> series = func_img.agg_data()
773+ >>> series.shape
774+ (642, 10)
775+
776+ In the case of a GIFTI file with unknown data arrays, it may be preferable
777+ to specify the intent code, so that a time series array is always returned:
778+
779+ >>> series_2 = func_img.agg_data('NIFTI_INTENT_TIME_SERIES')
780+ >>> series_3 = func_img.agg_data('time series')
781+ >>> series_4 = func_img.agg_data(2001)
782+ >>> np.array_equal(series, series_2)
783+ True
784+ >>> np.array_equal(series, series_3)
785+ True
786+ >>> np.array_equal(series, series_4)
787+ True
788+
789+ Requesting a data array from a GIFTI file with no matching intent codes
790+ will result in an empty tuple:
791+
792+ >>> surf_img.agg_data('time series')
793+ ()
794+ >>> func_img.agg_data('triangle')
795+ ()
796+ """
797+
798+ # Allow multiple intents to specify the order
799+ # e.g., agg_data(('pointset', 'triangle')) ensures consistent order
800+
801+ if isinstance (intent_code , tuple ):
802+ return tuple (self .agg_data (intent_code = code ) for code in intent_code )
803+
804+ darrays = self .darrays if intent_code is None else self .get_arrays_from_intent (intent_code )
805+ all_data = tuple (da .data for da in darrays )
806+ all_intent = {intent_codes .niistring [da .intent ] for da in darrays }
807+
808+ if all_intent == {'NIFTI_INTENT_TIME_SERIES' }: # stack when the gifti is a timeseries
809+ return np .column_stack (all_data )
810+
811+ if len (all_data ) == 1 :
812+ all_data = all_data [0 ]
813+
814+ return all_data
815+
681816 @deprecate_with_version (
682817 'getArraysFromIntent method deprecated. '
683818 "Use get_arrays_from_intent instead." ,
0 commit comments