The source file meshless_interpolation_function.h with comments:

        #ifndef LIBMESH_MESHLESS_INTERPOLATION_FUNCTION_H
        #define LIBMESH_MESHLESS_INTERPOLATION_FUNCTION_H
        
Local Includes
        #include "libmesh/function_base.h"
        #include "libmesh/meshfree_interpolation.h"
        #include "libmesh/threads.h"
        
        
C++ includes
        #include <cstddef>
        
        namespace libMesh
        {
        
        
        
Forward Declarations
        template <typename T>
        class DenseVector;
        
        
------------------------------------------------------------ MeshlessInterpolationFunction class definition
        class MeshlessInterpolationFunction : public FunctionBase<Number>
        {
        private:
          const MeshfreeInterpolation &_mfi;
          mutable std::vector<Point> _pts;
          mutable std::vector<Number> _vals;
          Threads::spin_mutex &_mutex;
          
        public:
        
          /**
           * Constructor.  Requires a \p \pMeshlessInterpolation object.
           */
          MeshlessInterpolationFunction (const MeshfreeInterpolation &mfi,
        				 Threads::spin_mutex &mutex) :
          _mfi (mfi),
          _mutex(mutex)
          {}
          
        
          /**
           * The actual initialization process.
           */
          void init ();
        
          /**
           * Clears the function.
           */
          void clear ();
        
          /**
           * Returns a new deep copy of the function.
           */
          virtual AutoPtr<FunctionBase<Number> > clone () const;
        
          /**
           * @returns the value at point \p p and time
           * \p time, which defaults to zero.
           */
          Number operator() (const Point& p,
        		     const Real time=0.);
        
          /**
           * Like before, but returns the values in a
           * writable reference.
           */
          void operator() (const Point& p,
        		   const Real time,
        		   DenseVector<Number>& output);
        
        };
        
        
        
------------------------------------------------------------ MeshlessInterpolationFunction inline methods
        inline
        Number MeshlessInterpolationFunction::operator() (const Point& p,
        						  const Real /* time */)
        {
          _pts.clear();
          _pts.push_back(p);
          _vals.resize(1);
        
          Threads::spin_mutex::scoped_lock lock(_mutex);
        
          _mfi.interpolate_field_data (_mfi.field_variables(),
        			       _pts, _vals);
        
          return _vals.front();
        }
        
        
        
        inline
        void MeshlessInterpolationFunction::operator() (const Point& p,
        						const Real time,
        						DenseVector<Number>& output)
        {
          output.resize(1);
          output(0) = (*this)(p,time);
          return;
        }
        
        
        
        inline
        void MeshlessInterpolationFunction::init ()
        {
        }
        
        
        
        inline
        void MeshlessInterpolationFunction::clear ()
        {
        }
        
        
        
        inline
        AutoPtr<FunctionBase<Number> >
        MeshlessInterpolationFunction::clone () const
        { 
          return AutoPtr<FunctionBase<Number> > (new MeshlessInterpolationFunction (_mfi, _mutex) );
        }
        
        
        } // namespace libMesh
        
        
        #endif // LIBMESH_MESHLESS_INTERPOLATION_FUNCTION_H
        



The source file miscellaneous_ex8.C with comments:

        #include "libmesh/libmesh.h"
        #include "libmesh/meshfree_interpolation.h"
        #include "libmesh/mesh.h"
        #include "libmesh/equation_systems.h"
        #include "libmesh/numeric_vector.h"
        #include "libmesh/tecplot_io.h"
        #include "libmesh/threads.h"
        #include "meshless_interpolation_function.h"
        
C++ includes
        #include <cstdlib>
        
        
Bring in everything from the libMesh namespace
        using namespace libMesh;
        
        
        void create_random_point_cloud (const unsigned int Npts,
        				std::vector<Point> &pts,
        				const Real max_range = 10)
        {
          std::cout << "Generating "<< Npts << " point cloud...";
          pts.resize(Npts);
        
          for (size_t i=0;i<Npts;i++)
            {
              pts[i](0) = max_range * (std::rand() % 1000) / Real(1000);
              pts[i](1) = max_range * (std::rand() % 1000) / Real(1000);
              pts[i](2) = max_range * (std::rand() % 1000) / Real(1000);
            }
          std::cout << "done\n";
        }
        
        
        
        Real exact_solution_u (const Point &p)
        {
          const Real
            x = p(0),
            y = p(1),
            z = p(2);
          
          return (x*x*x   +
        	  y*y*y*y +
        	  z*z*z*z*z);
        }
        
        
        
        Real exact_solution_v (const Point &p)
        {
          const Real
            x = p(0),
            y = p(1),
            z = p(2);
          
          return (x*x   +
        	  y*y +
        	  z*z*z);
        }
        
        Number exact_value (const Point& p,
                            const Parameters&,
                            const std::string&,
                            const std::string&)
        {
          return exact_solution_v(p);
        }
        
We now define the function which provides the initialization routines for the "Convection-Diffusion" system. This handles things like setting initial conditions and boundary conditions.
        void init_sys(EquationSystems& es,
                      const std::string& system_name)
        {
Get a reference to the Convection-Diffusion system object.
          System & system =
            es.get_system<System>(system_name);
          
          system.project_solution(exact_value, NULL, es.parameters);
        }
        
        
        
        
        int main(int argc, char** argv)
        {
Skip this example if we do not meet certain requirements
        #if LIBMESH_DIM != 3
          std::cout << "This example requires 3D support - skipping.\n";
          return 77;
        #endif
        #ifndef LIBMESH_HAVE_ZLIB_H
          std::cout << "This example requires zlib support - skipping.\n";
          return 77;
        #endif
        
Initialize libMesh.
          LibMeshInit init (argc, argv);
          {
Demonstration case 1
            {
              std::vector<Point>       tgt_pts;
              std::vector<Number>      tgt_data;
              std::vector<std::string> field_vars;
              
              field_vars.push_back("u");
              field_vars.push_back("v");
              
              InverseDistanceInterpolation<3> idi (/* n_interp_pts = */ 8,
        					   /* power =        */ 2);
              
              idi.set_field_variables (field_vars);
              
              create_random_point_cloud (1e5,
        				 idi.get_source_points());
              
Explicitly set the data values we will interpolate from
              {
        	const std::vector<Point> &src_pts  (idi.get_source_points());
        	std::vector<Number>      &src_vals (idi.get_source_vals());
        	
        	src_vals.clear(); src_vals.reserve(src_pts.size());
        	
        	for (std::vector<Point>::const_iterator pt_it=src_pts.begin();
        	     pt_it != src_pts.end(); ++pt_it)
        	  {
        	    src_vals.push_back (exact_solution_u (*pt_it));
        	    src_vals.push_back (exact_solution_v (*pt_it));
        	  }	  
              }
        
              std::cout << idi;
              
Interpolate to some other random points, and evaluate the result
              {
        	create_random_point_cloud (10,
        				   tgt_pts);
        	
        	idi.interpolate_field_data (field_vars,
        				    tgt_pts,
        				    tgt_data);
        	
              
        	std::vector<Number>::const_iterator v_it=tgt_data.begin();
        	
        	for (std::vector<Point>::const_iterator  p_it=tgt_pts.begin();
        	     p_it!=tgt_pts.end(); ++p_it)
        	  {
        	    std::cout << "\nAt target point " << *p_it
        		      << "\n u_interp=" << *v_it
        		      << ", u_exact="  << exact_solution_u(*p_it);
        	    ++v_it;
        	    std::cout << "\n v_interp=" << *v_it
        		      << ", v_exact="  << exact_solution_v(*p_it)
        		      << std::endl;
        	    ++v_it;
        	  }
              }
            }
        
        
Demonstration case 2
            {
              Mesh mesh_a, mesh_b;
        
              mesh_a.read("struct.ucd.gz"); mesh_b.read("unstruct.ucd.gz");
        
Create equation systems objects.
              EquationSystems
        	es_a(mesh_a), es_b(mesh_b);
        
              System
        	&sys_a = es_a.add_system<System>("src_system"),
        	&sys_b = es_b.add_system<System>("dest_system");
        
              sys_a.add_variable ("Cp", FIRST);
              sys_b.add_variable ("Cp", FIRST);
        
              sys_a.attach_init_function (init_sys);
              es_a.init();
              
Write out the initial conditions.
              TecplotIO(mesh_a).write_equation_systems ("src.dat",
        						es_a);
        
              InverseDistanceInterpolation<3> idi (/* n_interp_pts = */ 4,
        					   /* power =        */ 2);
        
              std::vector<Point>  &src_pts  (idi.get_source_points());
              std::vector<Number> &src_vals (idi.get_source_vals());
              std::vector<std::string> field_vars;      
              field_vars.push_back("Cp");
              idi.set_field_variables(field_vars);
        
We now will loop over every node in the source mesh and add it to a source point list, along with the solution
              {
        	MeshBase::const_node_iterator nd  = mesh_a.local_nodes_begin();
        	MeshBase::const_node_iterator end = mesh_a.local_nodes_end();
        
        	for (; nd!=end; ++nd)
        	  {
        	    const Node *node(*nd);
        	    src_pts.push_back(*node);
        	    src_vals.push_back(sys_a.current_solution(node->dof_number(0,0,0)));
        	  }	  	
              }
        
We have only set local values - prepare for use by gathering remote gata
              idi.prepare_for_use();
        
Create a MeshlessInterpolationFunction that uses our InverseDistanceInterpolation object. Since each MeshlessInterpolationFunction shares the same InverseDistanceInterpolation object in a threaded environment we must also provide a locking mechanism.
              Threads::spin_mutex mutex;
              MeshlessInterpolationFunction mif(idi, mutex);
        
project the solution onto system b
              es_b.init();
              sys_b.project_solution (&mif);
              
Write the result
              TecplotIO(mesh_b).write_equation_systems ("dest.dat",
        						es_b);
            }
        
        
            
          }
          return 0;
        }



The source file meshless_interpolation_function.h without comments:

 
  #ifndef LIBMESH_MESHLESS_INTERPOLATION_FUNCTION_H
  #define LIBMESH_MESHLESS_INTERPOLATION_FUNCTION_H
  
  #include "libmesh/function_base.h"
  #include "libmesh/meshfree_interpolation.h"
  #include "libmesh/threads.h"
  
  
  #include <cstddef>
  
  namespace libMesh
  {
  
  
  
  template <typename T>
  class DenseVector;
  
  
  class MeshlessInterpolationFunction : public FunctionBase<Number>
  {
  private:
    const MeshfreeInterpolation &_mfi;
    mutable std::vector<Point> _pts;
    mutable std::vector<Number> _vals;
    Threads::spin_mutex &_mutex;
    
  public:
  
    /**
     * Constructor.  Requires a \p \pMeshlessInterpolation object.
     */
    MeshlessInterpolationFunction (const MeshfreeInterpolation &mfi,
  				 Threads::spin_mutex &mutex) :
    _mfi (mfi),
    _mutex(mutex)
    {}
    
  
    /**
     * The actual initialization process.
     */
    void init ();
  
    /**
     * Clears the function.
     */
    void clear ();
  
    /**
     * Returns a new deep copy of the function.
     */
    virtual AutoPtr<FunctionBase<Number> > clone () const;
  
    /**
     * @returns the value at point \p p and time
     * \p time, which defaults to zero.
     */
    Number operator() (const Point& p,
  		     const Real time=0.);
  
    /**
     * Like before, but returns the values in a
     * writable reference.
     */
    void operator() (const Point& p,
  		   const Real time,
  		   DenseVector<Number>& output);
  
  };
  
  
  
  inline
  Number MeshlessInterpolationFunction::operator() (const Point& p,
  						  const Real /* time */)
  {
    _pts.clear();
    _pts.push_back(p);
    _vals.resize(1);
  
    Threads::spin_mutex::scoped_lock lock(_mutex);
  
    _mfi.interpolate_field_data (_mfi.field_variables(),
  			       _pts, _vals);
  
    return _vals.front();
  }
  
  
  
  inline
  void MeshlessInterpolationFunction::operator() (const Point& p,
  						const Real time,
  						DenseVector<Number>& output)
  {
    output.resize(1);
    output(0) = (*this)(p,time);
    return;
  }
  
  
  
  inline
  void MeshlessInterpolationFunction::init ()
  {
  }
  
  
  
  inline
  void MeshlessInterpolationFunction::clear ()
  {
  }
  
  
  
  inline
  AutoPtr<FunctionBase<Number> >
  MeshlessInterpolationFunction::clone () const
  { 
    return AutoPtr<FunctionBase<Number> > (new MeshlessInterpolationFunction (_mfi, _mutex) );
  }
  
  
  } // namespace libMesh
  
  
  #endif // LIBMESH_MESHLESS_INTERPOLATION_FUNCTION_H
  



The source file miscellaneous_ex8.C without comments:

 
  #include "libmesh/libmesh.h"
  #include "libmesh/meshfree_interpolation.h"
  #include "libmesh/mesh.h"
  #include "libmesh/equation_systems.h"
  #include "libmesh/numeric_vector.h"
  #include "libmesh/tecplot_io.h"
  #include "libmesh/threads.h"
  #include "meshless_interpolation_function.h"
  
  #include <cstdlib>
  
  
  using namespace libMesh;
  
  
  void create_random_point_cloud (const unsigned int Npts,
  				std::vector<Point> &pts,
  				const Real max_range = 10)
  {
    std::cout << "Generating "<< Npts << " point cloud...";
    pts.resize(Npts);
  
    for (size_t i=0;i<Npts;i++)
      {
        pts[i](0) = max_range * (std::rand() % 1000) / Real(1000);
        pts[i](1) = max_range * (std::rand() % 1000) / Real(1000);
        pts[i](2) = max_range * (std::rand() % 1000) / Real(1000);
      }
    std::cout << "done\n";
  }
  
  
  
  Real exact_solution_u (const Point &p)
  {
    const Real
      x = p(0),
      y = p(1),
      z = p(2);
    
    return (x*x*x   +
  	  y*y*y*y +
  	  z*z*z*z*z);
  }
  
  
  
  Real exact_solution_v (const Point &p)
  {
    const Real
      x = p(0),
      y = p(1),
      z = p(2);
    
    return (x*x   +
  	  y*y +
  	  z*z*z);
  }
  
  Number exact_value (const Point& p,
                      const Parameters&,
                      const std::string&,
                      const std::string&)
  {
    return exact_solution_v(p);
  }
  
  void init_sys(EquationSystems& es,
                const std::string& system_name)
  {
    System & system =
      es.get_system<System>(system_name);
    
    system.project_solution(exact_value, NULL, es.parameters);
  }
  
  
  
  
  int main(int argc, char** argv)
  {
  #if LIBMESH_DIM != 3
    std::cout << "This example requires 3D support - skipping.\n";
    return 77;
  #endif
  #ifndef LIBMESH_HAVE_ZLIB_H
    std::cout << "This example requires zlib support - skipping.\n";
    return 77;
  #endif
  
    LibMeshInit init (argc, argv);
    {
      {
        std::vector<Point>       tgt_pts;
        std::vector<Number>      tgt_data;
        std::vector<std::string> field_vars;
        
        field_vars.push_back("u");
        field_vars.push_back("v");
        
        InverseDistanceInterpolation<3> idi (/* n_interp_pts = */ 8,
  					   /* power =        */ 2);
        
        idi.set_field_variables (field_vars);
        
        create_random_point_cloud (1e5,
  				 idi.get_source_points());
        
        {
  	const std::vector<Point> &src_pts  (idi.get_source_points());
  	std::vector<Number>      &src_vals (idi.get_source_vals());
  	
  	src_vals.clear(); src_vals.reserve(src_pts.size());
  	
  	for (std::vector<Point>::const_iterator pt_it=src_pts.begin();
  	     pt_it != src_pts.end(); ++pt_it)
  	  {
  	    src_vals.push_back (exact_solution_u (*pt_it));
  	    src_vals.push_back (exact_solution_v (*pt_it));
  	  }	  
        }
  
        std::cout << idi;
        
        {
  	create_random_point_cloud (10,
  				   tgt_pts);
  	
  	idi.interpolate_field_data (field_vars,
  				    tgt_pts,
  				    tgt_data);
  	
        
  	std::vector<Number>::const_iterator v_it=tgt_data.begin();
  	
  	for (std::vector<Point>::const_iterator  p_it=tgt_pts.begin();
  	     p_it!=tgt_pts.end(); ++p_it)
  	  {
  	    std::cout << "\nAt target point " << *p_it
  		      << "\n u_interp=" << *v_it
  		      << ", u_exact="  << exact_solution_u(*p_it);
  	    ++v_it;
  	    std::cout << "\n v_interp=" << *v_it
  		      << ", v_exact="  << exact_solution_v(*p_it)
  		      << std::endl;
  	    ++v_it;
  	  }
        }
      }
  
  
      {
        Mesh mesh_a, mesh_b;
  
        mesh_a.read("struct.ucd.gz"); mesh_b.read("unstruct.ucd.gz");
  
        EquationSystems
  	es_a(mesh_a), es_b(mesh_b);
  
        System
  	&sys_a = es_a.add_system<System>("src_system"),
  	&sys_b = es_b.add_system<System>("dest_system");
  
        sys_a.add_variable ("Cp", FIRST);
        sys_b.add_variable ("Cp", FIRST);
  
        sys_a.attach_init_function (init_sys);
        es_a.init();
        
        TecplotIO(mesh_a).write_equation_systems ("src.dat",
  						es_a);
  
        InverseDistanceInterpolation<3> idi (/* n_interp_pts = */ 4,
  					   /* power =        */ 2);
  
        std::vector<Point>  &src_pts  (idi.get_source_points());
        std::vector<Number> &src_vals (idi.get_source_vals());
        std::vector<std::string> field_vars;      
        field_vars.push_back("Cp");
        idi.set_field_variables(field_vars);
  
        {
  	MeshBase::const_node_iterator nd  = mesh_a.local_nodes_begin();
  	MeshBase::const_node_iterator end = mesh_a.local_nodes_end();
  
  	for (; nd!=end; ++nd)
  	  {
  	    const Node *node(*nd);
  	    src_pts.push_back(*node);
  	    src_vals.push_back(sys_a.current_solution(node->dof_number(0,0,0)));
  	  }	  	
        }
  
        idi.prepare_for_use();
  
        Threads::spin_mutex mutex;
        MeshlessInterpolationFunction mif(idi, mutex);
  
        es_b.init();
        sys_b.project_solution (&mif);
        
        TecplotIO(mesh_b).write_equation_systems ("dest.dat",
  						es_b);
      }
  
  
      
    }
    return 0;
  }



The console output of the program:

***************************************************************
* Running Example miscellaneous_ex8:
*  mpirun -np 2 example-devel  -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
***************************************************************
 
Generating 100000 point cloud...done
MeshfreeInterpolation
 n_source_points()=100000
 n_field_variables()=2
  variables = u v 
Generating 10 point cloud...done
*** Warning, This code is untested, experimental, or likely to see future API changes: src/solution_transfer/meshfree_interpolation.C, line 283, compiled Feb  5 2013 at 13:19:19 ***

At target point (x,y,z)=(    9.61,     7.45,     5.24)
 u_interp=8241.22, u_exact=7918.57
 v_interp=298.061, v_exact=291.732

At target point (x,y,z)=(     3.5,     9.71,     8.11)
 u_interp=44209.4, u_exact=44016
 v_interp=641.42, v_exact=639.946

At target point (x,y,z)=(     5.3,      8.4,     4.36)
 u_interp=6706.14, u_exact=6703.14
 v_interp=183.748, v_exact=181.532

At target point (x,y,z)=(    3.57,     5.66,     5.31)
 u_interp=5324.74, u_exact=5293.34
 v_interp=194.632, v_exact=194.502

At target point (x,y,z)=(    7.66,     3.88,     8.31)
 u_interp=39033.6, u_exact=40304.4
 v_interp=635.944, v_exact=647.586

At target point (x,y,z)=(    9.54,     7.83,     1.32)
 u_interp=4706.06, u_exact=4631.04
 v_interp=155.44, v_exact=154.62

At target point (x,y,z)=(    8.32,     1.21,     6.09)
 u_interp=8806.23, u_exact=8955.03
 v_interp=293.686, v_exact=296.553

At target point (x,y,z)=(    0.32,     2.44,     1.59)
 u_interp=45.7344, u_exact=45.6403
 v_interp=10.1975, v_exact=10.0757

At target point (x,y,z)=(    2.14,     4.04,     8.72)
 u_interp=51856, u_exact=50693.8
 v_interp=692.691, v_exact=683.956

At target point (x,y,z)=(     9.4,     0.51,     2.22)
 u_interp=883.149, u_exact=884.574
 v_interp=99.3335, v_exact=99.5611
************************************************************************************************************************
***             WIDEN YOUR WINDOW TO 120 CHARACTERS.  Use 'enscript -r -fCourier9' to print this document            ***
************************************************************************************************************************

---------------------------------------------- PETSc Performance Summary: ----------------------------------------------

/workspace/libmesh/examples/miscellaneous/miscellaneous_ex8/.libs/lt-example-devel on a intel-12. named hbar.ices.utexas.edu with 2 processors, by benkirk Tue Feb  5 13:41:42 2013
Using Petsc Release Version 3.3.0, Patch 2, Fri Jul 13 15:42:00 CDT 2012 

                         Max       Max/Min        Avg      Total 
Time (sec):           1.899e+00      1.00000   1.899e+00
Objects:              1.500e+01      1.00000   1.500e+01
Flops:                0.000e+00      0.00000   0.000e+00  0.000e+00
Flops/sec:            0.000e+00      0.00000   0.000e+00  0.000e+00
MPI Messages:         6.000e+00      1.00000   6.000e+00  1.200e+01
MPI Message Lengths:  2.696e+03      1.00000   4.493e+02  5.392e+03
MPI Reductions:       3.100e+01      1.00000

Flop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)
                            e.g., VecAXPY() for real vectors of length N --> 2N flops
                            and VecAXPY() for complex vectors of length N --> 8N flops

Summary of Stages:   ----- Time ------  ----- Flops -----  --- Messages ---  -- Message Lengths --  -- Reductions --
                        Avg     %Total     Avg     %Total   counts   %Total     Avg         %Total   counts   %Total 
 0:      Main Stage: 1.8985e+00 100.0%  0.0000e+00   0.0%  1.200e+01 100.0%  4.493e+02      100.0%  3.000e+01  96.8% 

------------------------------------------------------------------------------------------------------------------------
See the 'Profiling' chapter of the users' manual for details on interpreting output.
Phase summary info:
   Count: number of times phase was executed
   Time and Flops: Max - maximum over all processors
                   Ratio - ratio of maximum to minimum over all processors
   Mess: number of messages sent
   Avg. len: average message length
   Reduct: number of global reductions
   Global: entire computation
   Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().
      %T - percent time in this phase         %f - percent flops in this phase
      %M - percent messages in this phase     %L - percent message lengths in this phase
      %R - percent reductions in this phase
   Total Mflop/s: 10e-6 * (sum of flops over all processors)/(max time over all processors)
------------------------------------------------------------------------------------------------------------------------
Event                Count      Time (sec)     Flops                             --- Global ---  --- Stage ---   Total
                   Max Ratio  Max     Ratio   Max  Ratio  Mess   Avg len Reduct  %T %f %M %L %R  %T %f %M %L %R Mflop/s
------------------------------------------------------------------------------------------------------------------------

--- Event Stage 0: Main Stage

VecCopy                2 1.0 4.0531e-06 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
VecSet                 4 1.0 4.7684e-06 1.2 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
VecAssemblyBegin       4 1.0 2.6240e-0317.4 0.00e+00 0.0 0.0e+00 0.0e+00 1.2e+01  0  0  0  0 39   0  0  0  0 40     0
VecAssemblyEnd         4 1.0 1.9789e-05 1.1 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
VecScatterBegin        2 1.0 4.0054e-05 1.0 0.00e+00 0.0 4.0e+00 9.0e+02 0.0e+00  0  0 33 66  0   0  0 33 66  0     0
VecScatterEnd          2 1.0 8.1062e-06 1.1 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
------------------------------------------------------------------------------------------------------------------------

Memory usage is given in bytes:

Object Type          Creations   Destructions     Memory  Descendants' Mem.
Reports information only for process 0.

--- Event Stage 0: Main Stage

              Vector     6              6        60144     0
      Vector Scatter     2              2         2072     0
           Index Set     4              4         3576     0
   IS L to G Mapping     2              2         1128     0
              Viewer     1              0            0     0
========================================================================================================================
Average time to get PetscTime(): 9.53674e-08
Average time for MPI_Barrier(): 1.19209e-06
Average time for zero size MPI_Send(): 1.54972e-05
#PETSc Option Table entries:
-ksp_right_pc
-log_summary
-pc_type bjacobi
-sub_pc_factor_levels 4
-sub_pc_factor_zeropivot 0
-sub_pc_type ilu
#End of PETSc Option Table entries
Compiled without FORTRAN kernels
Compiled with full precision matrices (default)
sizeof(short) 2 sizeof(int) 4 sizeof(long) 8 sizeof(void*) 8 sizeof(PetscScalar) 8 sizeof(PetscInt) 4
Configure run at: Thu Nov  8 11:21:02 2012
Configure options: --with-debugging=false --COPTFLAGS=-O3 --CXXOPTFLAGS=-O3 --FOPTFLAGS=-O3 --with-clanguage=C++ --with-shared-libraries=1 --with-mpi-dir=/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1 --with-mumps=true --download-mumps=1 --with-metis=true --download-metis=1 --with-parmetis=true --download-parmetis=1 --with-superlu=true --download-superlu=1 --with-superludir=true --download-superlu_dist=1 --with-blacs=true --download-blacs=1 --with-scalapack=true --download-scalapack=1 --with-hypre=true --download-hypre=1 --with-blas-lib="[/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_intel_lp64.so,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_sequential.so,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_core.so]" --with-lapack-lib="[/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_lapack95_lp64.a]"
-----------------------------------------
Libraries compiled on Thu Nov  8 11:21:02 2012 on daedalus.ices.utexas.edu 
Machine characteristics: Linux-2.6.32-279.1.1.el6.x86_64-x86_64-with-redhat-6.3-Carbon
Using PETSc directory: /opt/apps/ossw/libraries/petsc/petsc-3.3-p2
Using PETSc arch: intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt
-----------------------------------------

Using C compiler: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpicxx  -wd1572 -O3   -fPIC   ${COPTFLAGS} ${CFLAGS}
Using Fortran compiler: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpif90  -fPIC -O3   ${FOPTFLAGS} ${FFLAGS} 
-----------------------------------------

Using include paths: -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/include -I/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/include
-----------------------------------------

Using C linker: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpicxx
Using Fortran linker: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpif90
Using libraries: -Wl,-rpath,/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -L/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -lpetsc -lX11 -Wl,-rpath,/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -L/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -lcmumps -ldmumps -lsmumps -lzmumps -lmumps_common -lpord -lHYPRE -lpthread -lsuperlu_dist_3.0 -lparmetis -lmetis -lscalapack -lblacs -lsuperlu_4.3 -Wl,-rpath,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64 -L/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64 -lmkl_lapack95_lp64 -lmkl_intel_lp64 -lmkl_sequential -lmkl_core -Wl,-rpath,/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/lib -L/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/lib -Wl,-rpath,/opt/apps/sysnet/intel/12.1/composer_xe_2011_sp1.7.256/compiler/lib/intel64 -L/opt/apps/sysnet/intel/12.1/composer_xe_2011_sp1.7.256/compiler/lib/intel64 -Wl,-rpath,/usr/lib/gcc/x86_64-redhat-linux/4.4.6 -L/usr/lib/gcc/x86_64-redhat-linux/4.4.6 -lmpichf90 -lifport -lifcore -lm -lm -lmpichcxx -ldl -lmpich -lopa -lmpl -lrt -lpthread -limf -lsvml -lipgo -ldecimal -lcilkrts -lstdc++ -lgcc_s -lirc -lirc_s -ldl 
-----------------------------------------


 ----------------------------------------------------------------------------------------------------------------------
| Processor id:   0                                                                                                    |
| Num Processors: 2                                                                                                    |
| Time:           Tue Feb  5 13:41:42 2013                                                                             |
| OS:             Linux                                                                                                |
| HostName:       hbar.ices.utexas.edu                                                                                 |
| OS Release:     2.6.32-279.1.1.el6.x86_64                                                                            |
| OS Version:     #1 SMP Tue Jul 10 11:24:23 CDT 2012                                                                  |
| Machine:        x86_64                                                                                               |
| Username:       benkirk                                                                                              |
| Configuration:  ./configure  '--enable-everything'                                                                   |
|  '--prefix=/workspace/libmesh/install'                                                                               |
|  'CXX=mpicxx'                                                                                                        |
|  'CC=mpicc'                                                                                                          |
|  'F77=mpif77'                                                                                                        |
|  'FC=mpif90'                                                                                                         |
|  'PETSC_DIR=/opt/apps/ossw/libraries/petsc/petsc-3.3-p2'                                                             |
|  'PETSC_ARCH=intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt'                                                |
|  'SLEPC_DIR=/opt/apps/ossw/libraries/slepc/slepc-3.3-p2-petsc-3.3-p2-cxx-opt'                                        |
|  'TRILINOS_DIR=/opt/apps/ossw/libraries/trilinos/trilinos-10.12.2/sl6/intel-12.1/mpich2-1.4.1p1/mkl-intel-10.3.12.361'|
|  'VTK_DIR=/opt/apps/ossw/libraries/vtk/vtk-5.10.0/sl6/intel-12.1'                                                    |
 ----------------------------------------------------------------------------------------------------------------------
 ------------------------------------------------------------------------------------------------------------
| libMesh Performance: Alive time=1.90843, Active time=1.80016                                               |
 ------------------------------------------------------------------------------------------------------------
| Event                          nCalls    Total Time  Avg Time    Total Time  Avg Time    % of Active Time  |
|                                          w/o Sub     w/o Sub     With Sub    With Sub    w/o S    With S   |
|------------------------------------------------------------------------------------------------------------|
|                                                                                                            |
|                                                                                                            |
| DofMap                                                                                                     |
|   add_neighbors_to_send_list() 2         0.0824      0.041176    0.0898      0.044924    4.57     4.99     |
|   create_dof_constraints()     2         0.0076      0.003775    0.0076      0.003775    0.42     0.42     |
|   distribute_dofs()            2         0.0798      0.039892    0.2417      0.120840    4.43     13.43    |
|   dof_indices()                9157      0.2083      0.000023    0.2083      0.000023    11.57    11.57    |
|   prepare_send_list()          2         0.0001      0.000062    0.0001      0.000062    0.01     0.01     |
|   reinit()                     2         0.1601      0.080027    0.1601      0.080027    8.89     8.89     |
|                                                                                                            |
| EquationSystems                                                                                            |
|   build_solution_vector()      2         0.0164      0.008189    0.1165      0.058256    0.91     6.47     |
|                                                                                                            |
| InverseDistanceInterpolation<>                                                                             |
|   construct_kd_tree()          2         0.3420      0.171010    0.3420      0.171010    19.00    19.00    |
|   interpolate_field_data()     8437      0.0538      0.000006    0.0538      0.000006    2.99     2.99     |
|                                                                                                            |
| Mesh                                                                                                       |
|   find_neighbors()             2         0.1490      0.074516    0.1497      0.074828    8.28     8.31     |
|   read()                       2         0.0669      0.033452    0.0669      0.033452    3.72     3.72     |
|   renumber_nodes_and_elem()    4         0.0065      0.001629    0.0065      0.001629    0.36     0.36     |
|                                                                                                            |
| MeshCommunication                                                                                          |
|   broadcast()                  2         0.0729      0.036434    0.1337      0.066847    4.05     7.43     |
|   compute_hilbert_indices()    4         0.1365      0.034120    0.1365      0.034120    7.58     7.58     |
|   find_global_indices()        4         0.0509      0.012722    0.1993      0.049815    2.83     11.07    |
|   parallel_sort()              4         0.0099      0.002474    0.0101      0.002529    0.55     0.56     |
|                                                                                                            |
| MeshOutput                                                                                                 |
|   write_equation_systems()     2         0.0001      0.000046    0.1465      0.073248    0.01     8.14     |
|                                                                                                            |
| MeshfreeInterpolation                                                                                      |
|   gather_remote_data()         1         0.0005      0.000510    0.0006      0.000563    0.03     0.03     |
|                                                                                                            |
| MetisPartitioner                                                                                           |
|   partition()                  2         0.1809      0.090426    0.2808      0.140380    10.05    15.60    |
|                                                                                                            |
| Parallel                                                                                                   |
|   allgather()                  18        0.0012      0.000069    0.0013      0.000072    0.07     0.07     |
|   broadcast()                  12        0.0012      0.000103    0.0012      0.000103    0.07     0.07     |
|   max(scalar)                  204       0.0006      0.000003    0.0006      0.000003    0.04     0.04     |
|   max(vector)                  47        0.0003      0.000007    0.0008      0.000016    0.02     0.04     |
|   min(bool)                    241       0.0007      0.000003    0.0007      0.000003    0.04     0.04     |
|   min(scalar)                  198       0.0690      0.000349    0.0690      0.000349    3.83     3.83     |
|   min(vector)                  47        0.0004      0.000009    0.0009      0.000020    0.02     0.05     |
|   probe()                      21        0.0021      0.000099    0.0021      0.000099    0.12     0.12     |
|   receive()                    21        0.0003      0.000013    0.0024      0.000113    0.02     0.13     |
|   send()                       21        0.0001      0.000005    0.0001      0.000005    0.01     0.01     |
|   send_receive()               28        0.0004      0.000014    0.0029      0.000103    0.02     0.16     |
|   sum()                        20        0.0002      0.000011    0.0007      0.000033    0.01     0.04     |
|                                                                                                            |
| Parallel::Request                                                                                          |
|   wait()                       21        0.0001      0.000003    0.0001      0.000003    0.00     0.00     |
|                                                                                                            |
| Partitioner                                                                                                |
|   set_node_processor_ids()     2         0.0124      0.006222    0.0186      0.009324    0.69     1.04     |
|   set_parent_processor_ids()   2         0.0133      0.006660    0.0133      0.006660    0.74     0.74     |
|                                                                                                            |
| System                                                                                                     |
|   project_vector()             2         0.0435      0.021737    0.2063      0.103134    2.42     11.46    |
|                                                                                                            |
| TecplotIO                                                                                                  |
|   write_nodal_data()           2         0.0298      0.014897    0.0298      0.014897    1.66     1.66     |
 ------------------------------------------------------------------------------------------------------------
| Totals:                        18542     1.8002                                          100.00            |
 ------------------------------------------------------------------------------------------------------------

 
***************************************************************
* Done Running Example miscellaneous_ex8:
*  mpirun -np 2 example-devel  -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
***************************************************************

Site Created By: libMesh Developers
Last modified: February 05 2013 19:42:40 UTC

Hosted By:
SourceForge.net Logo