Browse Source

Bug fix in MPI_Node::Truncate()

- During truncate the parent vectors were not updated.

Other minor changes:
- Profile counts flop, memory even when not enabled.
- Profile::[Add_FLOP,Add_MEM,Enable] return original value.
- Move node_id to TreeNode from FMM_Node.
- In Tree, memgr initialized with zero size to reduce memory usage.
- Remove MPI_Comm_dup from MPI_Tree. We assume user will not free the
  communicator until after the tree is destroyed.
Dhairya Malhotra 10 years ago
parent
commit
94bce731d3
7 changed files with 39 additions and 25 deletions
  1. 0 1
      include/fmm_node.hpp
  2. 4 4
      include/mpi_node.txx
  3. 2 4
      include/mpi_tree.hpp
  4. 3 3
      include/profile.hpp
  5. 1 1
      include/tree.hpp
  6. 2 0
      include/tree_node.hpp
  7. 27 12
      src/profile.cpp

+ 0 - 1
include/fmm_node.hpp

@@ -171,7 +171,6 @@ class FMM_Node: public Node{
   Vector<size_t> trg_scatter;
 
   std::vector<std::vector<FMM_Node*> > interac_list;
-  size_t node_id; //For translating node pointer to index.
 
  private:
 

+ 4 - 4
include/mpi_node.txx

@@ -209,7 +209,7 @@ void MPI_Node<T>::Truncate(){
           Vector<Real_t>& chld_vec=*chld_pt_coord[i][j];
           vec_size+=chld_vec.Dim();
         }
-        Vector<Real_t> vec=*pt_coord[j];
+        Vector<Real_t>& vec=*pt_coord[j];
         vec.ReInit(vec_size);
 
         vec_size=0;
@@ -227,7 +227,7 @@ void MPI_Node<T>::Truncate(){
           Vector<Real_t>& chld_vec=*chld_pt_value[i][j];
           vec_size+=chld_vec.Dim();
         }
-        Vector<Real_t> vec=*pt_value[j];
+        Vector<Real_t>& vec=*pt_value[j];
         vec.ReInit(vec_size);
 
         vec_size=0;
@@ -245,7 +245,7 @@ void MPI_Node<T>::Truncate(){
           Vector<size_t>& chld_vec=*chld_pt_scatter[i][j];
           vec_size+=chld_vec.Dim();
         }
-        Vector<size_t> vec=*pt_scatter[j];
+        Vector<size_t>& vec=*pt_scatter[j];
         vec.ReInit(vec_size);
 
         vec_size=0;
@@ -519,7 +519,7 @@ void MPI_Node<T>::VTU_Data(VTUData_t& vtu_data, std::vector<Node_t*>& nodes, int
         value.push_back(n->pt_value[i]);
       }
     }
-    size_t value_dof=value.size()/point_cnt;
+    size_t value_dof=(value.size()?value.size()/point_cnt:0);
     assert(value_dof*point_cnt==value.size());
     for(size_t nid=0;nid<nodes.size();nid++){
       Node_t* n=nodes[nid];

+ 2 - 4
include/mpi_tree.hpp

@@ -38,14 +38,12 @@ class MPI_Tree: public Tree<TreeNode>{
   /**
    * \brief Constructor.
    */
-  MPI_Tree(MPI_Comm c): Tree<Node_t>() {MPI_Comm_dup(c,&comm);}
+  MPI_Tree(MPI_Comm c): Tree<Node_t>() {comm=c;}
 
   /**
    * \brief Virtual destructor.
    */
-  virtual ~MPI_Tree() {
-    MPI_Comm_free(&comm);
-  }
+  virtual ~MPI_Tree() {}
 
   /**
    * \brief Initialize the distributed MPI tree.

+ 3 - 3
include/profile.hpp

@@ -24,11 +24,11 @@ namespace pvfmm{
 class Profile{
   public:
 
-    static void Add_FLOP(long long inc);
+    static long long Add_FLOP(long long inc);
 
-    static void Add_MEM(long long inc);
+    static long long Add_MEM(long long inc);
 
-    static void Enable(bool state){enable_state=state;};
+    static bool Enable(bool state);
 
     static void Tic(const char* name_, const MPI_Comm* comm_=NULL,bool sync_=false, int level=0);
 

+ 1 - 1
include/tree.hpp

@@ -30,7 +30,7 @@ class Tree{
   /**
    * \brief Constructor.
    */
-  Tree(): dim(0), root_node(NULL), max_depth(MAX_DEPTH), memgr(DEVICE_BUFFER_SIZE*1024l*1024l) { };
+  Tree(): dim(0), root_node(NULL), max_depth(MAX_DEPTH), memgr(0) { };
 
   /**
    * \brief Virtual destructor.

+ 2 - 0
include/tree_node.hpp

@@ -128,6 +128,8 @@ class TreeNode{
    */
   void SetStatus(int flag);
 
+  size_t node_id; //For translating node pointer to index.
+
  protected:
 
   int dim;               //Dimension of the tree

+ 27 - 12
src/profile.cpp

@@ -20,24 +20,39 @@
 
 namespace pvfmm{
 
-void Profile::Add_FLOP(long long inc){
-#if __PROFILE__ >= 0
-  if(!enable_state) return;
+long long Profile::Add_FLOP(long long inc){
+  long long orig_val=FLOP;
+  #if __PROFILE__ >= 0
   #pragma omp critical (FLOP)
-  FLOP+=inc;
-#endif
+  {
+    orig_val=FLOP;
+    FLOP+=inc;
+  }
+  #endif
+  return orig_val;
 }
 
-void Profile::Add_MEM(long long inc){
-#if __PROFILE__ >= 0
-  if(!enable_state) return;
+long long Profile::Add_MEM(long long inc){
+  long long orig_val=MEM;
+  #if __PROFILE__ >= 0
   #pragma omp critical (MEM)
   {
-  MEM+=inc;
-  for(size_t i=0;i<max_mem.size();i++)
-    if(max_mem[i]<MEM) max_mem[i]=MEM;
+    orig_val=MEM;
+    MEM+=inc;
+    for(size_t i=0;i<max_mem.size();i++){
+      if(max_mem[i]<MEM) max_mem[i]=MEM;
+    }
   }
-#endif
+  #endif
+  return orig_val;
+}
+
+bool Profile::Enable(bool state){
+  bool orig_val=enable_state;
+  #if __PROFILE__ >= 0
+  enable_state=state;
+  #endif
+  return orig_val;
 }
 
 void Profile::Tic(const char* name_, const MPI_Comm* comm_,bool sync_, int verbose){