RNifti
Fast R and C++ Access to NIfTI Images
NiftiImage.h
1 #ifndef _NIFTI_IMAGE_H_
2 #define _NIFTI_IMAGE_H_
3 
4 
5 #ifdef USING_R
6 
7 #include <Rcpp.h>
8 
9 #else
10 
11 #define R_NegInf -INFINITY
12 
13 #include <stdint.h>
14 #include <cstddef>
15 #include <cmath>
16 #include <string>
17 #include <sstream>
18 #include <vector>
19 #include <stdexcept>
20 #include <algorithm>
21 #include <map>
22 #include <locale>
23 #include <limits>
24 
25 #endif
26 
27 
28 #include "niftilib/nifti1_io.h"
29 
38 namespace RNifti {
39 
40 namespace internal {
41 
42 struct vec3
43 {
44  float v[3];
45 
46  vec3 operator-() const
47  {
48  vec3 r;
49  r.v[0] = -v[0];
50  r.v[1] = -v[1];
51  r.v[2] = -v[2];
52  return r;
53  }
54 };
55 
56 inline mat33 topLeftCorner (const mat44 &matrix)
57 {
58  mat33 newMatrix;
59  for (int i=0; i<3; i++)
60  {
61  for (int j=0; j<3; j++)
62  newMatrix.m[i][j] = matrix.m[i][j];
63  }
64  return newMatrix;
65 }
66 
67 inline vec3 matrixVectorProduct (const mat33 &matrix, const vec3 &vector)
68 {
69  vec3 newVector;
70  for (int i=0; i<3; i++)
71  {
72  newVector.v[i] = 0.0;
73  for (int j=0; j<3; j++)
74  newVector.v[i] += matrix.m[i][j] * vector.v[j];
75  }
76  return newVector;
77 }
78 
79 } // internal namespace
80 
86 {
87 public:
92  struct Block
93  {
94  const NiftiImage &image;
95  const int dimension;
96  const int index;
105  Block (const NiftiImage &image, const int dimension, const int index)
107  {
108  if (dimension != image->ndim)
109  throw std::runtime_error("Blocks must be along the last dimension in the image");
110  }
111 
120  Block & operator= (const NiftiImage &source)
121  {
122  if (source->datatype != image->datatype)
123  throw std::runtime_error("New data does not have the same datatype as the target block");
124  if (source->scl_slope != image->scl_slope || source->scl_inter != image->scl_inter)
125  throw std::runtime_error("New data does not have the same scale parameters as the target block");
126 
127  size_t blockSize = 1;
128  for (int i=1; i<dimension; i++)
129  blockSize *= image->dim[i];
130 
131  if (blockSize != source->nvox)
132  throw std::runtime_error("New data does not have the same size as the target block");
133 
134  blockSize *= image->nbyper;
135  memcpy(static_cast<char*>(image->data) + blockSize*index, source->data, blockSize);
136  return *this;
137  }
138 
146  template <typename TargetType>
147  std::vector<TargetType> getData (const bool useSlope = true) const;
148  };
149 
150 #ifdef USING_R
151 
157  static short sexpTypeToNiftiType (const int sexpType)
158  {
159  if (sexpType == INTSXP || sexpType == LGLSXP)
160  return DT_INT32;
161  else if (sexpType == REALSXP)
162  return DT_FLOAT64;
163  else
164  throw std::runtime_error("Array elements must be numeric");
165  }
166 #endif
167 
173  static mat33 xformToRotation (const mat44 matrix)
174  {
175  float qb, qc, qd, qfac;
176  nifti_mat44_to_quatern(matrix, &qb, &qc, &qd, NULL, NULL, NULL, NULL, NULL, NULL, &qfac);
177  mat44 rotationMatrix = nifti_quatern_to_mat44(qb, qc, qd, 0.0, 0.0, 0.0, 1.0, 1.0, 1.0, qfac);
178  return internal::topLeftCorner(rotationMatrix);
179  }
180 
186  static std::string xformToString (const mat44 matrix)
187  {
188  int icode, jcode, kcode;
189  nifti_mat44_to_orientation(matrix, &icode, &jcode, &kcode);
190 
191  int codes[3] = { icode, jcode, kcode };
192  std::string result("---");
193  for (int i=0; i<3; i++)
194  {
195  switch (codes[i])
196  {
197  case NIFTI_L2R: result[i] = 'R'; break;
198  case NIFTI_R2L: result[i] = 'L'; break;
199  case NIFTI_P2A: result[i] = 'A'; break;
200  case NIFTI_A2P: result[i] = 'P'; break;
201  case NIFTI_I2S: result[i] = 'S'; break;
202  case NIFTI_S2I: result[i] = 'I'; break;
203  }
204  }
205  return result;
206  }
207 
214  static int fileVersion (const std::string &path)
215  {
216  nifti_1_header *header = nifti_read_header(path.c_str(), NULL, false);
217  if (header == NULL)
218  return -1;
219  else
220  {
221  int version = NIFTI_VERSION(*header);
222  if (version == 0)
223  {
224  // NIfTI-2 has a 540-byte header - check for this or its byte-swapped equivalent
225  if (header->sizeof_hdr == 540 || header->sizeof_hdr == 469893120)
226  {
227  // The magic number has moved in NIfTI-2, so find it by byte offset
228  const char *magic = (char *) header + 4;
229  if (strncmp(magic,"ni2",3) == 0 || strncmp(magic,"n+2",3) == 0)
230  version = 2;
231  }
232  else if (!nifti_hdr_looks_good(header))
233  {
234  // Not plausible as ANALYZE, so return -1
235  version = -1;
236  }
237  }
238  return version;
239  }
240  }
241 
242 
243 protected:
244  nifti_image *image;
245  bool persistent;
251  void copy (const nifti_image *source);
252 
257  void copy (const NiftiImage &source);
258 
263  void copy (const Block &source);
264 
265 
266 #ifdef USING_R
267 
273  void initFromNiftiS4 (const Rcpp::RObject &object, const bool copyData = true);
274 
280  void initFromMriImage (const Rcpp::RObject &object, const bool copyData = true);
281 
286  void initFromList (const Rcpp::RObject &object);
287 
293  void initFromArray (const Rcpp::RObject &object, const bool copyData = true);
294 
295 #endif
296 
301  void updatePixdim (const std::vector<float> &pixdim);
302 
307  void setPixunits (const std::vector<std::string> &pixunits);
308 
309 public:
314  : image(NULL), persistent(false) {}
315 
320  NiftiImage (const NiftiImage &source)
321  : image(NULL), persistent(false)
322  {
323  this->copy(source);
324 #ifndef NDEBUG
325  Rprintf("Creating NiftiImage with pointer %p (from NiftiImage)\n", this->image);
326 #endif
327  }
328 
333  NiftiImage (const Block &source)
334  : image(NULL), persistent(false)
335  {
336  this->copy(source);
337 #ifndef NDEBUG
338  Rprintf("Creating NiftiImage with pointer %p (from Block)\n", this->image);
339 #endif
340  }
341 
348  NiftiImage (nifti_image * const image, const bool copy = false)
349  : image(NULL), persistent(false)
350  {
351  if (copy)
352  this->copy(image);
353  else
354  this->image = image;
355 #ifndef NDEBUG
356  Rprintf("Creating NiftiImage with pointer %p (from pointer)\n", this->image);
357 #endif
358  }
359 
366  NiftiImage (const std::string &path, const bool readData = true)
367  : persistent(false)
368  {
369  this->image = nifti_image_read(path.c_str(), readData);
370  if (this->image == NULL)
371  throw std::runtime_error("Failed to read image from path " + path);
372 #ifndef NDEBUG
373  Rprintf("Creating NiftiImage with pointer %p (from string)\n", this->image);
374 #endif
375  }
376 
383  NiftiImage (const std::string &path, const std::vector<int> &volumes);
384 
385 #ifdef USING_R
386 
391  NiftiImage (const SEXP object, const bool readData = true);
392 #endif
393 
397  virtual ~NiftiImage ()
398  {
399  if (!persistent)
400  {
401 #ifndef NDEBUG
402  Rprintf("Freeing NiftiImage with pointer %p\n", this->image);
403 #endif
404  nifti_image_free(image);
405  }
406  }
407 
411  operator const nifti_image* () const { return image; }
412 
416  operator nifti_image* () { return image; }
417 
421  const nifti_image * operator-> () const { return image; }
422 
426  nifti_image * operator-> () { return image; }
427 
433  {
434  copy(source);
435 #ifndef NDEBUG
436  Rprintf("Creating NiftiImage with pointer %p (from NiftiImage)\n", this->image);
437 #endif
438  return *this;
439  }
440 
446  NiftiImage & operator= (const Block &source)
447  {
448  copy(source);
449 #ifndef NDEBUG
450  Rprintf("Creating NiftiImage with pointer %p (from Block)\n", this->image);
451 #endif
452  return *this;
453  }
454 
460  {
461  this->persistent = persistent;
462 #ifndef NDEBUG
463  if (persistent)
464  Rprintf("Setting NiftiImage with pointer %p to be persistent\n", this->image);
465 #endif
466  return *this;
467  }
468 
472  bool isNull () const { return (image == NULL); }
473 
477  bool isPersistent () const { return persistent; }
478 
482  bool isDataScaled () const { return (image != NULL && image->scl_slope != 0.0 && (image->scl_slope != 1.0 || image->scl_inter != 0.0)); }
483 
487  int nDims () const
488  {
489  if (image == NULL)
490  return 0;
491  else
492  return image->ndim;
493  }
494 
499  std::vector<int> dim () const
500  {
501  if (image == NULL)
502  return std::vector<int>();
503  else
504  return std::vector<int>(image->dim+1, image->dim+image->ndim+1);
505  }
506 
511  std::vector<float> pixdim () const
512  {
513  if (image == NULL)
514  return std::vector<float>();
515  else
516  return std::vector<float>(image->pixdim+1, image->pixdim+image->ndim+1);
517  }
518 
526  {
527  int ndim = image->ndim;
528  while (image->dim[ndim] < 2)
529  ndim--;
530  image->dim[0] = image->ndim = ndim;
531 
532  return *this;
533  }
534 
542  template <typename TargetType>
543  std::vector<TargetType> getData (const bool useSlope = true) const;
544 
551  NiftiImage & changeDatatype (const short datatype, const bool useSlope = false);
552 
559  NiftiImage & changeDatatype (const std::string &datatype, const bool useSlope = false);
560 
568  template <typename SourceType>
569  NiftiImage & replaceData (const std::vector<SourceType> &data, const short datatype = DT_NONE);
570 
575  {
576  nifti_image_unload(image);
577  return *this;
578  }
579 
585  NiftiImage & rescale (const std::vector<float> &scales);
586 
594  NiftiImage & reorient (const int i, const int j, const int k);
595 
603  NiftiImage & reorient (const std::string &orientation);
604 
605 #ifdef USING_R
606 
610  NiftiImage & update (const Rcpp::RObject &object);
611 #endif
612 
618  mat44 xform (const bool preferQuaternion = true) const;
619 
623  int nBlocks () const
624  {
625  if (image == NULL)
626  return 0;
627  else
628  return image->dim[image->ndim];
629  }
630 
638  const Block block (const int i) const { return Block(*this, nDims(), i); }
639 
647  Block block (const int i) { return Block(*this, nDims(), i); }
648 
654  const Block slice (const int i) const { return Block(*this, 3, i); }
655 
661  Block slice (const int i) { return Block(*this, 3, i); }
662 
668  const Block volume (const int i) const { return Block(*this, 4, i); }
669 
675  Block volume (const int i) { return Block(*this, 4, i); }
676 
682  void toFile (const std::string fileName, const short datatype = DT_NONE) const;
683 
689  void toFile (const std::string fileName, const std::string &datatype) const;
690 
691 #ifdef USING_R
692 
697  Rcpp::RObject toArray () const;
698 
704  Rcpp::RObject toPointer (const std::string label) const;
705 
712  Rcpp::RObject toArrayOrPointer (const bool internal, const std::string label) const;
713 
718  Rcpp::RObject headerToList () const;
719 
720 #endif
721 
722 };
723 
724 // Include helper functions
725 #include "lib/NiftiImage_internal.h"
726 
727 inline void NiftiImage::copy (const nifti_image *source)
728 {
729  if (image != NULL && !persistent)
730  nifti_image_free(image);
731 
732  if (source == NULL)
733  image = NULL;
734  else
735  {
736  image = nifti_copy_nim_info(source);
737  if (source->data != NULL)
738  {
739  size_t dataSize = nifti_get_volsize(source);
740  image->data = calloc(1, dataSize);
741  memcpy(image->data, source->data, dataSize);
742  }
743  }
744 
745  persistent = false;
746 }
747 
748 inline void NiftiImage::copy (const NiftiImage &source)
749 {
750  const nifti_image *sourceStruct = source;
751  copy(sourceStruct);
752 }
753 
754 inline void NiftiImage::copy (const Block &source)
755 {
756  if (image != NULL && !persistent)
757  nifti_image_free(image);
758 
759  const nifti_image *sourceStruct = source.image;
760  if (sourceStruct == NULL)
761  image = NULL;
762  else
763  {
764  image = nifti_copy_nim_info(sourceStruct);
765  image->dim[0] = source.image->dim[0] - 1;
766  image->dim[source.dimension] = 1;
767  image->pixdim[source.dimension] = 1.0;
768  nifti_update_dims_from_array(image);
769 
770  if (sourceStruct->data != NULL)
771  {
772  size_t blockSize = nifti_get_volsize(image);
773  image->data = calloc(1, blockSize);
774  memcpy(image->data, static_cast<char*>(source.image->data) + blockSize*source.index, blockSize);
775  }
776  }
777 
778  persistent = false;
779 }
780 
781 #ifdef USING_R
782 
783 // Convert an S4 "nifti" object, as defined in the oro.nifti package, to a "nifti_image" struct
784 inline void NiftiImage::initFromNiftiS4 (const Rcpp::RObject &object, const bool copyData)
785 {
786  nifti_1_header header;
787  header.sizeof_hdr = 348;
788 
789  const std::vector<short> dims = object.slot("dim_");
790  for (int i=0; i<8; i++)
791  header.dim[i] = dims[i];
792 
793  header.intent_p1 = object.slot("intent_p1");
794  header.intent_p2 = object.slot("intent_p2");
795  header.intent_p3 = object.slot("intent_p3");
796  header.intent_code = object.slot("intent_code");
797 
798  header.datatype = object.slot("datatype");
799  header.bitpix = object.slot("bitpix");
800 
801  header.slice_start = object.slot("slice_start");
802  header.slice_end = object.slot("slice_end");
803  header.slice_code = Rcpp::as<int>(object.slot("slice_code"));
804  header.slice_duration = object.slot("slice_duration");
805 
806  const std::vector<float> pixdims = object.slot("pixdim");
807  for (int i=0; i<8; i++)
808  header.pixdim[i] = pixdims[i];
809  header.xyzt_units = Rcpp::as<int>(object.slot("xyzt_units"));
810 
811  header.vox_offset = object.slot("vox_offset");
812 
813  header.scl_slope = object.slot("scl_slope");
814  header.scl_inter = object.slot("scl_inter");
815  header.toffset = object.slot("toffset");
816 
817  header.cal_max = object.slot("cal_max");
818  header.cal_min = object.slot("cal_min");
819  header.glmax = header.glmin = 0;
820 
821  strncpy(header.descrip, Rcpp::as<std::string>(object.slot("descrip")).c_str(), 79);
822  header.descrip[79] = '\0';
823  strncpy(header.aux_file, Rcpp::as<std::string>(object.slot("aux_file")).c_str(), 23);
824  header.aux_file[23] = '\0';
825  strncpy(header.intent_name, Rcpp::as<std::string>(object.slot("intent_name")).c_str(), 15);
826  header.intent_name[15] = '\0';
827  strncpy(header.magic, Rcpp::as<std::string>(object.slot("magic")).c_str(), 3);
828  header.magic[3] = '\0';
829 
830  header.qform_code = object.slot("qform_code");
831  header.sform_code = object.slot("sform_code");
832 
833  header.quatern_b = object.slot("quatern_b");
834  header.quatern_c = object.slot("quatern_c");
835  header.quatern_d = object.slot("quatern_d");
836  header.qoffset_x = object.slot("qoffset_x");
837  header.qoffset_y = object.slot("qoffset_y");
838  header.qoffset_z = object.slot("qoffset_z");
839 
840  const std::vector<float> srow_x = object.slot("srow_x");
841  const std::vector<float> srow_y = object.slot("srow_y");
842  const std::vector<float> srow_z = object.slot("srow_z");
843  for (int i=0; i<4; i++)
844  {
845  header.srow_x[i] = srow_x[i];
846  header.srow_y[i] = srow_y[i];
847  header.srow_z[i] = srow_z[i];
848  }
849 
850  if (header.datatype == DT_UINT8 || header.datatype == DT_INT16 || header.datatype == DT_INT32 || header.datatype == DT_INT8 || header.datatype == DT_UINT16 || header.datatype == DT_UINT32)
851  header.datatype = DT_INT32;
852  else if (header.datatype == DT_FLOAT32 || header.datatype == DT_FLOAT64)
853  header.datatype = DT_FLOAT64; // This assumes that sizeof(double) == 8
854  else
855  throw std::runtime_error("Data type is not supported");
856 
857  this->image = nifti_convert_nhdr2nim(header, NULL);
858 
859  const SEXP data = PROTECT(object.slot(".Data"));
860  if (!copyData || Rf_length(data) <= 1)
861  this->image->data = NULL;
862  else
863  {
864  const size_t dataSize = nifti_get_volsize(this->image);
865  this->image->data = calloc(1, dataSize);
866  if (header.datatype == DT_INT32)
867  {
868  Rcpp::IntegerVector intData(data);
869  std::copy(intData.begin(), intData.end(), static_cast<int32_t*>(this->image->data));
870  }
871  else
872  {
873  Rcpp::DoubleVector doubleData(data);
874  std::copy(doubleData.begin(), doubleData.end(), static_cast<double*>(this->image->data));
875  }
876  }
877  UNPROTECT(1);
878 }
879 
880 inline void NiftiImage::initFromMriImage (const Rcpp::RObject &object, const bool copyData)
881 {
882  Rcpp::Reference mriImage(object);
883  Rcpp::Function getXform = mriImage.field("getXform");
884  Rcpp::NumericMatrix xform = getXform();
885 
886  this->image = NULL;
887 
888  if (Rf_length(mriImage.field("tags")) > 0)
889  initFromList(mriImage.field("tags"));
890 
891  Rcpp::RObject data = mriImage.field("data");
892  if (data.inherits("SparseArray"))
893  {
894  Rcpp::Language call("as.array", data);
895  data = call.eval();
896  }
897 
898  const short datatype = (Rf_isNull(data) ? DT_INT32 : sexpTypeToNiftiType(data.sexp_type()));
899 
900  int dims[8] = { 0, 0, 0, 0, 0, 0, 0, 0 };
901  const std::vector<int> dimVector = mriImage.field("imageDims");
902  const int nDims = std::min(7, int(dimVector.size()));
903  dims[0] = nDims;
904  size_t nVoxels = 1;
905  for (int i=0; i<nDims; i++)
906  {
907  dims[i+1] = dimVector[i];
908  nVoxels *= dimVector[i];
909  }
910 
911  if (this->image == NULL)
912  this->image = nifti_make_new_nim(dims, datatype, FALSE);
913  else
914  {
915  std::copy(dims, dims+8, this->image->dim);
916  this->image->datatype = datatype;
917  nifti_datatype_sizes(image->datatype, &image->nbyper, NULL);
918  }
919 
920  if (copyData && !Rf_isNull(data))
921  {
922  // NB: nifti_get_volsize() will not be right here if there were tags
923  const size_t dataSize = nVoxels * image->nbyper;
924  this->image->data = calloc(1, dataSize);
925  if (datatype == DT_INT32)
926  memcpy(this->image->data, INTEGER(data), dataSize);
927  else
928  memcpy(this->image->data, REAL(data), dataSize);
929  }
930  else
931  this->image->data = NULL;
932 
933  const std::vector<float> pixdimVector = mriImage.field("voxelDims");
934  const int pixdimLength = pixdimVector.size();
935  for (int i=0; i<std::min(pixdimLength,nDims); i++)
936  this->image->pixdim[i+1] = std::abs(pixdimVector[i]);
937 
938  const std::vector<std::string> pixunitsVector = mriImage.field("voxelDimUnits");
939  setPixunits(pixunitsVector);
940 
941  if (xform.rows() != 4 || xform.cols() != 4)
942  this->image->qform_code = this->image->sform_code = 0;
943  else
944  {
945  mat44 matrix;
946  for (int i=0; i<4; i++)
947  {
948  for (int j=0; j<4; j++)
949  matrix.m[i][j] = static_cast<float>(xform(i,j));
950  }
951 
952  this->image->qto_xyz = matrix;
953  this->image->qto_ijk = nifti_mat44_inverse(image->qto_xyz);
954  nifti_mat44_to_quatern(image->qto_xyz, &image->quatern_b, &image->quatern_c, &image->quatern_d, &image->qoffset_x, &image->qoffset_y, &image->qoffset_z, NULL, NULL, NULL, &image->qfac);
955 
956  this->image->sto_xyz = matrix;
957  this->image->sto_ijk = nifti_mat44_inverse(image->sto_xyz);
958 
959  this->image->qform_code = this->image->sform_code = 2;
960  }
961 }
962 
963 inline void NiftiImage::initFromList (const Rcpp::RObject &object)
964 {
965  Rcpp::List list(object);
966  nifti_1_header *header = nifti_make_new_header(NULL, DT_FLOAT64);
967 
968  internal::updateHeader(header, list);
969 
970  this->image = nifti_convert_nhdr2nim(*header, NULL);
971  this->image->data = NULL;
972  free(header);
973 }
974 
975 inline void NiftiImage::initFromArray (const Rcpp::RObject &object, const bool copyData)
976 {
977  int dims[8] = { 0, 0, 0, 0, 0, 0, 0, 0 };
978  const std::vector<int> dimVector = object.attr("dim");
979 
980  const int nDims = std::min(7, int(dimVector.size()));
981  dims[0] = nDims;
982  for (int i=0; i<nDims; i++)
983  dims[i+1] = dimVector[i];
984 
985  const short datatype = sexpTypeToNiftiType(object.sexp_type());
986  this->image = nifti_make_new_nim(dims, datatype, int(copyData));
987 
988  if (copyData)
989  {
990  const size_t dataSize = nifti_get_volsize(image);
991  if (datatype == DT_INT32)
992  memcpy(this->image->data, INTEGER(object), dataSize);
993  else
994  memcpy(this->image->data, REAL(object), dataSize);
995  }
996  else
997  this->image->data = NULL;
998 
999  if (object.hasAttribute("pixdim"))
1000  {
1001  const std::vector<float> pixdimVector = object.attr("pixdim");
1002  const int pixdimLength = pixdimVector.size();
1003  for (int i=0; i<std::min(pixdimLength,nDims); i++)
1004  this->image->pixdim[i+1] = pixdimVector[i];
1005  }
1006 
1007  if (object.hasAttribute("pixunits"))
1008  {
1009  const std::vector<std::string> pixunitsVector = object.attr("pixunits");
1010  setPixunits(pixunitsVector);
1011  }
1012 }
1013 
1014 inline NiftiImage::NiftiImage (const SEXP object, const bool readData)
1015  : persistent(false)
1016 {
1017  Rcpp::RObject imageObject(object);
1018  bool resolved = false;
1019 
1020  if (imageObject.hasAttribute(".nifti_image_ptr"))
1021  {
1022  Rcpp::XPtr<NiftiImage> imagePtr(SEXP(imageObject.attr(".nifti_image_ptr")));
1023  NiftiImage *ptr = imagePtr;
1024  if (ptr != NULL)
1025  {
1026  this->image = ptr->image;
1027  this->persistent = true;
1028  resolved = true;
1029 
1030  if (imageObject.hasAttribute("dim"))
1031  update(imageObject);
1032  }
1033  else if (Rf_isString(object))
1034  throw std::runtime_error("Internal image is not valid");
1035  else
1036  Rf_warning("Ignoring invalid internal pointer");
1037  }
1038 
1039  if (!resolved)
1040  {
1041  if (Rf_isNull(object))
1042  this->image = NULL;
1043  else if (Rf_isString(object))
1044  {
1045  const std::string path = Rcpp::as<std::string>(object);
1046  this->image = nifti_image_read(path.c_str(), readData);
1047  if (this->image == NULL)
1048  throw std::runtime_error("Failed to read image from path " + path);
1049  }
1050  else if (imageObject.inherits("nifti"))
1051  initFromNiftiS4(imageObject, readData);
1052  else if (imageObject.inherits("anlz"))
1053  throw std::runtime_error("Cannot currently convert objects of class \"anlz\"");
1054  else if (imageObject.inherits("MriImage"))
1055  initFromMriImage(imageObject, readData);
1056  else if (Rf_isVectorList(object))
1057  initFromList(imageObject);
1058  else if (imageObject.hasAttribute("dim"))
1059  initFromArray(imageObject, readData);
1060  else
1061  throw std::runtime_error("Cannot convert object of class \"" + Rcpp::as<std::string>(imageObject.attr("class")) + "\" to a nifti_image");
1062  }
1063 
1064  if (this->image != NULL)
1065  nifti_update_dims_from_array(this->image);
1066 
1067 #ifndef NDEBUG
1068  Rprintf("Creating NiftiImage with pointer %p (from SEXP)\n", this->image);
1069 #endif
1070 }
1071 
1072 #endif // USING_R
1073 
1074 inline NiftiImage::NiftiImage (const std::string &path, const std::vector<int> &volumes)
1075  : persistent(false)
1076 {
1077  if (volumes.empty())
1078  throw std::runtime_error("The vector of volumes is empty");
1079 
1080  nifti_brick_list brickList;
1081  this->image = nifti_image_read_bricks(path.c_str(), volumes.size(), &volumes[0], &brickList);
1082  if (this->image == NULL)
1083  throw std::runtime_error("Failed to read image from path " + path);
1084 
1085  size_t brickSize = image->nbyper * image->nx * image->ny * image->nz;
1086  image->data = calloc(1, nifti_get_volsize(image));
1087  for (int i=0; i<brickList.nbricks; i++)
1088  memcpy((char *) image->data + i * brickSize, brickList.bricks[i], brickSize);
1089  nifti_free_NBL(&brickList);
1090 
1091 #ifndef NDEBUG
1092  Rprintf("Creating NiftiImage with pointer %p (from string and volume vector)\n", this->image);
1093 #endif
1094 }
1095 
1096 inline void NiftiImage::updatePixdim (const std::vector<float> &pixdim)
1097 {
1098  const int nDims = image->dim[0];
1099  const std::vector<float> origPixdim(image->pixdim+1, image->pixdim+4);
1100 
1101  for (int i=1; i<8; i++)
1102  image->pixdim[i] = 0.0;
1103 
1104  const int pixdimLength = pixdim.size();
1105  for (int i=0; i<std::min(pixdimLength,nDims); i++)
1106  image->pixdim[i+1] = pixdim[i];
1107 
1108  if (!std::equal(origPixdim.begin(), origPixdim.begin() + std::min(3,nDims), pixdim.begin()))
1109  {
1110  mat33 scaleMatrix;
1111  for (int i=0; i<3; i++)
1112  {
1113  for (int j=0; j<3; j++)
1114  {
1115  if (i != j)
1116  scaleMatrix.m[i][j] = 0.0;
1117  else if (i >= nDims)
1118  scaleMatrix.m[i][j] = 1.0;
1119  else
1120  scaleMatrix.m[i][j] = pixdim[i] / origPixdim[i];
1121  }
1122  }
1123 
1124  if (image->qform_code > 0)
1125  {
1126  mat33 prod = nifti_mat33_mul(scaleMatrix, internal::topLeftCorner(image->qto_xyz));
1127  for (int i=0; i<3; i++)
1128  {
1129  for (int j=0; j<3; j++)
1130  image->qto_xyz.m[i][j] = prod.m[i][j];
1131  }
1132  image->qto_ijk = nifti_mat44_inverse(image->qto_xyz);
1133  nifti_mat44_to_quatern(image->qto_xyz, &image->quatern_b, &image->quatern_c, &image->quatern_d, &image->qoffset_x, &image->qoffset_y, &image->qoffset_z, NULL, NULL, NULL, &image->qfac);
1134  }
1135 
1136  if (image->sform_code > 0)
1137  {
1138  mat33 prod = nifti_mat33_mul(scaleMatrix, internal::topLeftCorner(image->sto_xyz));
1139  for (int i=0; i<3; i++)
1140  {
1141  for (int j=0; j<3; j++)
1142  image->sto_xyz.m[i][j] = prod.m[i][j];
1143  }
1144  image->sto_ijk = nifti_mat44_inverse(image->sto_xyz);
1145  }
1146  }
1147 }
1148 
1149 inline void NiftiImage::setPixunits (const std::vector<std::string> &pixunits)
1150 {
1151  for (size_t i=0; i<pixunits.size(); i++)
1152  {
1153  if (pixunits[i] == "m")
1154  image->xyz_units = NIFTI_UNITS_METER;
1155  else if (pixunits[i] == "mm")
1156  image->xyz_units = NIFTI_UNITS_MM;
1157  else if (pixunits[i] == "um")
1158  image->xyz_units = NIFTI_UNITS_MICRON;
1159  else if (pixunits[i] == "s")
1160  image->time_units = NIFTI_UNITS_SEC;
1161  else if (pixunits[i] == "ms")
1162  image->time_units = NIFTI_UNITS_MSEC;
1163  else if (pixunits[i] == "us")
1164  image->time_units = NIFTI_UNITS_USEC;
1165  else if (pixunits[i] == "Hz")
1166  image->time_units = NIFTI_UNITS_HZ;
1167  else if (pixunits[i] == "ppm")
1168  image->time_units = NIFTI_UNITS_PPM;
1169  else if (pixunits[i] == "rad/s")
1170  image->time_units = NIFTI_UNITS_RADS;
1171  }
1172 }
1173 
1174 inline NiftiImage & NiftiImage::rescale (const std::vector<float> &scales)
1175 {
1176  std::vector<float> pixdim(image->pixdim+1, image->pixdim+4);
1177 
1178  for (int i=0; i<std::min(3, int(scales.size())); i++)
1179  {
1180  if (scales[i] != 1.0)
1181  {
1182  pixdim[i] /= scales[i];
1183  image->dim[i+1] = static_cast<int>(std::floor(image->dim[i+1] * scales[i]));
1184  }
1185  }
1186 
1188  nifti_update_dims_from_array(image);
1189 
1190  // Data vector is now the wrong size, so drop it
1191  if (!persistent)
1192  nifti_image_unload(image);
1193 
1194  image->scl_slope = 0.0;
1195  image->scl_inter = 0.0;
1196 
1197  return *this;
1198 }
1199 
1200 inline NiftiImage & NiftiImage::reorient (const int icode, const int jcode, const int kcode)
1201 {
1202  if (this->isNull())
1203  return *this;
1204  if (image->qform_code == 0 && image->sform_code == 0)
1205  {
1206  Rf_warning("Image qform and sform codes are both zero, so it cannot be reoriented");
1207  return *this;
1208  }
1209 
1210  int used[6] = { 0, 0, 0, 0, 0, 0 };
1211  used[icode-1] = 1;
1212  used[jcode-1] = 1;
1213  used[kcode-1] = 1;
1214  if (used[0]+used[1] != 1 || used[2]+used[3] != 1 || used[4]+used[5] != 1)
1215  throw std::runtime_error("Each canonical axis should be used exactly once");
1216 
1217  const int codes[3] = { icode, jcode, kcode };
1218  const mat44 native = this->xform();
1219 
1220  // Calculate the origin, which requires inverting the current xform
1221  // Here we use a simplified formula that exploits blockwise inversion and the nature of xforms
1222  internal::vec3 origin;
1223  for (int i=0; i<3; i++)
1224  origin.v[i] = native.m[i][3];
1225  origin = -internal::matrixVectorProduct(nifti_mat33_inverse(internal::topLeftCorner(native)), origin);
1226 
1227  // Create a target xform (rotation matrix only)
1228  mat33 target;
1229  for (int j=0; j<3; j++)
1230  {
1231  for (int i=0; i<3; i++)
1232  target.m[i][j] = 0.0;
1233 
1234  switch (codes[j])
1235  {
1236  case NIFTI_L2R: target.m[0][j] = 1.0; break;
1237  case NIFTI_R2L: target.m[0][j] = -1.0; break;
1238  case NIFTI_P2A: target.m[1][j] = 1.0; break;
1239  case NIFTI_A2P: target.m[1][j] = -1.0; break;
1240  case NIFTI_I2S: target.m[2][j] = 1.0; break;
1241  case NIFTI_S2I: target.m[2][j] = -1.0; break;
1242  }
1243  }
1244 
1245  // Extract (inverse of) canonical axis matrix from native xform
1246  int nicode, njcode, nkcode;
1247  nifti_mat44_to_orientation(native, &nicode, &njcode, &nkcode);
1248  int ncodes[3] = { nicode, njcode, nkcode };
1249  mat33 nativeAxesTransposed;
1250  for (int i=0; i<3; i++)
1251  {
1252  for (int j=0; j<3; j++)
1253  nativeAxesTransposed.m[i][j] = 0.0;
1254 
1255  switch (ncodes[i])
1256  {
1257  case NIFTI_L2R: nativeAxesTransposed.m[i][0] = 1.0; break;
1258  case NIFTI_R2L: nativeAxesTransposed.m[i][0] = -1.0; break;
1259  case NIFTI_P2A: nativeAxesTransposed.m[i][1] = 1.0; break;
1260  case NIFTI_A2P: nativeAxesTransposed.m[i][1] = -1.0; break;
1261  case NIFTI_I2S: nativeAxesTransposed.m[i][2] = 1.0; break;
1262  case NIFTI_S2I: nativeAxesTransposed.m[i][2] = -1.0; break;
1263  }
1264  }
1265 
1266  // Check for no-op case
1267  if (icode == nicode && jcode == njcode && kcode == nkcode)
1268  return *this;
1269 
1270  // The transform is t(approx_old_xform) %*% target_xform
1271  // The new xform is old_xform %*% transform
1272  // NB: "transform" is really 4x4, but the last row is simple and the last column is filled below
1273  mat33 transform = nifti_mat33_mul(nativeAxesTransposed, target);
1274  mat44 result;
1275  for (int i=0; i<4; i++)
1276  {
1277  for (int j=0; j<3; j++)
1278  result.m[i][j] = native.m[i][0] * transform.m[0][j] + native.m[i][1] * transform.m[1][j] + native.m[i][2] * transform.m[2][j];
1279 
1280  result.m[3][i] = i == 3 ? 1.0 : 0.0;
1281  }
1282 
1283  // Extract the mapping between dimensions and the signs
1284  // These vectors are all indexed in the target space, except "revsigns"
1285  int locs[3], signs[3], newdim[3], revsigns[3];
1286  float newpixdim[3];
1287  double maxes[3] = { R_NegInf, R_NegInf, R_NegInf };
1288  internal::vec3 offset;
1289  for (int j=0; j<3; j++)
1290  {
1291  // Find the largest absolute value in each column, which gives the old dimension corresponding to each new dimension
1292  for (int i=0; i<3; i++)
1293  {
1294  const double value = static_cast<double>(transform.m[i][j]);
1295  if (fabs(value) > maxes[j])
1296  {
1297  maxes[j] = fabs(value);
1298  signs[j] = value > 0.0 ? 1 : -1;
1299  locs[j] = i;
1300  }
1301  }
1302 
1303  // Obtain the sign for the reverse mapping
1304  revsigns[locs[j]] = signs[j];
1305 
1306  // Permute dim and pixdim
1307  newdim[j] = image->dim[locs[j]+1];
1308  newpixdim[j] = image->pixdim[locs[j]+1];
1309 
1310  // Flip and/or permute the origin
1311  if (signs[j] < 0)
1312  offset.v[j] = image->dim[locs[j]+1] - origin.v[locs[j]] - 1.0;
1313  else
1314  offset.v[j] = origin.v[locs[j]];
1315  }
1316 
1317  // Convert the origin back to an xform offset and insert it
1318  offset = -internal::matrixVectorProduct(internal::topLeftCorner(result), offset);
1319  for (int i=0; i<3; i++)
1320  result.m[i][3] = offset.v[i];
1321 
1322  // Update the xforms with nonzero codes
1323  if (image->qform_code > 0)
1324  {
1325  image->qto_xyz = result;
1326  image->qto_ijk = nifti_mat44_inverse(image->qto_xyz);
1327  nifti_mat44_to_quatern(image->qto_xyz, &image->quatern_b, &image->quatern_c, &image->quatern_d, &image->qoffset_x, &image->qoffset_y, &image->qoffset_z, NULL, NULL, NULL, &image->qfac);
1328  }
1329  if (image->sform_code > 0)
1330  {
1331  image->sto_xyz = result;
1332  image->sto_ijk = nifti_mat44_inverse(image->sto_xyz);
1333  }
1334 
1335  // Calculate strides: the step in target space associated with each dimension in source space
1336  ptrdiff_t strides[3];
1337  strides[locs[0]] = 1;
1338  strides[locs[1]] = strides[locs[0]] * image->dim[locs[0]+1];
1339  strides[locs[2]] = strides[locs[1]] * image->dim[locs[1]+1];
1340 
1341  // Permute the data (if present)
1342  if (image->data != NULL)
1343  {
1344  size_t volSize = size_t(image->nx * image->ny * image->nz);
1345  size_t nVolumes = std::max(size_t(1), image->nvox / volSize);
1346 
1347  const std::vector<double> oldData = this->getData<double>();
1348  std::vector<double> newData(image->nvox);
1349 
1350  // Where the sign is negative we need to start at the end of the dimension
1351  size_t volStart = 0;
1352  for (int i=0; i<3; i++)
1353  {
1354  if (revsigns[i] < 0)
1355  volStart += (image->dim[i+1] - 1) * strides[i];
1356  }
1357 
1358  // Iterate over the data and place it into a new vector
1359  std::vector<double>::const_iterator it = oldData.begin();
1360  for (size_t v=0; v<nVolumes; v++)
1361  {
1362  for (int k=0; k<image->nz; k++)
1363  {
1364  ptrdiff_t offset = k * strides[2] * revsigns[2];
1365  for (int j=0; j<image->ny; j++)
1366  {
1367  for (int i=0; i<image->nx; i++)
1368  {
1369  newData[volStart + offset] = *it++;
1370  offset += strides[0] * revsigns[0];
1371  }
1372  offset += strides[1] * revsigns[1] - image->nx * strides[0] * revsigns[0];
1373  }
1374  }
1375  volStart += volSize;
1376  }
1377 
1378  // Replace the existing data in the image
1379  this->replaceData(newData);
1380  }
1381 
1382  // Copy new dims and pixdims in
1383  // NB: Old dims are used above, so this must happen last
1384  std::copy(newdim, newdim+3, image->dim+1);
1385  std::copy(newpixdim, newpixdim+3, image->pixdim+1);
1386  nifti_update_dims_from_array(image);
1387 
1388  return *this;
1389 }
1390 
1391 inline NiftiImage & NiftiImage::reorient (const std::string &orientation)
1392 {
1393  if (orientation.length() != 3)
1394  throw std::runtime_error("Orientation string should have exactly three characters");
1395 
1396  int codes[3];
1397  for (int i=0; i<3; i++)
1398  {
1399  switch(orientation[i])
1400  {
1401  case 'r': case 'R': codes[i] = NIFTI_L2R; break;
1402  case 'l': case 'L': codes[i] = NIFTI_R2L; break;
1403  case 'a': case 'A': codes[i] = NIFTI_P2A; break;
1404  case 'p': case 'P': codes[i] = NIFTI_A2P; break;
1405  case 's': case 'S': codes[i] = NIFTI_I2S; break;
1406  case 'i': case 'I': codes[i] = NIFTI_S2I; break;
1407 
1408  default:
1409  throw std::runtime_error("Orientation string is invalid");
1410  }
1411  }
1412 
1413  return reorient(codes[0], codes[1], codes[2]);
1414 }
1415 
1416 #ifdef USING_R
1417 
1418 inline NiftiImage & NiftiImage::update (const Rcpp::RObject &object)
1419 {
1420  if (Rf_isVectorList(object))
1421  {
1422  Rcpp::List list(object);
1423  nifti_1_header *header = NULL;
1424  if (this->isNull())
1425  {
1426  header = nifti_make_new_header(NULL, DT_FLOAT64);
1427  internal::updateHeader(header, list, true);
1428  }
1429  else
1430  {
1431  header = (nifti_1_header *) calloc(1, sizeof(nifti_1_header));
1432  *header = nifti_convert_nim2nhdr(image);
1433  internal::updateHeader(header, list, true);
1434  }
1435 
1436  if (header != NULL)
1437  {
1438  nifti_image *newImage = nifti_convert_nhdr2nim(*header, NULL);
1439  if (this->image->data != NULL)
1440  {
1441  size_t dataSize = nifti_get_volsize(image);
1442  newImage->data = calloc(1, dataSize);
1443  memcpy(newImage->data, image->data, dataSize);
1444  }
1445  if (!persistent)
1446  nifti_image_free(image);
1447  this->image = newImage;
1448  }
1449  }
1450  else if (object.hasAttribute("dim"))
1451  {
1452  for (int i=0; i<8; i++)
1453  image->dim[i] = 0;
1454  const std::vector<int> dimVector = object.attr("dim");
1455 
1456  const int nDims = std::min(7, int(dimVector.size()));
1457  image->dim[0] = nDims;
1458  for (int i=0; i<nDims; i++)
1459  image->dim[i+1] = dimVector[i];
1460 
1461  if (object.hasAttribute("pixdim"))
1462  {
1463  const std::vector<float> pixdimVector = object.attr("pixdim");
1464  updatePixdim(pixdimVector);
1465  }
1466 
1467  if (object.hasAttribute("pixunits"))
1468  {
1469  const std::vector<std::string> pixunitsVector = object.attr("pixunits");
1470  setPixunits(pixunitsVector);
1471  }
1472 
1473  // This NIfTI-1 library function clobbers dim[0] if the last dimension is unitary; we undo that here
1474  nifti_update_dims_from_array(image);
1475  image->dim[0] = image->ndim = nDims;
1476 
1477  image->datatype = NiftiImage::sexpTypeToNiftiType(object.sexp_type());
1478  nifti_datatype_sizes(image->datatype, &image->nbyper, NULL);
1479 
1480  if (!persistent)
1481  nifti_image_unload(image);
1482 
1483  const size_t dataSize = nifti_get_volsize(image);
1484  image->data = calloc(1, dataSize);
1485  if (image->datatype == DT_INT32)
1486  memcpy(image->data, INTEGER(object), dataSize);
1487  else
1488  memcpy(image->data, REAL(object), dataSize);
1489 
1490  image->scl_slope = 0.0;
1491  image->scl_inter = 0.0;
1492  }
1493 
1494  return *this;
1495 }
1496 
1497 #endif // USING_R
1498 
1499 inline mat44 NiftiImage::xform (const bool preferQuaternion) const
1500 {
1501  if (image == NULL)
1502  {
1503  mat44 matrix;
1504  for (int i=0; i<4; i++)
1505  {
1506  for (int j=0; j<4; j++)
1507  matrix.m[i][j] = 0.0;
1508  }
1509  return matrix;
1510  }
1511  else if (image->qform_code <= 0 && image->sform_code <= 0)
1512  {
1513  // No qform or sform so use pixdim (NB: other software may assume differently)
1514  mat44 matrix;
1515  for (int i=0; i<4; i++)
1516  {
1517  for (int j=0; j<4; j++)
1518  {
1519  if (i != j)
1520  matrix.m[i][j] = 0.0;
1521  else if (i == 3)
1522  matrix.m[3][3] = 1.0;
1523  else
1524  matrix.m[i][j] = (image->pixdim[i+1]==0.0 ? 1.0 : image->pixdim[i+1]);
1525  }
1526  }
1527  return matrix;
1528  }
1529  else if ((preferQuaternion && image->qform_code > 0) || image->sform_code <= 0)
1530  return image->qto_xyz;
1531  else
1532  return image->sto_xyz;
1533 }
1534 
1535 template <typename TargetType>
1536 inline std::vector<TargetType> NiftiImage::Block::getData (const bool useSlope) const
1537 {
1538  if (image.isNull())
1539  return std::vector<TargetType>();
1540 
1541  size_t blockSize = 1;
1542  for (int i=1; i<dimension; i++)
1543  blockSize *= image->dim[i];
1544 
1545  std::vector<TargetType> data(blockSize);
1546  internal::DataConverter<TargetType> *converter = NULL;
1547  if (useSlope && image.isDataScaled())
1548  converter = new internal::DataConverter<TargetType>(image->scl_slope, image->scl_inter);
1549 
1550  internal::convertData<TargetType>(image->data, image->datatype, blockSize, data.begin(), blockSize*index, converter);
1551 
1552  delete converter;
1553  return data;
1554 }
1555 
1556 template <typename TargetType>
1557 inline std::vector<TargetType> NiftiImage::getData (const bool useSlope) const
1558 {
1559  if (this->isNull())
1560  return std::vector<TargetType>();
1561 
1562  std::vector<TargetType> data(image->nvox);
1563  internal::DataConverter<TargetType> *converter = NULL;
1564  if (useSlope && this->isDataScaled())
1565  converter = new internal::DataConverter<TargetType>(image->scl_slope, image->scl_inter);
1566 
1567  internal::convertData<TargetType>(image->data, image->datatype, image->nvox, data.begin(), 0, converter);
1568 
1569  delete converter;
1570  return data;
1571 }
1572 
1573 inline NiftiImage & NiftiImage::changeDatatype (const short datatype, const bool useSlope)
1574 {
1575  if (this->isNull() || image->datatype == datatype)
1576  return *this;
1577 
1578  if (useSlope && this->isDataScaled())
1579  throw std::runtime_error("Resetting the slope and intercept for an image with them already set is not supported");
1580 
1581  if (image->data != NULL)
1582  {
1583  int bytesPerPixel;
1584  nifti_datatype_sizes(datatype, &bytesPerPixel, NULL);
1585  void *data = calloc(image->nvox, bytesPerPixel);
1586 
1587  switch (datatype)
1588  {
1589  case DT_UINT8:
1590  {
1591  internal::DataConverter<uint8_t> converter(useSlope ? internal::DataConverter<uint8_t>::IndexMode : internal::DataConverter<uint8_t>::CastMode);
1592  internal::convertData(image->data, image->datatype, image->nvox, static_cast<uint8_t *>(data), 0, &converter);
1593  image->scl_slope = static_cast<float>(converter.getSlope());
1594  image->scl_inter = static_cast<float>(converter.getIntercept());
1595  break;
1596  }
1597 
1598  case DT_INT16:
1599  {
1600  internal::DataConverter<int16_t> converter(useSlope ? internal::DataConverter<int16_t>::IndexMode : internal::DataConverter<int16_t>::CastMode);
1601  internal::convertData(image->data, image->datatype, image->nvox, static_cast<int16_t *>(data), 0, &converter);
1602  image->scl_slope = static_cast<float>(converter.getSlope());
1603  image->scl_inter = static_cast<float>(converter.getIntercept());
1604  break;
1605  }
1606 
1607  case DT_INT32:
1608  {
1609  internal::DataConverter<int32_t> converter(useSlope ? internal::DataConverter<int32_t>::IndexMode : internal::DataConverter<int32_t>::CastMode);
1610  internal::convertData(image->data, image->datatype, image->nvox, static_cast<int32_t *>(data), 0, &converter);
1611  image->scl_slope = static_cast<float>(converter.getSlope());
1612  image->scl_inter = static_cast<float>(converter.getIntercept());
1613  break;
1614  }
1615 
1616  case DT_FLOAT32:
1617  {
1618  internal::DataConverter<float> converter(useSlope ? internal::DataConverter<float>::IndexMode : internal::DataConverter<float>::CastMode);
1619  internal::convertData(image->data, image->datatype, image->nvox, static_cast<float *>(data), 0, &converter);
1620  image->scl_slope = static_cast<float>(converter.getSlope());
1621  image->scl_inter = static_cast<float>(converter.getIntercept());
1622  break;
1623  }
1624 
1625  case DT_FLOAT64:
1626  {
1627  internal::DataConverter<double> converter(useSlope ? internal::DataConverter<double>::IndexMode : internal::DataConverter<double>::CastMode);
1628  internal::convertData(image->data, image->datatype, image->nvox, static_cast<double *>(data), 0, &converter);
1629  image->scl_slope = static_cast<float>(converter.getSlope());
1630  image->scl_inter = static_cast<float>(converter.getIntercept());
1631  break;
1632  }
1633 
1634  case DT_INT8:
1635  {
1636  internal::DataConverter<int8_t> converter(useSlope ? internal::DataConverter<int8_t>::IndexMode : internal::DataConverter<int8_t>::CastMode);
1637  internal::convertData(image->data, image->datatype, image->nvox, static_cast<int8_t *>(data), 0, &converter);
1638  image->scl_slope = static_cast<float>(converter.getSlope());
1639  image->scl_inter = static_cast<float>(converter.getIntercept());
1640  break;
1641  }
1642 
1643  case DT_UINT16:
1644  {
1645  internal::DataConverter<uint16_t> converter(useSlope ? internal::DataConverter<uint16_t>::IndexMode : internal::DataConverter<uint16_t>::CastMode);
1646  internal::convertData(image->data, image->datatype, image->nvox, static_cast<uint16_t *>(data), 0, &converter);
1647  image->scl_slope = static_cast<float>(converter.getSlope());
1648  image->scl_inter = static_cast<float>(converter.getIntercept());
1649  break;
1650  }
1651 
1652  case DT_UINT32:
1653  {
1654  internal::DataConverter<uint32_t> converter(useSlope ? internal::DataConverter<uint32_t>::IndexMode : internal::DataConverter<uint32_t>::CastMode);
1655  internal::convertData(image->data, image->datatype, image->nvox, static_cast<uint32_t *>(data), 0, &converter);
1656  image->scl_slope = static_cast<float>(converter.getSlope());
1657  image->scl_inter = static_cast<float>(converter.getIntercept());
1658  break;
1659  }
1660 
1661  case DT_INT64:
1662  {
1663  internal::DataConverter<int64_t> converter(useSlope ? internal::DataConverter<int64_t>::IndexMode : internal::DataConverter<int64_t>::CastMode);
1664  internal::convertData(image->data, image->datatype, image->nvox, static_cast<int64_t *>(data), 0, &converter);
1665  image->scl_slope = static_cast<float>(converter.getSlope());
1666  image->scl_inter = static_cast<float>(converter.getIntercept());
1667  break;
1668  }
1669 
1670  case DT_UINT64:
1671  {
1672  internal::DataConverter<uint64_t> converter(useSlope ? internal::DataConverter<uint64_t>::IndexMode : internal::DataConverter<uint64_t>::CastMode);
1673  internal::convertData(image->data, image->datatype, image->nvox, static_cast<uint64_t *>(data), 0, &converter);
1674  image->scl_slope = static_cast<float>(converter.getSlope());
1675  image->scl_inter = static_cast<float>(converter.getIntercept());
1676  break;
1677  }
1678 
1679  default:
1680  throw std::runtime_error("Unsupported data type (" + std::string(nifti_datatype_string(datatype)) + ")");
1681  }
1682 
1683  nifti_image_unload(image);
1684  image->data = data;
1685  }
1686 
1687  image->datatype = datatype;
1688  nifti_datatype_sizes(datatype, &image->nbyper, &image->swapsize);
1689 
1690  return *this;
1691 }
1692 
1693 inline NiftiImage & NiftiImage::changeDatatype (const std::string &datatype, const bool useSlope)
1694 {
1695  return changeDatatype(internal::stringToDatatype(datatype), useSlope);
1696 }
1697 
1698 template <typename SourceType>
1699 inline NiftiImage & NiftiImage::replaceData (const std::vector<SourceType> &data, const short datatype)
1700 {
1701  if (this->isNull())
1702  return *this;
1703  else if (data.size() != image->nvox)
1704  throw std::runtime_error("New data length does not match the number of voxels in the image");
1705 
1706  if (datatype != DT_NONE)
1707  {
1708  nifti_image_unload(image);
1709  image->datatype = datatype;
1710  nifti_datatype_sizes(datatype, &image->nbyper, &image->swapsize);
1711  }
1712 
1713  if (image->data == NULL)
1714  image->data = calloc(image->nvox, image->nbyper);
1715  internal::replaceData<SourceType>(data.begin(), data.end(), image->data, image->datatype);
1716 
1717  image->scl_slope = 0.0;
1718  image->scl_inter = 0.0;
1719  image->cal_min = static_cast<float>(*std::min_element(data.begin(), data.end()));
1720  image->cal_max = static_cast<float>(*std::max_element(data.begin(), data.end()));
1721 
1722  return *this;
1723 }
1724 
1725 inline void NiftiImage::toFile (const std::string fileName, const short datatype) const
1726 {
1727  const bool changingDatatype = (datatype != DT_NONE && !this->isNull() && datatype != image->datatype);
1728 
1729  // Copy the source image only if the datatype will be changed
1730  NiftiImage imageToWrite(image, changingDatatype);
1731 
1732  if (changingDatatype)
1733  imageToWrite.changeDatatype(datatype, true);
1734  else
1735  imageToWrite.setPersistence(true);
1736 
1737  const int status = nifti_set_filenames(imageToWrite, fileName.c_str(), false, true);
1738  if (status != 0)
1739  throw std::runtime_error("Failed to set filenames for NIfTI object");
1740  nifti_image_write(imageToWrite);
1741 }
1742 
1743 inline void NiftiImage::toFile (const std::string fileName, const std::string &datatype) const
1744 {
1745  toFile(fileName, internal::stringToDatatype(datatype));
1746 }
1747 
1748 #ifdef USING_R
1749 
1750 inline Rcpp::RObject NiftiImage::toArray () const
1751 {
1752  Rcpp::RObject array;
1753 
1754  if (this->isNull())
1755  return array;
1756  else if (this->isDataScaled())
1757  {
1758  array = internal::imageDataToArray<REALSXP>(image);
1759  std::transform(REAL(array), REAL(array)+Rf_length(array), REAL(array), internal::DataConverter<double>(image->scl_slope,image->scl_inter));
1760  }
1761  else
1762  {
1763  switch (image->datatype)
1764  {
1765  case DT_UINT8:
1766  case DT_INT16:
1767  case DT_INT32:
1768  case DT_INT8:
1769  case DT_UINT16:
1770  case DT_UINT32:
1771  case DT_INT64:
1772  case DT_UINT64:
1773  array = internal::imageDataToArray<INTSXP>(image);
1774  break;
1775 
1776  case DT_FLOAT32:
1777  case DT_FLOAT64:
1778  array = internal::imageDataToArray<REALSXP>(image);
1779  break;
1780 
1781  default:
1782  throw std::runtime_error("Unsupported data type (" + std::string(nifti_datatype_string(image->datatype)) + ")");
1783  }
1784  }
1785 
1786  internal::addAttributes(array, image);
1787  array.attr("class") = Rcpp::CharacterVector::create("niftiImage", "array");
1788  return array;
1789 }
1790 
1791 inline Rcpp::RObject NiftiImage::toPointer (const std::string label) const
1792 {
1793  if (this->isNull())
1794  return Rcpp::RObject();
1795  else
1796  {
1797  Rcpp::RObject string = Rcpp::wrap(label);
1798  internal::addAttributes(string, image, false);
1799  string.attr("class") = Rcpp::CharacterVector::create("internalImage", "niftiImage");
1800  return string;
1801  }
1802 }
1803 
1804 inline Rcpp::RObject NiftiImage::toArrayOrPointer (const bool internal, const std::string label) const
1805 {
1806  return (internal ? toPointer(label) : toArray());
1807 }
1808 
1809 inline Rcpp::RObject NiftiImage::headerToList () const
1810 {
1811  if (this->image == NULL)
1812  return Rcpp::RObject();
1813 
1814  nifti_1_header header = nifti_convert_nim2nhdr(this->image);
1815  Rcpp::List result;
1816 
1817  result["sizeof_hdr"] = header.sizeof_hdr;
1818 
1819  result["dim_info"] = int(header.dim_info);
1820  result["dim"] = std::vector<short>(header.dim, header.dim+8);
1821 
1822  result["intent_p1"] = header.intent_p1;
1823  result["intent_p2"] = header.intent_p2;
1824  result["intent_p3"] = header.intent_p3;
1825  result["intent_code"] = header.intent_code;
1826 
1827  result["datatype"] = header.datatype;
1828  result["bitpix"] = header.bitpix;
1829 
1830  result["slice_start"] = header.slice_start;
1831  result["pixdim"] = std::vector<float>(header.pixdim, header.pixdim+8);
1832  result["vox_offset"] = header.vox_offset;
1833  result["scl_slope"] = header.scl_slope;
1834  result["scl_inter"] = header.scl_inter;
1835  result["slice_end"] = header.slice_end;
1836  result["slice_code"] = int(header.slice_code);
1837  result["xyzt_units"] = int(header.xyzt_units);
1838  result["cal_max"] = header.cal_max;
1839  result["cal_min"] = header.cal_min;
1840  result["slice_duration"] = header.slice_duration;
1841  result["toffset"] = header.toffset;
1842  result["descrip"] = std::string(header.descrip, 80);
1843  result["aux_file"] = std::string(header.aux_file, 24);
1844 
1845  result["qform_code"] = header.qform_code;
1846  result["sform_code"] = header.sform_code;
1847  result["quatern_b"] = header.quatern_b;
1848  result["quatern_c"] = header.quatern_c;
1849  result["quatern_d"] = header.quatern_d;
1850  result["qoffset_x"] = header.qoffset_x;
1851  result["qoffset_y"] = header.qoffset_y;
1852  result["qoffset_z"] = header.qoffset_z;
1853  result["srow_x"] = std::vector<float>(header.srow_x, header.srow_x+4);
1854  result["srow_y"] = std::vector<float>(header.srow_y, header.srow_y+4);
1855  result["srow_z"] = std::vector<float>(header.srow_z, header.srow_z+4);
1856 
1857  result["intent_name"] = std::string(header.intent_name, 16);
1858  result["magic"] = std::string(header.magic, 4);
1859 
1860  result.attr("class") = Rcpp::CharacterVector::create("niftiHeader");
1861 
1862  return result;
1863 }
1864 
1865 #endif // USING_R
1866 
1867 } // namespace
1868 
1869 #endif
NiftiImage & reorient(const int i, const int j, const int k)
Reorient the image by permuting dimensions and potentially reversing some.
Definition: NiftiImage.h:1200
NiftiImage()
Default constructor.
Definition: NiftiImage.h:313
const int index
The location along dimension.
Definition: NiftiImage.h:96
static std::string xformToString(const mat44 matrix)
Convert a 4x4 xform matrix to a string describing its canonical axes.
Definition: NiftiImage.h:186
std::vector< int > dim() const
Return the dimensions of the image.
Definition: NiftiImage.h:499
const nifti_image * operator->() const
Allows a NiftiImage object to be treated as a pointer to a const nifti_image.
Definition: NiftiImage.h:421
bool persistent
Marker of persistence, which determines whether the nifti_image should be freed on destruction...
Definition: NiftiImage.h:245
std::vector< float > pixdim() const
Return the dimensions of the pixels or voxels in the image.
Definition: NiftiImage.h:511
NiftiImage(const std::string &path, const bool readData=true)
Initialise using a path string.
Definition: NiftiImage.h:366
NiftiImage & rescale(const std::vector< float > &scales)
Rescale the image, changing its image dimensions and pixel dimensions.
Definition: NiftiImage.h:1174
NiftiImage(nifti_image *const image, const bool copy=false)
Initialise using an existing nifti_image pointer.
Definition: NiftiImage.h:348
Definition: NiftiImage.h:38
Block & operator=(const NiftiImage &source)
Copy assignment operator, which allows a block in one image to be replaced with the contents of anoth...
Definition: NiftiImage.h:120
bool isNull() const
Determine whether or not the internal pointer is NULL.
Definition: NiftiImage.h:472
Inner class referring to a subset of an image.
Definition: NiftiImage.h:92
Block slice(const int i)
Extract a slice block from a 3D image.
Definition: NiftiImage.h:661
NiftiImage & replaceData(const std::vector< SourceType > &data, const short datatype=DT_NONE)
Replace the pixel data in the image with the contents of a vector.
Definition: NiftiImage.h:1699
bool isDataScaled() const
Determine whether nontrivial scale and slope parameters are set.
Definition: NiftiImage.h:482
const Block slice(const int i) const
Extract a slice block from a 3D image.
Definition: NiftiImage.h:654
const Block volume(const int i) const
Extract a volume block from a 4D image.
Definition: NiftiImage.h:668
bool isPersistent() const
Determine whether or not the image is marked as persistent.
Definition: NiftiImage.h:477
std::vector< TargetType > getData(const bool useSlope=true) const
Extract a vector of data from the image, casting it to any required element type. ...
Definition: NiftiImage.h:1557
std::vector< TargetType > getData(const bool useSlope=true) const
Extract a vector of data from a block, casting it to any required element type.
Definition: NiftiImage.h:1536
void setPixunits(const std::vector< std::string > &pixunits)
Modify the pixel dimension units.
Definition: NiftiImage.h:1149
NiftiImage(const Block &source)
Initialise from a block, copying in the data.
Definition: NiftiImage.h:333
virtual ~NiftiImage()
Destructor which frees the wrapped pointer, unless the object is marked as persistent.
Definition: NiftiImage.h:397
Thin wrapper around a C-style nifti_image struct that allows C++-style destruction.
Definition: NiftiImage.h:85
NiftiImage(const NiftiImage &source)
Copy constructor.
Definition: NiftiImage.h:320
int nBlocks() const
Return the number of blocks in the image.
Definition: NiftiImage.h:623
const int dimension
The dimension along which the block applies (which should be the last)
Definition: NiftiImage.h:95
const NiftiImage & image
The parent image.
Definition: NiftiImage.h:94
NiftiImage & operator=(const NiftiImage &source)
Copy assignment operator, which copies from its argument.
Definition: NiftiImage.h:432
static int fileVersion(const std::string &path)
Get the NIfTI format version used by the file at the specified path.
Definition: NiftiImage.h:214
static mat33 xformToRotation(const mat44 matrix)
Extract the pure rotation part of a 4x4 xform matrix.
Definition: NiftiImage.h:173
void toFile(const std::string fileName, const short datatype=DT_NONE) const
Write the image to a NIfTI-1 file.
Definition: NiftiImage.h:1725
void updatePixdim(const std::vector< float > &pixdim)
Modify the pixel dimensions, and potentially the xform matrices to match.
Definition: NiftiImage.h:1096
const Block block(const int i) const
Extract a block from the image.
Definition: NiftiImage.h:638
NiftiImage & setPersistence(const bool persistent)
Marked the image as persistent, so that it can be passed back to R.
Definition: NiftiImage.h:459
Definition: NiftiImage.h:42
NiftiImage & changeDatatype(const short datatype, const bool useSlope=false)
Change the datatype of the image, casting the pixel data if present.
Definition: NiftiImage.h:1573
void copy(const nifti_image *source)
Copy the contents of a nifti_image to create a new image.
Definition: NiftiImage.h:727
mat44 xform(const bool preferQuaternion=true) const
Obtain an xform matrix, indicating the orientation of the image.
Definition: NiftiImage.h:1499
nifti_image * image
The wrapped nifti_image pointer.
Definition: NiftiImage.h:244
int nDims() const
Return the number of dimensions in the image.
Definition: NiftiImage.h:487
Block volume(const int i)
Extract a volume block from a 4D image.
Definition: NiftiImage.h:675
NiftiImage & drop()
Drop unitary dimensions.
Definition: NiftiImage.h:525
Block(const NiftiImage &image, const int dimension, const int index)
Standard constructor for this class.
Definition: NiftiImage.h:105
NiftiImage & dropData()
Drop the data from the image, retaining only the metadata.
Definition: NiftiImage.h:574
Block block(const int i)
Extract a block from the image.
Definition: NiftiImage.h:647