bool extract_boundary::run(viennamesh::algorithm_handle &)
  {
    mesh_handle input_mesh = get_required_input<mesh_handle>("mesh");
    mesh_handle output_mesh = make_data<mesh_handle>();

    point_container hole_points;
    seed_point_container seed_points;

    std::cerr << "Extract boundary of " << input_mesh.size() << " partitions" << std::endl;
    //MY IMPLEMENTATION
    output_mesh.resize(input_mesh.size());
    for (size_t i = 0; i < input_mesh.size(); ++i)
    {

      hole_points.clear();
      seed_points.clear();

      input_mesh(i).set_full_layout();
      viennagrid::extract_boundary( input_mesh(i), output_mesh(i), viennagrid::facet_dimension(input_mesh(i)) );
      viennagrid::extract_seed_points( input_mesh(i), seed_points );
  //     viennagrid::extract_hole_points( input_mesh(), hole_points );

     // set_output( "mesh", output_mesh ); //MY IMPLEMENTATION

  //     if (!hole_points.empty())
  //     {
  //       info(1) << "Extracted " << hole_points.size() << " hole points" << std::endl;
  //       for (point_container::const_iterator it = hole_points.begin(); it != hole_points.end(); ++it)
  //         info(1) << "   " << *it << std::endl;
  //
  //       point_handle output_hole_points = make_data<point>();
  //       output_hole_points.set( hole_points );
  //       set_output( "hole_points", output_hole_points );
  //     }

      if (!seed_points.empty())
      {
        info(1) << "Extracted " << seed_points.size() << " seed points" << std::endl;
        for (seed_point_container::const_iterator it = seed_points.begin(); it != seed_points.end(); ++it)
          info(1) << "   " << (*it).first << " -> " << (*it).second << std::endl;


        seed_point_handle output_seed_points = make_data<seed_point>();
        output_seed_points.set( seed_points );
        set_output( "seed_points", output_seed_points );
      }
    }//end of for over input_mesh.size()


    set_output( "mesh", output_mesh );
    return true;
  }
  bool extract_symmetric_slice_3d::run(viennamesh::algorithm_handle &)
  {
    mesh_handle input_mesh = get_required_input<mesh_handle>("mesh");
    int geometric_dimension = viennagrid::geometric_dimension( input_mesh() );

    if (geometric_dimension != 3)
      return false;

    double tol = 1e-6;
    if (get_input<double>("tolerance").valid())
      tol = get_input<double>("tolerance")();

    typedef viennagrid::mesh                                                MeshType;
    typedef viennagrid::result_of::point<MeshType>::type                    PointType;


    PointType axis = get_required_input<point>("axis")();
    axis.normalize();

    int rotational_frequency = get_required_input<int>("rotational_frequency")();
    double angle = 2*M_PI/rotational_frequency;

    PointType N[2];
    N[0] = viennagrid::make_point(0,1,0);
    if ( std::abs(viennagrid::inner_prod(axis,N[0])) > 1.0-tol )
      N[0] = viennagrid::make_point(-1,0,0);
    N[0] -= axis * viennagrid::inner_prod(axis,N[0]);
    N[0].normalize();

    N[1] = -rotate(N[0], point(geometric_dimension), axis, angle);


    info(1) << "Using rotational frequency " << rotational_frequency << std::endl;
    info(1) << "Angle = " << angle << std::endl;
    info(1) << "Axis = " << axis << std::endl;
    info(1) << "Normal[0] = " << N[0] << std::endl;
    info(1) << "Normal[1] = " << N[1] << std::endl;


    mesh_handle tmp = cut( input_mesh(), N[0], tol );
    mesh_handle output_mesh = cut( tmp(), N[1], tol );

    set_output( "mesh", output_mesh );

    return true;
  }
Example #3
0
static PyArrayObject *
_get_transform_mesh(PyObject *py_affine, npy_intp *dims)
{
    /* TODO: Could we get away with float, rather than double, arrays here? */

    /* Given a non-affine transform object, create a mesh that maps
    every pixel in the output image to the input image.  This is used
    as a lookup table during the actual resampling. */

    PyObject *py_inverse = NULL;
    npy_intp out_dims[3];

    out_dims[0] = dims[0] * dims[1];
    out_dims[1] = 2;

    py_inverse = PyObject_CallMethod(
        py_affine, (char *)"inverted", (char *)"", NULL);
    if (py_inverse == NULL) {
        return NULL;
    }

    numpy::array_view<double, 2> input_mesh(out_dims);
    double *p = (double *)input_mesh.data();

    for (npy_intp y = 0; y < dims[0]; ++y) {
        for (npy_intp x = 0; x < dims[1]; ++x) {
            *p++ = (double)x;
            *p++ = (double)y;
        }
    }

    PyObject *output_mesh =
        PyObject_CallMethod(
            py_inverse, (char *)"transform", (char *)"O",
            (char *)input_mesh.pyobj(), NULL);

    Py_DECREF(py_inverse);

    if (output_mesh == NULL) {
        return NULL;
    }

    PyArrayObject *output_mesh_array =
        (PyArrayObject *)PyArray_ContiguousFromAny(
            output_mesh, NPY_DOUBLE, 2, 2);

    Py_DECREF(output_mesh);

    if (output_mesh_array == NULL) {
        return NULL;
    }

    return output_mesh_array;
}
    bool make_mesh::run_impl()
    {
      output_parameter_proxy<netgen::mesh> omp(output_mesh);

      StdCaptureHandle capture_handle;

      if (omp != input_mesh)
        omp() = input_mesh();


      ::netgen::MeshingParameters mesh_parameters;

      if (cell_size.valid())
      {
        ::netgen::Point3d pmin;
        ::netgen::Point3d pmax;
        omp()().GetBox(pmin, pmax);

        double bb_volume = std::abs( (pmax.X() - pmin.X()) * (pmax.Y() - pmin.Y()) * (pmax.Z() - pmin.Z()) );
        double cell_volume = cell_size()*cell_size()*cell_size();
        double max_cell_count = 100000000;

        if ( cell_volume * max_cell_count < bb_volume )
        {
          warning(1) << "Cell size is too small and might result in too much elements" << std::endl;
          warning(1) << "Cell size                = " << cell_size() << std::endl;
          warning(1) << "Mesh max cell count      = " << max_cell_count << std::endl;
          warning(1) << "Mesh bounding box        = " << pmin << "," << pmax << std::endl;
          warning(1) << "Mesh bounding box volume = " << bb_volume << std::endl;
          warning(1) << "Mesh max cell volume     = " << cell_volume * max_cell_count << std::endl;
        }

        mesh_parameters.maxh = cell_size();
      }

      try
      {
        omp()().CalcLocalH(mesh_parameters.grading);
        MeshVolume (mesh_parameters, omp()());
        RemoveIllegalElements (omp()());
        OptimizeVolume (mesh_parameters, omp()());
      }
      catch (::netgen::NgException const & ex)
      {
        error(1) << "Netgen Error: " << ex.What() << std::endl;
        return false;
      }

      return true;
    }
  bool check_hull_topology::run(viennamesh::algorithm_handle &)
  {
    mesh_handle input_mesh = get_required_input<mesh_handle>("mesh");
    mesh_handle output_mesh = make_data<mesh_handle>();

    typedef viennagrid::mesh                                                            MeshType;

    typedef viennagrid::result_of::const_element_range<MeshType>::type                  ConstElementRangeType;
    typedef viennagrid::result_of::iterator<ConstElementRangeType>::type                ConstElementIteratorType;

    typedef viennagrid::result_of::const_coboundary_range<MeshType>::type               ConstCoboundaryElementRangeType;
    typedef viennagrid::result_of::iterator<ConstCoboundaryElementRangeType>::type      ConstCoboundaryElementIteratorType;

    viennagrid::result_of::element_copy_map<>::type copy_map(output_mesh(), true);


    ConstElementRangeType lines( input_mesh(), 1 );
    for (ConstElementIteratorType lit = lines.begin(); lit != lines.end(); ++lit)
    {
      ConstCoboundaryElementRangeType coboundary_triangles( input_mesh(), *lit, 2 );
      if (coboundary_triangles.size() < 2)
      {
        std::cout << "Line " << *lit << " has less than 2 co-boundary triangles" << std::endl;

        for (ConstCoboundaryElementIteratorType ctit = coboundary_triangles.begin(); ctit != coboundary_triangles.end(); ++ctit)
        {
          copy_map(*ctit);
        }
      }

    }

    set_output( "mesh", output_mesh );

    return true;
  }
  bool uniform_refine::run(viennamesh::algorithm_handle &)
  {
    mesh_handle input_mesh = get_required_input<mesh_handle>("mesh");
    if (!input_mesh.valid())
      return false;

    mesh_handle output_mesh = make_data<mesh_handle>();

    if (output_mesh == input_mesh)
      return false;

    viennagrid::cell_refine_uniformly(input_mesh(), output_mesh());

    set_output( "mesh", output_mesh );

    return true;
  }
  bool map_regions::run(viennamesh::algorithm_handle &)
  {
    mesh_handle input_mesh = get_required_input<mesh_handle>("mesh");
    string_handle input_region_mapping = get_required_input<string_handle>("region_mapping");

    std::map<std::string, std::string> region_mapping;
    {
      std::string tmp = input_region_mapping();

      std::list<std::string> mappings;
      boost::algorithm::split( mappings, tmp, boost::is_any_of(";") );
//       std::list<std::string> mappings = stringtools::split_string( input_region_mapping(), ";" );

      for (std::list<std::string>::const_iterator sit = mappings.begin(); sit != mappings.end(); ++sit)
      {
        std::list<std::string> from_to;
        boost::algorithm::split( from_to, *sit, boost::is_any_of(",") );

//         std::list<std::string> from_to = stringtools::split_string( *sit, "," );
        if (from_to.size() != 2)
        {
          return false;
        }

        std::list<std::string>::const_iterator it = from_to.begin();
        std::string from = *it;
        ++it;
        std::string to = *it;

        region_mapping[from] = to;
      }
    }

    for (std::map<std::string, std::string>::const_iterator it = region_mapping.begin(); it != region_mapping.end(); ++it)
    {
      info(1) << "Mapping region " << it->first << " to " << it->second << std::endl;
    }

    mesh_handle output_mesh = make_data<mesh_handle>();

    map_regions_impl( input_mesh(), output_mesh(), region_mapping );

    set_output( "mesh", output_mesh );
    return true;
  }
  bool stretch_middle::run(viennamesh::algorithm_handle &)
  {
    mesh_handle input_mesh = get_required_input<mesh_handle>("mesh");

    data_handle<double> stretch = get_required_input<double>("stretch");
    data_handle<double> middle_tolerance = get_required_input<double>("middle_tolerance");

    mesh_handle output_mesh = make_data<mesh_handle>();

    if (output_mesh != input_mesh)
      viennagrid::copy( input_mesh(), output_mesh() );

    typedef viennagrid::mesh MeshType;
    typedef viennagrid::result_of::point<MeshType>::type PointType;
    typedef viennagrid::result_of::const_element_range<MeshType>::type ConstElementRangeType;
    typedef viennagrid::result_of::iterator<ConstElementRangeType>::type ConstElementIteratorType;

    ConstElementRangeType vertices( output_mesh(), 0 );
    for (ConstElementIteratorType vit = vertices.begin(); vit != vertices.end(); ++vit)
    {
      PointType pt = viennagrid::get_point(*vit);

      if (std::abs(pt[0]) >= middle_tolerance())
      {
        if (pt[0] < 0)
          pt[0] -= stretch();
        else
          pt[0] += stretch();
      }

      viennagrid::set_point(*vit, pt);
    }

    set_output( "mesh", output_mesh );

    return true;
  }
  bool symmetry_detection_3d::run(viennamesh::algorithm_handle &)
  {


//
//     std::cout << "dynamic: " << jacobi_polynom<double>(4,2,2) << std::endl;
//     std::cout << "static: " << static_jacobi_polynom<double,4>(2,2) << std::endl;
//
//     return true;






    mesh_handle input_mesh = get_required_input<mesh_handle>("mesh");
    int geometric_dimension = viennagrid::geometric_dimension( input_mesh() );
    int cell_dimension = viennagrid::cell_dimension( input_mesh() );

    data_handle<int> p = get_required_input<int>("p");
    data_handle<double> relative_integrate_tolerance = get_required_input<double>("relative_integrate_tolerance");
    data_handle<double> absolute_integrate_tolerance = get_required_input<double>("absolute_integrate_tolerance");
//     data_handle<int> max_iteration_count = get_required_input<int>("max_iteration_count");
    data_handle<double> mirror_symmetry_tolerance = get_required_input<double>("mirror_symmetry_tolerance");
    data_handle<double> rotational_symmetry_tolerance = get_required_input<double>("rotational_symmetry_tolerance");

    if (geometric_dimension != 3)
      return false;

    if (cell_dimension != 2)
      return false;



    typedef viennagrid::mesh    MeshType;
    typedef point               PointType;

    typedef viennagrid::result_of::const_vertex_range<MeshType>::type       ConstVertexRangeType;
    typedef viennagrid::result_of::iterator<ConstVertexRangeType>::type     ConstVertexRangeIterator;


    double max_size = 0.0;
    {
      ConstVertexRangeType vertices(input_mesh());
      for (ConstVertexRangeIterator vit = vertices.begin(); vit != vertices.end(); ++vit)
      {
        double cur_size = viennagrid::norm_2( viennagrid::get_point(*vit) );
        if (cur_size > max_size)
          max_size = cur_size;
      }
    }

    info(1) << "Before start" << std::endl;

    MeshType mesh;
    viennagrid::copy( input_mesh(), mesh );
    viennagrid::scale( mesh, 1.0/max_size );

    info(1) << "After copy/scale" << std::endl;


    viennautils::Timer timer;
    timer.start();
    RealGeneralizedMoment m_real(2*p(), mesh);
//     , relative_integrate_tolerance(), absolute_integrate_tolerance(), max_iteration_count());

    info(1) << "After calculating generalized moment (!!! took " << timer.get() << "sec !!!)" << std::endl;

    double sphere_radius = 1.0;
    if (get_input<double>("sphere_radius").valid())
      sphere_radius = get_input<double>("sphere_radius")();

    MeshType sphere;
    viennagrid::make_sphere_hull( sphere, viennagrid::make_point(0,0,0), sphere_radius, 4 );

    viennagrid::quantity_field gradient_field_real(0, 1);
    gradient_field_real.set_name("gradient_real");

    ConstVertexRangeType vertices(sphere);
    for (ConstVertexRangeIterator vit = vertices.begin(); vit != vertices.end(); ++vit)
    {
      PointType const & pt = viennagrid::get_point(*vit);

      double theta;
      double phi;
      double r;
      to_spherical(pt, theta, phi, r);

      double grad_real = m_real.grad(theta, phi, 1e-2);
      gradient_field_real.set(*vit, grad_real);
    }

//     {
//       int bench_count = 100000;
//       std::vector<double> v(bench_count);
//       viennamesh::LoggingStack s("bench");
//
//       for (int i = 0; i != bench_count; ++i)
//         v[i] = m_real.grad(i*0.1, i*0.2, 1e-2);
//     }

    info(1) << "After calculating sphere" << std::endl;

    set_output("sphere", sphere);
    set_output("mesh", mesh);

    quantity_field_handle quantities = make_data<viennagrid::quantity_field>();
    quantities.set(gradient_field_real);
    set_output("sphere_quantities", quantities);


//     m_real.print();
//     std::cout << std::endl;
//     std::cout << "m_real hast mirror symmetry: " << std::boolalpha << m_real.z_mirror_symmetry( mirror_symmetry_tolerance() ) << std::endl;
//     m_real.rotation_symmetry_angles();
//     // rotational_symmetry_tolerance() );
//     std::cout << std::endl;



    data_handle<viennamesh_point> rotation_vector = get_input<viennamesh_point>("rotation_vector");
    data_handle<int> rotational_frequencies = get_input<int>("rotational_frequencies");

    if (rotation_vector.valid())
    {
      for (int i = 0; i != rotation_vector.size(); ++i)
      {
        point new_z = rotation_vector(i);
        info(1) << "Using rotation vector " << new_z << std::endl;
        RealGeneralizedMoment rotated_m = m_real.get_rotated(new_z);

//         rotated_m.print();
//         std::cout << std::endl;
        info(1) << "rotated_m (z = "<< new_z << ") hast mirror symmetry: " << std::boolalpha << rotated_m.z_mirror_symmetry( mirror_symmetry_tolerance() ) << std::endl;
//         rotated_m.rotation_symmetry_angles();
//         rotated_m.rotation_symmetry_angles( rotational_symmetry_tolerance() );
        rotated_m.check_rotation_symmetry(M_PI);

        if (rotational_frequencies.valid())
        {
          for (int i = 0; i != rotational_frequencies.size(); ++i)
          {
            int rotational_frequency = rotational_frequencies(i);
            double angle = 2*M_PI/rotational_frequency;
            info(1) << "Using rotational frequency " << rotational_frequency << " (angle = " << angle << ") error = " << rotated_m.check_rotation_symmetry(angle) << std::endl;
          }
        }
      }
    }


    return true;
  }
Example #10
0
    bool make_statistic::run(viennamesh::algorithm_handle &)
    {
        info(1) << "make_statistic started:" << std::endl;

        mesh_handle input_mesh = get_required_input<mesh_handle>("mesh");

        /*original mesh for comparison statistics (hausdorff distance, curvature difference...)*/
        mesh_handle original_mesh = get_input<mesh_handle>("original_mesh");

        /*All provided cell quality metric types are implemented as header files, which are located in statistics/metrics.
         * Note that most metrics are only implemented for triangles*/
        data_handle<viennamesh_string> metric_type = get_required_input<viennamesh_string>("metric_type");

        /*Decision threshold for triangle qualitity classification
         * E.g., radius_ratio < 1.5
         * Whether '<' or '>' is used for comparison is deduced from metric_ordering_tag*/
        data_handle<viennagrid_numeric> good_element_threshold = get_input<viennagrid_numeric>("good_element_threshold");

        /* comprehensive mesh quality metric is defined as
             alpha * (1 - good_elements_counted/count) + beta * min_dist_rms + gamma * mean_curvature + delta * volume_deviation;
        */
        data_handle<viennagrid_numeric> alpha = get_input<viennagrid_numeric>("alpha");
        data_handle<viennagrid_numeric> beta = get_input<viennagrid_numeric>("beta");
        data_handle<viennagrid_numeric> gamma = get_input<viennagrid_numeric>("gamma");
        data_handle<viennagrid_numeric> delta = get_input<viennagrid_numeric>("delta");

#if HIST
        data_handle<viennagrid_numeric> histogram_bins = get_input<viennagrid_numeric>("histogram_bin");
        data_handle<viennagrid_numeric> histogram_min = get_input<viennagrid_numeric>("histogram_min");
        data_handle<viennagrid_numeric> histogram_max = get_input<viennagrid_numeric>("histogram_max");
        data_handle<int> histogram_bin_count = get_input<int>("histogram_bin_count");
#endif


        typedef viennagrid::mesh                                  MeshType;
        typedef viennagrid::result_of::element<MeshType>::type    ElementType; //=Triangle or Tetrahedron

        typedef viennamesh::statistic<viennagrid_numeric>         StatisticType;
        StatisticType statistic;

#if HIST
        typedef StatisticType::histogram_type                     HistogramType;


        if (histogram_bins.valid())
        {
            std::vector<viennagrid_numeric> bins;
            for (int i = 0; i != histogram_bins.size(); ++i)
                bins.push_back( histogram_bins(i) );

            statistic.set_histogram( StatisticType::histogram_type::make(bins.begin(), bins.end()) );
        }
        else if (histogram_min.valid() && histogram_max.valid() && histogram_bin_count.valid())
        {
            statistic.set_histogram( StatisticType::histogram_type::make_uniform(histogram_min(), histogram_max(), histogram_bin_count()) );
        }
        else
        {
            //If histograms are about to be added again, returing false here will prevent the basic statistics features from working!
            error(1) << "No histogram configuration provided" << std::endl;
            return false;
        }
#endif





        //Good element threshold input was set, call statistics with the appropriate metric
        if(good_element_threshold.valid())
        {
            viennamesh::LoggingStack stack( std::string("High quality cell counter with metric type \"") + metric_type() + "\"" );

            if (metric_type() == "aspect_ratio")
                statistic.cell_quality_count<viennamesh::aspect_ratio_tag>( input_mesh(), viennamesh::aspect_ratio<ElementType>, good_element_threshold());
            else if (metric_type() == "min_angle")
                statistic.cell_quality_count<viennamesh::min_angle_tag>( input_mesh(), viennamesh::min_angle<ElementType>, good_element_threshold());
            else if (metric_type() == "max_angle")
                statistic.cell_quality_count<viennamesh::max_angle_tag>( input_mesh(), viennamesh::max_angle<ElementType>, good_element_threshold());
            else if (metric_type() == "min_dihedral_angle")
                statistic.cell_quality_count<viennamesh::min_dihedral_angle_tag>( input_mesh(), viennamesh::min_dihedral_angle<ElementType>, good_element_threshold());
            else if (metric_type() == "radius_edge_ratio")
                statistic.cell_quality_count<viennamesh::radius_edge_ratio_tag>( input_mesh(), viennamesh::radius_edge_ratio<ElementType>, good_element_threshold());
            else if (metric_type() == "radius_ratio")
                statistic.cell_quality_count<viennamesh::radius_ratio_tag>( input_mesh(), viennamesh::radius_ratio<ElementType>, good_element_threshold());
            else if (metric_type() == "perimeter_inradius_ratio")
                statistic.cell_quality_count<viennamesh::perimeter_inradius_ratio_tag>( input_mesh(), viennamesh::perimeter_inradius_ratio<ElementType>, good_element_threshold());
            else if (metric_type() == "edge_ratio")
                statistic.cell_quality_count<viennamesh::edge_ratio_tag>( input_mesh(), viennamesh::edge_ratio<ElementType>, good_element_threshold());
            else if (metric_type() == "circum_perimeter_ratio")
                statistic.cell_quality_count<viennamesh::circum_perimeter_ratio_tag>( input_mesh(), viennamesh::circum_perimeter_ratio<ElementType>, good_element_threshold());
            else if (metric_type() == "stretch")
                statistic.cell_quality_count<viennamesh::stretch_tag>( input_mesh(), viennamesh::stretch<ElementType>, good_element_threshold());
            else if (metric_type() == "skewness")
                statistic.cell_quality_count<viennamesh::skewness_tag>( input_mesh(), viennamesh::skewness<ElementType>, good_element_threshold());
            else
            {
                error(1) << "Metric type \"" << metric_type() << "\" is not supported for cell quality classifaction" << std::endl;
                return false;
            }
        }
        else //no triangle classifiction takes place
        {
            viennamesh::LoggingStack stack( std::string("Cell statistics with metric type \"") + metric_type() + "\"" );

            if (metric_type() == "aspect_ratio")
                statistic.cell_stats( input_mesh(), viennamesh::aspect_ratio<ElementType> );
            else if (metric_type() == "min_angle")
                statistic.cell_stats( input_mesh(), viennamesh::min_angle<ElementType> );
            else if (metric_type() == "max_angle")
                statistic.cell_stats( input_mesh(), viennamesh::max_angle<ElementType> );
            else if (metric_type() == "min_dihedral_angle")
                statistic.cell_stats( input_mesh(), viennamesh::min_dihedral_angle<ElementType> );
            else if (metric_type() == "radius_edge_ratio")
                statistic.cell_stats( input_mesh(), viennamesh::radius_edge_ratio<ElementType> );
            else if (metric_type() == "radius_ratio")
                statistic.cell_stats( input_mesh(), viennamesh::radius_ratio<ElementType> );
            else if (metric_type() == "perimeter_inradius_ratio")
                statistic.cell_stats( input_mesh(), viennamesh::perimeter_inradius_ratio<ElementType> );
            else if (metric_type() == "edge_ratio")
                statistic.cell_stats( input_mesh(), viennamesh::edge_ratio<ElementType> );
            else if (metric_type() == "circum_perimeter_ratio")
                statistic.cell_stats( input_mesh(), viennamesh::circum_perimeter_ratio<ElementType> );
            else if (metric_type() == "stretch")
                statistic.cell_stats( input_mesh(), viennamesh::stretch<ElementType> );
            else if (metric_type() == "skewness")
                statistic.cell_stats( input_mesh(), viennamesh::skewness<ElementType> );
            else
            {
                error(1) << "Metric type \"" << metric_type() << "\" is not supported" << std::endl;
                return false;
            }
        }

        if(original_mesh.valid())//a second mesh is set, calculate mesh comparison measures
        {
            viennamesh::LoggingStack stack( std::string("Calculation of Mesh Comparison Measures") );

            statistic.mesh_comparison_quality(input_mesh(), original_mesh());



            ConstTriangleRange tr_orig(original_mesh());
            ConstTriangleRange tr(input_mesh());

            StatisticType statistic_orig;
            statistic_orig.cell_stats( original_mesh(), viennamesh::aspect_ratio<ElementType> );


            //use median of orig mesh for input mesh triangle shape characterization
            statistic.cell_quality_count<viennamesh::aspect_ratio_tag>( input_mesh(), viennamesh::aspect_ratio<ElementType>, statistic_orig.median());

            if(alpha.valid() && beta.valid() && gamma.valid() && delta.valid() )
            {
                statistic.set_mesh_quality_weights(alpha(), beta(), gamma(), delta());
                info(5) << "values for comprehensive mesh quality metric: alpha = " << alpha()
                        << ", beta = " << beta() << ", gamma = " << gamma() << ", delta = "<< delta() <<  std::endl;

            }
            else
            {
                info(5) << "default values for comprehensive mesh quality metric used: alpha = 0.25, beta = 20, gamma = 1.0, delta = 1.3" << std::endl;
            }

            set_output("minimum_distance_rms", statistic.min_dist_rms());
            set_output("mean_curvature_difference", statistic.mean_curvature());
            set_output("area_deviation", statistic.volume_deviation());
            set_output("triangle_shape", statistic.good_elements()/statistic.count());
            set_output("mesh_quality_metric", statistic.mesh_quality_metric());
        }


        info(5) << statistic << "\n";

#if HIST
        statistic.normalize();
        std::vector<viennagrid_numeric> bins;
        for (HistogramType::const_iterator bit = statistic.histogram().begin(); bit != statistic.histogram().end(); ++bit)
            bins.push_back( (*bit).second );
        bins.push_back( statistic.histogram().overflow_bin() );

        data_handle<viennagrid_numeric> output_bins = make_data<viennagrid_numeric>();
        output_bins.set( bins );
        set_output( "bins", output_bins );
#endif

        set_output( "min", statistic.min() );
        set_output( "max", statistic.max() );
        set_output( "mean", statistic.mean() );
        set_output( "median", statistic.median());

        return true;
    }
    bool chessboard_coloring::run(viennamesh::algorithm_handle &)
    {
        viennamesh::info(1) << name() << std::endl;

        //get input mesh and create output mesh
        mesh_handle input_mesh = get_required_input<mesh_handle>("mesh");
        mesh_handle output_mesh = make_data<mesh_handle>();

        //Typedefs
        typedef viennagrid::mesh                                                                  MeshType;
        // typedef viennagrid::result_of::element<MeshType>::type                                    VertexType;
        typedef viennagrid::result_of::element<MeshType>::type                                    EdgeType;
        typedef viennagrid::result_of::element<MeshType>::type                                    TriangleType;

        //typedef viennagrid::result_of::element_range<MeshType, 0>::type                           VertexRange;

        typedef viennagrid::result_of::element_range<MeshType, 2>::type                           TriangleRange; 
        typedef viennagrid::result_of::iterator<TriangleRange>::type                              TriangleIterator; 

        typedef viennagrid::result_of::neighbor_range<MeshType, 1, 2>::type                       TriangleNeighborRangeType;
        typedef viennagrid::result_of::iterator<TriangleNeighborRangeType>::type                  TriangleNeighborIterator;
                        
        //get number of vertices and triangles
        int num_vertices = viennagrid::vertex_count(input_mesh());
        int num_triangles = viennagrid::element_count(input_mesh(), 2);
/*
        viennamesh::info(2) << "Number of vertices in mesh : " << num_vertices << std::endl;
        viennamesh::info(2) << "Number of triangles in mesh: " << num_triangles << std::endl;
*/
        //set up accessor
        std::vector<bool> color_triangles(num_triangles, false);
        std::vector<bool> touched (num_triangles, false);

        viennagrid::result_of::accessor<std::vector<bool>, TriangleType>::type color_accessor(color_triangles);

        //set first element to color "black" (BOOL = TRUE)
        color_accessor.set(viennagrid::cells(input_mesh())[0], true);

        //iterate triangles in mesh
        TriangleRange triangles(input_mesh());

        //Iterate over all triangles in the mesh
        viennagrid_element_id * triangle_ids_begin;
        viennagrid_element_id * triangle_ids_end;
        viennagrid_dimension topological_dimension = viennagrid::cell_dimension( input_mesh() );	

        viennagrid_element_id * neighbor_begin;
        viennagrid_element_id * neighbor_end;

        viennagrid_dimension connector = 1;

        viennagrid_mesh_elements_get(input_mesh().internal(), topological_dimension, &triangle_ids_begin, &triangle_ids_end);	

        viennagrid::quantity_field color_field(2,1);
       // viennagrid_quantity_field_create(&color_field);
        color_field.set_name("color");

        //viennagrid_quantity_field_init(color_field, 2, VIENNAGRID_QUANTITY_FIELD_TYPE_NUMERIC, 1 , VIENNAGRID_QUANTITY_FIELD_STORAGE_DENSE);

        //get triangles from mesh

        for (viennagrid_element_id *tri = triangle_ids_begin; tri != triangle_ids_end; ++tri)
        {
            int tri_index = viennagrid_index_from_element_id( *tri );

            //std::cout << tri_index << std::endl;

            if ( !touched[tri_index])
            { 
                //color_accessor.set(viennagrid::cells(input_mesh())[tri_index], true);
                color_triangles[tri_index] = true;
                touched[tri_index] = true;

                double clr = 1;

                color_field.set(*tri, clr);

                viennagrid_element_neighbor_elements(input_mesh().internal(), *tri, 0, 2, &neighbor_begin, &neighbor_end);

                for (viennagrid_element_id *n_tri = neighbor_begin; n_tri != neighbor_end; ++n_tri)
                {
                    int n_tri_index = viennagrid_index_from_element_id( *n_tri );

                    //std::cout << "  " << n_tri_index << std::endl;

                    viennagrid_element_id * n_neighbor_begin;
                    viennagrid_element_id * n_neighbor_end;

                    viennagrid_element_neighbor_elements(input_mesh().internal(), *n_tri, 0, 2, &n_neighbor_begin, &n_neighbor_end);

                    for (viennagrid_element_id *n_n_tri = n_neighbor_begin; n_n_tri != n_neighbor_end; ++n_n_tri)
                    {
                        int n_n_tri_index = viennagrid_index_from_element_id( *n_n_tri );

                        if ( !touched[n_n_tri_index])
                        { 
                            //color_accessor.set(viennagrid::cells(input_mesh())[n_tri_index], false);
                            touched[n_n_tri_index] = true;

                            clr = 0;
                            color_field.set(*n_n_tri, clr);
                        }
                    }
                }
            }
        }

        output_mesh = input_mesh;
        quantity_field_handle quantities = make_data<viennagrid::quantity_field>();

        quantities.set(color_field);
   
        set_output("color_field", quantities);
        set_output("mesh", output_mesh());

        return true;
    } //end of bool chessboard_coloring::run(viennamesh::algorithm_handle &)
  bool merge_meshes::run(viennamesh::algorithm_handle &)
  {
    mesh_handle output_mesh = make_data<mesh_handle>();
    mesh_handle input_mesh = get_input<mesh_handle>("mesh");

    bool region_offset = true;
    if ( get_input<bool>("region_offset").valid() )
      region_offset = get_input<bool>("region_offset")();

    double tolerance = 1e-6;
    if ( get_input<double>("tolerance").valid() )
      tolerance = get_input<double>("tolerance")();

    info(1) << "Using region offset: " << std::boolalpha << region_offset << std::endl;

    int merged_count = 0;
    if (input_mesh.valid())
    {
      int mesh_count = input_mesh.size();
      for (int i = 0; i != mesh_count; ++i)
      {
        merge_meshes_impl( input_mesh(i), output_mesh(), (merged_count == 0 ? -1 : tolerance), region_offset );
        ++merged_count;
      }
    }

    int mesh_index = 0;
    while (true)
    {
      std::string input_name = "mesh" + boost::lexical_cast<std::string>(mesh_index);
      mesh_handle another_input_mesh = get_input<mesh_handle>(input_name);

      if (!another_input_mesh.valid())
        break;

      int mesh_count = another_input_mesh.size();
      for (int i = 0; i != mesh_count; ++i)
      {
        merge_meshes_impl( another_input_mesh(i), output_mesh(), (merged_count == 0 ? -1 : tolerance), region_offset );
        ++merged_count;
      }

      ++mesh_index;
    }

//



//     if (input_mesh)
//       merge_meshes_impl( input_mesh(), output_mesh() );
//
//     int index = 0;
//     input_mesh = get_input<mesh_handle>("mesh[" + lexical_cast<std::string>(index++) + "]");
//     while (input_mesh)
//     {
//       merge_meshes_impl( input_mesh(), output_mesh() );
//       input_mesh = get_input<mesh_handle>("mesh[" + lexical_cast<std::string>(index++) + "]");
//     }

    info(1) << "Merged " << merged_count << " meshes" << std::endl;

    set_output( "mesh", output_mesh );

    return true;
  }
  bool center_mesh::run(viennamesh::algorithm_handle &)
  {

    mesh_handle input_mesh = get_input<mesh_handle>("mesh");
    if (input_mesh.valid())
    {
      mesh_handle output_mesh = make_data<mesh_handle>();

      typedef viennagrid::mesh                                                MeshType;
      typedef viennagrid::result_of::point<MeshType>::type                    PointType;

      typedef viennagrid::result_of::const_vertex_range<MeshType>::type       ConstVertexRangeType;
      typedef viennagrid::result_of::iterator<ConstVertexRangeType>::type     ConstVertexRangeIterator;


      viennagrid::copy( input_mesh(), output_mesh() );
      PointType mesh_centroid = viennagrid::centroid( output_mesh() );

      info(1) << "Mesh centroid " << mesh_centroid << std::endl;

      ConstVertexRangeType vertices( output_mesh() );
      for (ConstVertexRangeIterator vit = vertices.begin(); vit != vertices.end(); ++vit)
        viennagrid::set_point( *vit, viennagrid::get_point(*vit) - mesh_centroid );


      set_output( "mesh", output_mesh );
    }


    data_handle<viennagrid_plc> input_geometry = get_input<viennagrid_plc>("geometry");
    if (input_geometry.valid())
    {
      data_handle<viennagrid_plc> output_geometry = make_data<viennagrid_plc>();

      viennagrid_dimension geometric_dimension;
      viennagrid_plc_geometric_dimension_get( input_geometry(), &geometric_dimension );

      viennagrid_numeric * src_coords;
      viennagrid_plc_vertex_coords_pointer( input_geometry(), &src_coords );

      viennagrid_int vertex_count;
      viennagrid_plc_element_count( input_geometry(), 0, &vertex_count );

      std::vector<viennagrid_numeric> center( geometric_dimension, 0 );

      for (viennagrid_int i = 0; i != vertex_count; ++i)
      {
        for (viennagrid_dimension d = 0; d != geometric_dimension; ++d)
          center[d] += src_coords[ geometric_dimension*i + d ];
      }

      for (viennagrid_dimension d = 0; d != geometric_dimension; ++d)
        center[d] /= vertex_count;


      viennagrid_plc_copy( input_geometry(), output_geometry() );

      viennagrid_numeric * dst_coords;
      viennagrid_plc_vertex_coords_pointer( output_geometry(), &dst_coords );

      for (viennagrid_int i = 0; i != vertex_count; ++i)
      {
        for (viennagrid_dimension d = 0; d != geometric_dimension; ++d)
          dst_coords[ geometric_dimension*i + d ] -= center[d];
      }


      set_output( "geometry", output_geometry );
    }

    return true;
  }