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