network.h 11 KB
#ifndef STIM_NETWORK_H
#define STIM_NETWORK_H

#include <stdlib.h>
#include <assert.h>
#include <sstream>
#include <fstream>
#include <algorithm>
#include <string.h>
#include <math.h>
#include <stim/math/vec3.h>
#include <stim/visualization/obj.h>
#include <stim/visualization/cylinder.h>
#include <ANN/ANN.h>
#include <boost/tuple/tuple.hpp>


namespace stim{
/** This is the a class that interfaces with gl_spider in order to store the currently
 *   segmented network. The following data is stored and can be extracted:
 *   1)Network geometry and centerline.
 *   2)Network connectivity (a graph of nodes and edges), reconstructed using ANN library.
*/


template<typename T>
class network{

	///Each edge is a fiber with two nodes.
	///Each node is an in index to the endpoint of the fiber in the nodes array.
	class edge : public cylinder<T>
	{
		public:
		unsigned v[2];		//unique id's designating the starting and ending
		// default constructor
		 edge() : cylinder<T>(){v[1] = -1; v[0] = -1;}
		/// Constructor - creates an edge from a list of points by calling the stim::fiber constructor

		///@param p is an array of positions in space
		edge(std::vector< stim::vec3<T> > p) : cylinder<T>(p){}

		/// Copy constructor creates an edge from a fiber
		edge(stim::cylinder<T> f) : cylinder<T>(f) {}

		/// Resamples an edge by calling the fiber resampling function
		edge resample(T spacing){
			edge e(cylinder<T>::resample(spacing));	//call the fiber->edge constructor
			e.v[0] = v[0];					//copy the vertex data
			e.v[1] = v[1];

			return e;						//return the new edge
		}

		/// Output the edge information as a string
		std::string str(){
			std::stringstream ss;
			ss<<"("<<cylinder<T>::size()<<")\tl = "<<length()<<"\t"<<v[0]<<"----"<<v[1];
			return ss.str();
		}

	};

	///Node class that stores the physical position of the node as well as the edges it is connected to (edges that connect to it), As well as any additional data necessary.
	class vertex : public stim::vec3<T>
	{
		public:
			//std::vector<unsigned int> edges;		//indices of edges connected to this node.
			std::vector<unsigned int> e[2];			//indices of edges going out (e[0]) and coming in (e[1])
			//stim::vec3<T> p;						//position of this node in physical space.

			//constructor takes a stim::vec
			vertex(stim::vec3<T> p) : stim::vec3<T>(p){}

			/// Output the vertex information as a string
			std::string	str(){
				std::stringstream ss;
				ss<<"\t(x, y, z) = "<<stim::vec3<T>::str();

				if(e[0].size() > 0){
					ss<<"\t> ";
					for(unsigned int o = 0; o < e[0].size(); o++)
						ss<<e[0][o]<<" ";
				}
				if(e[1].size() > 0){
					ss<<"\t< ";
					for(unsigned int i = 0; i < e[1].size(); i++)
						ss<<e[1][i]<<" ";
				}

				return ss.str();
			}

	};

protected:

	std::vector<edge> E;       //list of edges
	std::vector<vertex> V;	    //list of vertices.

public:

	///Returns the number of edges in the network.
	unsigned int edges(){
		return E.size();
	}

	///Returns the number of nodes in the network.
	unsigned int vertices(){
		return V.size();
	}

	stim::cylinder<T> get_cylinder(unsigned f){
		return E[f];					//return the specified edge (casting it to a fiber)
	}

	//load a network from an OBJ file
	void load_obj(std::string filename){

		stim::obj<T> O;			//create an OBJ object
		O.load(filename);		//load the OBJ file as an object

		std::vector<unsigned> id2vert;	//this list stores the OBJ vertex ID associated with each network vertex

		unsigned i[2];					//temporary, IDs associated with the first and last points in an OBJ line

		//for each line in the OBJ object
		for(unsigned int l = 1; l <= O.numL(); l++){

			std::vector< stim::vec<T> > c;				//allocate an array of points for the vessel centerline
			O.getLine(l, c);							//get the fiber centerline

			std::vector< stim::vec3<T> > c3(c.size());
			for(size_t j = 0; j < c.size(); j++)
				c3[j] = c[j];

			edge new_edge = c3;							//create an edge from the given centerline
			unsigned int I = new_edge.size();			//calculate the number of points on the centerline

			//get the first and last vertex IDs for the line
			std::vector< unsigned > id;					//create an array to store the centerline point IDs
			O.getLinei(l, id);							//get the list of point IDs for the line
			i[0] = id.front();							//get the OBJ ID for the first element of the line
			i[1] = id.back();							//get the OBJ ID for the last element of the line

			std::vector<unsigned>::iterator it;			//create an iterator for searching the id2vert array
			unsigned it_idx;							//create an integer for the id2vert entry index

			//find out if the nodes for this fiber have already been created
			it = find(id2vert.begin(), id2vert.end(), i[0]);	//look for the first node
			it_idx = std::distance(id2vert.begin(), it);
			if(it == id2vert.end()){							//if i[0] hasn't already been used
				vertex new_vertex = new_edge[0];				//create a new vertex, assign it a position
				new_vertex.e[0].push_back(E.size());			//add the current edge as outgoing
				new_edge.v[0] = V.size();						//add the new vertex to the edge
				V.push_back(new_vertex);						//add the new vertex to the vertex list
				id2vert.push_back(i[0]);						//add the ID to the ID->vertex conversion list
			}
			else{												//if the vertex already exists
				V[it_idx].e[0].push_back(E.size());				//add the current edge as outgoing
				new_edge.v[0] = it_idx;
			}

			it = find(id2vert.begin(), id2vert.end(), i[1]);	//look for the second ID
			it_idx = std::distance(id2vert.begin(), it);
			if(it == id2vert.end()){							//if i[1] hasn't already been used
				vertex new_vertex = new_edge[I-1];							//create a new vertex, assign it a position
				new_vertex.e[1].push_back(E.size());						//add the current edge as incoming
				new_edge.v[1] = V.size();
				V.push_back(new_vertex);									//add the new vertex to the vertex list
				id2vert.push_back(i[1]);						//add the ID to the ID->vertex conversion list
			}
			else{												//if the vertex already exists
				V[it_idx].e[1].push_back(E.size());				//add the current edge as incoming
				new_edge.v[1] = it_idx;
			}

			E.push_back(new_edge);								//push the edge to the list

		}
	}

	/// Output the network as a string
	std::string str(){

		std::stringstream ss;
		ss<<"Nodes ("<<V.size()<<")--------"<<std::endl;
		for(unsigned int v = 0; v < V.size(); v++){
			ss<<"\t"<<v<<V[v].str()<<std::endl;
		}

		ss<<"Edges ("<<E.size()<<")--------"<<std::endl;
		for(unsigned e = 0; e < E.size(); e++){
			ss<<"\t"<<e<<E[e].str()<<std::endl;
		}

		return ss.str();
	}

	/// This function resamples all fibers in a network given a desired minimum spacing
	/// @param spacing is the minimum distance between two points on the network
	stim::network<T> resample(T spacing){
		stim::network<T> n;					//create a new network that will be an exact copy, with resampled fibers
		n.V = V;							//copy all vertices

		n.E.resize(edges());				//allocate space for the edge list

		//copy all fibers, resampling them in the process
		for(unsigned e = 0; e < edges(); e++){				//for each edge in the edge list
			n.E[e] = E[e].resample(spacing);				//resample the edge and copy it to the new network
		}

		return n;							              //return the resampled network
	}



	/// Calculate the total number of points on all edges.
	unsigned total_points(){
		unsigned n = 0;
		for(unsigned e = 0; e < E.size(); e++)
			n += E[e].size();
		return n;
	}

	// gaussian function
	float gaussianFunction(float x, float std=25){ return exp(-x/(2*std*std));} // by default std = 25

    // stim 3d vector to annpoint of 3 dimensions
	void stim2ann(ANNpoint &a, stim::vec3<T> b){
		a[0] = b[0];
		a[1] = b[1];
		a[2] = b[2];
	}

	/// Calculate the average magnitude across the entire network.
	/// @param m is the magnitude value to use. The default is 0 (usually radius).
	T average(unsigned m = 0){

		T M, L;												//allocate space for the total magnitude and length
		M = L = 0;											//initialize both the initial magnitude and length to zero
		for(unsigned e = 0; e < E.size(); e++){				//for each edge in the network
			M += E[e].integrate(m);							//get the integrated magnitude
			L += E[e].length();								//get the edge length
		}

		return M / L;					//return the average magnitude
	}

	/// This function compares two networks and returns the percentage of the current network that is missing from A.

	/// @param A is the network to compare to - the field is generated for A
	/// @param sigma is the user-defined tolerance value - smaller values provide a stricter comparison
	stim::network<T> compare(stim::network<T> A, float sigma){

		stim::network<T> R;								//generate a network storing the result of the comparison
		R = (*this);									//initialize the result with the current network

		//generate a KD-tree for network A
		float metric = 0.0;                               // initialize metric to be returned after comparing the networks
		ANNkd_tree* kdt;                                 // initialize a pointer to a kd tree
		double **c;						                 // centerline (array of double pointers) - points on kdtree must be double
		unsigned int n_data = A.total_points();          // set the number of points
		c = (double**) malloc(sizeof(double*) * n_data); // allocate the array pointer
		for(unsigned int i = 0; i < n_data; i++)		 // allocate space for each point of 3 dimensions
			c[i] = (double*) malloc(sizeof(double) * 3);

		unsigned t = 0;
		for(unsigned e = 0; e < A.E.size(); e++){					//for each edge in the network
			for(unsigned p = 0; p < A.E[e].size(); p++){			//for each point in the edge
				for(unsigned d = 0; d < 3; d++){				//for each coordinate

					c[t][d] = A.E[e][p][d];
				}
				t++;
			}
		}

		//compare each point in the current network to the field produced by A
		ANNpointArray pts = (ANNpointArray)c;           // create an array of data points of type double
		kdt = new ANNkd_tree(pts, n_data, 3);			// build a KD tree using the annpointarray
		double eps = 0; // error bound
		ANNdistArray dists = new ANNdist[1];     // near neighbor distances
		ANNidxArray nnIdx = new ANNidx[1];				// near neighbor indices // allocate near neigh indices

		stim::vec3<T> p0, p1;
		float m1;
		float M = 0;											//stores the total metric value
		float L = 0;											//stores the total network length
		ANNpoint queryPt = annAllocPt(3);
		for(unsigned e = 0; e < R.E.size(); e++){					//for each edge in A
			R.E[e].add_mag(0);									//add a new magnitude for the metric

			for(unsigned p = 0; p < R.E[e].size(); p++){			//for each point in the edge

				p1 = R.E[e][p];									//get the next point in the edge
				stim2ann(queryPt, p1);
				kdt->annkSearch( queryPt, 1, nnIdx, dists, eps);	//find the distance between A and the current network
				m1 = 1.0f - gaussianFunction((float)dists[0], sigma);		//calculate the metric value based on the distance
				R.E[e].set_mag(m1, p, 1);						//set the error for the second point in the segment

			}
		}

		return R;		//return the resulting network
	}

	/// Returns the number of magnitude values stored in each edge. This should be uniform across the network.
	unsigned nmags(){
		return E[0].nmags();
	}


};		//end stim::network class
};		//end stim namespace
#endif