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 #define ISNAN(x) (x != x)
13 
14 #include <stdint.h>
15 #include <cstddef>
16 #include <cmath>
17 #include <string>
18 #include <sstream>
19 #include <vector>
20 #include <stdexcept>
21 #include <algorithm>
22 #include <map>
23 #include <locale>
24 #include <limits>
25 
26 #endif
27 
28 
29 #include "niftilib/nifti1_io.h"
30 
39 namespace RNifti {
40 
46 {
47 public:
52  struct Block
53  {
54  const NiftiImage &image;
55  const int dimension;
56  const int index;
65  Block (const NiftiImage &image, const int dimension, const int index)
67  {
68  if (dimension != image->ndim)
69  throw std::runtime_error("Blocks must be along the last dimension in the image");
70  }
71 
80  Block & operator= (const NiftiImage &source)
81  {
82  if (source->datatype != image->datatype)
83  throw std::runtime_error("New data does not have the same datatype as the target block");
84  if (source->scl_slope != image->scl_slope || source->scl_inter != image->scl_inter)
85  throw std::runtime_error("New data does not have the same scale parameters as the target block");
86 
87  size_t blockSize = 1;
88  for (int i=1; i<dimension; i++)
89  blockSize *= image->dim[i];
90 
91  if (blockSize != source->nvox)
92  throw std::runtime_error("New data does not have the same size as the target block");
93 
94  blockSize *= image->nbyper;
95  memcpy(static_cast<char*>(image->data) + blockSize*index, source->data, blockSize);
96  return *this;
97  }
98 
102  template <typename TargetType>
103  std::vector<TargetType> getData () const;
104  };
105 
106 #ifdef USING_R
107 
113  static short sexpTypeToNiftiType (const int sexpType)
114  {
115  if (sexpType == INTSXP || sexpType == LGLSXP)
116  return DT_INT32;
117  else if (sexpType == REALSXP)
118  return DT_FLOAT64;
119  else
120  throw std::runtime_error("Array elements must be numeric");
121  }
122 #endif
123 
129  static std::string xformToString (const mat44 matrix)
130  {
131  int icode, jcode, kcode;
132  nifti_mat44_to_orientation(matrix, &icode, &jcode, &kcode);
133 
134  int codes[3] = { icode, jcode, kcode };
135  std::string result("---");
136  for (int i=0; i<3; i++)
137  {
138  switch (codes[i])
139  {
140  case NIFTI_L2R: result[i] = 'R'; break;
141  case NIFTI_R2L: result[i] = 'L'; break;
142  case NIFTI_P2A: result[i] = 'A'; break;
143  case NIFTI_A2P: result[i] = 'P'; break;
144  case NIFTI_I2S: result[i] = 'S'; break;
145  case NIFTI_S2I: result[i] = 'I'; break;
146  }
147  }
148  return result;
149  }
150 
157  static int fileVersion (const std::string &path)
158  {
159  nifti_1_header *header = nifti_read_header(path.c_str(), NULL, false);
160  if (header == NULL)
161  return -1;
162  else
163  {
164  int version = NIFTI_VERSION(*header);
165  if (version == 0)
166  {
167  // NIfTI-2 has a 540-byte header - check for this or its byte-swapped equivalent
168  if (header->sizeof_hdr == 540 || header->sizeof_hdr == 469893120)
169  {
170  // The magic number has moved in NIfTI-2, so find it by byte offset
171  const char *magic = (char *) header + 4;
172  if (strncmp(magic,"ni2",3) == 0 || strncmp(magic,"n+2",3) == 0)
173  version = 2;
174  }
175  else if (!nifti_hdr_looks_good(header))
176  {
177  // Not plausible as ANALYZE, so return -1
178  version = -1;
179  }
180  }
181  return version;
182  }
183  }
184 
185 
186 protected:
187  nifti_image *image;
188  bool persistent;
194  void copy (const nifti_image *source);
195 
200  void copy (const NiftiImage &source);
201 
206  void copy (const Block &source);
207 
208 
209 #ifdef USING_R
210 
216  void initFromNiftiS4 (const Rcpp::RObject &object, const bool copyData = true);
217 
223  void initFromMriImage (const Rcpp::RObject &object, const bool copyData = true);
224 
229  void initFromList (const Rcpp::RObject &object);
230 
236  void initFromArray (const Rcpp::RObject &object, const bool copyData = true);
237 
238 #endif
239 
244  void updatePixdim (const std::vector<float> &pixdim);
245 
250  void setPixunits (const std::vector<std::string> &pixunits);
251 
252 public:
257  : image(NULL), persistent(false) {}
258 
263  NiftiImage (const NiftiImage &source)
264  : image(NULL), persistent(false)
265  {
266  this->copy(source);
267 #ifndef NDEBUG
268  Rprintf("Creating NiftiImage with pointer %p (from NiftiImage)\n", this->image);
269 #endif
270  }
271 
276  NiftiImage (const Block &source)
277  : image(NULL), persistent(false)
278  {
279  this->copy(source);
280 #ifndef NDEBUG
281  Rprintf("Creating NiftiImage with pointer %p (from Block)\n", this->image);
282 #endif
283  }
284 
291  NiftiImage (nifti_image * const image, const bool copy = false)
292  : image(NULL), persistent(false)
293  {
294  if (copy)
295  this->copy(image);
296  else
297  this->image = image;
298 #ifndef NDEBUG
299  Rprintf("Creating NiftiImage with pointer %p (from pointer)\n", this->image);
300 #endif
301  }
302 
309  NiftiImage (const std::string &path, const bool readData = true)
310  : persistent(false)
311  {
312  this->image = nifti_image_read(path.c_str(), readData);
313  if (this->image == NULL)
314  throw std::runtime_error("Failed to read image from path " + path);
315 #ifndef NDEBUG
316  Rprintf("Creating NiftiImage with pointer %p (from string)\n", this->image);
317 #endif
318  }
319 
326  NiftiImage (const std::string &path, const std::vector<int> &volumes);
327 
328 #ifdef USING_R
329 
334  NiftiImage (const SEXP object, const bool readData = true);
335 #endif
336 
340  virtual ~NiftiImage ()
341  {
342  if (!persistent)
343  {
344 #ifndef NDEBUG
345  Rprintf("Freeing NiftiImage with pointer %p\n", this->image);
346 #endif
347  nifti_image_free(image);
348  }
349  }
350 
354  operator const nifti_image* () const { return image; }
355 
359  operator nifti_image* () { return image; }
360 
364  const nifti_image * operator-> () const { return image; }
365 
369  nifti_image * operator-> () { return image; }
370 
376  {
377  copy(source);
378 #ifndef NDEBUG
379  Rprintf("Creating NiftiImage with pointer %p (from NiftiImage)\n", this->image);
380 #endif
381  return *this;
382  }
383 
389  NiftiImage & operator= (const Block &source)
390  {
391  copy(source);
392 #ifndef NDEBUG
393  Rprintf("Creating NiftiImage with pointer %p (from Block)\n", this->image);
394 #endif
395  return *this;
396  }
397 
403  {
404  this->persistent = persistent;
405 #ifndef NDEBUG
406  if (persistent)
407  Rprintf("Setting NiftiImage with pointer %p to be persistent\n", this->image);
408 #endif
409  return *this;
410  }
411 
415  bool isNull () const { return (image == NULL); }
416 
420  bool isPersistent () const { return persistent; }
421 
425  bool isDataScaled () const { return (image != NULL && image->scl_slope != 0.0 && (image->scl_slope != 1.0 || image->scl_inter != 0.0)); }
426 
430  int nDims () const
431  {
432  if (image == NULL)
433  return 0;
434  else
435  return image->ndim;
436  }
437 
442  std::vector<int> dim () const
443  {
444  if (image == NULL)
445  return std::vector<int>();
446  else
447  return std::vector<int>(image->dim+1, image->dim+image->ndim+1);
448  }
449 
454  std::vector<float> pixdim () const
455  {
456  if (image == NULL)
457  return std::vector<float>();
458  else
459  return std::vector<float>(image->pixdim+1, image->pixdim+image->ndim+1);
460  }
461 
469  {
470  int ndim = image->ndim;
471  while (image->dim[ndim] < 2)
472  ndim--;
473  image->dim[0] = image->ndim = ndim;
474 
475  return *this;
476  }
477 
481  template <typename TargetType>
482  std::vector<TargetType> getData () const;
483 
488  NiftiImage & changeDatatype (const short datatype);
489 
494  NiftiImage & changeDatatype (const std::string &datatype);
495 
503  template <typename SourceType>
504  NiftiImage & replaceData (const std::vector<SourceType> &data, const short datatype = DT_NONE);
505 
510  {
511  nifti_image_unload(image);
512  return *this;
513  }
514 
520  NiftiImage & rescale (const std::vector<float> &scales);
521 
529  NiftiImage & reorient (const int i, const int j, const int k);
530 
538  NiftiImage & reorient (const std::string &orientation);
539 
540 #ifdef USING_R
541 
545  NiftiImage & update (const Rcpp::RObject &object);
546 #endif
547 
553  mat44 xform (const bool preferQuaternion = true) const;
554 
558  int nBlocks () const
559  {
560  if (image == NULL)
561  return 0;
562  else
563  return image->dim[image->ndim];
564  }
565 
573  const Block block (const int i) const { return Block(*this, nDims(), i); }
574 
582  Block block (const int i) { return Block(*this, nDims(), i); }
583 
589  const Block slice (const int i) const { return Block(*this, 3, i); }
590 
596  Block slice (const int i) { return Block(*this, 3, i); }
597 
603  const Block volume (const int i) const { return Block(*this, 4, i); }
604 
610  Block volume (const int i) { return Block(*this, 4, i); }
611 
617  void toFile (const std::string fileName, const short datatype = DT_NONE) const;
618 
624  void toFile (const std::string fileName, const std::string &datatype) const;
625 
626 #ifdef USING_R
627 
632  Rcpp::RObject toArray () const;
633 
639  Rcpp::RObject toPointer (const std::string label) const;
640 
647  Rcpp::RObject toArrayOrPointer (const bool internal, const std::string label) const;
648 
653  Rcpp::RObject headerToList () const;
654 
655 #endif
656 
657 };
658 
659 // Include helper functions
660 #include "lib/NiftiImage_internal.h"
661 
662 inline void NiftiImage::copy (const nifti_image *source)
663 {
664  if (image != NULL && !persistent)
665  nifti_image_free(image);
666 
667  if (source == NULL)
668  image = NULL;
669  else
670  {
671  image = nifti_copy_nim_info(source);
672  if (source->data != NULL)
673  {
674  size_t dataSize = nifti_get_volsize(source);
675  image->data = calloc(1, dataSize);
676  memcpy(image->data, source->data, dataSize);
677  }
678  }
679 
680  persistent = false;
681 }
682 
683 inline void NiftiImage::copy (const NiftiImage &source)
684 {
685  const nifti_image *sourceStruct = source;
686  copy(sourceStruct);
687 }
688 
689 inline void NiftiImage::copy (const Block &source)
690 {
691  if (image != NULL && !persistent)
692  nifti_image_free(image);
693 
694  const nifti_image *sourceStruct = source.image;
695  if (sourceStruct == NULL)
696  image = NULL;
697  else
698  {
699  image = nifti_copy_nim_info(sourceStruct);
700  image->dim[0] = source.image->dim[0] - 1;
701  image->dim[source.dimension] = 1;
702  image->pixdim[source.dimension] = 1.0;
703  nifti_update_dims_from_array(image);
704 
705  if (sourceStruct->data != NULL)
706  {
707  size_t blockSize = nifti_get_volsize(image);
708  image->data = calloc(1, blockSize);
709  memcpy(image->data, static_cast<char*>(source.image->data) + blockSize*source.index, blockSize);
710  }
711  }
712 
713  persistent = false;
714 }
715 
716 #ifdef USING_R
717 
718 // Convert an S4 "nifti" object, as defined in the oro.nifti package, to a "nifti_image" struct
719 inline void NiftiImage::initFromNiftiS4 (const Rcpp::RObject &object, const bool copyData)
720 {
721  nifti_1_header header;
722  header.sizeof_hdr = 348;
723 
724  const std::vector<short> dims = object.slot("dim_");
725  for (int i=0; i<8; i++)
726  header.dim[i] = dims[i];
727 
728  header.intent_p1 = object.slot("intent_p1");
729  header.intent_p2 = object.slot("intent_p2");
730  header.intent_p3 = object.slot("intent_p3");
731  header.intent_code = object.slot("intent_code");
732 
733  header.datatype = object.slot("datatype");
734  header.bitpix = object.slot("bitpix");
735 
736  header.slice_start = object.slot("slice_start");
737  header.slice_end = object.slot("slice_end");
738  header.slice_code = Rcpp::as<int>(object.slot("slice_code"));
739  header.slice_duration = object.slot("slice_duration");
740 
741  const std::vector<float> pixdims = object.slot("pixdim");
742  for (int i=0; i<8; i++)
743  header.pixdim[i] = pixdims[i];
744  header.xyzt_units = Rcpp::as<int>(object.slot("xyzt_units"));
745 
746  header.vox_offset = object.slot("vox_offset");
747 
748  header.scl_slope = object.slot("scl_slope");
749  header.scl_inter = object.slot("scl_inter");
750  header.toffset = object.slot("toffset");
751 
752  header.cal_max = object.slot("cal_max");
753  header.cal_min = object.slot("cal_min");
754  header.glmax = header.glmin = 0;
755 
756  strncpy(header.descrip, Rcpp::as<std::string>(object.slot("descrip")).c_str(), 79);
757  header.descrip[79] = '\0';
758  strncpy(header.aux_file, Rcpp::as<std::string>(object.slot("aux_file")).c_str(), 23);
759  header.aux_file[23] = '\0';
760  strncpy(header.intent_name, Rcpp::as<std::string>(object.slot("intent_name")).c_str(), 15);
761  header.intent_name[15] = '\0';
762  strncpy(header.magic, Rcpp::as<std::string>(object.slot("magic")).c_str(), 3);
763  header.magic[3] = '\0';
764 
765  header.qform_code = object.slot("qform_code");
766  header.sform_code = object.slot("sform_code");
767 
768  header.quatern_b = object.slot("quatern_b");
769  header.quatern_c = object.slot("quatern_c");
770  header.quatern_d = object.slot("quatern_d");
771  header.qoffset_x = object.slot("qoffset_x");
772  header.qoffset_y = object.slot("qoffset_y");
773  header.qoffset_z = object.slot("qoffset_z");
774 
775  const std::vector<float> srow_x = object.slot("srow_x");
776  const std::vector<float> srow_y = object.slot("srow_y");
777  const std::vector<float> srow_z = object.slot("srow_z");
778  for (int i=0; i<4; i++)
779  {
780  header.srow_x[i] = srow_x[i];
781  header.srow_y[i] = srow_y[i];
782  header.srow_z[i] = srow_z[i];
783  }
784 
785  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)
786  header.datatype = DT_INT32;
787  else if (header.datatype == DT_FLOAT32 || header.datatype == DT_FLOAT64)
788  header.datatype = DT_FLOAT64; // This assumes that sizeof(double) == 8
789  else
790  throw std::runtime_error("Data type is not supported");
791 
792  this->image = nifti_convert_nhdr2nim(header, NULL);
793 
794  const SEXP data = PROTECT(object.slot(".Data"));
795  if (!copyData || Rf_length(data) <= 1)
796  this->image->data = NULL;
797  else
798  {
799  const size_t dataSize = nifti_get_volsize(this->image);
800  this->image->data = calloc(1, dataSize);
801  if (header.datatype == DT_INT32)
802  {
803  Rcpp::IntegerVector intData(data);
804  std::copy(intData.begin(), intData.end(), static_cast<int32_t*>(this->image->data));
805  }
806  else
807  {
808  Rcpp::DoubleVector doubleData(data);
809  std::copy(doubleData.begin(), doubleData.end(), static_cast<double*>(this->image->data));
810  }
811  }
812  UNPROTECT(1);
813 }
814 
815 inline void NiftiImage::initFromMriImage (const Rcpp::RObject &object, const bool copyData)
816 {
817  Rcpp::Reference mriImage(object);
818  Rcpp::Function getXform = mriImage.field("getXform");
819  Rcpp::NumericMatrix xform = getXform();
820 
821  this->image = NULL;
822 
823  if (Rf_length(mriImage.field("tags")) > 0)
824  initFromList(mriImage.field("tags"));
825 
826  Rcpp::RObject data = mriImage.field("data");
827  if (data.inherits("SparseArray"))
828  {
829  Rcpp::Language call("as.array", data);
830  data = call.eval();
831  }
832 
833  const short datatype = (Rf_isNull(data) ? DT_INT32 : sexpTypeToNiftiType(data.sexp_type()));
834 
835  int dims[8] = { 0, 0, 0, 0, 0, 0, 0, 0 };
836  const std::vector<int> dimVector = mriImage.field("imageDims");
837  const int nDims = std::min(7, int(dimVector.size()));
838  dims[0] = nDims;
839  size_t nVoxels = 1;
840  for (int i=0; i<nDims; i++)
841  {
842  dims[i+1] = dimVector[i];
843  nVoxels *= dimVector[i];
844  }
845 
846  if (this->image == NULL)
847  this->image = nifti_make_new_nim(dims, datatype, FALSE);
848  else
849  {
850  std::copy(dims, dims+8, this->image->dim);
851  this->image->datatype = datatype;
852  nifti_datatype_sizes(image->datatype, &image->nbyper, NULL);
853  }
854 
855  if (copyData && !Rf_isNull(data))
856  {
857  // NB: nifti_get_volsize() will not be right here if there were tags
858  const size_t dataSize = nVoxels * image->nbyper;
859  this->image->data = calloc(1, dataSize);
860  if (datatype == DT_INT32)
861  memcpy(this->image->data, INTEGER(data), dataSize);
862  else
863  memcpy(this->image->data, REAL(data), dataSize);
864  }
865  else
866  this->image->data = NULL;
867 
868  const std::vector<float> pixdimVector = mriImage.field("voxelDims");
869  const int pixdimLength = pixdimVector.size();
870  for (int i=0; i<std::min(pixdimLength,nDims); i++)
871  this->image->pixdim[i+1] = std::abs(pixdimVector[i]);
872 
873  const std::vector<std::string> pixunitsVector = mriImage.field("voxelDimUnits");
874  setPixunits(pixunitsVector);
875 
876  if (xform.rows() != 4 || xform.cols() != 4)
877  this->image->qform_code = this->image->sform_code = 0;
878  else
879  {
880  mat44 matrix;
881  for (int i=0; i<4; i++)
882  {
883  for (int j=0; j<4; j++)
884  matrix.m[i][j] = static_cast<float>(xform(i,j));
885  }
886 
887  this->image->qto_xyz = matrix;
888  this->image->qto_ijk = nifti_mat44_inverse(image->qto_xyz);
889  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);
890 
891  this->image->sto_xyz = matrix;
892  this->image->sto_ijk = nifti_mat44_inverse(image->sto_xyz);
893 
894  this->image->qform_code = this->image->sform_code = 2;
895  }
896 }
897 
898 inline void NiftiImage::initFromList (const Rcpp::RObject &object)
899 {
900  Rcpp::List list(object);
901  nifti_1_header *header = nifti_make_new_header(NULL, DT_FLOAT64);
902 
903  internal::updateHeader(header, list);
904 
905  this->image = nifti_convert_nhdr2nim(*header, NULL);
906  this->image->data = NULL;
907  free(header);
908 }
909 
910 inline void NiftiImage::initFromArray (const Rcpp::RObject &object, const bool copyData)
911 {
912  int dims[8] = { 0, 0, 0, 0, 0, 0, 0, 0 };
913  const std::vector<int> dimVector = object.attr("dim");
914 
915  const int nDims = std::min(7, int(dimVector.size()));
916  dims[0] = nDims;
917  for (int i=0; i<nDims; i++)
918  dims[i+1] = dimVector[i];
919 
920  const short datatype = sexpTypeToNiftiType(object.sexp_type());
921  this->image = nifti_make_new_nim(dims, datatype, int(copyData));
922 
923  if (copyData)
924  {
925  const size_t dataSize = nifti_get_volsize(image);
926  if (datatype == DT_INT32)
927  memcpy(this->image->data, INTEGER(object), dataSize);
928  else
929  memcpy(this->image->data, REAL(object), dataSize);
930  }
931  else
932  this->image->data = NULL;
933 
934  if (object.hasAttribute("pixdim"))
935  {
936  const std::vector<float> pixdimVector = object.attr("pixdim");
937  const int pixdimLength = pixdimVector.size();
938  for (int i=0; i<std::min(pixdimLength,nDims); i++)
939  this->image->pixdim[i+1] = pixdimVector[i];
940  }
941 
942  if (object.hasAttribute("pixunits"))
943  {
944  const std::vector<std::string> pixunitsVector = object.attr("pixunits");
945  setPixunits(pixunitsVector);
946  }
947 }
948 
949 inline NiftiImage::NiftiImage (const SEXP object, const bool readData)
950  : persistent(false)
951 {
952  Rcpp::RObject imageObject(object);
953  bool resolved = false;
954 
955  if (imageObject.hasAttribute(".nifti_image_ptr"))
956  {
957  Rcpp::XPtr<NiftiImage> imagePtr(SEXP(imageObject.attr(".nifti_image_ptr")));
958  NiftiImage *ptr = imagePtr;
959  if (ptr != NULL)
960  {
961  this->image = ptr->image;
962  this->persistent = true;
963  resolved = true;
964 
965  if (imageObject.hasAttribute("dim"))
966  update(imageObject);
967  }
968  else if (Rf_isString(object))
969  throw std::runtime_error("Internal image is not valid");
970  else
971  Rf_warning("Ignoring invalid internal pointer");
972  }
973 
974  if (!resolved)
975  {
976  if (Rf_isNull(object))
977  this->image = NULL;
978  else if (Rf_isString(object))
979  {
980  const std::string path = Rcpp::as<std::string>(object);
981  this->image = nifti_image_read(path.c_str(), readData);
982  if (this->image == NULL)
983  throw std::runtime_error("Failed to read image from path " + path);
984  }
985  else if (imageObject.inherits("nifti"))
986  initFromNiftiS4(imageObject, readData);
987  else if (imageObject.inherits("anlz"))
988  throw std::runtime_error("Cannot currently convert objects of class \"anlz\"");
989  else if (imageObject.inherits("MriImage"))
990  initFromMriImage(imageObject, readData);
991  else if (Rf_isVectorList(object))
992  initFromList(imageObject);
993  else if (imageObject.hasAttribute("dim"))
994  initFromArray(imageObject, readData);
995  else
996  throw std::runtime_error("Cannot convert object of class \"" + Rcpp::as<std::string>(imageObject.attr("class")) + "\" to a nifti_image");
997  }
998 
999  if (this->image != NULL)
1000  nifti_update_dims_from_array(this->image);
1001 
1002 #ifndef NDEBUG
1003  Rprintf("Creating NiftiImage with pointer %p (from SEXP)\n", this->image);
1004 #endif
1005 }
1006 
1007 #endif // USING_R
1008 
1009 inline NiftiImage::NiftiImage (const std::string &path, const std::vector<int> &volumes)
1010  : persistent(false)
1011 {
1012  if (volumes.empty())
1013  throw std::runtime_error("The vector of volumes is empty");
1014 
1015  nifti_brick_list brickList;
1016  this->image = nifti_image_read_bricks(path.c_str(), volumes.size(), &volumes[0], &brickList);
1017  if (this->image == NULL)
1018  throw std::runtime_error("Failed to read image from path " + path);
1019 
1020  size_t brickSize = image->nbyper * image->nx * image->ny * image->nz;
1021  image->data = calloc(1, nifti_get_volsize(image));
1022  for (int i=0; i<brickList.nbricks; i++)
1023  memcpy((char *) image->data + i * brickSize, brickList.bricks[i], brickSize);
1024  nifti_free_NBL(&brickList);
1025 
1026 #ifndef NDEBUG
1027  Rprintf("Creating NiftiImage with pointer %p (from string and volume vector)\n", this->image);
1028 #endif
1029 }
1030 
1031 inline void NiftiImage::updatePixdim (const std::vector<float> &pixdim)
1032 {
1033  const int nDims = image->dim[0];
1034  const std::vector<float> origPixdim(image->pixdim+1, image->pixdim+4);
1035 
1036  for (int i=1; i<8; i++)
1037  image->pixdim[i] = 0.0;
1038 
1039  const int pixdimLength = pixdim.size();
1040  for (int i=0; i<std::min(pixdimLength,nDims); i++)
1041  image->pixdim[i+1] = pixdim[i];
1042 
1043  if (!std::equal(origPixdim.begin(), origPixdim.begin() + std::min(3,nDims), pixdim.begin()))
1044  {
1045  mat33 scaleMatrix;
1046  for (int i=0; i<3; i++)
1047  {
1048  for (int j=0; j<3; j++)
1049  {
1050  if (i != j)
1051  scaleMatrix.m[i][j] = 0.0;
1052  else if (i >= nDims)
1053  scaleMatrix.m[i][j] = 1.0;
1054  else
1055  scaleMatrix.m[i][j] = pixdim[i] / origPixdim[i];
1056  }
1057  }
1058 
1059  if (image->qform_code > 0)
1060  {
1061  mat33 prod = nifti_mat33_mul(scaleMatrix, internal::topLeftCorner(image->qto_xyz));
1062  for (int i=0; i<3; i++)
1063  {
1064  for (int j=0; j<3; j++)
1065  image->qto_xyz.m[i][j] = prod.m[i][j];
1066  }
1067  image->qto_ijk = nifti_mat44_inverse(image->qto_xyz);
1068  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);
1069  }
1070 
1071  if (image->sform_code > 0)
1072  {
1073  mat33 prod = nifti_mat33_mul(scaleMatrix, internal::topLeftCorner(image->sto_xyz));
1074  for (int i=0; i<3; i++)
1075  {
1076  for (int j=0; j<3; j++)
1077  image->sto_xyz.m[i][j] = prod.m[i][j];
1078  }
1079  image->sto_ijk = nifti_mat44_inverse(image->sto_xyz);
1080  }
1081  }
1082 }
1083 
1084 inline void NiftiImage::setPixunits (const std::vector<std::string> &pixunits)
1085 {
1086  for (size_t i=0; i<pixunits.size(); i++)
1087  {
1088  if (pixunits[i] == "m")
1089  image->xyz_units = NIFTI_UNITS_METER;
1090  else if (pixunits[i] == "mm")
1091  image->xyz_units = NIFTI_UNITS_MM;
1092  else if (pixunits[i] == "um")
1093  image->xyz_units = NIFTI_UNITS_MICRON;
1094  else if (pixunits[i] == "s")
1095  image->time_units = NIFTI_UNITS_SEC;
1096  else if (pixunits[i] == "ms")
1097  image->time_units = NIFTI_UNITS_MSEC;
1098  else if (pixunits[i] == "us")
1099  image->time_units = NIFTI_UNITS_USEC;
1100  else if (pixunits[i] == "Hz")
1101  image->time_units = NIFTI_UNITS_HZ;
1102  else if (pixunits[i] == "ppm")
1103  image->time_units = NIFTI_UNITS_PPM;
1104  else if (pixunits[i] == "rad/s")
1105  image->time_units = NIFTI_UNITS_RADS;
1106  }
1107 }
1108 
1109 inline NiftiImage & NiftiImage::rescale (const std::vector<float> &scales)
1110 {
1111  std::vector<float> pixdim(image->pixdim+1, image->pixdim+4);
1112 
1113  for (int i=0; i<std::min(3, int(scales.size())); i++)
1114  {
1115  if (scales[i] != 1.0)
1116  {
1117  pixdim[i] /= scales[i];
1118  image->dim[i+1] = static_cast<int>(std::floor(image->dim[i+1] * scales[i]));
1119  }
1120  }
1121 
1123  nifti_update_dims_from_array(image);
1124 
1125  // Data vector is now the wrong size, so drop it
1126  if (!persistent)
1127  nifti_image_unload(image);
1128 
1129  image->scl_slope = 0.0;
1130  image->scl_inter = 0.0;
1131 
1132  return *this;
1133 }
1134 
1135 inline NiftiImage & NiftiImage::reorient (const int icode, const int jcode, const int kcode)
1136 {
1137  if (this->isNull())
1138  return *this;
1139  if (image->qform_code == 0 && image->sform_code == 0)
1140  {
1141  Rf_warning("Image qform and sform codes are both zero, so it cannot be reoriented");
1142  return *this;
1143  }
1144 
1145  int used[6] = { 0, 0, 0, 0, 0, 0 };
1146  used[icode-1] = 1;
1147  used[jcode-1] = 1;
1148  used[kcode-1] = 1;
1149  if (used[0]+used[1] != 1 || used[2]+used[3] != 1 || used[4]+used[5] != 1)
1150  throw std::runtime_error("Each canonical axis should be used exactly once");
1151 
1152  const int codes[3] = { icode, jcode, kcode };
1153  const mat44 native = this->xform();
1154 
1155  // Create a target xform (rotation matrix only)
1156  mat33 target;
1157  for (int j=0; j<3; j++)
1158  {
1159  for (int i=0; i<3; i++)
1160  target.m[i][j] = 0.0;
1161 
1162  switch (codes[j])
1163  {
1164  case NIFTI_L2R: target.m[0][j] = 1.0; break;
1165  case NIFTI_R2L: target.m[0][j] = -1.0; break;
1166  case NIFTI_P2A: target.m[1][j] = 1.0; break;
1167  case NIFTI_A2P: target.m[1][j] = -1.0; break;
1168  case NIFTI_I2S: target.m[2][j] = 1.0; break;
1169  case NIFTI_S2I: target.m[2][j] = -1.0; break;
1170  }
1171  }
1172 
1173  // Extract (inverse of) canonical axis matrix from native xform
1174  int nicode, njcode, nkcode;
1175  nifti_mat44_to_orientation(native, &nicode, &njcode, &nkcode);
1176  int ncodes[3] = { nicode, njcode, nkcode };
1177  mat33 nativeAxesTransposed;
1178  for (int i=0; i<3; i++)
1179  {
1180  for (int j=0; j<3; j++)
1181  nativeAxesTransposed.m[i][j] = 0.0;
1182 
1183  switch (ncodes[i])
1184  {
1185  case NIFTI_L2R: nativeAxesTransposed.m[i][0] = 1.0; break;
1186  case NIFTI_R2L: nativeAxesTransposed.m[i][0] = -1.0; break;
1187  case NIFTI_P2A: nativeAxesTransposed.m[i][1] = 1.0; break;
1188  case NIFTI_A2P: nativeAxesTransposed.m[i][1] = -1.0; break;
1189  case NIFTI_I2S: nativeAxesTransposed.m[i][2] = 1.0; break;
1190  case NIFTI_S2I: nativeAxesTransposed.m[i][2] = -1.0; break;
1191  }
1192  }
1193 
1194  // Check for no-op case
1195  if (icode == nicode && jcode == njcode && kcode == nkcode)
1196  return *this;
1197 
1198  // The transform is t(approx_old_xform) %*% target_xform
1199  // The new xform is old_xform %*% transform
1200  // NB: "transform" is really 4x4, but the last row and column are filled implicitly during the multiplication loop
1201  mat33 transform = nifti_mat33_mul(nativeAxesTransposed, target);
1202  mat44 result;
1203  for (int i=0; i<4; i++)
1204  {
1205  for (int j=0; j<3; j++)
1206  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];
1207 
1208  result.m[i][3] = native.m[i][3];
1209  }
1210 
1211  // Update the xforms with nonzero codes
1212  if (image->qform_code > 0)
1213  {
1214  image->qto_xyz = result;
1215  image->qto_ijk = nifti_mat44_inverse(image->qto_xyz);
1216  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);
1217  }
1218  if (image->sform_code > 0)
1219  {
1220  image->sto_xyz = result;
1221  image->sto_ijk = nifti_mat44_inverse(image->sto_xyz);
1222  }
1223 
1224  // Extract the mapping between dimensions and the signs
1225  int locs[3], signs[3], newdim[3];
1226  float newpixdim[3];
1227  double maxes[3] = { R_NegInf, R_NegInf, R_NegInf };
1228  for (int j=0; j<3; j++)
1229  {
1230  for (int i=0; i<3; i++)
1231  {
1232  const double value = static_cast<double>(transform.m[i][j]);
1233  if (fabs(value) > maxes[j])
1234  {
1235  maxes[j] = fabs(value);
1236  signs[j] = value > 0.0 ? 1 : -1;
1237  locs[j] = i;
1238  }
1239  }
1240 
1241  // Permute dim and pixdim
1242  newdim[j] = image->dim[locs[j]+1];
1243  newpixdim[j] = image->pixdim[locs[j]+1];
1244  }
1245 
1246  // Calculate strides in target space
1247  ptrdiff_t strides[3];
1248  strides[locs[0]] = 1;
1249  for (int n=1; n<3; n++)
1250  strides[locs[n]] = strides[locs[n-1]] * image->dim[locs[n-1]+1];
1251 
1252  if (image->data != NULL)
1253  {
1254  size_t volSize = size_t(image->nx * image->ny * image->nz);
1255  size_t nVolumes = std::max(size_t(1), image->nvox / volSize);
1256 
1257  const std::vector<double> oldData = this->getData<double>();
1258  std::vector<double> newData(image->nvox);
1259 
1260  // Where the sign is negative we need to start at the end of the dimension
1261  size_t volStart = 0;
1262  for (int i=0; i<3; i++)
1263  {
1264  if (signs[i] < 0)
1265  volStart += (image->dim[i+1] - 1) * strides[i];
1266  }
1267 
1268  // Iterate over the data and place it into a new vector
1269  std::vector<double>::const_iterator it = oldData.begin();
1270  for (size_t v=0; v<nVolumes; v++)
1271  {
1272  for (int k=0; k<image->nz; k++)
1273  {
1274  ptrdiff_t offset = k * strides[2] * signs[2];
1275  for (int j=0; j<image->ny; j++)
1276  {
1277  for (int i=0; i<image->nx; i++)
1278  {
1279  newData[volStart + offset] = *it++;
1280  offset += strides[0] * signs[0];
1281  }
1282  offset += strides[1] * signs[1] - image->nx * strides[0] * signs[0];
1283  }
1284  }
1285  volStart += volSize;
1286  }
1287 
1288  // Replace the existing data in the image
1289  this->replaceData(newData);
1290  }
1291 
1292  // Copy new dims and pixdims in
1293  // NB: Old dims are used above, so this must happen last
1294  std::copy(newdim, newdim+3, image->dim+1);
1295  std::copy(newpixdim, newpixdim+3, image->pixdim+1);
1296  nifti_update_dims_from_array(image);
1297 
1298  return *this;
1299 }
1300 
1301 inline NiftiImage & NiftiImage::reorient (const std::string &orientation)
1302 {
1303  if (orientation.length() != 3)
1304  throw std::runtime_error("Orientation string should have exactly three characters");
1305 
1306  int codes[3];
1307  for (int i=0; i<3; i++)
1308  {
1309  switch(orientation[i])
1310  {
1311  case 'r': case 'R': codes[i] = NIFTI_L2R; break;
1312  case 'l': case 'L': codes[i] = NIFTI_R2L; break;
1313  case 'a': case 'A': codes[i] = NIFTI_P2A; break;
1314  case 'p': case 'P': codes[i] = NIFTI_A2P; break;
1315  case 's': case 'S': codes[i] = NIFTI_I2S; break;
1316  case 'i': case 'I': codes[i] = NIFTI_S2I; break;
1317 
1318  default:
1319  throw std::runtime_error("Orientation string is invalid");
1320  }
1321  }
1322 
1323  return reorient(codes[0], codes[1], codes[2]);
1324 }
1325 
1326 #ifdef USING_R
1327 
1328 inline NiftiImage & NiftiImage::update (const Rcpp::RObject &object)
1329 {
1330  if (Rf_isVectorList(object))
1331  {
1332  Rcpp::List list(object);
1333  nifti_1_header *header = NULL;
1334  if (this->isNull())
1335  {
1336  header = nifti_make_new_header(NULL, DT_FLOAT64);
1337  internal::updateHeader(header, list, true);
1338  }
1339  else
1340  {
1341  header = (nifti_1_header *) calloc(1, sizeof(nifti_1_header));
1342  *header = nifti_convert_nim2nhdr(image);
1343  internal::updateHeader(header, list, true);
1344  }
1345 
1346  if (header != NULL)
1347  {
1348  nifti_image *newImage = nifti_convert_nhdr2nim(*header, NULL);
1349  if (this->image->data != NULL)
1350  {
1351  size_t dataSize = nifti_get_volsize(image);
1352  newImage->data = calloc(1, dataSize);
1353  memcpy(newImage->data, image->data, dataSize);
1354  }
1355  if (!persistent)
1356  nifti_image_free(image);
1357  this->image = newImage;
1358  }
1359  }
1360  else if (object.hasAttribute("dim"))
1361  {
1362  for (int i=0; i<8; i++)
1363  image->dim[i] = 0;
1364  const std::vector<int> dimVector = object.attr("dim");
1365 
1366  const int nDims = std::min(7, int(dimVector.size()));
1367  image->dim[0] = nDims;
1368  for (int i=0; i<nDims; i++)
1369  image->dim[i+1] = dimVector[i];
1370 
1371  if (object.hasAttribute("pixdim"))
1372  {
1373  const std::vector<float> pixdimVector = object.attr("pixdim");
1374  updatePixdim(pixdimVector);
1375  }
1376 
1377  if (object.hasAttribute("pixunits"))
1378  {
1379  const std::vector<std::string> pixunitsVector = object.attr("pixunits");
1380  setPixunits(pixunitsVector);
1381  }
1382 
1383  // This NIfTI-1 library function clobbers dim[0] if the last dimension is unitary; we undo that here
1384  nifti_update_dims_from_array(image);
1385  image->dim[0] = image->ndim = nDims;
1386 
1387  image->datatype = NiftiImage::sexpTypeToNiftiType(object.sexp_type());
1388  nifti_datatype_sizes(image->datatype, &image->nbyper, NULL);
1389 
1390  if (!persistent)
1391  nifti_image_unload(image);
1392 
1393  const size_t dataSize = nifti_get_volsize(image);
1394  image->data = calloc(1, dataSize);
1395  if (image->datatype == DT_INT32)
1396  memcpy(image->data, INTEGER(object), dataSize);
1397  else
1398  memcpy(image->data, REAL(object), dataSize);
1399 
1400  image->scl_slope = 0.0;
1401  image->scl_inter = 0.0;
1402  }
1403 
1404  return *this;
1405 }
1406 
1407 #endif // USING_R
1408 
1409 inline mat44 NiftiImage::xform (const bool preferQuaternion) const
1410 {
1411  if (image == NULL)
1412  {
1413  mat44 matrix;
1414  for (int i=0; i<4; i++)
1415  {
1416  for (int j=0; j<4; j++)
1417  matrix.m[i][j] = 0.0;
1418  }
1419  return matrix;
1420  }
1421  else if (image->qform_code <= 0 && image->sform_code <= 0)
1422  {
1423  // No qform or sform so use pixdim (NB: other software may assume differently)
1424  mat44 matrix;
1425  for (int i=0; i<4; i++)
1426  {
1427  for (int j=0; j<4; j++)
1428  {
1429  if (i != j)
1430  matrix.m[i][j] = 0.0;
1431  else if (i == 3)
1432  matrix.m[3][3] = 1.0;
1433  else
1434  matrix.m[i][j] = (image->pixdim[i+1]==0.0 ? 1.0 : image->pixdim[i+1]);
1435  }
1436  }
1437  return matrix;
1438  }
1439  else if ((preferQuaternion && image->qform_code > 0) || image->sform_code <= 0)
1440  return image->qto_xyz;
1441  else
1442  return image->sto_xyz;
1443 }
1444 
1445 template <typename TargetType>
1446 inline std::vector<TargetType> NiftiImage::Block::getData () const
1447 {
1448  if (image.isNull())
1449  return std::vector<TargetType>();
1450 
1451  size_t blockSize = 1;
1452  for (int i=1; i<dimension; i++)
1453  blockSize *= image->dim[i];
1454 
1455  std::vector<TargetType> data(blockSize);
1456  internal::convertData<TargetType>(image->data, image->datatype, blockSize, data.begin(), blockSize*index);
1457 
1458  if (image.isDataScaled())
1459  std::transform(data.begin(), data.end(), data.begin(), internal::DataRescaler<TargetType>(image->scl_slope,image->scl_inter));
1460 
1461  return data;
1462 }
1463 
1464 template <typename TargetType>
1465 inline std::vector<TargetType> NiftiImage::getData () const
1466 {
1467  if (this->isNull())
1468  return std::vector<TargetType>();
1469 
1470  std::vector<TargetType> data(image->nvox);
1471  internal::convertData<TargetType>(image->data, image->datatype, image->nvox, data.begin());
1472 
1473  if (this->isDataScaled())
1474  std::transform(data.begin(), data.end(), data.begin(), internal::DataRescaler<TargetType>(image->scl_slope,image->scl_inter));
1475 
1476  return data;
1477 }
1478 
1479 inline NiftiImage & NiftiImage::changeDatatype (const short datatype)
1480 {
1481  if (this->isNull() || image->datatype == datatype)
1482  return *this;
1483 
1484  if (image->data != NULL)
1485  {
1486  int bytesPerPixel;
1487  nifti_datatype_sizes(datatype, &bytesPerPixel, NULL);
1488  void *data = calloc(image->nvox, bytesPerPixel);
1489 
1490  switch (datatype)
1491  {
1492  case DT_UINT8:
1493  internal::convertData<uint8_t>(image->data, image->datatype, image->nvox, static_cast<uint8_t *>(data));
1494  break;
1495 
1496  case DT_INT16:
1497  internal::convertData<int16_t>(image->data, image->datatype, image->nvox, static_cast<int16_t *>(data));
1498  break;
1499 
1500  case DT_INT32:
1501  internal::convertData<int32_t>(image->data, image->datatype, image->nvox, static_cast<int32_t *>(data));
1502  break;
1503 
1504  case DT_FLOAT32:
1505  internal::convertData<float>(image->data, image->datatype, image->nvox, static_cast<float *>(data));
1506  break;
1507 
1508  case DT_FLOAT64:
1509  internal::convertData<double>(image->data, image->datatype, image->nvox, static_cast<double *>(data));
1510  break;
1511 
1512  case DT_INT8:
1513  internal::convertData<int8_t>(image->data, image->datatype, image->nvox, static_cast<int8_t *>(data));
1514  break;
1515 
1516  case DT_UINT16:
1517  internal::convertData<uint16_t>(image->data, image->datatype, image->nvox, static_cast<uint16_t *>(data));
1518  break;
1519 
1520  case DT_UINT32:
1521  internal::convertData<uint32_t>(image->data, image->datatype, image->nvox, static_cast<uint32_t *>(data));
1522  break;
1523 
1524  case DT_INT64:
1525  internal::convertData<int64_t>(image->data, image->datatype, image->nvox, static_cast<int64_t *>(data));
1526  break;
1527 
1528  case DT_UINT64:
1529  internal::convertData<uint64_t>(image->data, image->datatype, image->nvox, static_cast<uint64_t *>(data));
1530  break;
1531 
1532  default:
1533  throw std::runtime_error("Unsupported data type (" + std::string(nifti_datatype_string(datatype)) + ")");
1534  }
1535 
1536  nifti_image_unload(image);
1537  image->data = data;
1538  }
1539 
1540  image->datatype = datatype;
1541  nifti_datatype_sizes(datatype, &image->nbyper, &image->swapsize);
1542 
1543  return *this;
1544 }
1545 
1546 inline NiftiImage & NiftiImage::changeDatatype (const std::string &datatype)
1547 {
1548  return changeDatatype(internal::stringToDatatype(datatype));
1549 }
1550 
1551 template <typename SourceType>
1552 inline NiftiImage & NiftiImage::replaceData (const std::vector<SourceType> &data, const short datatype)
1553 {
1554  if (this->isNull())
1555  return *this;
1556  else if (data.size() != image->nvox)
1557  throw std::runtime_error("New data length does not match the number of voxels in the image");
1558 
1559  if (datatype != DT_NONE)
1560  {
1561  nifti_image_unload(image);
1562  image->datatype = datatype;
1563  nifti_datatype_sizes(datatype, &image->nbyper, &image->swapsize);
1564  }
1565 
1566  if (image->data == NULL)
1567  image->data = calloc(image->nvox, image->nbyper);
1568  internal::replaceData<SourceType>(data.begin(), data.end(), image->data, image->datatype);
1569 
1570  image->scl_slope = 0.0;
1571  image->scl_inter = 0.0;
1572  image->cal_min = static_cast<float>(*std::min_element(data.begin(), data.end()));
1573  image->cal_max = static_cast<float>(*std::max_element(data.begin(), data.end()));
1574 
1575  return *this;
1576 }
1577 
1578 inline void NiftiImage::toFile (const std::string fileName, const short datatype) const
1579 {
1580  // Copy the source image only if the datatype will be changed
1581  NiftiImage imageToWrite(image, datatype != DT_NONE);
1582 
1583  if (datatype == DT_NONE)
1584  imageToWrite.setPersistence(true);
1585  else
1586  imageToWrite.changeDatatype(datatype);
1587 
1588  const int status = nifti_set_filenames(imageToWrite, fileName.c_str(), false, true);
1589  if (status != 0)
1590  throw std::runtime_error("Failed to set filenames for NIfTI object");
1591  nifti_image_write(imageToWrite);
1592 }
1593 
1594 inline void NiftiImage::toFile (const std::string fileName, const std::string &datatype) const
1595 {
1596  toFile(fileName, internal::stringToDatatype(datatype));
1597 }
1598 
1599 #ifdef USING_R
1600 
1601 inline Rcpp::RObject NiftiImage::toArray () const
1602 {
1603  Rcpp::RObject array;
1604 
1605  if (this->isNull())
1606  return array;
1607  else if (this->isDataScaled())
1608  {
1609  array = internal::imageDataToArray<REALSXP>(image);
1610  std::transform(REAL(array), REAL(array)+Rf_length(array), REAL(array), internal::DataRescaler<double>(image->scl_slope,image->scl_inter));
1611  }
1612  else
1613  {
1614  switch (image->datatype)
1615  {
1616  case DT_UINT8:
1617  case DT_INT16:
1618  case DT_INT32:
1619  case DT_INT8:
1620  case DT_UINT16:
1621  case DT_UINT32:
1622  case DT_INT64:
1623  case DT_UINT64:
1624  array = internal::imageDataToArray<INTSXP>(image);
1625  break;
1626 
1627  case DT_FLOAT32:
1628  case DT_FLOAT64:
1629  array = internal::imageDataToArray<REALSXP>(image);
1630  break;
1631 
1632  default:
1633  throw std::runtime_error("Unsupported data type (" + std::string(nifti_datatype_string(image->datatype)) + ")");
1634  }
1635  }
1636 
1637  internal::addAttributes(array, image);
1638  array.attr("class") = Rcpp::CharacterVector::create("niftiImage", "array");
1639  return array;
1640 }
1641 
1642 inline Rcpp::RObject NiftiImage::toPointer (const std::string label) const
1643 {
1644  if (this->isNull())
1645  return Rcpp::RObject();
1646  else
1647  {
1648  Rcpp::RObject string = Rcpp::wrap(label);
1649  internal::addAttributes(string, image, false);
1650  string.attr("class") = Rcpp::CharacterVector::create("internalImage", "niftiImage");
1651  return string;
1652  }
1653 }
1654 
1655 inline Rcpp::RObject NiftiImage::toArrayOrPointer (const bool internal, const std::string label) const
1656 {
1657  return (internal ? toPointer(label) : toArray());
1658 }
1659 
1660 inline Rcpp::RObject NiftiImage::headerToList () const
1661 {
1662  if (this->image == NULL)
1663  return Rcpp::RObject();
1664 
1665  nifti_1_header header = nifti_convert_nim2nhdr(this->image);
1666  Rcpp::List result;
1667 
1668  result["sizeof_hdr"] = header.sizeof_hdr;
1669 
1670  result["dim_info"] = int(header.dim_info);
1671  result["dim"] = std::vector<short>(header.dim, header.dim+8);
1672 
1673  result["intent_p1"] = header.intent_p1;
1674  result["intent_p2"] = header.intent_p2;
1675  result["intent_p3"] = header.intent_p3;
1676  result["intent_code"] = header.intent_code;
1677 
1678  result["datatype"] = header.datatype;
1679  result["bitpix"] = header.bitpix;
1680 
1681  result["slice_start"] = header.slice_start;
1682  result["pixdim"] = std::vector<float>(header.pixdim, header.pixdim+8);
1683  result["vox_offset"] = header.vox_offset;
1684  result["scl_slope"] = header.scl_slope;
1685  result["scl_inter"] = header.scl_inter;
1686  result["slice_end"] = header.slice_end;
1687  result["slice_code"] = int(header.slice_code);
1688  result["xyzt_units"] = int(header.xyzt_units);
1689  result["cal_max"] = header.cal_max;
1690  result["cal_min"] = header.cal_min;
1691  result["slice_duration"] = header.slice_duration;
1692  result["toffset"] = header.toffset;
1693  result["descrip"] = std::string(header.descrip, 80);
1694  result["aux_file"] = std::string(header.aux_file, 24);
1695 
1696  result["qform_code"] = header.qform_code;
1697  result["sform_code"] = header.sform_code;
1698  result["quatern_b"] = header.quatern_b;
1699  result["quatern_c"] = header.quatern_c;
1700  result["quatern_d"] = header.quatern_d;
1701  result["qoffset_x"] = header.qoffset_x;
1702  result["qoffset_y"] = header.qoffset_y;
1703  result["qoffset_z"] = header.qoffset_z;
1704  result["srow_x"] = std::vector<float>(header.srow_x, header.srow_x+4);
1705  result["srow_y"] = std::vector<float>(header.srow_y, header.srow_y+4);
1706  result["srow_z"] = std::vector<float>(header.srow_z, header.srow_z+4);
1707 
1708  result["intent_name"] = std::string(header.intent_name, 16);
1709  result["magic"] = std::string(header.magic, 4);
1710 
1711  result.attr("class") = Rcpp::CharacterVector::create("niftiHeader");
1712 
1713  return result;
1714 }
1715 
1716 #endif // USING_R
1717 
1718 } // namespace
1719 
1720 #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:1135
NiftiImage()
Default constructor.
Definition: NiftiImage.h:256
const int index
The location along dimension.
Definition: NiftiImage.h:56
static std::string xformToString(const mat44 matrix)
Convert a 4x4 xform matrix to a string describing its canonical axes.
Definition: NiftiImage.h:129
std::vector< int > dim() const
Return the dimensions of the image.
Definition: NiftiImage.h:442
const nifti_image * operator->() const
Allows a NiftiImage object to be treated as a pointer to a const nifti_image.
Definition: NiftiImage.h:364
bool persistent
Marker of persistence, which determines whether the nifti_image should be freed on destruction...
Definition: NiftiImage.h:188
std::vector< float > pixdim() const
Return the dimensions of the pixels or voxels in the image.
Definition: NiftiImage.h:454
NiftiImage(const std::string &path, const bool readData=true)
Initialise using a path string.
Definition: NiftiImage.h:309
NiftiImage & rescale(const std::vector< float > &scales)
Rescale the image, changing its image dimensions and pixel dimensions.
Definition: NiftiImage.h:1109
NiftiImage(nifti_image *const image, const bool copy=false)
Initialise using an existing nifti_image pointer.
Definition: NiftiImage.h:291
Definition: NiftiImage.h:39
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:80
bool isNull() const
Determine whether or not the internal pointer is NULL.
Definition: NiftiImage.h:415
Inner class referring to a subset of an image.
Definition: NiftiImage.h:52
Block slice(const int i)
Extract a slice block from a 3D image.
Definition: NiftiImage.h:596
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:1552
std::vector< TargetType > getData() const
Extract a vector of data from the image, casting it to any required element type. ...
Definition: NiftiImage.h:1465
bool isDataScaled() const
Determine whether nontrivial scale and slope parameters are set.
Definition: NiftiImage.h:425
const Block slice(const int i) const
Extract a slice block from a 3D image.
Definition: NiftiImage.h:589
const Block volume(const int i) const
Extract a volume block from a 4D image.
Definition: NiftiImage.h:603
bool isPersistent() const
Determine whether or not the image is marked as persistent.
Definition: NiftiImage.h:420
void setPixunits(const std::vector< std::string > &pixunits)
Modify the pixel dimension units.
Definition: NiftiImage.h:1084
NiftiImage(const Block &source)
Initialise from a block, copying in the data.
Definition: NiftiImage.h:276
virtual ~NiftiImage()
Destructor which frees the wrapped pointer, unless the object is marked as persistent.
Definition: NiftiImage.h:340
Thin wrapper around a C-style nifti_image struct that allows C++-style destruction.
Definition: NiftiImage.h:45
NiftiImage(const NiftiImage &source)
Copy constructor.
Definition: NiftiImage.h:263
int nBlocks() const
Return the number of blocks in the image.
Definition: NiftiImage.h:558
const int dimension
The dimension along which the block applies (which should be the last)
Definition: NiftiImage.h:55
const NiftiImage & image
The parent image.
Definition: NiftiImage.h:54
NiftiImage & operator=(const NiftiImage &source)
Copy assignment operator, which copies from its argument.
Definition: NiftiImage.h:375
static int fileVersion(const std::string &path)
Get the NIfTI format version used by the file at the specified path.
Definition: NiftiImage.h:157
void toFile(const std::string fileName, const short datatype=DT_NONE) const
Write the image to a NIfTI-1 file.
Definition: NiftiImage.h:1578
void updatePixdim(const std::vector< float > &pixdim)
Modify the pixel dimensions, and potentially the xform matrices to match.
Definition: NiftiImage.h:1031
const Block block(const int i) const
Extract a block from the image.
Definition: NiftiImage.h:573
NiftiImage & setPersistence(const bool persistent)
Marked the image as persistent, so that it can be passed back to R.
Definition: NiftiImage.h:402
void copy(const nifti_image *source)
Copy the contents of a nifti_image to create a new image.
Definition: NiftiImage.h:662
mat44 xform(const bool preferQuaternion=true) const
Obtain an xform matrix, indicating the orientation of the image.
Definition: NiftiImage.h:1409
nifti_image * image
The wrapped nifti_image pointer.
Definition: NiftiImage.h:187
NiftiImage & changeDatatype(const short datatype)
Change the datatype of the image, casting the pixel data if present.
Definition: NiftiImage.h:1479
std::vector< TargetType > getData() const
Extract a vector of data from a block, casting it to any required element type.
Definition: NiftiImage.h:1446
int nDims() const
Return the number of dimensions in the image.
Definition: NiftiImage.h:430
Block volume(const int i)
Extract a volume block from a 4D image.
Definition: NiftiImage.h:610
NiftiImage & drop()
Drop unitary dimensions.
Definition: NiftiImage.h:468
Block(const NiftiImage &image, const int dimension, const int index)
Standard constructor for this class.
Definition: NiftiImage.h:65
NiftiImage & dropData()
Drop the data from the image, retaining only the metadata.
Definition: NiftiImage.h:509
Block block(const int i)
Extract a block from the image.
Definition: NiftiImage.h:582