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