Parcourir la source

minor changes

Dhairya Malhotra il y a 11 ans
Parent
commit
d802dd28bb
4 fichiers modifiés avec 11 ajouts et 13 suppressions
  1. 1 1
      INSTALL
  2. 2 2
      include/mem_utils.txx
  3. 8 8
      include/mpi_tree.txx
  4. 0 2
      include/pvfmm_common.hpp

+ 1 - 1
INSTALL

@@ -88,7 +88,7 @@ are given below:
 
 `Stampede (TACC)'
      module load fftw3
-     ./configure CXXFLAGS="-mavx -wd3218 -wd2571 -wd2570 -wd2568 -wd161" --with-fftw="$TACC_FFTW3_DIR" FLIBS=" "
+     ./configure CXXFLAGS="-mavx -wd3218 -wd2570" --with-fftw="$TACC_FFTW3_DIR" FLIBS=" "
 
 `Ronaldo (ICES)'
      ./configure CXXFLAGS="-msse4" --with-fftw="$FFTW_DIR"

+ 2 - 2
include/mem_utils.txx

@@ -25,9 +25,9 @@ namespace mem{
     assert(alignment <= 0x8000);
     size_t size=size_*sizeof(T);
     uintptr_t r = (uintptr_t)malloc(size + --alignment + 2);
-    assert(r!=0);
+    //if (!r) return NULL;
+    ASSERT_WITH_MSG(r!=0, "malloc failed.");
     uintptr_t o = (uintptr_t)(r + 2 + alignment) & ~(uintptr_t)alignment;
-    if (!r) return NULL;
     ((uint16_t*)o)[-1] = (uint16_t)(o-r);
     return (T*)o;
     //return (T*)fftw_malloc(size);

+ 8 - 8
include/mpi_tree.txx

@@ -1121,6 +1121,11 @@ bool MPI_Tree<TreeNode>::CheckTree(){
   MPI_Comm_rank(*Comm(),&myrank);
   MPI_Comm_size(*Comm(),&np);
   std::vector<MortonId> mins=GetMins();
+
+  std::stringstream st;
+  st<<"PID_"<<myrank<<" : ";
+  std::string str;
+
   Node_t* n=this->PostorderFirst();
   while(n!=NULL){
     if(myrank<np-1) if(n->GetMortonId().getDFD()>=mins[myrank+1])break;
@@ -1140,7 +1145,8 @@ bool MPI_Tree<TreeNode>::CheckTree(){
   }
   while(n!=NULL){
     if(n->IsLeaf() && !n->IsGhost()){
-      assert(false);
+      st<<"non-ghost leaf node "<<n->GetMortonId()<<"; after last node.";
+      str=st.str(); ASSERT_WITH_MSG(false,str.c_str());
     }
     n=this->PostorderNxt(n);
   }
@@ -1259,17 +1265,11 @@ void MPI_Tree<TreeNode>::ConstructLET(BoundaryType bndry){
   //std::cout<<"Shared = "<<shared_nodes.size()<<'\n';
 
   // Pack shared nodes.
-#ifdef NDEBUG
   static std::vector<char> shrd_buff_vec0(omp_p*64l*1024l*1024l);
   static std::vector<char> shrd_buff_vec1(omp_p*128l*1024l*1024l);
   static std::vector<char> send_buff_vec(omp_p*64l*1024l*1024l); char* send_buff;
   static std::vector<char> recv_buff_vec(omp_p*64l*1024l*1024l); char* recv_buff;
-#else
-  std::vector<char> shrd_buff_vec0(omp_p*64l*1024l*1024l);
-  std::vector<char> shrd_buff_vec1(omp_p*128l*1024l*1024l);
-  std::vector<char> send_buff_vec(omp_p*64l*1024l*1024l); char* send_buff;
-  std::vector<char> recv_buff_vec(omp_p*64l*1024l*1024l); char* recv_buff;
-#endif
+
   std::vector<PackedData> shrd_data;
   size_t max_data_size=0;
   {

+ 0 - 2
include/pvfmm_common.hpp

@@ -15,8 +15,6 @@
 #define NULL 0
 #endif
 
-//#define LAZY_MAT_COMP //TODO: Implement this.
-
 //Disable assert checks.
 //#define NDEBUG