ISMRM Raw Data Format (ISMRMRD)

Contents:

Preamble

A prerequisite for sharing magnetic resonance (imaging) reconstruction algorithms and code is a common raw data format. This document describes such a common raw data format and attempts to capture the data fields that are require to describe enough details about the magnetic resonance experiment to reconstruct images from the data. This standard was developed by a subcommittee of the ISMRM Sedona 2013 workshop. Comments and requests for additions/modifications can be sent to:

  • Michael S. Hansen (michael DOT hansen AT nih DOT gov)
  • Wally Block (wblock AT cae DOT wisc DOT edu)
  • Mark Griswold (mag46 AT case DOT edu)
  • Brian Hargreaves (bah AT stanford DOT edu)
  • Peter Boernert (peter DOT boernert AT philips DOT com)
  • Sebastian Kozerke (kozerke AT biomed DOT ee DOT ethz DOT ch)
  • Craig Meyer (cmeyer AT virginia DOT edu)
  • Doug Noll (dnoll AT umich DOT edu)
  • Jim Pipe (Jim.Pipe AT DignityHealth DOT org)

Developers/Contributors

  • Michael S. Hansen, National Institutes of Health, USA
  • Nick Zwart, Barrow Neurological Institute, Phoenix, Arizona
  • Souheil Inati, National Institutes of Health, USA
  • Joe Naegele, National Institutes of Health, USA
  • Kaveh Vahedipour, Juelich Research Centre, Juelich, Germany

Obtaining and Installing

The source code, examples, and example datasets can be found on the ISMRM Raw Data Github website.

To download the source code, clone the git archive:

git clone https://github.com/ismrmrd/ismrmrd

API Documentation can be found at https://ismrmrd.github.io/ismrmrd/api/.

Dependencies

The ISMRM Raw Data format is described by an XML schema and some C-style structs with fixed memory layout and as such does not have dependencies. However, it uses HDF5 files for storage and a C++ library for reading and writing the ISMRMRD files is included in this distribution. Furthermore, since the XML header is defined with an XML schema, we encourage using XML data binding when writing software using the format. To compile all components of this distribution you need:

Note

It is only necessary to install the dependencies if you wish to develop compiled C/C++ software, which uses the ISMRMRD format. The format can be read in Matlab without installing any additional software.

Linux installation

The dependencies mentioned above should be included in most linux distributions. On Ubuntu you can install all required dependencies with:

sudo apt-get install libhdf5-serial-dev h5utils cmake cmake-curses-gui libboost-all-dev doxygen git

After installation of dependencies, the library can be installed with:

git clone https://github.com/ismrmrd/ismrmrd
cd ismrmrd/
mkdir build
cd build
cmake ../
make
sudo make install

This will install the library in /usr/local/ by default. To specify an alternative installation directory, pass -D CMAKE_INSTALL_PREFIX=<install dir> to cmake.

Mac OSX Installation

There are numerous different package management systems for Mac. In this example, we have used Homebrew (http://brew.sh/). First install the dependencies:

brew install wget hdf5 boost cmake doxygen fftw

Then download and compile:

git clone https://github.com/ismrmrd/ismrmrd
cd ismrmrd
mkdir build
cd build/
cmake ../
make
make install

This will install the library in /usr/local/ by default. To specify an alternative installation directory, pass -D CMAKE_INSTALL_PREFIX=<install dir> to cmake.

Windows Installation

Setting up a Windows development environment is usually a bit more challenging than working on Unix platforms where most library dependencies are easily installed with package management systems (see above). The general Windows installation instructions (you may have to make adjustments for your setup) is as follows:

This can seem a bit daunting, we have included a Windows powershell script, which you can use to guide you through the installation process.

After installing all dependencies, download the code, e.g. from a git bash shell:

git clone https://github.com/ismrmrd/ismrmrd
cd ismrmrd/
mkdir build
cd build/
cmake-gui.exe

Last command will open CMake’s graphical user interface. Hit the configure button and deal with the dependencies that CMake is unable to find. Hit configure again and repeat the process until CMake has enough information to configure. Once the configuration is complete, you can hit generate to generate a Visual Studio project, which you can open and use to build ISMRMRD. There are step-by-step commands included in the powershell script below to guide you through the CMake configuration and build process from the command line. The command line CMake configuration line (assuming you have installed with the paths above), would look something like (backslashes are just there to break the command over multiple lines):

cmake -G"Visual Studio 10 Win64" \
    -DBOOST_ROOT=C:/Program Files/boost/boost_1_51 \
    -DFFTW3_INCLUDE_DIR=C:/MRILibraries/fftw3 \
    -DFFTW3F_LIBRARY=C:/MRILibraries/fftw3/libfftw3f-3.lib ../

Again, you may have to adjust for your specific installation paths. After generating the Visual Studio project, you can build from a Visual Studio Command Prompt with:

msbuild .\ISMRMRD.sln /p:Configuration=Release

Overview

The raw data format combines a mix of flexible data structures (XML header) and fixed structures (equivalent to C-structs). A raw data set consist mainly of 2 sections:

  1. A flexible XML format document that can contain an arbitrary number of fields and accommodate everything from simple values (b-values, etc.) to entire vendor protocols, etc. This purpose of this XML document is to provide parameters that may be meaningful for some experiments but not for others. This XML format is defined by an XML Schema Definition file (ismrmrd.xsd).
  2. Raw data section. This section contains all the acquired data in the experiment. Each data item is preceded by a C-struct with encoding numbers, etc. Following this data header is a channel header and data for each acquired channel. The raw data headers are defined in a C/C++ header file (ismrmrd.h)

In addition to these sections, the ISMRMRD format also specifies an image header for storing reconstructed images and the accompanying C++ library provides a convenient way of writing such images into HDF5 files along with generic arrays for storing less well defined data structures, e.g. coil sensitivity maps or other calibration data.

Flexible Data Header

The flexible data structure is defined by the xml schema definition in schema/ismrmrd.xsd (schema is included in appendix below).

An example of an XML file for a Cartesian 3D acquisition could look like:

<?xml version="1.0"?>
<ismrmrdHeader xmlns="http://www.ismrm.org/ISMRMRD" xmlns:xsi="http://www.w3.org/2001/XMLSchema-instance" xmlns:xs="http://www.w3.org/2001/XMLSchema" xsi:schemaLocation="http://www.ismrm.org/ISMRMRD ismrmrd.xsd">
  <subjectInformation>
    <patientName>phantom</patientName>
    <patientWeight_kg>70.3068</patientWeight_kg>
  </subjectInformation>
  <acquisitionSystemInformation>
    <systemVendor>SIEMENS</systemVendor>
    <systemModel>Avanto</systemModel>
    <systemFieldStrength_T>1.494</systemFieldStrength_T>
    <receiverChannels>32</receiverChannels>
    <relativeReceiverNoiseBandwidth>0.79</relativeReceiverNoiseBandwidth>
  </acquisitionSystemInformation>
  <experimentalConditions>
    <H1resonanceFrequency_Hz>63642459</H1resonanceFrequency_Hz>
  </experimentalConditions>
  <encoding>
    <trajectory>cartesian</trajectory>
    <encodedSpace>
      <matrixSize>
        <x>256</x>
        <y>140</y>
        <z>80</z>
      </matrixSize>
      <fieldOfView_mm>
        <x>600</x>
        <y>328.153125</y>
        <z>160</z>
      </fieldOfView_mm>
    </encodedSpace>
    <reconSpace>
      <matrixSize>
        <x>128</x>
        <y>116</y>
        <z>64</z>
      </matrixSize>
      <fieldOfView_mm>
        <x>300</x>
        <y>271.875</y>
        <z>128</z>
      </fieldOfView_mm>
    </reconSpace>
    <encodingLimits>
      <kspace_encoding_step_1>
        <minimum>0</minimum>
        <maximum>83</maximum>
        <center>28</center>
      </kspace_encoding_step_1>
      <kspace_encoding_step_2>
        <minimum>0</minimum>
        <maximum>45</maximum>
        <center>20</center>
      </kspace_encoding_step_2>
      <slice>
        <minimum>0</minimum>
        <maximum>0</maximum>
        <center>0</center>
      </slice>
      <set>
        <minimum>0</minimum>
        <maximum>0</maximum>
        <center>0</center>
      </set>
    </encodingLimits>
  </encoding>
  <parallelImaging>
    <accelerationFactor>
      <kspace_encoding_step_1>1</kspace_encoding_step_1>
      <kspace_encoding_step_2>1</kspace_encoding_step_2>
    </accelerationFactor>
    <calibrationMode>other</calibrationMode>
  </parallelImaging>
  <sequenceParameters>
    <TR>4.6</TR>
    <TE>2.35</TE>
    <TI>300</TI>
  </sequenceParameters>
</ismrmrdHeader>

The most critical elements for image reconstruction are contained in the <encoding> section of the document, which describes the encoded spaced and also the target reconstructed space. Along with the <encodingLimits> this section allows the reconstruction program to determine matrix sizes, oversampling factors, partial Fourier, etc. In the example above, data is acquired with two-fold oversampling in the read-out (x) direction, which is reflected in the larger matrix size in the encoded space compared to the reconstruction space. The field of view is also twice as large in the encoded space. For the first phase encoding dimension (y), we have a combination of oversampling (20%), reduced phase resolution (only 83 lines of k-space acquired, and partial Fourier sampling, which is reflected in the asymmetric center of the encoding limits of the <kspace_encoding_step_1>. Specifically, the data lines would be placed into the encoding space like this:

0                                     70                                         139
|-------------------------------------|-------------------------------------------|
                      ****************************************************
                      ^               ^                                  ^
                      0              28                                  83

After FFT, only the central 116 lines are kept, i.e. there is a reduced field of view in the phase encoding direction. Center and encoding limits for the readout dimension is not given in the XML header. This is to accommodate sequences where the center of the readout may change from readout to readout (alternating directions of readout). There is a field on the individual data headers (see below) to indicate the center of the readout.

An experiment can have multiple encoding spaces and it is possible to indicate on each acquired data readout, which encoding space the data belongs to (see below).

In addition to the defined field in the xml header, it is possible to add an arbitrary number of user defined parameters to accommodate special sequence parameters. Please consult the xml schema to see how user parameters are defined. Briefly, the XML header can have a section at the end which looks like:

<userParameters>
  <userParameterLong>
    <name>MyVar1</name><value>1003</value>
  </userParameterLong>
  <userParameterLong>
    <name>MyVar2</name><value>1999</value>
  </userParameterLong>
  <userParameterDouble>
    <name>MyDoubleVar</name><value>87.6676</value>
  </userParameterDouble>
</userParameters>

Fixed Data structures

Each raw data acquisition is preceded by the following fixed layout structure:

    uint16_t version;                                    /**< First unsigned int indicates the version */
    uint64_t flags;                                      /**< bit field with flags */
    uint32_t measurement_uid;                            /**< Unique ID for the measurement */
    uint32_t scan_counter;                               /**< Current acquisition number in the measurement */
    uint32_t acquisition_time_stamp;                     /**< Acquisition clock */
    uint32_t physiology_time_stamp[ISMRMRD_PHYS_STAMPS]; /**< Physiology time stamps, e.g. ecg, breating, etc. */
    uint16_t number_of_samples;                          /**< Number of samples acquired */
    uint16_t available_channels;                         /**< Available coils */
    uint16_t active_channels;                            /**< Active coils on current acquisiton */
    uint64_t channel_mask[ISMRMRD_CHANNEL_MASKS];        /**< Mask to indicate which channels are active. Support for 1024 channels */
    uint16_t discard_pre;                                /**< Samples to be discarded at the beginning of  acquisition */
    uint16_t discard_post;                               /**< Samples to be discarded at the end of acquisition */
    uint16_t center_sample;                              /**< Sample at the center of k-space */
    uint16_t encoding_space_ref;                         /**< Reference to an encoding space, typically only one per acquisition */
    uint16_t trajectory_dimensions;                      /**< Indicates the dimensionality of the trajectory vector (0 means no trajectory) */
    float sample_time_us;                                /**< Time between samples in micro seconds, sampling BW */
    float position[3];                                   /**< Three-dimensional spatial offsets from isocenter */
    float read_dir[3];                                   /**< Directional cosines of the readout/frequency encoding */
    float phase_dir[3];                                  /**< Directional cosines of the phase */
    float slice_dir[3];                                  /**< Directional cosines of the slice direction */
    float patient_table_position[3];                     /**< Patient table off-center */
    ISMRMRD_EncodingCounters idx;                        /**< Encoding loop counters, see above */
    int32_t user_int[ISMRMRD_USER_INTS];                 /**< Free user parameters */
    float user_float[ISMRMRD_USER_FLOATS];               /**< Free user parameters */

Where EncodingCounters are defined as:

    uint16_t kspace_encode_step_1;    /**< e.g. phase encoding line number */
    uint16_t kspace_encode_step_2;    /**< e.g. partition encoding number */
    uint16_t average;                 /**< e.g. signal average number */
    uint16_t slice;                   /**< e.g. imaging slice number */
    uint16_t contrast;                /**< e.g. echo number in multi-echo */
    uint16_t phase;                   /**< e.g. cardiac phase number */
    uint16_t repetition;              /**< e.g. dynamic number for dynamic scanning */
    uint16_t set;                     /**< e.g. flow encoding set */
    uint16_t segment;                 /**< e.g. segment number for segmented acquisition */
    uint16_t user[ISMRMRD_USER_INTS]; /**< Free user parameters */

The interpretation of some of these fields may vary from sequence to sequence, i.e. for a Cartesian sequence, kspace_encode_step_1 would be the phase encoding step, for a spiral sequence where phase encoding direction does not make sense, it would be the spiral interleave number. The encoding_space_ref enables the user to tie an acquisition to a specific encoding space (see above) in case there are multiple, e.g. in situations where a calibration scan may be integrated in the acquisition.

The flags field is a bit mask, which in principle can be used freely by the user, but suggested flag values are given in ismrmrd.h, it is recommended not to use already designated flag bits for custom purposes. There are a set of bits reserved for prototyping (bits 57-64), please see ismrmrd.h for details.

The header contains a trajectory_dimensions field. If the value of this field is larger than 0, it means that trajectories are stored with each individual acquisition. For a 2D acquisition, the trajectory_dimensions would typically be 2 and the convention (for memory layout) is that the header is followed immediately by the trajectory before the complex data. There is an example of how this memory layout could be implemented with a C++ class in the ismrmrd.h file:

class Acquisition
{

//....

AcquisitionHeader head_; //Header, see above

float* traj_;            //Trajectory, elements = head_.trajectory_dimensions*head_.number_of_samples
                         //   [kx,ky,kx,ky.....]        (for head_.trajectory_dimensions = 2)
                         //   [kx,ky,kz,kx,ky,kz,.....] (for head_.trajectory_dimensions = 3)

float* data_;            //Actual data, elements = head_.number_of_samples*head_.active_channels*2
                         //   [re,im,re,im,.....,re,im,re,im,.....,re,im,re,im,.....]
                         //    ---channel 1-------channel 2---------channel 3-----

};

This suggested memory layout is only a suggestion. The HDF5 interface (see below) can be used to read the data into many different data structures. In fact, the user can choose to read only part of the header or not read the data, etc.

As mentioned above, the ISMRMRD format also suggests a way to store reconstructed images (or maybe image data used for calibration). An ImageHeader structure is defined in ismrmrd.h:

    uint16_t version;                                    /**< First unsigned int indicates the version */
    uint16_t data_type;                                  /**< e.g. unsigned short, float, complex float, etc. */
    uint64_t flags;                                      /**< bit field with flags */
    uint32_t measurement_uid;                            /**< Unique ID for the measurement  */
    uint16_t matrix_size[3];                             /**< Pixels in the 3 spatial dimensions */
    float field_of_view[3];                              /**< Size (in mm) of the 3 spatial dimensions */
    uint16_t channels;                                   /**< Number of receive channels */
    float position[3];                                   /**< Three-dimensional spatial offsets from isocenter */
    float read_dir[3];                                   /**< Directional cosines of the readout/frequency encoding */
    float phase_dir[3];                                  /**< Directional cosines of the phase */
    float slice_dir[3];                                  /**< Directional cosines of the slice direction */
    float patient_table_position[3];                     /**< Patient table off-center */
    uint16_t average;                                    /**< e.g. signal average number */
    uint16_t slice;                                      /**< e.g. imaging slice number */
    uint16_t contrast;                                   /**< e.g. echo number in multi-echo */
    uint16_t phase;                                      /**< e.g. cardiac phase number */
    uint16_t repetition;                                 /**< e.g. dynamic number for dynamic scanning */
    uint16_t set;                                        /**< e.g. flow encodning set */
    uint32_t acquisition_time_stamp;                     /**< Acquisition clock */
    uint32_t physiology_time_stamp[ISMRMRD_PHYS_STAMPS]; /**< Physiology time stamps, e.g. ecg, breathing, etc. */
    uint16_t image_type;                                 /**< e.g. magnitude, phase, complex, real, imag, etc. */
    uint16_t image_index;                                /**< e.g. image number in series of images  */
    uint16_t image_series_index;                         /**< e.g. series number */
    int32_t user_int[ISMRMRD_USER_INTS];                 /**< Free user parameters */
    float user_float[ISMRMRD_USER_FLOATS];               /**< Free user parameters */
    uint32_t attribute_string_len;                       /**< Length of attributes string */

In a similar fashion to the raw data acquisition data, the intention is to store a header followed by the image data. Since image data can be in several different format (e.g. float, complex, etc.), the memory layout is less well defined but can be described as:

template <typename T> class Image {

ImageHeader head_;     //ImageHeader as defined above
T* data_;              //Data, array of size (matrix_size[0]*matrix_size[1]*matrix_size[2]*channels),
                       //first spatial dimension is fastest changing array index, channels outer most (slowest changing).
};

File Storage

The ISMRM Raw Data format is stored in HDF5 format. Details on this format can be found at the HDF5 website. Briefly it is a hierarchical data format (much like a file system), which can contain multiple variable organized in groups (like folders in a file system). The variables can contain arrays of data values, custom defined structs, or simple text fields. It is the convention (but not a requirement) that the ISMRMRD datasets are stored in a group called /dataset. The XML configuration is stored in the variable /dataset/xml and the data is stored in /dataset/data. HDF5 files can be viewed with the HDFView application which is available on the HDF5 website for multiple platforms. Files can also be read directly in Matlab, in fact Matlab uses (since file format v7.3) HDF5 as its internal data format in the .mat files. As an example the data from an ISMRMRD file with name myfile.h5 can be read in matlab with a command like:

>> data = h5read('simple_gre.h5', '/dataset/data');
>> data

data =

head: [1x1 struct]
traj: {1x1281 cell}
data: {1x1281 cell}

 >> data.head

 ans =

                version: [1x1281 uint16]
                  flags: [1x1281 uint64]
        measurement_uid: [1x1281 uint32]
           scan_counter: [1x1281 uint32]
 acquisition_time_stamp: [1x1281 uint32]
  physiology_time_stamp: [3x1281 uint32]
      number_of_samples: [1x1281 uint16]
     available_channels: [1x1281 uint16]
        active_channels: [1x1281 uint16]
           channel_mask: [16x1281 uint64]
            discard_pre: [1x1281 uint16]
           discard_post: [1x1281 uint16]
          center_sample: [1x1281 uint16]
     encoding_space_ref: [1x1281 uint16]
  trajectory_dimensions: [1x1281 uint16]
         sample_time_us: [1x1281 single]
               position: [3x1281 single]
               read_dir: [3x1281 single]
              phase_dir: [3x1281 single]
              slice_dir: [3x1281 single]
 patient_table_position: [3x1281 single]
                    idx: [1x1 struct]
               user_int: [8x1281 int32]
             user_float: [8x1281 single]

 >>

The HDF5 file format can be access from C, C++, and java using the libraries provided on the HDF5 website. The ISMRMRD distribution also comes with some C++ wrappers that can be used for easy access (read and write) from C++ programs. See below.

In addition to storing acquisition data and images as defined by the headers above, the HDF5 format also enables storage of generic multi-dimensional arrays. The ISMRMRD format does not explicitly define how such data should be stored, but leaves it open for the user to add variables and data as dictated by a given application.

C++ Support Library

To enable easy prototyping of C++ software using the ISMRMRD data format, a simple C++ wrapper class is provided (defined in dataset.h):

class EXPORTISMRMRD Dataset {
public:
    // Constructor and destructor
    Dataset(const char* filename, const char* groupname, bool create_file_if_needed = true);
    ~Dataset();
    
    // Methods
    // XML Header
    void writeHeader(const std::string &xmlstring);
    void readHeader(std::string& xmlstring);
    // Acquisitions
    void appendAcquisition(const Acquisition &acq);
    void readAcquisition(uint32_t index, Acquisition &acq);
    uint32_t getNumberOfAcquisitions();
    // Images
    template <typename T> void appendImage(const std::string &var, const Image<T> &im);
    void appendImage(const std::string &var, const ISMRMRD_Image *im);
    template <typename T> void readImage(const std::string &var, uint32_t index, Image<T> &im);
    uint32_t getNumberOfImages(const std::string &var);
    // NDArrays
    template <typename T> void appendNDArray(const std::string &var, const NDArray<T> &arr);
    void appendNDArray(const std::string &var, const ISMRMRD_NDArray *arr);
    template <typename T> void readNDArray(const std::string &var, uint32_t index, NDArray<T> &arr);
    uint32_t getNumberOfNDArrays(const std::string &var);

protected:
    ISMRMRD_Dataset dset_;

Using this wrapper, C++ applications can be programmed as:

// Open dataset
ISMRMRD::Dataset d(datafile.c_str(), "dataset", false);

std::string xml;
d.readHeader(xml);
ISMRMRD::IsmrmrdHeader hdr;
ISMRMRD::deserialize(xml.c_str(),hdr);

// Do something with the header

unsigned int number_of_acquisitions = d.getNumberOfAcquisitions();
ISMRMRD::Acquisition acq;
for (unsigned int i = 0; i < number_of_acquisitions; i++) {
    // Read one acquisition at a time
    d.readAcquisition(i, acq);

    // Do something with the data
}

Since the XML header is defined in the schema/ismrmrd.xsd file, it can be parsed with numerous xml parsing libraries. The ISMRMRD library includes an API that allows for programmatically deserializing, manipulating, and serializing the XML header. See the code in the utilities directory for examples of how to use the XML API.

C++ Example Applications

The distribution includes two example applications, one that creates a simple 2D single-channel dataset from scratch and one that reconstructs this dataset (you need FFTW installed to compile these test applications). The data generation application looks like this (generate_cartesian_shepp_logan.cpp):

int main(int argc, char** argv)
{
	/** TODO
	 *
	 *  Noise samples
	 *  Acceleration
	 *  k-space coordinates
	 *
	 */

	unsigned int matrix_size; //Matrix size
	unsigned int ncoils;      //Number of coils
	unsigned int ros;           //Readout ovesampling
	unsigned int repetitions;
	unsigned int acc_factor;
	float noise_level;
	std::string outfile;
	std::string dataset;
	bool store_coordinates;
	bool noise_calibration;

	po::options_description desc("Allowed options");
	desc.add_options()
	    ("help,h", "produce help message")
	    ("matrix,m", po::value<unsigned int>(&matrix_size)->default_value(256), "Matrix Size")
	    ("coils,c", po::value<unsigned int>(&ncoils)->default_value(8), "Number of Coils")
	    ("oversampling,O", po::value<unsigned int>(&ros)->default_value(2), "Readout oversampling")
	    ("repetitions,r", po::value<unsigned int>(&repetitions)->default_value(1), "Repetitions")
	    ("acceleration,a", po::value<unsigned int>(&acc_factor)->default_value(1), "Acceleration factor")
	    ("noise-level,n", po::value<float>(&noise_level)->default_value(0.05f,"0.05"), "Noise Level")
	    ("output,o", po::value<std::string>(&outfile)->default_value("testdata.h5"), "Output File Name")
	    ("dataset,d", po::value<std::string>(&dataset)->default_value("dataset"), "Output Dataset Name")
	    ("noise-calibration,C", po::value<bool>(&noise_calibration)->zero_tokens(), "Add noise calibration")
	    ("k-coordinates,k",  po::value<bool>(&store_coordinates)->zero_tokens(), "Store k-space coordinates")
	;

	po::variables_map vm;
	po::store(po::parse_command_line(argc, argv, desc), vm);
	po::notify(vm);

	if (vm.count("help")) {
	    std::cout << desc << "\n";
	    return 1;
	}

	std::cout << "Generating Cartesian Shepp Logan Phantom!!!" << std::endl;
	std::cout << "Acceleration: " << acc_factor << std::endl;

	boost::shared_ptr<NDArray<complex_float_t> > phantom = shepp_logan_phantom(matrix_size);
	boost::shared_ptr<NDArray<complex_float_t> > coils = generate_birdcage_sensititivies(matrix_size, ncoils, 1.5);

	std::vector<size_t> dims;
	dims.push_back(matrix_size*ros); //oversampling in the readout direction
	dims.push_back(matrix_size);
	dims.push_back(ncoils);

	NDArray<complex_float_t> coil_images(dims);
	memset(coil_images.getDataPtr(), 0, coil_images.getDataSize());

	for (unsigned int c = 0; c < ncoils; c++) {
            for (unsigned int y = 0; y < matrix_size; y++) {
                for (unsigned int x = 0; x < matrix_size; x++) {
                    uint16_t xout = x + (matrix_size*ros-matrix_size)/2;
                    coil_images(xout,y,c) = (*phantom)(x,y) * (*coils)(x,y,c);
                }
            }
	}

        //Let's append the data to the file
        //Create if needed
	Dataset d(outfile.c_str(),dataset.c_str(), true);
	Acquisition acq;
        size_t readout = matrix_size*ros;
        
	if (noise_calibration)
        {
            acq.resize(readout, ncoils);
            memset((void *)acq.getDataPtr(), 0, acq.getDataSize());
            acq.setFlag(ISMRMRD_ACQ_IS_NOISE_MEASUREMENT);
            add_noise(acq,noise_level);
            acq.sample_time_us() = 5.0;
            d.appendAcquisition(acq);
	}
        
        if (store_coordinates) {
            acq.resize(readout, ncoils, 2);
        }
        else {
            acq.resize(readout, ncoils);
        }
        memset((void*)acq.getDataPtr(), 0, acq.getDataSize());
        
        acq.available_channels() = ncoils;
	acq.center_sample() = (readout>>1);


        for (unsigned int r = 0; r < repetitions; r++) {
            for (unsigned int a = 0; a < acc_factor; a++) {
                NDArray<complex_float_t> cm = coil_images;
                fft2c(cm);

                add_noise(cm,noise_level);
                for (size_t i = a; i < matrix_size; i+=acc_factor) {
                    acq.clearAllFlags();
                    
                    //Set some flags
                    if (i == a) {
                        acq.setFlag(ISMRMRD_ACQ_FIRST_IN_SLICE);
                    }
                    if (i >= (matrix_size-acc_factor)) {
                        acq.setFlag(ISMRMRD_ACQ_LAST_IN_SLICE);
                    }
                    acq.idx().kspace_encode_step_1 = i;
                    acq.idx().repetition = r*acc_factor + a;
                    acq.sample_time_us() = 5.0;
                    for (size_t c = 0; c < ncoils; c++) {
                        for (size_t s = 0; s < readout; s++) {
                            acq.data(s,c) = cm(s,i,c);
                        }
                    }
                    
                    if (store_coordinates) {
                        float ky = (1.0*i-(matrix_size>>1))/(1.0*matrix_size);
                        for (size_t x = 0; x < readout; x++) {
                            float kx = (1.0*x-(readout>>1))/(1.0*readout);
                            acq.traj(0,x) = kx;
                            acq.traj(1,x) = ky;
                        }
                    }
                    d.appendAcquisition(acq);
                }
            }
	}

	//Let's create a header, we will use the C++ classes in ismrmrd/xml.h
	IsmrmrdHeader h;
        h.version = ISMRMRD_XMLHDR_VERSION;
	h.experimentalConditions.H1resonanceFrequency_Hz = 63500000; //~1.5T        

	AcquisitionSystemInformation sys;
	sys.institutionName = "ISMRM Synthetic Imaging Lab";
	sys.receiverChannels = ncoils;
	h.acquisitionSystemInformation = sys;

	//Create an encoding section
        Encoding e;
        e.encodedSpace.matrixSize.x = readout;
        e.encodedSpace.matrixSize.y = matrix_size;
        e.encodedSpace.matrixSize.z = 1;
        e.encodedSpace.fieldOfView_mm.x = 600;
        e.encodedSpace.fieldOfView_mm.y = 300;
        e.encodedSpace.fieldOfView_mm.z = 6;
        e.reconSpace.matrixSize.x = readout/2;
        e.reconSpace.matrixSize.y = matrix_size;
        e.reconSpace.matrixSize.z = 1;
        e.reconSpace.fieldOfView_mm.x = 300;
        e.reconSpace.fieldOfView_mm.y = 300;
        e.reconSpace.fieldOfView_mm.z = 6;
        e.trajectory = "cartesian";
        e.encodingLimits.kspace_encoding_step_1 = Limit(0, matrix_size-1,(matrix_size>>1));
        e.encodingLimits.repetition = Limit(0, repetitions*acc_factor,0);
        
	//e.g. parallel imaging
	if (acc_factor > 1) {
            ParallelImaging parallel;
            parallel.accelerationFactor.kspace_encoding_step_1 = acc_factor;
            parallel.accelerationFactor.kspace_encoding_step_2 = 1;
            parallel.calibrationMode = "interleaved";
            e.parallelImaging = parallel;
	}

	//Add the encoding section to the header
	h.encoding.push_back(e);

	//Add any additional fields that you may want would go here....

	//Serialize the header
        std::stringstream str;
        ISMRMRD::serialize( h, str);
        std::string xml_header = str.str();
        //std::cout << xml_header << std::endl;
        
	//Write the header to the data file.
	d.writeHeader(xml_header);

        //Write out some arrays for convenience
        d.appendNDArray("phantom", *phantom);
        d.appendNDArray("csm", *coils);
        d.appendNDArray("coil_images", coil_images);
        
	return 0;
}

To reconstruct this synthetic dataset, you can use the test reconstruction application (recon_cartesian_2d.cpp):

int main(int argc, char** argv)
{
    if (argc < 2) {
        print_usage(argv[0]);
        return -1;
    }

    std::string datafile(argv[1]);

    std::cout << "Simple ISMRMRD Reconstruction program" << std::endl;
    std::cout << "   - filename: " << datafile << std::endl;

    //Let's open the existing dataset
    ISMRMRD::Dataset d(datafile.c_str(),"dataset", false);

    std::string xml;
    d.readHeader(xml);
    ISMRMRD::IsmrmrdHeader hdr;
    ISMRMRD::deserialize(xml.c_str(),hdr);

    //Let's print some information from the header
    if (hdr.version) {
        std::cout << "XML Header version: " << hdr.version << std::endl;
    }
    else {
        std::cout << "XML Header unspecified version." << std::endl;
    }
    
    if (hdr.encoding.size() != 1) {
        std::cout << "Number of encoding spaces: " << hdr.encoding.size() << std::endl;
        std::cout << "This simple reconstruction application only supports one encoding space" << std::endl;
        return -1;
    }

    ISMRMRD::EncodingSpace e_space = hdr.encoding[0].encodedSpace;
    ISMRMRD::EncodingSpace r_space = hdr.encoding[0].reconSpace;

    std::cout << "Encoding Matrix Size        : [" << e_space.matrixSize.x << ", " << e_space.matrixSize.y << ", " << e_space.matrixSize.z << "]" << std::endl;
    std::cout << "Reconstruction Matrix Size  : [" << r_space.matrixSize.x << ", " << r_space.matrixSize.y << ", " << r_space.matrixSize.z << "]" << std::endl;
    std::cout << "Number of acquisitions      : " << d.getNumberOfAcquisitions() << std::endl;

    if (e_space.matrixSize.z != 1) {
        std::cout << "This simple reconstruction application only supports 2D encoding spaces" << std::endl;
        return -1;
    }

    //Allocate a buffer for the data
    std::vector<size_t> dims;
    dims.push_back(e_space.matrixSize.x);
    dims.push_back(e_space.matrixSize.y);
    ISMRMRD::NDArray<complex_float_t> buffer(dims);
    
    //Now loop through and copy data
    unsigned int number_of_acquisitions = d.getNumberOfAcquisitions();
    ISMRMRD::Acquisition acq;
    for (unsigned int i = 0; i < number_of_acquisitions; i++) {
        //Read one acquisition at a time
        d.readAcquisition(i, acq);

        //Copy data, we should probably be more careful here and do more tests....
        //We are not considering multiple channels here.
        unsigned int offset = acq.idx().kspace_encode_step_1*dims[0];
        memcpy(&buffer.getDataPtr()[offset], acq.getDataPtr(),sizeof(complex_float_t)*dims[0]);
    }

    //Let's FFT the k-space to image
    fftwf_complex* tmp = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex)*buffer.getNumberOfElements());

    if (!tmp) {
        std::cout << "Error allocating temporary storage for FFTW" << std::endl;
        return -1;
    }

    //FFTSHIFT
    fftshift(reinterpret_cast<complex_float_t*>(tmp), buffer.getDataPtr(), dims[0], dims[1]);

    //Create the FFTW plan
    fftwf_plan p = fftwf_plan_dft_2d(dims[1], dims[0], tmp ,tmp, FFTW_BACKWARD, FFTW_ESTIMATE);

    //Execute the FFT
    fftwf_execute(p);

    //FFTSHIFT
    fftshift( buffer.getDataPtr(), reinterpret_cast<std::complex<float>*>(tmp), dims[0], dims[1]);

    //Clean up.
    fftwf_destroy_plan(p);
    fftwf_free(tmp);

    //Allocate an image
    ISMRMRD::Image<float> img_out(r_space.matrixSize.x, r_space.matrixSize.y, 1, 1);

    //f there is oversampling in the readout direction remove it
    //Take the magnitude
    size_t offset = ((e_space.matrixSize.x - r_space.matrixSize.x)>>1);
    for (unsigned int y = 0; y < r_space.matrixSize.y; y++) {
        for (unsigned int x = 0; x < r_space.matrixSize.x; x++) {
            img_out(x,y) = std::abs(buffer(x+offset, y));
        }
    }
    
    // The following are extra guidance we can put in the image header
    img_out.setImageType(ISMRMRD::ISMRMRD_IMTYPE_MAGNITUDE);
    img_out.setSlice(0);
    img_out.setFieldOfView(r_space.fieldOfView_mm.x, r_space.fieldOfView_mm.y, r_space.fieldOfView_mm.z);
    //And so on
    
    //Let's write the reconstructed image into the same data file
    d.appendImage("myimage", img_out);

    return 0;
}

Matlab Example Code and Datasets

The examples folder contains some matlab code to illustrate simple interaction with the ISMRMRD data format. The examples use test data sets, wich can be downloaded from the Github website. Go to the examples/data folder and type the following to download the data:

wget https://sourceforge.net/projects/ismrmrd/files/data/3D_partial_fourier.h5
wget https://sourceforge.net/projects/ismrmrd/files/data/simple_gre.h5
wget https://sourceforge.net/projects/ismrmrd/files/data/simple_spiral.h5

For instance, to reconstruct a 2D Cartesian acquisition (10 image repetitions), type (from the examples/matlab folder):

>> images = simple_cartesian_recon('../data/simple_gre.h5');
Reconstructing image 1....done
Reconstructing image 2....done
Reconstructing image 3....done
Reconstructing image 4....done
Reconstructing image 5....done
Reconstructing image 6....done
Reconstructing image 7....done
Reconstructing image 8....done
Reconstructing image 9....done
Reconstructing image 10....done
>>

You should see one of the reconstructed images display. An example is also given of a 3D acquisition with partial Fourier, phase and slice oversampling, etc. Reconstruct this dataset with:

>> images = simple_cartesian_recon('../data/3D_partial_fourier.h5');
Reconstructing image 1....done

The center slice of the volume should be displayed at the end of the reconstruction.

Finally, there is also a spiral dataset. This dataset illustrates how the flexible section of the <trajectoryDescription> can be used to add user defined parameters and an identifier to describe the trajectory. This dataset is also an example of storing the trajectory with the data for direct reconstruction. Reconstruct this dataset with:

>> images = simple_spiral_recon('../data/simple_spiral.h5');
Reconstructing image 1....done
Reconstructing image 2....done
Reconstructing image 3....done
Reconstructing image 4....done
Reconstructing image 5....done
Reconstructing image 6....done
Reconstructing image 7....done
Reconstructing image 8....done
Reconstructing image 9....done
Reconstructing image 10....done
>>

Appendix

XML Schema Definition

<?xml version="1.0" encoding="UTF-8" standalone="no"?>
<xs:schema xmlns="http://www.ismrm.org/ISMRMRD" xmlns:xs="http://www.w3.org/2001/XMLSchema" elementFormDefault="qualified" targetNamespace="http://www.ismrm.org/ISMRMRD" version="1">

    <xs:element name="ismrmrdHeader" type="ismrmrdHeader"/>

    <xs:complexType name="ismrmrdHeader">
      <xs:sequence>
        <xs:element maxOccurs="1" minOccurs="0" name="version" type="xs:long"/>
        <xs:element maxOccurs="1" minOccurs="0" name="subjectInformation" type="subjectInformationType"/>
        <xs:element maxOccurs="1" minOccurs="0" name="studyInformation" type="studyInformationType"/>
        <xs:element maxOccurs="1" minOccurs="0" name="measurementInformation" type="measurementInformationType"/>
        <xs:element maxOccurs="1" minOccurs="0" name="acquisitionSystemInformation" type="acquisitionSystemInformationType"/>
        <xs:element maxOccurs="1" minOccurs="1" name="experimentalConditions" type="experimentalConditionsType"/>
        <xs:element maxOccurs="unbounded" minOccurs="1" name="encoding" type="encoding"/>
        <xs:element maxOccurs="1" minOccurs="0" name="sequenceParameters" type="sequenceParametersType"/>
        <xs:element maxOccurs="1" minOccurs="0" name="userParameters" type="userParameters"/>
      </xs:sequence>
    </xs:complexType>

  <xs:complexType name="subjectInformationType">
    <xs:all>
      <xs:element minOccurs="0" name="patientName" type="xs:string"/>
      <xs:element minOccurs="0" name="patientWeight_kg" type="xs:float"/>
      <xs:element minOccurs="0" name="patientID" type="xs:string"/>
      <xs:element minOccurs="0" name="patientBirthdate" type="xs:date"/>
      <xs:element minOccurs="0" name="patientGender">
        <xs:simpleType>
          <xs:restriction base="xs:string">
            <xs:pattern value="[MFO]"/>
          </xs:restriction>
        </xs:simpleType>
      </xs:element>
    </xs:all>
  </xs:complexType>

  <xs:complexType name="studyInformationType">
    <xs:all>
      <xs:element minOccurs="0" maxOccurs="1" name="studyDate" type="xs:date"/>
      <xs:element minOccurs="0" maxOccurs="1" name="studyTime" type="xs:time"/>
      <xs:element minOccurs="0" maxOccurs="1" name="studyID" type="xs:string"/>
      <xs:element minOccurs="0" maxOccurs="1" name="accessionNumber" type="xs:long"/>
      <xs:element minOccurs="0" maxOccurs="1" name="referringPhysicianName" type="xs:string"/>
      <xs:element minOccurs="0" maxOccurs="1" name="studyDescription" type="xs:string"/>
      <xs:element minOccurs="0" maxOccurs="1" name="studyInstanceUID" type="xs:string"/>
    </xs:all>
  </xs:complexType>

  <xs:complexType name="measurementInformationType">
    <xs:sequence>
      <xs:element minOccurs="0" name="measurementID" type="xs:string"/>
      <xs:element minOccurs="0" name="seriesDate" type="xs:date"/>
      <xs:element minOccurs="0" name="seriesTime" type="xs:time"/>
      <xs:element minOccurs="1" name="patientPosition">
        <xs:simpleType>
          <xs:restriction base="xs:string">
            <xs:enumeration value="HFP"/>
            <xs:enumeration value="HFS"/>
            <xs:enumeration value="HFDR"/>
            <xs:enumeration value="HFDL"/>
            <xs:enumeration value="FFP"/>
            <xs:enumeration value="FFS"/>
            <xs:enumeration value="FFDR"/>
            <xs:enumeration value="FFDL"/>
          </xs:restriction>
        </xs:simpleType>
      </xs:element>
      <xs:element minOccurs="0" name="initialSeriesNumber" type="xs:long"/>
      <xs:element minOccurs="0" name="protocolName" type="xs:string"/>
      <xs:element minOccurs="0" name="seriesDescription" type="xs:string"/>
      <xs:element maxOccurs="unbounded" minOccurs="0" name="measurementDependency" type="measurementDependencyType"/>
      <xs:element minOccurs="0" name="seriesInstanceUIDRoot" type="xs:string"/>
      <xs:element minOccurs="0" name="frameOfReferenceUID" type="xs:string"/>
      <xs:element minOccurs="0" name="referencedImageSequence" type="referencedImageSequence"/>
    </xs:sequence>
  </xs:complexType>

  <xs:complexType name="measurementDependencyType">
    <xs:sequence>
        <xs:element maxOccurs="1" minOccurs="1" name="dependencyType" type="xs:string"/>
        <xs:element maxOccurs="1" minOccurs="1" name="measurementID" type="xs:string"/>
    </xs:sequence>
  </xs:complexType>

  <xs:complexType name="acquisitionSystemInformationType">
    <xs:all>
      <xs:element minOccurs="0" name="systemVendor" type="xs:string"/>
      <xs:element minOccurs="0" name="systemModel" type="xs:string"/>
      <xs:element minOccurs="0" name="systemFieldStrength_T" type="xs:float"/>
      <xs:element minOccurs="0" name="relativeReceiverNoiseBandwidth" type="xs:float"/>
      <xs:element minOccurs="0" name="receiverChannels" type="xs:unsignedShort"/>
      <xs:element minOccurs="0" name="institutionName" type="xs:string"/>
      <xs:element minOccurs="0" name="stationName" type="xs:string"/>
    </xs:all>
  </xs:complexType>

  <xs:complexType name="experimentalConditionsType">
    <xs:all>
      <xs:element name="H1resonanceFrequency_Hz" type="xs:long"/>
    </xs:all>
  </xs:complexType>

  <xs:complexType name="encoding">
    <xs:all>
      <xs:element maxOccurs="1" minOccurs="1" name="encodedSpace" type="encodingSpaceType"/>
      <xs:element maxOccurs="1" minOccurs="1" name="reconSpace" type="encodingSpaceType"/>
      <xs:element maxOccurs="1" minOccurs="1" name="encodingLimits" type="encodingLimitsType"/>
      <xs:element maxOccurs="1" minOccurs="1" name="trajectory" type="trajectoryType"/>
      <xs:element maxOccurs="1" minOccurs="0" name="trajectoryDescription" type="trajectoryDescriptionType"/>
      <xs:element maxOccurs="1" minOccurs="0" name="parallelImaging" type="parallelImagingType"/>
      <xs:element maxOccurs="1" minOccurs="0" name="echoTrainLength" type="xs:long"/>
    </xs:all>
  </xs:complexType>

  <xs:complexType name="encodingSpaceType">
    <xs:all>
      <xs:element maxOccurs="1" minOccurs="1" name="matrixSize" type="matrixSize"/>
      <xs:element maxOccurs="1" minOccurs="1" name="fieldOfView_mm" type="fieldOfView_mm"/>
    </xs:all>
  </xs:complexType>

  <xs:complexType name="matrixSize">
    <xs:sequence>
      <xs:element default="1" maxOccurs="1" minOccurs="1" name="x" type="xs:unsignedShort"/>
      <xs:element default="1" maxOccurs="1" minOccurs="1" name="y" type="xs:unsignedShort"/>
      <xs:element default="1" maxOccurs="1" minOccurs="1" name="z" type="xs:unsignedShort"/>
    </xs:sequence>
  </xs:complexType>

  <xs:complexType name="fieldOfView_mm">
    <xs:sequence>
      <xs:element maxOccurs="1" minOccurs="1" name="x" type="xs:float"/>
      <xs:element maxOccurs="1" minOccurs="1" name="y" type="xs:float"/>
      <xs:element maxOccurs="1" minOccurs="1" name="z" type="xs:float"/>
    </xs:sequence>
  </xs:complexType>

  <xs:complexType name="limitType">
    <xs:all>
      <xs:element default="0" name="minimum" type="xs:unsignedShort"/>
      <xs:element default="0" name="maximum" type="xs:unsignedShort"/>
      <xs:element default="0" name="center" type="xs:unsignedShort"/>
    </xs:all>
  </xs:complexType>

  <xs:complexType name="encodingLimitsType">
    <xs:all>
      <xs:element maxOccurs="1" minOccurs="0" name="kspace_encoding_step_0" type="limitType"/>
      <xs:element maxOccurs="1" minOccurs="0" name="kspace_encoding_step_1" type="limitType"/>
      <xs:element maxOccurs="1" minOccurs="0" name="kspace_encoding_step_2" type="limitType"/>
      <xs:element maxOccurs="1" minOccurs="0" name="average" type="limitType"/>
      <xs:element maxOccurs="1" minOccurs="0" name="slice" type="limitType"/>
      <xs:element maxOccurs="1" minOccurs="0" name="contrast" type="limitType"/>
      <xs:element maxOccurs="1" minOccurs="0" name="phase" type="limitType"/>
      <xs:element maxOccurs="1" minOccurs="0" name="repetition" type="limitType"/>
      <xs:element maxOccurs="1" minOccurs="0" name="set" type="limitType"/>
      <xs:element maxOccurs="1" minOccurs="0" name="segment" type="limitType"/>
    </xs:all>
  </xs:complexType>

  <xs:simpleType name="trajectoryType">
    <xs:restriction base="xs:string">
      <xs:enumeration value="cartesian"/>
      <xs:enumeration value="epi"/>
      <xs:enumeration value="radial"/>
      <xs:enumeration value="goldenangle"/>
      <xs:enumeration value="spiral"/>
      <xs:enumeration value="other"/>
    </xs:restriction>
  </xs:simpleType>

  <xs:complexType name="trajectoryDescriptionType">
    <xs:sequence>
      <xs:element maxOccurs="1" minOccurs="1" name="identifier" type="xs:string"/>
      <xs:element maxOccurs="unbounded" minOccurs="0" name="userParameterLong" type="userParameterLongType"/>
      <xs:element maxOccurs="unbounded" minOccurs="0" name="userParameterDouble" type="userParameterDoubleType"/>
      <xs:element maxOccurs="1" minOccurs="0" name="comment" type="xs:string"/>
    </xs:sequence>
  </xs:complexType>

  <xs:complexType name="sequenceParametersType">
    <xs:sequence>
      <xs:element minOccurs="1" maxOccurs="unbounded" type="xs:float" name="TR"/>
      <xs:element minOccurs="1" maxOccurs="unbounded" type="xs:float" name="TE"/>
      <xs:element minOccurs="0" maxOccurs="unbounded" type="xs:float" name="TI"/>
      <xs:element minOccurs="0" maxOccurs="unbounded" type="xs:float" name="flipAngle_deg"/>

    </xs:sequence>
  </xs:complexType>

  <xs:complexType name="userParameterLongType">
    <xs:all>
      <xs:element name="name" type="xs:string"/>
      <xs:element name="value" type="xs:long"/>
    </xs:all>
  </xs:complexType>

  <xs:complexType name="userParameterDoubleType">
    <xs:all>
      <xs:element name="name" type="xs:string"/>
      <xs:element name="value" type="xs:double"/>
    </xs:all>
  </xs:complexType>

  <xs:complexType name="userParameterStringType">
    <xs:all>
      <xs:element name="name" type="xs:string"/>
      <xs:element name="value" type="xs:string"/>
    </xs:all>
  </xs:complexType>

  <xs:complexType name="userParameterBase64Type">
    <xs:all>
      <xs:element name="name" type="xs:string"/>
      <xs:element name="value" type="xs:base64Binary"/>
    </xs:all>
  </xs:complexType>

  <xs:complexType name="referencedImageSequence">
    <xs:sequence>
      <xs:element minOccurs="0" maxOccurs="unbounded" name="referencedSOPInstanceUID" type="xs:string"/>
    </xs:sequence>
  </xs:complexType>
      
  <xs:complexType name="userParameters">
    <xs:sequence>
      <xs:element maxOccurs="unbounded" minOccurs="0" name="userParameterLong" type="userParameterLongType"/>
      <xs:element maxOccurs="unbounded" minOccurs="0" name="userParameterDouble" type="userParameterDoubleType"/>
      <xs:element maxOccurs="unbounded" minOccurs="0" name="userParameterString" type="userParameterStringType"/>
      <xs:element maxOccurs="unbounded" minOccurs="0" name="userParameterBase64" type="userParameterBase64Type"/>
    </xs:sequence>
  </xs:complexType>

  <xs:complexType name="accelerationFactorType">
    <xs:all>
      <xs:element name="kspace_encoding_step_1" type="xs:unsignedShort"/>
      <xs:element name="kspace_encoding_step_2" type="xs:unsignedShort"/>
    </xs:all>
  </xs:complexType>

  <xs:simpleType name="calibrationModeType">
    <xs:restriction base="xs:string">
      <xs:enumeration value="embedded"/>
      <xs:enumeration value="interleaved"/>
      <xs:enumeration value="separate"/>
      <xs:enumeration value="external"/>
      <xs:enumeration value="other"/>
    </xs:restriction>
  </xs:simpleType>

  <xs:simpleType name="interleavingDimensionType">
    <xs:restriction base="xs:string">
      <xs:enumeration value="phase"/>
      <xs:enumeration value="repetition"/>
      <xs:enumeration value="contrast"/>
      <xs:enumeration value="average"/>
      <xs:enumeration value="other"/>
    </xs:restriction>
  </xs:simpleType>

  <xs:complexType name="parallelImagingType">
  	<xs:sequence>
  	 <xs:element type="accelerationFactorType" name="accelerationFactor"/>
  	 <xs:element maxOccurs="1" minOccurs="0" type="calibrationModeType" name="calibrationMode"/>
  	 <xs:element maxOccurs="1" minOccurs="0" type="interleavingDimensionType" name="interleavingDimension"/>
  	</xs:sequence>
  </xs:complexType>
</xs:schema>