#pragma once#include<fstream>#include<vtkContourFilter.h>#include<vtkMassProperties.h>#include<itkImage.h>#include<ParticleRegionDomain.h>#include<itkImageToVTKImageFilter.h>#include<itkZeroCrossingImageFilter.h>#include<itkImageRegionConstIteratorWithIndex.h>#include<chrono>// we have to undef foreach here because both Qt and OpenVDB define foreach#undef foreach#ifndef Q_MOC_RUN#include<openvdb/openvdb.h>#include<openvdb/tools/SignedFloodFill.h>#include<openvdb/tools/Interpolation.h>#include<openvdb/tools/GridOperators.h>#include<openvdb/math/Math.h>#include<openvdb/math/Transform.h>#endifnamespaceshapeworks{template<classT>classParticleImageDomain:publicParticleRegionDomain{public:usingPointer=std::shared_ptr<ParticleImageDomain>;usingImageType=itk::Image<T,DIMENSION>;usingPointType=ParticleRegionDomain::PointType;voidSetImage(ImageType*I,doublenarrow_band){this->m_FixedDomain=false;//this->Modified();openvdb::initialize();// It is safe to initialize multiple times.// Set a large background value, so that we quickly catch particles outside or on the edge the narrow band.// (Downside: its more difficult to display the correct location of the point of failure.)m_VDBImage=openvdb::FloatGrid::create(1e8);m_VDBImage->setGridClass(openvdb::GRID_LEVEL_SET);autovdbAccessor=m_VDBImage->getAccessor();// Save properties of the Image needed for the optimizerm_Size=I->GetRequestedRegion().GetSize();m_Spacing=I->GetSpacing();m_Origin=I->GetOrigin();m_Index=I->GetRequestedRegion().GetIndex();// Transformation from index space to world spaceopenvdb::math::Mat4fmat;mat.setIdentity();mat.postScale(openvdb::Vec3f(m_Spacing[0],m_Spacing[1],m_Spacing[2]));mat.postTranslate(openvdb::Vec3f(m_Origin[0],m_Origin[1],m_Origin[2]));constautoxform=openvdb::math::Transform::createLinearTransform(mat);m_VDBImage->setTransform(xform);itk::ImageRegionIterator<ImageType>it(I,I->GetRequestedRegion());it.GoToBegin();while(!it.IsAtEnd()){constautoidx=it.GetIndex();constautopixel=it.Get();if(abs(pixel)>narrow_band){++it;continue;}constautocoord=openvdb::Coord(idx[0],idx[1],idx[2]);vdbAccessor.setValue(coord,pixel);++it;}typenameImageType::PointTypel0;I->TransformIndexToPhysicalPoint(m_Index,l0);for(unsignedinti=0;i<DIMENSION;i++)m_Index[i]+=m_Size[i]-1;typenameImageType::PointTypeu0;I->TransformIndexToPhysicalPoint(m_Index,u0);// Cast points to higher precision if needed. Parent class uses doubles// because they are compared directly with points in the particle system,// which are always double precision.PointTypel;PointTypeu;for(unsignedinti=0;i<DIMENSION;i++){l[i]=static_cast<double>(l0[i]);u[i]=static_cast<double>(u0[i]);}this->SetLowerBound(l);this->SetUpperBound(u);// Precompute and save values that are used in parts of the optimizerthis->SetupImageForCrossingPointUpdate(I);// this->UpdateZeroCrossingPoint(I);this->UpdateSurfaceArea(I);}inlinedoubleGetSurfaceArea()constoverride{throwstd::runtime_error("Surface area is not computed currently.");returnm_SurfaceArea;}inlinePointTypeGetOrigin()const{returnm_Origin;}inlinetypenameImageType::SizeTypeGetSize()const{returnm_Size;}inlinetypenameImageType::SpacingTypeGetSpacing()const{returnm_Spacing;}inlinetypenameImageType::RegionType::IndexTypeGetIndex()const{returnm_Index;}inlinePointTypeGetValidLocationNear(PointTypep)constoverride{// todo why is this function ignoring the argument? Also see Optimize::AddSinglePointreturnm_ZeroCrossingPoint;}inlineTSample(constPointType&p)const{if(this->IsInsideBuffer(p)){constautocoord=this->ToVDBCoord(p);returnopenvdb::tools::BoxSampler::sample(m_VDBImage->tree(),coord);}else{std::ostringstreammessage;\
message<<"Domain "<<m_DomainID<<": "<<m_DomainName<<" : Distance transform queried for a Point, "<<p<<", outside the given image domain. Consider increasing the narrow band";throwstd::runtime_error(message.str());}}inlinedoubleGetMaxDiameter()constoverride{doublebestRadius=0;doublemaxdim=0;for(unsignedinti=0;i<ImageType::ImageDimension;i++){if(GetSize()[i]>maxdim){maxdim=GetSize()[i];bestRadius=maxdim*GetSpacing()[i];}}returnbestRadius;}voidDeleteImages()override{m_VDBImage=0;}// Updates zero crossing points. Raster scans candidate zero crossing points, and finds one that does not violate any constraints.voidUpdateZeroCrossingPoint()override{for(size_ti=0;i<m_possible_zero_crossings.size();i++){this->m_ZeroCrossingPoint=m_possible_zero_crossings[i];if(!this->GetConstraints()->isAnyViolated(this->m_ZeroCrossingPoint)){//std::cout << "Chosen initial point " << this->m_ZeroCrossingPoint << std::endl;break;}}if(this->GetConstraints()->isAnyViolated(this->m_ZeroCrossingPoint)){std::cerr<<"A particle initialization violates at least one constraint. Make sure at least one point satisfies all constraints"<<std::endl;}}protected:openvdb::FloatGrid::PtrGetVDBImage()const{returnm_VDBImage;}ParticleImageDomain(){}virtual~ParticleImageDomain(){};voidPrintSelf(std::ostream&os,itk::Indentindent)const{ParticleRegionDomain::PrintSelf(os,indent);os<<indent<<"VDB Active Voxels = "<<m_VDBImage->activeVoxelCount()<<std::endl;}inlineopenvdb::math::Transform::Ptrtransform()const{returnthis->m_VDBImage->transformPtr();}// Converts a coordinate from an ITK Image point in world space to the corresponding// coordinate in OpenVDB Index space. Raises an exception if the narrow band is not// sufficiently large to sample the point.inlineopenvdb::Vec3RToVDBCoord(constPointType&p)const{constautoworldCoord=openvdb::Vec3R(p[0],p[1],p[2]);constautoidxCoord=this->transform()->worldToIndex(worldCoord);// Make sure the coordinate is part of the narrow bandif(m_VDBImage->tree().isValueOff(openvdb::Coord::round(idxCoord))){// `isValueOff` requires an integer coordinate// If multiple threads crash here at the same time, the error message displayed is just "terminate called recursively",// which isn't helpful. So we std::cerr the error to make sure its printed to the console.std::cerr<<"Sampled point outside the narrow band: "<<p<<std::endl;std::ostringstreammessage;\
message<<"Attempt to sample at a point outside the narrow band: "<<p<<". Consider increasing the narrow band";throwstd::runtime_error(message.str());}returnidxCoord;}private:openvdb::FloatGrid::Ptrm_VDBImage;typenameImageType::SizeTypem_Size;typenameImageType::SpacingTypem_Spacing;PointTypem_Origin;PointTypem_ZeroCrossingPoint;typenameImageType::RegionType::IndexTypem_Index;// Index defining the corner of the regiondoublem_SurfaceArea;std::vector<PointType>m_possible_zero_crossings;// Computes possible zero crossing points. Later on, one can find the ones that do not violate constraints.voidSetupImageForCrossingPointUpdate(ImageType*I){typenameitk::ZeroCrossingImageFilter<ImageType,ImageType>::Pointerzc=itk::ZeroCrossingImageFilter<ImageType,ImageType>::New();zc->SetInput(I);zc->Update();typenameitk::ImageRegionConstIteratorWithIndex<ImageType>zcIt(zc->GetOutput(),zc->GetOutput()->GetRequestedRegion());for(zcIt.GoToReverseBegin();!zcIt.IsAtReverseEnd();--zcIt){if(zcIt.Get()==1.0){PointTypepos;I->TransformIndexToPhysicalPoint(zcIt.GetIndex(),pos);this->m_ZeroCrossingPoint=pos;m_possible_zero_crossings.push_back(pos);}}}voidUpdateSurfaceArea(ImageType*I){//TODO: This code has been copied from Optimize.cpp. It does not work/* typename itk::ImageToVTKImageFilter < ImageType > ::Pointer itk2vtkConnector; itk2vtkConnector = itk::ImageToVTKImageFilter < ImageType > ::New(); itk2vtkConnector->SetInput(I); vtkSmartPointer < vtkContourFilter > ls = vtkSmartPointer < vtkContourFilter > ::New(); ls->SetInputData(itk2vtkConnector->GetOutput()); ls->SetValue(0, 0.0); ls->Update(); vtkSmartPointer < vtkMassProperties > mp = vtkSmartPointer < vtkMassProperties > ::New(); mp->SetInputData(ls->GetOutput()); mp->Update(); m_SurfaceArea = mp->GetSurfaceArea(); */}};}// end namespace itk