RNifti
Fast R and C++ Access to NIfTI Images
NiftiImage.h
1 #ifndef _NIFTI_IMAGE_H_
2 #define _NIFTI_IMAGE_H_
3 
4 #include <Rcpp.h>
5 
6 #include "niftilib/nifti1_io.h"
7 
21 {
22 public:
27  struct Block
28  {
29  const NiftiImage &image;
30  const int dimension;
31  const int index;
40  Block (const NiftiImage &image, const int dimension, const int index)
41  : image(image), dimension(dimension), index(index)
42  {
43  if (dimension != image->ndim)
44  throw std::runtime_error("Blocks must be along the last dimension in the image");
45  }
46 
55  Block & operator= (const NiftiImage &source)
56  {
57  if (source->datatype != image->datatype)
58  throw std::runtime_error("New data does not have the same datatype as the target block");
59 
60  size_t blockSize = 1;
61  for (int i=1; i<dimension; i++)
62  blockSize *= image->dim[i];
63 
64  if (blockSize != source->nvox)
65  throw std::runtime_error("New data does not have the same size as the target block");
66 
67  blockSize *= image->nbyper;
68  memcpy(static_cast<char*>(image->data) + blockSize*index, source->data, blockSize);
69  return *this;
70  }
71  };
72 
79  static short sexpTypeToNiftiType (const int sexpType)
80  {
81  if (sexpType == INTSXP || sexpType == LGLSXP)
82  return DT_INT32;
83  else if (sexpType == REALSXP)
84  return DT_FLOAT64;
85  else
86  throw std::runtime_error("Array elements must be numeric");
87  }
88 
89 protected:
90  nifti_image *image;
91  bool persistent;
97  void copy (nifti_image * const source);
98 
103  void copy (const NiftiImage &source);
104 
109  void copy (const Block &source);
110 
116  void initFromNiftiS4 (const Rcpp::RObject &object, const bool copyData = true);
117 
123  void initFromMriImage (const Rcpp::RObject &object, const bool copyData = true);
124 
129  void initFromList (const Rcpp::RObject &object);
130 
136  void initFromArray (const Rcpp::RObject &object, const bool copyData = true);
137 
142  void updatePixdim (const std::vector<float> &pixdim);
143 
148  void setPixunits (const std::vector<std::string> &pixunits);
149 
150 public:
155  : image(NULL), persistent(false) {}
156 
161  NiftiImage (const NiftiImage &source)
162  : image(NULL), persistent(false)
163  {
164  this->copy(source);
165 #ifndef NDEBUG
166  Rprintf("Creating NiftiImage with pointer %p (from NiftiImage)\n", this->image);
167 #endif
168  }
169 
176  NiftiImage (nifti_image * const image, const bool copy = false)
177  : image(NULL), persistent(false)
178  {
179  if (copy)
180  this->copy(image);
181  else
182  this->image = image;
183 #ifndef NDEBUG
184  Rprintf("Creating NiftiImage with pointer %p (from pointer)\n", this->image);
185 #endif
186  }
187 
194  NiftiImage (const std::string &path, const bool readData = true)
195  : persistent(false)
196  {
197  this->image = nifti_image_read(path.c_str(), readData);
198  if (this->image == NULL)
199  throw std::runtime_error("Failed to read image from path " + path);
200 #ifndef NDEBUG
201  Rprintf("Creating NiftiImage with pointer %p (from string)\n", this->image);
202 #endif
203  }
204 
210  NiftiImage (const SEXP object, const bool readData = true);
211 
216  {
217  if (!persistent)
218  {
219 #ifndef NDEBUG
220  Rprintf("Freeing NiftiImage with pointer %p\n", this->image);
221 #endif
222  nifti_image_free(image);
223  }
224  }
225 
229  operator nifti_image* () const { return image; }
230 
234  nifti_image * operator-> () const { return image; }
235 
241  {
242  copy(source);
243 #ifndef NDEBUG
244  Rprintf("Creating NiftiImage with pointer %p (from NiftiImage)\n", this->image);
245 #endif
246  return *this;
247  }
248 
254  NiftiImage & operator= (const Block &source)
255  {
256  copy(source);
257 #ifndef NDEBUG
258  Rprintf("Creating NiftiImage with pointer %p (from Block)\n", this->image);
259 #endif
260  return *this;
261  }
262 
267  void setPersistence (const bool persistent)
268  {
269  this->persistent = persistent;
270 #ifndef NDEBUG
271  if (persistent)
272  Rprintf("Setting NiftiImage with pointer %p to be persistent\n", this->image);
273 #endif
274  }
275 
279  bool isNull () const { return (image == NULL); }
280 
284  bool isPersistent () const { return persistent; }
285 
289  int nDims () const
290  {
291  if (image == NULL)
292  return 0;
293  else
294  return image->ndim;
295  }
296 
304  {
305  int ndim = image->ndim;
306  while (image->dim[ndim] < 2)
307  ndim--;
308  image->dim[0] = image->ndim = ndim;
309 
310  return *this;
311  }
312 
318  void rescale (const std::vector<float> &scales);
319 
324  void update (const SEXP array);
325 
331  mat44 xform (const bool preferQuaternion = true) const;
332 
338  const Block slice (const int i) const { return Block(*this, 3, i); }
339 
345  Block slice (const int i) { return Block(*this, 3, i); }
346 
352  const Block volume (const int i) const { return Block(*this, 4, i); }
353 
359  Block volume (const int i) { return Block(*this, 4, i); }
360 
366  void toFile (const std::string fileName, const short datatype) const;
367 
373  void toFile (const std::string fileName, const std::string &datatype) const;
374 
379  Rcpp::RObject toArray () const;
380 
386  Rcpp::RObject toPointer (const std::string label) const;
387 
394  Rcpp::RObject toArrayOrPointer (const bool internal, const std::string label) const;
395 
400  Rcpp::RObject headerToList () const;
401 };
402 
403 inline void NiftiImage::copy (nifti_image * const source)
404 {
405  if (image != NULL)
406  nifti_image_free(image);
407 
408  if (source == NULL)
409  image = NULL;
410  else
411  {
412  image = nifti_copy_nim_info(source);
413  if (source->data != NULL)
414  {
415  size_t dataSize = nifti_get_volsize(source);
416  image->data = calloc(1, dataSize);
417  memcpy(image->data, source->data, dataSize);
418  }
419  }
420 }
421 
422 inline void NiftiImage::copy (const NiftiImage &source)
423 {
424  nifti_image *sourceStruct = source;
425  copy(sourceStruct);
426 }
427 
428 inline void NiftiImage::copy (const Block &source)
429 {
430  if (image != NULL)
431  nifti_image_free(image);
432 
433  nifti_image *sourceStruct = source.image;
434  if (sourceStruct == NULL)
435  image = NULL;
436  else
437  {
438  image = nifti_copy_nim_info(sourceStruct);
439  image->dim[0] = source.image->dim[0] - 1;
440  image->dim[source.dimension] = 1;
441  image->pixdim[source.dimension] = 1.0;
442  nifti_update_dims_from_array(image);
443 
444  if (sourceStruct->data != NULL)
445  {
446  size_t blockSize = nifti_get_volsize(image);
447  image->data = calloc(1, blockSize);
448  memcpy(image->data, static_cast<char*>(source.image->data) + blockSize*source.index, blockSize);
449  }
450  }
451 }
452 
453 // Convert an S4 "nifti" object, as defined in the oro.nifti package, to a "nifti_image" struct
454 inline void NiftiImage::initFromNiftiS4 (const Rcpp::RObject &object, const bool copyData)
455 {
456  nifti_1_header header;
457  header.sizeof_hdr = 348;
458 
459  const std::vector<short> dims = object.slot("dim_");
460  for (int i=0; i<8; i++)
461  header.dim[i] = dims[i];
462 
463  header.intent_p1 = object.slot("intent_p1");
464  header.intent_p2 = object.slot("intent_p2");
465  header.intent_p3 = object.slot("intent_p3");
466  header.intent_code = object.slot("intent_code");
467 
468  header.datatype = object.slot("datatype");
469  header.bitpix = object.slot("bitpix");
470 
471  header.slice_start = object.slot("slice_start");
472  header.slice_end = object.slot("slice_end");
473  header.slice_code = Rcpp::as<int>(object.slot("slice_code"));
474  header.slice_duration = object.slot("slice_duration");
475 
476  const std::vector<float> pixdims = object.slot("pixdim");
477  for (int i=0; i<8; i++)
478  header.pixdim[i] = pixdims[i];
479  header.xyzt_units = Rcpp::as<int>(object.slot("xyzt_units"));
480 
481  header.vox_offset = object.slot("vox_offset");
482 
483  header.scl_slope = object.slot("scl_slope");
484  header.scl_inter = object.slot("scl_inter");
485  header.toffset = object.slot("toffset");
486 
487  header.cal_max = object.slot("cal_max");
488  header.cal_min = object.slot("cal_min");
489  header.glmax = header.glmin = 0;
490 
491  strncpy(header.descrip, Rcpp::as<std::string>(object.slot("descrip")).c_str(), 79);
492  header.descrip[79] = '\0';
493  strncpy(header.aux_file, Rcpp::as<std::string>(object.slot("aux_file")).c_str(), 23);
494  header.aux_file[23] = '\0';
495  strncpy(header.intent_name, Rcpp::as<std::string>(object.slot("intent_name")).c_str(), 15);
496  header.intent_name[15] = '\0';
497  strncpy(header.magic, Rcpp::as<std::string>(object.slot("magic")).c_str(), 3);
498  header.magic[3] = '\0';
499 
500  header.qform_code = object.slot("qform_code");
501  header.sform_code = object.slot("sform_code");
502 
503  header.quatern_b = object.slot("quatern_b");
504  header.quatern_c = object.slot("quatern_c");
505  header.quatern_d = object.slot("quatern_d");
506  header.qoffset_x = object.slot("qoffset_x");
507  header.qoffset_y = object.slot("qoffset_y");
508  header.qoffset_z = object.slot("qoffset_z");
509 
510  const std::vector<float> srow_x = object.slot("srow_x");
511  const std::vector<float> srow_y = object.slot("srow_y");
512  const std::vector<float> srow_z = object.slot("srow_z");
513  for (int i=0; i<4; i++)
514  {
515  header.srow_x[i] = srow_x[i];
516  header.srow_y[i] = srow_y[i];
517  header.srow_z[i] = srow_z[i];
518  }
519 
520  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)
521  header.datatype = DT_INT32;
522  else if (header.datatype == DT_FLOAT32 || header.datatype == DT_FLOAT64)
523  header.datatype = DT_FLOAT64; // This assumes that sizeof(double) == 8
524  else
525  throw std::runtime_error("Data type is not supported");
526 
527  this->image = nifti_convert_nhdr2nim(header, NULL);
528 
529  const SEXP data = PROTECT(object.slot(".Data"));
530  if (!copyData || Rf_length(data) <= 1)
531  this->image->data = NULL;
532  else
533  {
534  const size_t dataSize = nifti_get_volsize(this->image);
535  this->image->data = calloc(1, dataSize);
536  if (header.datatype == DT_INT32)
537  {
538  Rcpp::IntegerVector intData(data);
539  std::copy(intData.begin(), intData.end(), static_cast<int32_t*>(this->image->data));
540  }
541  else
542  {
543  Rcpp::DoubleVector doubleData(data);
544  std::copy(doubleData.begin(), doubleData.end(), static_cast<double*>(this->image->data));
545  }
546  }
547  UNPROTECT(1);
548 }
549 
550 inline void NiftiImage::initFromMriImage (const Rcpp::RObject &object, const bool copyData)
551 {
552  Rcpp::Reference mriImage(object);
553  Rcpp::Function getXform = mriImage.field("getXform");
554  Rcpp::NumericMatrix xform = getXform();
555 
556  this->image = NULL;
557 
558  if (Rf_length(mriImage.field("tags")) > 0)
559  initFromList(mriImage.field("tags"));
560 
561  Rcpp::RObject data = mriImage.field("data");
562  if (data.inherits("SparseArray"))
563  {
564  Rcpp::Language call("as.array", data);
565  data = call.eval();
566  }
567 
568  const short datatype = (Rf_isNull(data) ? DT_INT32 : sexpTypeToNiftiType(data.sexp_type()));
569 
570  int dims[8] = { 0, 0, 0, 0, 0, 0, 0, 0 };
571  const std::vector<int> dimVector = mriImage.field("imageDims");
572  const int nDims = std::min(7, int(dimVector.size()));
573  dims[0] = nDims;
574  size_t nVoxels = 1;
575  for (int i=0; i<nDims; i++)
576  {
577  dims[i+1] = dimVector[i];
578  nVoxels *= dimVector[i];
579  }
580 
581  if (this->image == NULL)
582  this->image = nifti_make_new_nim(dims, datatype, FALSE);
583  else
584  {
585  std::copy(dims, dims+8, this->image->dim);
586  this->image->datatype = datatype;
587  nifti_datatype_sizes(image->datatype, &image->nbyper, NULL);
588  }
589 
590  if (copyData && !Rf_isNull(data))
591  {
592  // NB: nifti_get_volsize() will not be right here if there were tags
593  const size_t dataSize = nVoxels * image->nbyper;
594  this->image->data = calloc(1, dataSize);
595  if (datatype == DT_INT32)
596  memcpy(this->image->data, INTEGER(data), dataSize);
597  else
598  memcpy(this->image->data, REAL(data), dataSize);
599  }
600  else
601  this->image->data = NULL;
602 
603  const std::vector<float> pixdimVector = mriImage.field("voxelDims");
604  const int pixdimLength = pixdimVector.size();
605  for (int i=0; i<std::min(pixdimLength,nDims); i++)
606  this->image->pixdim[i+1] = std::abs(pixdimVector[i]);
607 
608  const std::vector<std::string> pixunitsVector = mriImage.field("voxelDimUnits");
609  setPixunits(pixunitsVector);
610 
611  if (xform.rows() != 4 || xform.cols() != 4)
612  this->image->qform_code = this->image->sform_code = 0;
613  else
614  {
615  mat44 matrix;
616  for (int i=0; i<4; i++)
617  {
618  for (int j=0; j<4; j++)
619  matrix.m[i][j] = static_cast<float>(xform(i,j));
620  }
621 
622  this->image->qto_xyz = matrix;
623  this->image->qto_ijk = nifti_mat44_inverse(image->qto_xyz);
624  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);
625 
626  this->image->sto_xyz = matrix;
627  this->image->sto_ijk = nifti_mat44_inverse(image->sto_xyz);
628 
629  this->image->qform_code = this->image->sform_code = 2;
630  }
631 }
632 
633 template <typename TargetType>
634 inline void copyIfPresent (const Rcpp::List &list, const std::set<std::string> names, const std::string &name, TargetType &target)
635 {
636  if (names.count(name) == 1)
637  target = Rcpp::as<TargetType>(list[name]);
638 }
639 
640 // Special case for char, because Rcpp tries to be too clever and convert it to a string
641 template <>
642 inline void copyIfPresent (const Rcpp::List &list, const std::set<std::string> names, const std::string &name, char &target)
643 {
644  if (names.count(name) == 1)
645  target = static_cast<char>(Rcpp::as<int>(list[name]));
646 }
647 
648 inline void NiftiImage::initFromList (const Rcpp::RObject &object)
649 {
650  Rcpp::List list(object);
651  nifti_1_header *header = nifti_make_new_header(NULL, DT_FLOAT64);
652 
653  const Rcpp::CharacterVector _names = list.names();
654  std::set<std::string> names;
655  for (Rcpp::CharacterVector::const_iterator it=_names.begin(); it!=_names.end(); it++)
656  names.insert(Rcpp::as<std::string>(*it));
657 
658  copyIfPresent(list, names, "sizeof_hdr", header->sizeof_hdr);
659 
660  copyIfPresent(list, names, "dim_info", header->dim_info);
661  if (names.count("dim") == 1)
662  {
663  std::vector<short> dim = list["dim"];
664  for (int i=0; i<std::min(dim.size(),size_t(8)); i++)
665  header->dim[i] = dim[i];
666  }
667 
668  copyIfPresent(list, names, "intent_p1", header->intent_p1);
669  copyIfPresent(list, names, "intent_p2", header->intent_p2);
670  copyIfPresent(list, names, "intent_p3", header->intent_p3);
671  copyIfPresent(list, names, "intent_code", header->intent_code);
672 
673  copyIfPresent(list, names, "datatype", header->datatype);
674  copyIfPresent(list, names, "bitpix", header->bitpix);
675 
676  copyIfPresent(list, names, "slice_start", header->slice_start);
677  if (names.count("pixdim") == 1)
678  {
679  std::vector<float> pixdim = list["pixdim"];
680  for (int i=0; i<std::min(pixdim.size(),size_t(8)); i++)
681  header->pixdim[i] = pixdim[i];
682  }
683  copyIfPresent(list, names, "vox_offset", header->vox_offset);
684  copyIfPresent(list, names, "scl_slope", header->scl_slope);
685  copyIfPresent(list, names, "scl_inter", header->scl_inter);
686  copyIfPresent(list, names, "slice_end", header->slice_end);
687  copyIfPresent(list, names, "slice_code", header->slice_code);
688  copyIfPresent(list, names, "xyzt_units", header->xyzt_units);
689  copyIfPresent(list, names, "cal_max", header->cal_max);
690  copyIfPresent(list, names, "cal_min", header->cal_min);
691  copyIfPresent(list, names, "slice_duration", header->slice_duration);
692  copyIfPresent(list, names, "toffset", header->toffset);
693 
694  if (names.count("descrip") == 1)
695  strcpy(header->descrip, Rcpp::as<std::string>(list["descrip"]).substr(0,79).c_str());
696  if (names.count("aux_file") == 1)
697  strcpy(header->aux_file, Rcpp::as<std::string>(list["aux_file"]).substr(0,23).c_str());
698 
699  copyIfPresent(list, names, "qform_code", header->qform_code);
700  copyIfPresent(list, names, "sform_code", header->sform_code);
701  copyIfPresent(list, names, "quatern_b", header->quatern_b);
702  copyIfPresent(list, names, "quatern_c", header->quatern_c);
703  copyIfPresent(list, names, "quatern_d", header->quatern_d);
704  copyIfPresent(list, names, "qoffset_x", header->qoffset_x);
705  copyIfPresent(list, names, "qoffset_y", header->qoffset_y);
706  copyIfPresent(list, names, "qoffset_z", header->qoffset_z);
707 
708  if (names.count("srow_x") == 1)
709  {
710  std::vector<float> srow_x = list["srow_x"];
711  for (int i=0; i<std::min(srow_x.size(),size_t(4)); i++)
712  header->srow_x[i] = srow_x[i];
713  }
714  if (names.count("srow_y") == 1)
715  {
716  std::vector<float> srow_y = list["srow_y"];
717  for (int i=0; i<std::min(srow_y.size(),size_t(4)); i++)
718  header->srow_y[i] = srow_y[i];
719  }
720  if (names.count("srow_z") == 1)
721  {
722  std::vector<float> srow_z = list["srow_z"];
723  for (int i=0; i<std::min(srow_z.size(),size_t(4)); i++)
724  header->srow_z[i] = srow_z[i];
725  }
726 
727  if (names.count("intent_name") == 1)
728  strcpy(header->intent_name, Rcpp::as<std::string>(list["intent_name"]).substr(0,15).c_str());
729  if (names.count("magic") == 1)
730  strcpy(header->magic, Rcpp::as<std::string>(list["magic"]).substr(0,3).c_str());
731 
732  this->image = nifti_convert_nhdr2nim(*header, NULL);
733  this->image->data = NULL;
734  free(header);
735 }
736 
737 inline void NiftiImage::initFromArray (const Rcpp::RObject &object, const bool copyData)
738 {
739  int dims[8] = { 0, 0, 0, 0, 0, 0, 0, 0 };
740  const std::vector<int> dimVector = object.attr("dim");
741 
742  const int nDims = std::min(7, int(dimVector.size()));
743  dims[0] = nDims;
744  for (int i=0; i<nDims; i++)
745  dims[i+1] = dimVector[i];
746 
747  const short datatype = sexpTypeToNiftiType(object.sexp_type());
748  this->image = nifti_make_new_nim(dims, datatype, int(copyData));
749 
750  if (copyData)
751  {
752  const size_t dataSize = nifti_get_volsize(image);
753  if (datatype == DT_INT32)
754  memcpy(this->image->data, INTEGER(object), dataSize);
755  else
756  memcpy(this->image->data, REAL(object), dataSize);
757  }
758  else
759  this->image->data = NULL;
760 
761  if (object.hasAttribute("pixdim"))
762  {
763  const std::vector<float> pixdimVector = object.attr("pixdim");
764  const int pixdimLength = pixdimVector.size();
765  for (int i=0; i<std::min(pixdimLength,nDims); i++)
766  this->image->pixdim[i+1] = pixdimVector[i];
767  }
768 
769  if (object.hasAttribute("pixunits"))
770  {
771  const std::vector<std::string> pixunitsVector = object.attr("pixunits");
772  setPixunits(pixunitsVector);
773  }
774 }
775 
776 inline NiftiImage::NiftiImage (const SEXP object, const bool readData)
777  : persistent(false)
778 {
779  Rcpp::RObject imageObject(object);
780  bool resolved = false;
781 
782  if (imageObject.hasAttribute(".nifti_image_ptr"))
783  {
784  Rcpp::XPtr<NiftiImage> imagePtr(SEXP(imageObject.attr(".nifti_image_ptr")));
785  if (imagePtr.get() != NULL)
786  {
787  this->image = *imagePtr;
788  this->persistent = true;
789  resolved = true;
790 
791  if (imageObject.hasAttribute("dim"))
792  update(object);
793  }
794  else if (Rf_isString(object))
795  throw std::runtime_error("Internal image is not valid");
796  else
797  Rf_warning("Ignoring invalid internal pointer");
798  }
799 
800  if (!resolved)
801  {
802  if (Rf_isNull(object))
803  this->image = NULL;
804  else if (Rf_isString(object))
805  {
806  const std::string path = Rcpp::as<std::string>(object);
807  this->image = nifti_image_read(path.c_str(), readData);
808  if (this->image == NULL)
809  throw std::runtime_error("Failed to read image from path " + path);
810  }
811  else if (imageObject.inherits("nifti"))
812  initFromNiftiS4(imageObject, readData);
813  else if (imageObject.inherits("MriImage"))
814  initFromMriImage(imageObject, readData);
815  else if (Rf_isVectorList(object))
816  initFromList(imageObject);
817  else if (imageObject.hasAttribute("dim"))
818  initFromArray(imageObject, readData);
819  else
820  throw std::runtime_error("Cannot convert object of class \"" + Rcpp::as<std::string>(imageObject.attr("class")) + "\" to a nifti_image");
821  }
822 
823  if (this->image != NULL)
824  nifti_update_dims_from_array(this->image);
825 
826 #ifndef NDEBUG
827  Rprintf("Creating NiftiImage with pointer %p (from SEXP)\n", this->image);
828 #endif
829 }
830 
831 inline mat33 topLeftCorner (const mat44 &matrix)
832 {
833  mat33 newMatrix;
834  for (int i=0; i<3; i++)
835  {
836  for (int j=0; j<3; j++)
837  newMatrix.m[i][j] = matrix.m[i][j];
838  }
839 
840  return newMatrix;
841 }
842 
843 inline void NiftiImage::updatePixdim (const std::vector<float> &pixdim)
844 {
845  const int nDims = image->dim[0];
846  const std::vector<float> origPixdim(image->pixdim+1, image->pixdim+4);
847 
848  for (int i=1; i<8; i++)
849  image->pixdim[i] = 0.0;
850 
851  const int pixdimLength = pixdim.size();
852  for (int i=0; i<std::min(pixdimLength,nDims); i++)
853  image->pixdim[i+1] = pixdim[i];
854 
855  if (!std::equal(origPixdim.begin(), origPixdim.begin() + std::min(3,nDims), pixdim.begin()))
856  {
857  mat33 scaleMatrix;
858  for (int i=0; i<3; i++)
859  {
860  for (int j=0; j<3; j++)
861  {
862  if (i != j)
863  scaleMatrix.m[i][j] = 0.0;
864  else if (i >= nDims)
865  scaleMatrix.m[i][j] = 1.0;
866  else
867  scaleMatrix.m[i][j] = pixdim[i] / origPixdim[i];
868  }
869  }
870 
871  if (image->qform_code > 0)
872  {
873  mat33 prod = nifti_mat33_mul(scaleMatrix, topLeftCorner(image->qto_xyz));
874  for (int i=0; i<3; i++)
875  {
876  for (int j=0; j<3; j++)
877  image->qto_xyz.m[i][j] = prod.m[i][j];
878  }
879  image->qto_ijk = nifti_mat44_inverse(image->qto_xyz);
880  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);
881  }
882 
883  if (image->sform_code > 0)
884  {
885  mat33 prod = nifti_mat33_mul(scaleMatrix, topLeftCorner(image->sto_xyz));
886  for (int i=0; i<3; i++)
887  {
888  for (int j=0; j<3; j++)
889  image->sto_xyz.m[i][j] = prod.m[i][j];
890  }
891  image->sto_ijk = nifti_mat44_inverse(image->sto_xyz);
892  }
893  }
894 }
895 
896 inline void NiftiImage::setPixunits (const std::vector<std::string> &pixunits)
897 {
898  for (int i=0; i<pixunits.size(); i++)
899  {
900  if (pixunits[i] == "m")
901  image->xyz_units = NIFTI_UNITS_METER;
902  else if (pixunits[i] == "mm")
903  image->xyz_units = NIFTI_UNITS_MM;
904  else if (pixunits[i] == "um")
905  image->xyz_units = NIFTI_UNITS_MICRON;
906  else if (pixunits[i] == "s")
907  image->time_units = NIFTI_UNITS_SEC;
908  else if (pixunits[i] == "ms")
909  image->time_units = NIFTI_UNITS_MSEC;
910  else if (pixunits[i] == "us")
911  image->time_units = NIFTI_UNITS_USEC;
912  else if (pixunits[i] == "Hz")
913  image->time_units = NIFTI_UNITS_HZ;
914  else if (pixunits[i] == "ppm")
915  image->time_units = NIFTI_UNITS_PPM;
916  else if (pixunits[i] == "rad/s")
917  image->time_units = NIFTI_UNITS_RADS;
918  }
919 }
920 
921 inline void NiftiImage::rescale (const std::vector<float> &scales)
922 {
923  std::vector<float> pixdim(image->pixdim+1, image->pixdim+4);
924 
925  for (int i=0; i<std::min(3, int(scales.size())); i++)
926  {
927  if (scales[i] != 1.0)
928  {
929  pixdim[i] /= scales[i];
930  image->dim[i+1] = static_cast<int>(std::floor(image->dim[i+1] * scales[i]));
931  }
932  }
933 
934  updatePixdim(pixdim);
935  nifti_update_dims_from_array(image);
936 
937  // Data vector is now the wrong size, so drop it
938  free(image->data);
939 }
940 
941 inline void NiftiImage::update (const SEXP array)
942 {
943  Rcpp::RObject object(array);
944  if (!object.hasAttribute("dim"))
945  return;
946 
947  for (int i=0; i<8; i++)
948  image->dim[i] = 0;
949  const std::vector<int> dimVector = object.attr("dim");
950 
951  const int nDims = std::min(7, int(dimVector.size()));
952  image->dim[0] = nDims;
953  for (int i=0; i<nDims; i++)
954  image->dim[i+1] = dimVector[i];
955 
956  if (object.hasAttribute("pixdim"))
957  {
958  const std::vector<float> pixdimVector = object.attr("pixdim");
959  updatePixdim(pixdimVector);
960  }
961 
962  if (object.hasAttribute("pixunits"))
963  {
964  const std::vector<std::string> pixunitsVector = object.attr("pixunits");
965  setPixunits(pixunitsVector);
966  }
967 
968  // This NIfTI-1 library function clobbers dim[0] if the last dimension is unitary; we undo that here
969  nifti_update_dims_from_array(image);
970  image->dim[0] = image->ndim = nDims;
971 
972  image->datatype = NiftiImage::sexpTypeToNiftiType(object.sexp_type());
973  nifti_datatype_sizes(image->datatype, &image->nbyper, NULL);
974 
975  free(image->data);
976 
977  const size_t dataSize = nifti_get_volsize(image);
978  image->data = calloc(1, dataSize);
979  if (image->datatype == DT_INT32)
980  memcpy(image->data, INTEGER(object), dataSize);
981  else
982  memcpy(image->data, REAL(object), dataSize);
983 }
984 
985 inline mat44 NiftiImage::xform (const bool preferQuaternion) const
986 {
987  if (image == NULL)
988  {
989  mat44 matrix;
990  for (int i=0; i<4; i++)
991  {
992  for (int j=0; j<4; j++)
993  matrix.m[i][j] = 0.0;
994  }
995  return matrix;
996  }
997  else if (image->qform_code <= 0 && image->sform_code <= 0)
998  {
999  // No qform or sform so use pixdim (NB: other software may assume differently)
1000  mat44 matrix;
1001  for (int i=0; i<4; i++)
1002  {
1003  for (int j=0; j<4; j++)
1004  {
1005  if (i != j)
1006  matrix.m[i][j] = 0.0;
1007  else if (i == 3)
1008  matrix.m[3][3] = 1.0;
1009  else
1010  matrix.m[i][j] = (image->pixdim[i+1]==0.0 ? 1.0 : image->pixdim[i+1]);
1011  }
1012  }
1013  return matrix;
1014  }
1015  else if ((preferQuaternion && image->qform_code > 0) || image->sform_code <= 0)
1016  return image->qto_xyz;
1017  else
1018  return image->sto_xyz;
1019 }
1020 
1021 template <typename SourceType, typename TargetType>
1022 inline TargetType convertValue (SourceType value)
1023 {
1024  return static_cast<TargetType>(value);
1025 }
1026 
1027 template <typename SourceType, typename TargetType>
1028 inline void convertArray (const SourceType *source, const size_t length, TargetType *target)
1029 {
1030  std::transform(source, source + length, target, convertValue<SourceType,TargetType>);
1031 }
1032 
1033 template <typename TargetType>
1034 inline void changeDatatype (nifti_image *image, const short datatype)
1035 {
1036  TargetType *data;
1037  const size_t dataSize = image->nvox * sizeof(TargetType);
1038  data = static_cast<TargetType *>(calloc(1, dataSize));
1039 
1040  switch (image->datatype)
1041  {
1042  case DT_UINT8:
1043  convertArray(static_cast<uint8_t *>(image->data), image->nvox, data);
1044  break;
1045 
1046  case DT_INT16:
1047  convertArray(static_cast<int16_t *>(image->data), image->nvox, data);
1048  break;
1049 
1050  case DT_INT32:
1051  convertArray(static_cast<int32_t *>(image->data), image->nvox, data);
1052  break;
1053 
1054  case DT_FLOAT32:
1055  convertArray(static_cast<float *>(image->data), image->nvox, data);
1056  break;
1057 
1058  case DT_FLOAT64:
1059  convertArray(static_cast<double *>(image->data), image->nvox, data);
1060  break;
1061 
1062  case DT_INT8:
1063  convertArray(static_cast<int8_t *>(image->data), image->nvox, data);
1064  break;
1065 
1066  case DT_UINT16:
1067  convertArray(static_cast<uint16_t *>(image->data), image->nvox, data);
1068  break;
1069 
1070  case DT_UINT32:
1071  convertArray(static_cast<uint32_t *>(image->data), image->nvox, data);
1072  break;
1073 
1074  case DT_INT64:
1075  convertArray(static_cast<int64_t *>(image->data), image->nvox, data);
1076  break;
1077 
1078  case DT_UINT64:
1079  convertArray(static_cast<uint64_t *>(image->data), image->nvox, data);
1080  break;
1081 
1082  default:
1083  throw std::runtime_error("Unsupported data type (" + std::string(nifti_datatype_string(image->datatype)) + ")");
1084  }
1085 
1086  free(image->data);
1087  image->data = data;
1088  image->datatype = datatype;
1089  nifti_datatype_sizes(datatype, &image->nbyper, &image->swapsize);
1090 }
1091 
1092 inline void NiftiImage::toFile (const std::string fileName, const short datatype) const
1093 {
1094  // Copy the source image only if the datatype will be changed
1095  NiftiImage imageToWrite(image, datatype != DT_NONE);
1096 
1097  switch (datatype)
1098  {
1099  case DT_NONE:
1100  imageToWrite.setPersistence(true);
1101  break;
1102 
1103  case DT_UINT8:
1104  changeDatatype<uint8_t>(imageToWrite, datatype);
1105  break;
1106 
1107  case DT_INT16:
1108  changeDatatype<int16_t>(imageToWrite, datatype);
1109  break;
1110 
1111  case DT_INT32:
1112  changeDatatype<int32_t>(imageToWrite, datatype);
1113  break;
1114 
1115  case DT_FLOAT32:
1116  changeDatatype<float>(imageToWrite, datatype);
1117  break;
1118 
1119  case DT_FLOAT64:
1120  changeDatatype<double>(imageToWrite, datatype);
1121  break;
1122 
1123  case DT_INT8:
1124  changeDatatype<int8_t>(imageToWrite, datatype);
1125  break;
1126 
1127  case DT_UINT16:
1128  changeDatatype<uint16_t>(imageToWrite, datatype);
1129  break;
1130 
1131  case DT_UINT32:
1132  changeDatatype<uint32_t>(imageToWrite, datatype);
1133  break;
1134 
1135  case DT_INT64:
1136  changeDatatype<int64_t>(imageToWrite, datatype);
1137  break;
1138 
1139  case DT_UINT64:
1140  changeDatatype<uint64_t>(imageToWrite, datatype);
1141  break;
1142 
1143  default:
1144  throw std::runtime_error("Unsupported data type (" + std::string(nifti_datatype_string(datatype)) + ")");
1145  }
1146 
1147  const int status = nifti_set_filenames(imageToWrite, fileName.c_str(), false, true);
1148  if (status != 0)
1149  throw std::runtime_error("Failed to set filenames for NIfTI object");
1150  nifti_image_write(imageToWrite);
1151 }
1152 
1153 inline void NiftiImage::toFile (const std::string fileName, const std::string &datatype) const
1154 {
1155  static std::map<std::string,short> datatypeCodes;
1156  if (datatypeCodes.empty())
1157  {
1158  datatypeCodes["auto"] = DT_NONE;
1159  datatypeCodes["none"] = DT_NONE;
1160  datatypeCodes["unknown"] = DT_NONE;
1161  datatypeCodes["uint8"] = DT_UINT8;
1162  datatypeCodes["char"] = DT_UINT8;
1163  datatypeCodes["int16"] = DT_INT16;
1164  datatypeCodes["short"] = DT_INT16;
1165  datatypeCodes["int32"] = DT_INT32;
1166  datatypeCodes["int"] = DT_INT32;
1167  datatypeCodes["float32"] = DT_FLOAT32;
1168  datatypeCodes["float"] = DT_FLOAT32;
1169  datatypeCodes["float64"] = DT_FLOAT64;
1170  datatypeCodes["double"] = DT_FLOAT64;
1171  datatypeCodes["int8"] = DT_INT8;
1172  datatypeCodes["uint16"] = DT_UINT16;
1173  datatypeCodes["uint32"] = DT_UINT32;
1174  datatypeCodes["int64"] = DT_INT64;
1175  datatypeCodes["uint64"] = DT_UINT64;
1176  }
1177 
1178  if (datatypeCodes.count(datatype) == 0)
1179  {
1180  std::ostringstream message;
1181  message << "Datatype \"" << datatype << "\" is not valid";
1182  Rf_warning(message.str().c_str());
1183  toFile(fileName, DT_NONE);
1184  }
1185  else
1186  toFile(fileName, datatypeCodes[datatype]);
1187 }
1188 
1189 template <typename SourceType, int SexpType>
1190 inline Rcpp::RObject imageDataToArray (const nifti_image *source)
1191 {
1192  if (source == NULL)
1193  return Rcpp::RObject();
1194  else if (source->data == NULL)
1195  {
1196  Rf_warning("Internal image contains no data - filling array with NAs");
1197  Rcpp::Vector<SexpType> array(static_cast<int>(source->nvox));
1198  // Rcpp's proxy infrastructure should handle converting NA_REAL to the appropriate NA
1199  std::fill(array.begin(), array.end(), NA_REAL);
1200  return array;
1201  }
1202  else
1203  {
1204  SourceType *original = static_cast<SourceType *>(source->data);
1205  Rcpp::Vector<SexpType> array(static_cast<int>(source->nvox));
1206 
1207  if (SexpType == INTSXP || SexpType == LGLSXP)
1208  std::transform(original, original + source->nvox, array.begin(), convertValue<SourceType,int>);
1209  else if (SexpType == REALSXP)
1210  std::transform(original, original + source->nvox, array.begin(), convertValue<SourceType,double>);
1211  else
1212  throw std::runtime_error("Only numeric arrays can be created");
1213 
1214  return array;
1215  }
1216 }
1217 
1218 inline void finaliseNiftiImage (SEXP xptr)
1219 {
1220  NiftiImage *object = (NiftiImage *) R_ExternalPtrAddr(xptr);
1221  object->setPersistence(false);
1222  delete object;
1223  R_ClearExternalPtr(xptr);
1224 }
1225 
1226 inline void addAttributes (Rcpp::RObject &object, nifti_image *source, const bool realDim = true)
1227 {
1228  const int nDims = source->dim[0];
1229  Rcpp::IntegerVector dim(source->dim+1, source->dim+1+nDims);
1230 
1231  if (realDim)
1232  object.attr("dim") = dim;
1233  else
1234  object.attr("imagedim") = dim;
1235 
1236  Rcpp::DoubleVector pixdim(source->pixdim+1, source->pixdim+1+nDims);
1237  object.attr("pixdim") = pixdim;
1238 
1239  if (source->xyz_units == NIFTI_UNITS_UNKNOWN && source->time_units == NIFTI_UNITS_UNKNOWN)
1240  object.attr("pixunits") = "Unknown";
1241  else
1242  {
1243  Rcpp::CharacterVector pixunits(2);
1244  pixunits[0] = nifti_units_string(source->xyz_units);
1245  pixunits[1] = nifti_units_string(source->time_units);
1246  object.attr("pixunits") = pixunits;
1247  }
1248 
1249  NiftiImage *wrappedSource = new NiftiImage(source, true);
1250  wrappedSource->setPersistence(true);
1251  Rcpp::XPtr<NiftiImage> xptr(wrappedSource);
1252  R_RegisterCFinalizerEx(SEXP(xptr), &finaliseNiftiImage, FALSE);
1253  object.attr(".nifti_image_ptr") = xptr;
1254 }
1255 
1256 inline Rcpp::RObject NiftiImage::toArray () const
1257 {
1258  Rcpp::RObject array;
1259 
1260  if (this->isNull())
1261  return array;
1262 
1263  switch (image->datatype)
1264  {
1265  case DT_UINT8:
1266  array = imageDataToArray<uint8_t,INTSXP>(image);
1267  break;
1268 
1269  case DT_INT16:
1270  array = imageDataToArray<int16_t,INTSXP>(image);
1271  break;
1272 
1273  case DT_INT32:
1274  array = imageDataToArray<int32_t,INTSXP>(image);
1275  break;
1276 
1277  case DT_FLOAT32:
1278  array = imageDataToArray<float,REALSXP>(image);
1279  break;
1280 
1281  case DT_FLOAT64:
1282  array = imageDataToArray<double,REALSXP>(image);
1283  break;
1284 
1285  case DT_INT8:
1286  array = imageDataToArray<int8_t,INTSXP>(image);
1287  break;
1288 
1289  case DT_UINT16:
1290  array = imageDataToArray<uint16_t,INTSXP>(image);
1291  break;
1292 
1293  case DT_UINT32:
1294  array = imageDataToArray<uint32_t,INTSXP>(image);
1295  break;
1296 
1297  case DT_INT64:
1298  array = imageDataToArray<int64_t,INTSXP>(image);
1299  break;
1300 
1301  case DT_UINT64:
1302  array = imageDataToArray<uint64_t,INTSXP>(image);
1303  break;
1304 
1305  default:
1306  throw std::runtime_error("Unsupported data type (" + std::string(nifti_datatype_string(image->datatype)) + ")");
1307  }
1308 
1309  addAttributes(array, image);
1310  const Rcpp::IntegerVector dim = array.attr("dim");
1311  array.attr("class") = Rcpp::CharacterVector::create("niftiImage");
1312  return array;
1313 }
1314 
1315 inline Rcpp::RObject NiftiImage::toPointer (const std::string label) const
1316 {
1317  if (this->isNull())
1318  return Rcpp::RObject();
1319  else
1320  {
1321  Rcpp::RObject string = Rcpp::wrap(label);
1322  addAttributes(string, image, false);
1323  string.attr("class") = Rcpp::CharacterVector::create("internalImage", "niftiImage");
1324  return string;
1325  }
1326 }
1327 
1328 inline Rcpp::RObject NiftiImage::toArrayOrPointer (const bool internal, const std::string label) const
1329 {
1330  return (internal ? toPointer(label) : toArray());
1331 }
1332 
1333 inline Rcpp::RObject NiftiImage::headerToList () const
1334 {
1335  if (this->image == NULL)
1336  return Rcpp::RObject();
1337 
1338  nifti_1_header header = nifti_convert_nim2nhdr(this->image);
1339  Rcpp::List result;
1340 
1341  result["sizeof_hdr"] = header.sizeof_hdr;
1342 
1343  result["dim_info"] = int(header.dim_info);
1344  result["dim"] = std::vector<short>(header.dim, header.dim+8);
1345 
1346  result["intent_p1"] = header.intent_p1;
1347  result["intent_p2"] = header.intent_p2;
1348  result["intent_p3"] = header.intent_p3;
1349  result["intent_code"] = header.intent_code;
1350 
1351  result["datatype"] = header.datatype;
1352  result["bitpix"] = header.bitpix;
1353 
1354  result["slice_start"] = header.slice_start;
1355  result["pixdim"] = std::vector<float>(header.pixdim, header.pixdim+8);
1356  result["vox_offset"] = header.vox_offset;
1357  result["scl_slope"] = header.scl_slope;
1358  result["scl_inter"] = header.scl_inter;
1359  result["slice_end"] = header.slice_end;
1360  result["slice_code"] = int(header.slice_code);
1361  result["xyzt_units"] = int(header.xyzt_units);
1362  result["cal_max"] = header.cal_max;
1363  result["cal_min"] = header.cal_min;
1364  result["slice_duration"] = header.slice_duration;
1365  result["toffset"] = header.toffset;
1366  result["descrip"] = std::string(header.descrip, 80);
1367  result["aux_file"] = std::string(header.aux_file, 24);
1368 
1369  result["qform_code"] = header.qform_code;
1370  result["sform_code"] = header.sform_code;
1371  result["quatern_b"] = header.quatern_b;
1372  result["quatern_c"] = header.quatern_c;
1373  result["quatern_d"] = header.quatern_d;
1374  result["qoffset_x"] = header.qoffset_x;
1375  result["qoffset_y"] = header.qoffset_y;
1376  result["qoffset_z"] = header.qoffset_z;
1377  result["srow_x"] = std::vector<float>(header.srow_x, header.srow_x+4);
1378  result["srow_y"] = std::vector<float>(header.srow_y, header.srow_y+4);
1379  result["srow_z"] = std::vector<float>(header.srow_z, header.srow_z+4);
1380 
1381  result["intent_name"] = std::string(header.intent_name, 16);
1382  result["magic"] = std::string(header.magic, 4);
1383 
1384  result.attr("class") = Rcpp::CharacterVector::create("niftiHeader");
1385 
1386  return result;
1387 }
1388 
1389 #endif
Block(const NiftiImage &image, const int dimension, const int index)
Standard constructor for this class.
Definition: NiftiImage.h:40
NiftiImage(const std::string &path, const bool readData=true)
Initialise using a path string.
Definition: NiftiImage.h:194
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:55
NiftiImage(const NiftiImage &source)
Copy constructor.
Definition: NiftiImage.h:161
Block slice(const int i)
Extract a slice block from a 3D image.
Definition: NiftiImage.h:345
Thin wrapper around a C-style nifti_image struct that allows C++-style destruction.
Definition: NiftiImage.h:20
const int index
The location along dimension.
Definition: NiftiImage.h:31
void setPersistence(const bool persistent)
Allows the image to be marked as persistent, so that it can be passed back to R.
Definition: NiftiImage.h:267
bool isPersistent() const
Determines whether or not the image is marked as persistent.
Definition: NiftiImage.h:284
nifti_image * image
The wrapped nifti_image pointer.
Definition: NiftiImage.h:90
void initFromNiftiS4(const Rcpp::RObject &object, const bool copyData=true)
Initialise the object from an S4 object of class "nifti".
Definition: NiftiImage.h:454
Inner class referring to a subset of an image.
Definition: NiftiImage.h:27
void copy(nifti_image *const source)
Copy the contents of a nifti_image to create a new image.
Definition: NiftiImage.h:403
NiftiImage()
Default constructor.
Definition: NiftiImage.h:154
Block volume(const int i)
Extract a volume block from a 4D image.
Definition: NiftiImage.h:359
bool persistent
Marker of persistence, which determines whether the nifti_image should be freed on destruction...
Definition: NiftiImage.h:91
const int dimension
The dimension along which the block applies (which should be the last)
Definition: NiftiImage.h:30
const NiftiImage & image
The parent image.
Definition: NiftiImage.h:29
void toFile(const std::string fileName, const short datatype) const
Write the image to a NIfTI-1 file.
Definition: NiftiImage.h:1092
void update(const SEXP array)
Update the image from an R array.
Definition: NiftiImage.h:941
NiftiImage(nifti_image *const image, const bool copy=false)
Initialise using an existing nifti_image pointer.
Definition: NiftiImage.h:176
static short sexpTypeToNiftiType(const int sexpType)
Convert between R SEXP object type and nifti_image datatype codes.
Definition: NiftiImage.h:79
void initFromArray(const Rcpp::RObject &object, const bool copyData=true)
Initialise the object from an R array.
Definition: NiftiImage.h:737
Rcpp::RObject headerToList() const
Create an R list containing raw image metadata.
Definition: NiftiImage.h:1333
~NiftiImage()
Destructor which frees the wrapped pointer, unless the object is marked as persistent.
Definition: NiftiImage.h:215
mat44 xform(const bool preferQuaternion=true) const
Obtain an xform matrix, indicating the orientation of the image.
Definition: NiftiImage.h:985
void rescale(const std::vector< float > &scales)
Rescales the image, changing its image dimensions and pixel dimensions.
Definition: NiftiImage.h:921
nifti_image * operator->() const
Allows a NiftiImage object to be treated as a pointer to a nifti_image.
Definition: NiftiImage.h:234
bool isNull() const
Determines whether or not the internal pointer is NULL.
Definition: NiftiImage.h:279
void updatePixdim(const std::vector< float > &pixdim)
Modify the pixel dimensions, and potentially the xform matrices to match.
Definition: NiftiImage.h:843
const Block slice(const int i) const
Extract a slice block from a 3D image.
Definition: NiftiImage.h:338
int nDims() const
Returns the number of dimensions in the image.
Definition: NiftiImage.h:289
Rcpp::RObject toArray() const
Create an R array from the image.
Definition: NiftiImage.h:1256
Rcpp::RObject toArrayOrPointer(const bool internal, const std::string label) const
A conditional method that calls either toArray or toPointer.
Definition: NiftiImage.h:1328
void initFromList(const Rcpp::RObject &object)
Initialise the object from an R list with named elements, which can only contain metadata.
Definition: NiftiImage.h:648
void initFromMriImage(const Rcpp::RObject &object, const bool copyData=true)
Initialise the object from a reference object of class "MriImage".
Definition: NiftiImage.h:550
NiftiImage & drop()
Drop unitary dimensions.
Definition: NiftiImage.h:303
const Block volume(const int i) const
Extract a volume block from a 4D image.
Definition: NiftiImage.h:352
Rcpp::RObject toPointer(const std::string label) const
Create an internal image to pass back to R.
Definition: NiftiImage.h:1315
void setPixunits(const std::vector< std::string > &pixunits)
Modify the pixel dimension units.
Definition: NiftiImage.h:896