19 #include <tclap/CmdLine.h> 22 #include <vtkSmartPointer.h> 23 #include <vtkPolyData.h> 24 #include <vtkMultiThreader.h> 25 #include <vtkPolyDataCollection.h> 26 #include <vtkTransformCollection.h> 29 #include <milxDeformableModel.h> 31 using namespace TCLAP;
33 typedef std::vector< std::string >::iterator stringiterator;
35 typedef unsigned char charPixelType;
36 typedef itk::Image<charPixelType, milx::imgDimension> charImageType;
39 enum operations {none = 0, convert, duplicate, cat, split, scale, decimate, smooth, laplacian, thresholdscalars, flip, diffscalars, copyscalars, diffscalarspairs, statscalars, removescalars, mse, procrustes, icp, clip,
40 voxelise, orient, flat};
75 int main(
int argc,
char *argv[])
79 milx::PrintInfo(
"--------------------------------------------------------");
80 milx::PrintInfo(
"SMILI Deformable Model Tool for Models/Surfaces/Meshes.");
84 milx::PrintInfo(
"Australian e-Health Research Centre, CSIRO, Australia.");
85 milx::PrintInfo(
"--------------------------------------------------------\n");
89 CmdLine cmd(
"A diagnostic tool for deformable models/surface operations",
' ',
milx::NumberToString(milx::Version));
92 ValueArg<size_t> threadsArg(
"",
"threads",
"Set he number of global threads to use.",
false,
milx::NumberOfProcessors(),
"Threads");
93 ValueArg<std::string> outputArg(
"o",
"output",
"Output Surface",
false,
"result.vtk",
"Output");
94 ValueArg<std::string> prefixArg(
"p",
"prefix",
"Output prefix for multiple output",
false,
"surface_",
"Output Prefix");
95 ValueArg<std::string> outputFormatArg(
"",
"outputformat",
"Specify the default output format for multiple outputs (vtk, vtp, ply, stl, default: same as input)",
false,
"vtk",
"Output format");
97 ValueArg<float> decimateArg(
"d",
"decimate",
"Decimate all the meshes provided using the Quadric Decimate algorithm with decimation factor.",
false, 0.5,
"Decimate");
98 ValueArg<float> smoothArg(
"",
"smooth",
"Smooth all the meshes provided using the Windowed Sinc Algorithm with iterations.",
false, 18,
"Smooth");
99 ValueArg<float> laplacianArg(
"",
"laplacian",
"Smooth all the meshes provided using the Laplacian Algorithm with iterations.",
false, 18,
"Laplacian");
100 ValueArg<float> scaleArg(
"s",
"scale",
"Scale the coordinates of the points by scale.",
false, 0.9,
"Scale");
101 ValueArg<float> thresholdAboveArg(
"",
"thresholdabove",
"Thresold scalars above value.",
false, 0.0,
"Above");
102 ValueArg<float> thresholdBelowArg(
"",
"thresholdbelow",
"Thresold scalars below value.",
false, 0.0,
"Below");
103 ValueArg<float> clipArg(
"",
"clip",
"Clip model based on scalars value (keeping only parts with value).",
false, 1.0,
"Clip");
104 MultiArg<std::string> componentsArg(
"",
"component",
"Surface is a component of the surfaces.",
false,
"Component");
106 std::vector<size_t> axesAllowed;
107 axesAllowed.push_back(0);
108 axesAllowed.push_back(1);
109 axesAllowed.push_back(2);
110 ValuesConstraint<size_t> allowedAxesVals( axesAllowed );
111 ValueArg<size_t> flipArg(
"f",
"flip",
"Flip the meshes in the axis provided (0: x-axis, 1: y-axis, 2: z-axis).",
false, 0, &allowedAxesVals);
113 SwitchArg partitionArg(
"",
"partitionList",
"Partition the list of surfaces provided into two for operating on one vs the other",
false);
115 SwitchArg duplicateArg(
"",
"duplicate",
"Simply open and save the input file(s). Supports single (-o) or multiple (-p) input.",
false);
116 SwitchArg convertArg(
"c",
"convert",
"Convert the model from the current format to the one given at output (from file extension). Supports single (-o) or multiple (-p) input.",
false);
117 SwitchArg concatenateArg(
"",
"cat",
"Concatenate N surfaces into single mesh with filename provided.",
false);
118 SwitchArg colourConcatenateArg(
"",
"colourcat",
"Concatenate N surfaces into single mesh with filename provided and colour them.",
false);
119 SwitchArg splitArg(
"",
"split",
"Split each surface given components.",
false);
120 SwitchArg diffScalarArg(
"",
"scalardiff",
"Compute the differences in Scalars.",
false);
121 SwitchArg statsScalarArg(
"",
"scalarstats",
"Compute statistics of scalars (mean, variance etc. per point) output mesh with stats as arrays.",
false);
122 SwitchArg removeScalarArg(
"",
"scalarremove",
"Remove the scalars.",
false);
123 SwitchArg copyScalarArg(
"",
"scalarcopy",
"Copy the Scalars from first mesh to all others while removing existing ones.",
false);
124 SwitchArg mseArg(
"",
"mse",
"Mean Squared Error of Points in models.",
false);
125 SwitchArg procrustesArg(
"",
"procrustes",
"Similarity alignment of surfaces assuming points have correspondence. Use --rigid if rigid alignment is required.",
false);
126 SwitchArg rigidArg(
"",
"rigid",
"Rigid alignment of surfaces assuming points have correspondence. Use with procrustes or other registration arguments.",
false);
127 SwitchArg icpArg(
"",
"icp",
"Iterative Closest Points alignment of surfaces assuming points don't have correspondence. Last surface in list is used as the reference 'fixed' surface.",
false);
128 SwitchArg saveTransformsArg(
"",
"savetransforms",
"Save the transformation matrix after surface alignment. For use with --icp",
false);
131 UnlabeledMultiArg<std::string> multinames(
"surfaces",
"Surfaces to operate on",
true,
"Surfaces");
134 ValueArg<float> voxeliseArg(
"",
"voxelise",
"Voxelise a model into image space with given spacing. Output will be in Nifty format (*.nii.gz).",
false, 0.5,
"Voxelise");
135 ValueArg<std::string> orientArg(
"",
"orient",
"Apply direction/orientation to a model from the given image",
false,
"image.nii.gz",
"Orient");
136 ValueArg<size_t> flatArg(
"",
"flatten",
"Flatten a model to either a sphere (0) or plane (1). Needs to be a genus zero (no topological hole) mesh.",
false, 0,
"Flatten");
139 cmd.add( threadsArg );
140 cmd.add( multinames );
142 cmd.add( outputArg );
143 cmd.add( prefixArg );
144 cmd.add( outputFormatArg );
145 cmd.add( componentsArg );
147 cmd.add( partitionArg );
148 cmd.add( saveTransformsArg );
152 std::vector<Arg*> xorlist;
153 xorlist.push_back(&duplicateArg);
154 xorlist.push_back(&convertArg);
155 xorlist.push_back(&concatenateArg);
156 xorlist.push_back(&colourConcatenateArg);
157 xorlist.push_back(&splitArg);
158 xorlist.push_back(&scaleArg);
159 xorlist.push_back(&decimateArg);
160 xorlist.push_back(&smoothArg);
161 xorlist.push_back(&laplacianArg);
162 xorlist.push_back(&diffScalarArg);
163 xorlist.push_back(&statsScalarArg);
164 xorlist.push_back(&removeScalarArg);
165 xorlist.push_back(©ScalarArg);
166 xorlist.push_back(&mseArg);
167 xorlist.push_back(&procrustesArg);
168 xorlist.push_back(&icpArg);
169 xorlist.push_back(&thresholdAboveArg);
170 xorlist.push_back(&thresholdBelowArg);
171 xorlist.push_back(&clipArg);
172 xorlist.push_back(&flipArg);
174 xorlist.push_back(&voxeliseArg);
175 xorlist.push_back(&orientArg);
176 #ifdef ITK_USE_REVIEW //Review only members 177 xorlist.push_back(&flatArg);
182 cmd.parse( argc, argv );
185 const size_t threads = threadsArg.getValue();
187 std::vector<std::string> filenames = multinames.getValue();
189 std::string outputName = outputArg.getValue();
190 std::string outputFormatName = outputFormatArg.getValue();
191 const std::string prefixName = prefixArg.getValue();
192 const float decimateFactor = decimateArg.getValue();
193 const float scaleFactor = scaleArg.getValue();
194 const float smoothIterations = smoothArg.getValue();
195 const float laplacianIterations = laplacianArg.getValue();
196 float thresholdAbove = thresholdAboveArg.getValue();
197 float thresholdBelow = thresholdBelowArg.getValue();
198 float clipValue = clipArg.getValue();
199 size_t flipAxis = flipArg.getValue();
200 std::vector<std::string> componentNames = componentsArg.getValue();
203 vtkMultiThreader::SetGlobalDefaultNumberOfThreads(threads);
207 const float voxelSpacing = voxeliseArg.getValue();
208 std::string orientName = orientArg.getValue();
209 const size_t flatMode = flatArg.getValue();
212 operations operation = none;
213 bool componentsValid =
false;
214 if(duplicateArg.isSet())
218 operation = duplicate;
220 if(convertArg.isSet())
228 if(concatenateArg.isSet() || colourConcatenateArg.isSet())
240 if(componentsArg.isSet())
242 if(!splitArg.isSet())
244 milx::PrintError(
"Components Argument Error: One or more of the operations");
250 for(stringiterator name = componentNames.begin(); name != componentNames.end(); name ++)
251 std::cout << *name <<
", ";
252 std::cout << std::endl;
257 if(!componentsArg.isSet())
259 milx::PrintError(
"Split Argument Error: Components of the given surfaces must be set.");
261 milx::PrintError(
"Set components via the component argument as many times as needed.");
266 componentsValid =
true;
271 milx::PrintInfo(
"Scaling coordinates of each point in the surfaces: ");
274 if(smoothArg.isSet())
280 if(laplacianArg.isSet())
284 operation = laplacian;
286 if(decimateArg.isSet())
290 operation = decimate;
292 if(thresholdAboveArg.isSet() || thresholdBelowArg.isSet())
294 if(!prefixArg.isSet())
296 milx::PrintError(
"Threshold Argument Error: Output Prefix (-p) must be provided.");
300 if(thresholdAboveArg.isSet() && !thresholdBelowArg.isSet())
302 thresholdBelow = -std::numeric_limits<float>::max();
304 else if(!thresholdAboveArg.isSet() && thresholdBelowArg.isSet())
306 thresholdAbove = std::numeric_limits<float>::max();
308 std::cout <<
"Threshold Above Value: " << thresholdAbove << std::endl;
309 std::cout <<
"Threshold Below Value: " << thresholdBelow << std::endl;
312 operation = thresholdscalars;
320 if(diffScalarArg.isSet())
324 operation = diffscalars;
332 if(statsScalarArg.isSet())
336 operation = statscalars;
338 if(removeScalarArg.isSet())
342 operation = removescalars;
344 if(copyScalarArg.isSet())
346 if(filenames.size() == 1)
348 milx::PrintError(
"Argument Error: Need another filename to copy from mesh A to mesh B.");
351 else if(filenames.size() == 2 && !outputArg.isSet())
353 milx::PrintError(
"Argument Error: Use Output (-o) for copying from mesh A to mesh B.");
357 else if(filenames.size() > 2 && !prefixArg.isSet())
359 milx::PrintError(
"Argument Error: Output Prefix (-p) must be provided for more than 2 meshes.");
366 operation = copyscalars;
368 if(rigidArg.isSet() && (!procrustesArg.isSet() && !icpArg.isSet()))
371 milx::PrintError(
"Rigid option can only be used with the Procrustes/ICP option");
374 if(procrustesArg.isSet())
376 if(!prefixArg.isSet())
378 milx::PrintError(
"Procrustes Argument Error: Output Prefix (-p) must be provided.");
384 operation = procrustes;
388 if(!prefixArg.isSet())
390 milx::PrintError(
"ICP Argument Error: Output Prefix (-p) must be provided.");
400 if(voxeliseArg.isSet())
402 if(!prefixArg.isSet())
404 milx::PrintError(
"Voxelise Argument Error: Output Prefix (-p) must be provided.");
408 if(!outputFormatArg.isSet())
411 outputFormatName =
"nii.gz";
415 operation = voxelise;
417 if(orientArg.isSet())
419 if(!prefixArg.isSet())
421 milx::PrintError(
"Orient Argument Error: Output Prefix (-p) must be provided.");
431 if(!prefixArg.isSet())
433 milx::PrintError(
"Orient Argument Error: Output Prefix (-p) must be provided.");
443 std::vector<std::string> altfilenames;
444 if(partitionArg.isSet())
446 for(
size_t j = filenames.size()/2; j < filenames.size(); j ++)
447 altfilenames.push_back( filenames[j] );
449 for(
size_t j = 0; j < altfilenames.size(); j ++)
450 filenames.pop_back();
453 std::cout <<
"Total Surfaces: " << filenames.size() << std::endl;
454 for(stringiterator name = filenames.begin(); name != filenames.end(); name ++)
455 std::cout << *name <<
", ";
456 std::cout << std::endl;
458 if(partitionArg.isSet())
460 if(!diffScalarArg.isSet())
462 milx::PrintError(
"Partitioning list of surfaces is only supported with the following operations:");
466 if(!prefixArg.isSet())
468 milx::PrintError(
"Partition Surfaces List Argument Error: Output Prefix (-p) must be provided.");
474 milx::PrintInfo(
"Operation will be applied in conjuction with surfaces: ");
475 std::cout <<
"Total Surfaces: " << altfilenames.size() << std::endl;
476 for(stringiterator name = altfilenames.begin(); name != altfilenames.end(); name ++)
477 std::cout << *name <<
", ";
478 std::cout << std::endl;
480 if(diffScalarArg.isSet())
481 operation = diffscalarspairs;
486 std::cerr <<
"Reading Surfaces... ";
487 vtkSmartPointer<vtkPolyDataCollection> collection;
490 std::cerr <<
"Done" << std::endl;
491 std::cout <<
"Read " << collection->GetNumberOfItems() <<
" models" << std::endl;
493 vtkSmartPointer<vtkPolyDataCollection> altcollection;
494 if(partitionArg.isSet())
496 std::cerr <<
"Reading Alterate Surfaces... ";
499 std::cerr <<
"Done" << std::endl;
502 vtkSmartPointer<vtkPolyDataCollection> components;
503 if(componentsArg.isSet() && componentsValid)
505 std::cerr <<
"Reading Components... ";
508 std::cerr <<
"Done" << std::endl;
512 itk::SmartPointer<charImageType> refImage;
513 if(orientArg.isSet())
515 std::cerr <<
"Reading Reference Image... ";
516 if( !milx::File::OpenImage<charImageType>(orientName, refImage) )
518 std::cerr <<
"Done" << std::endl;
521 if(colourConcatenateArg.isSet())
524 collection->InitTraversal();
525 for(
int j = 0; j < collection->GetNumberOfItems(); j ++)
527 vtkSmartPointer<vtkPolyData> mesh = collection->GetNextItem();
528 const size_t numPoints = mesh->GetNumberOfPoints();
530 vtkSmartPointer<vtkFloatArray> weights = vtkSmartPointer<vtkFloatArray>::New();
531 weights->SetNumberOfComponents(1);
532 weights->SetNumberOfTuples(numPoints);
533 weights->FillComponent(0, count);
535 mesh->GetPointData()->SetScalars(weights);
542 vtkSmartPointer<vtkTransformCollection> transformCollection = vtkSmartPointer<vtkTransformCollection>::New();
547 bool outputRequired =
true;
548 bool multiOutputRequired =
false;
549 bool standardOperation =
false;
550 bool imageOutputRequired =
false;
551 vtkSmartPointer<vtkPolyDataCollection> resultCollections;
552 std::vector< itk::SmartPointer<charImageType> > imgCollection;
554 if (collection->GetNumberOfItems() == 1)
556 outputRequired =
true;
557 multiOutputRequired =
false;
561 outputRequired =
false;
562 multiOutputRequired =
true;
565 std::cerr <<
"Applying... ";
569 standardOperation =
true;
573 standardOperation =
true;
578 outputRequired =
true;
579 multiOutputRequired =
false;
584 outputRequired =
false;
589 std::vector< vtkSmartPointer<vtkPolyDataCollection> > splitCollections;
595 typedef std::vector< vtkSmartPointer<vtkPolyDataCollection> >::iterator iterator;
596 stringiterator name = filenames.begin();
597 for(iterator surfaces = splitCollections.begin();
598 surfaces != splitCollections.end();
599 surfaces ++, name ++)
601 std::vector<std::string> splitNames;
603 std::string splitName = prefixName + milx::File::GetBaseName(*name) +
"_";
605 for(stringiterator componentname = componentNames.begin(); componentname != componentNames.end(); componentname ++)
607 splitNames.push_back(splitName + milx::File::GetBaseName(*componentname) +
"." + ext);
608 milx::PrintInfo(
"Wrote " + splitName + milx::File::GetBaseName(*componentname) +
"." + ext);
615 outputRequired =
false;
616 multiOutputRequired =
false;
622 standardOperation =
true;
628 standardOperation =
true;
634 standardOperation =
true;
640 standardOperation =
true;
643 case thresholdscalars:
646 standardOperation =
true;
652 standardOperation =
true;
657 bool xAxis =
false, yAxis =
false, zAxis =
false;
660 else if(flipAxis == 1)
666 standardOperation =
true;
674 case diffscalarspairs:
675 resultCollections = vtkSmartPointer<vtkPolyDataCollection>::New();
678 standardOperation =
true;
679 collection = resultCollections;
689 standardOperation =
true;
695 if(filenames.size() == 2)
697 outputRequired =
true;
698 multiOutputRequired =
false;
699 collection->InitTraversal();
700 collection->GetNextItem();
701 Model.
SetInput(collection->GetNextItem());
704 standardOperation =
true;
712 outputName = prefixName +
"_mean." + ext;
713 outputRequired =
true;
714 multiOutputRequired =
true;
721 standardOperation =
true;
728 outputRequired =
false;
729 multiOutputRequired =
true;
730 imageOutputRequired =
true;
736 outputRequired =
false;
737 multiOutputRequired =
true;
739 #ifdef ITK_USE_REVIEW //Review only members 741 Model.FlattenCollection(collection, flatMode);
743 outputRequired =
false;
744 multiOutputRequired =
true;
751 std::cerr <<
"Done" << std::endl;
753 if (collection->GetNumberOfItems() == 1 && standardOperation)
755 collection->InitTraversal();
756 Model.
Result() = collection->GetNextItem();
766 if(multiOutputRequired)
768 if(collection->GetNumberOfItems() > 0)
770 const int N = collection->GetNumberOfItems();
773 std::vector<std::string> names;
776 for(stringiterator filename = filenames.begin(); filename != filenames.end() && n < N; filename ++, n ++)
779 if (outputFormatArg.isSet() || imageOutputRequired) {
780 ext = outputFormatName;
785 names.push_back(prefixName + milx::File::GetBaseName(*filename) +
"." + ext);
786 milx::PrintInfo(
"Will be writing " + prefixName + milx::File::GetBaseName(*filename) +
"." + ext);
789 if(imageOutputRequired)
790 milx::File::SaveImages<charImageType>(names, imgCollection);
795 if (saveTransformsArg.isSet() && transformCollection->GetNumberOfItems() == collection->GetNumberOfItems())
797 const std::string ext =
"trsf.gz";
798 std::vector<std::string> tnames;
800 for(stringiterator name = names.begin(); name != names.end() && n < N; name ++, n++)
802 tnames.push_back(milx::File::StripFileExtension(*name) +
"." + ext);
void ScalarStatisticsCollection(vtkSmartPointer< vtkPolyDataCollection > collection)
Computes the statistics (min, max, variance etc.) of the scalars over the collection assuming the mes...
void ConcatenateCollection(vtkSmartPointer< vtkPolyDataCollection > collection)
Concatenates all the models in the collection together into one model.
void ScaleCollection(vtkSmartPointer< vtkPolyDataCollection > collection, const coordinateType scale)
Scales the coordinates of each point in all the models in the collection.
void PrintError(const std::string msg)
Displays a generic msg to standard error with carriage return.
void DecimateCollection(vtkSmartPointer< vtkPolyDataCollection > collection, const coordinateType factor)
Decimates all points in all the models in the collection using the Quadric algorithm algorithm...
This program Animates image slices and animates surfaces in the model view for creating movies...
static bool OpenModelCollection(std::vector< std::string > filenames, vtkSmartPointer< vtkPolyDataCollection > &collection)
Opens model files, which can be a VTK XML, Legacy VTK PolyData File (i.e. either a *...
void SplitCollection(vtkSmartPointer< vtkPolyDataCollection > collection, vtkSmartPointer< vtkPolyDataCollection > components, std::vector< vtkSmartPointer< vtkPolyDataCollection > > &splitCollections)
Splits each of the models in the collection into a collection of models.
vtkSmartPointer< vtkPolyData > & Result()
Returns the current model, i.e. the result of the latest operation.
vtkSmartPointer< vtkPolyDataCollection > ProcrustesAlignCollection(vtkSmartPointer< vtkPolyDataCollection > collection, bool rigid=false)
Aligns the collection to the mean mesh (computed internally) of the collection assuming the meshes ha...
void SmoothCollection(vtkSmartPointer< vtkPolyDataCollection > collection, const size_t iterations)
Smoothes all points in all the models in the collection using the Windowed Sinc algorithm.
static bool SaveTransformCollection(std::vector< std::string > filenames, vtkSmartPointer< vtkTransformCollection > collection, const bool ITK=false)
Saves transform files, which can be VTK .trsf or ITK from a collection.
static bool SaveModelCollection(std::vector< std::string > filenames, vtkSmartPointer< vtkPolyDataCollection > collection, const bool binary=false)
Saves model files, which can be a VTK XML, Legacy VTK PolyData File (i.e. either a *...
void SetInput(vtkSmartPointer< vtkPolyData > model)
Assigns the input model to class.
void ScalarThresholdCollection(vtkSmartPointer< vtkPolyDataCollection > collection, const coordinateType aboveVal, const coordinateType belowVal)
Computes the Threshold of scalar values of meshes in place over the collection wrt the first mesh ass...
void PrintWarning(const std::string msg)
Displays a generic msg to standard output with carriage return.
static bool SaveModel(const std::string filename, vtkSmartPointer< vtkPolyData > data, const bool binary=false)
Saves a model as a file, which can be a VTK XML or Legacy VTK PolyData File (i.e. either a *...
void PrintInfo(const std::string msg)
Displays a generic msg to standard output with carriage return.
std::string NumberToString(double num, unsigned zeroPad=0)
Number to string converter.
double MeanSquaredErrorCollection(vtkSmartPointer< vtkPolyDataCollection > collection)
Computes the MSE of all models in collection accummulatively.
void ScalarCopyCollection(vtkSmartPointer< vtkPolyDataCollection > collection)
Copies the scalars over the collection (from the first mesh) assuming the meshes have the same number...
void LaplacianCollection(vtkSmartPointer< vtkPolyDataCollection > collection, const size_t iterations)
Smoothes all points in all the models in the collection using the Laplacian algorithm.
void ScalarRemoveCollection(vtkSmartPointer< vtkPolyDataCollection > collection)
Removes the scalars over the collection assuming the meshes have the same number of points...
void FlipCollection(vtkSmartPointer< vtkPolyDataCollection > collection, const bool xAxis, const bool yAxis, const bool zAxis)
Flips all points in all the models in the collection using the axis provided.
static std::string GetFileExtension(const std::string &filename)
Returns the file extension (in lower case) of the given filename.
unsigned NumberOfProcessors()
Number of processors or cores on the machine.
void ClipCollection(vtkSmartPointer< vtkPolyDataCollection > collection, const coordinateType aboveVal, const coordinateType belowVal)
Computes the Clipping of models via the Threshold of scalar values in place over the collection...
void ScalarDifferenceCollection(vtkSmartPointer< vtkPolyDataCollection > collection)
Computes the difference in scalars cummulatively over the collection wrt the first mesh assuming the ...
vtkSmartPointer< vtkPolyDataCollection > IterativeClosestPointsAlignCollection(vtkSmartPointer< vtkPolyDataCollection > collection, bool rigid=false, vtkSmartPointer< vtkTransformCollection > tCollection=0)
Aligns the collection to the mean mesh (computed internally) of the collection without needing corres...