Browse Source

Fix error in balanceOctree(...)

- include/mpi_tree.txx: In balanceOctree(...), use new_size and new_rank
  after MPI_Comm_split(...)

- include/interac_list.txx: Fix error: more than one instance of
  overloaded function "pow" matches the argument list.

- include/kernel.hpp: include files vector.hpp and matrix.hpp

- include/mortonif.txx: declare max_int as static in contructor, so that
  it is not computed each time.
Dhairya Malhotra 10 năm trước cách đây
mục cha
commit
226e47b9ae
4 tập tin đã thay đổi với 19 bổ sung14 xóa
  1. 3 3
      include/interac_list.txx
  2. 2 0
      include/kernel.hpp
  3. 4 4
      include/mortonid.txx
  4. 10 7
      include/mpi_tree.txx

+ 3 - 3
include/interac_list.txx

@@ -95,8 +95,8 @@ std::vector<Perm_Type>& InteracList<Node_t>::PermutList(Mat_Type t, size_t i){
 template <class Node_t>
 std::vector<Node_t*> InteracList<Node_t>::BuildList(Node_t* n, Mat_Type t){
   std::vector<Node_t*> interac_list(ListCount(t),NULL);
-  int n_collg=(int)pow(3,dim);
-  int n_child=(int)pow(2,dim);
+  int n_collg=(int)pow(3.0,(int)dim);
+  int n_child=(int)pow(2.0,(int)dim);
   int rel_coord[3];
 
   switch (t){
@@ -367,7 +367,7 @@ int InteracList<Node_t>::class_hash(int* c_){
  */
 template <class Node_t>
 void InteracList<Node_t>::InitList(int max_r, int min_r, int step, Mat_Type t){
-  size_t count=(size_t)(pow((max_r*2)/step+1,dim)-(min_r>0?pow((min_r*2)/step-1,dim):0));
+  size_t count=(size_t)(pow((max_r*2)/step+1.0,(int)dim)-(min_r>0?pow((min_r*2)/step-1.0,(int)dim):0));
   Matrix<int>& M=rel_coord[t];
   M.Resize(count,dim);
   hash_lut[t].assign(PVFMM_MAX_COORD_HASH, -1);

+ 2 - 0
include/kernel.hpp

@@ -11,6 +11,8 @@
 
 #include <pvfmm_common.hpp>
 #include <mem_mgr.hpp>
+#include <vector.hpp>
+#include <matrix.hpp>
 
 #ifndef _PVFMM_FMM_KERNEL_HPP_
 #define _PVFMM_FMM_KERNEL_HPP_

+ 4 - 4
include/mortonid.txx

@@ -15,7 +15,7 @@ inline MortonId::MortonId(MortonId m, unsigned char depth_):x(m.x), y(m.y), z(m.
 
 template <class T>
 inline MortonId::MortonId(T x_f,T y_f, T z_f, unsigned char depth_): depth(depth_){
-  UINT_T max_int=((UINT_T)1)<<(MAX_DEPTH);
+  static UINT_T max_int=((UINT_T)1)<<(MAX_DEPTH);
   x=(UINT_T)floor(x_f*max_int);
   y=(UINT_T)floor(y_f*max_int);
   z=(UINT_T)floor(z_f*max_int);
@@ -23,7 +23,7 @@ inline MortonId::MortonId(T x_f,T y_f, T z_f, unsigned char depth_): depth(depth
 
 template <class T>
 inline MortonId::MortonId(T* coord, unsigned char depth_): depth(depth_){
-  UINT_T max_int=((UINT_T)1)<<(MAX_DEPTH);
+  static UINT_T max_int=((UINT_T)1)<<(MAX_DEPTH);
   x=(UINT_T)floor(coord[0]*max_int);
   y=(UINT_T)floor(coord[1]*max_int);
   z=(UINT_T)floor(coord[2]*max_int);
@@ -31,8 +31,8 @@ inline MortonId::MortonId(T* coord, unsigned char depth_): depth(depth_){
 
 template <class T>
 inline void MortonId::GetCoord(T* coord){
-  UINT_T max_int=((UINT_T)1)<<(MAX_DEPTH);
-  T s=1.0/((T)max_int);
+  static UINT_T max_int=((UINT_T)1)<<(MAX_DEPTH);
+  static T s=1.0/((T)max_int);
   coord[0]=x*s;
   coord[1]=y*s;
   coord[2]=z*s;

+ 10 - 7
include/mpi_tree.txx

@@ -866,20 +866,20 @@ inline int balanceOctree (std::vector<MortonId > &in, std::vector<MortonId > &ou
         MPI_Request recvRequest;
         MPI_Request sendRequest;
 
-        if(rank > 0) {
-          MPI_Irecv(&lastOnPrev, 1, par::Mpi_datatype<MortonId>::value(), rank-1, 1, comm, &recvRequest);
+        if(new_rank > 0) {
+          MPI_Irecv(&lastOnPrev, 1, par::Mpi_datatype<MortonId>::value(), new_rank-1, 1, new_comm, &recvRequest);
         }
-        if(rank < (size-1)) {
-          MPI_Issend( &lastOctant, 1, par::Mpi_datatype<MortonId>::value(), rank+1, 1, comm,  &sendRequest);
+        if(new_rank < (new_size-1)) {
+          MPI_Issend( &lastOctant, 1, par::Mpi_datatype<MortonId>::value(), new_rank+1, 1, new_comm,  &sendRequest);
         }
 
-        if(rank > 0) {
+        if(new_rank > 0) {
           MPI_Status statusWait;
           MPI_Wait(&recvRequest, &statusWait);
           nxt_mid = lastOnPrev.NextId();
         }
 
-        if(rank < (size-1)) {
+        if(new_rank < (new_size-1)) {
           MPI_Status statusWait;
           MPI_Wait(&sendRequest, &statusWait);
         }
@@ -902,7 +902,7 @@ inline int balanceOctree (std::vector<MortonId > &in, std::vector<MortonId > &ou
         }
         nxt_mid=out[i].NextId();
       }
-      if(rank==size-1){
+      if(new_rank==new_size-1){
         while(nxt_mid.GetDepth()>0){
           out1.push_back(nxt_mid);
           nxt_mid=nxt_mid.NextId();
@@ -910,6 +910,9 @@ inline int balanceOctree (std::vector<MortonId > &in, std::vector<MortonId > &ou
       }
       out.swap(out1);
     }
+    if(new_size<size){
+      par::partitionW<MortonId>(out, NULL , comm);
+    }
   }
 
   //////////////////////////////////////////////////////////////////////////////////////////////////