Parcourir la source

fix GCC warnings, null-terminated strings in fortran

Dhairya Malhotra il y a 5 ans
Parent
commit
1e03c0ab5f
5 fichiers modifiés avec 85 ajouts et 92 suppressions
  1. 2 1
      include/sctl/mem_mgr.txx
  2. 61 68
      include/sctl/tree.hpp
  3. 7 7
      src/libtree.cpp
  4. 6 6
      src/test-cpp.cpp
  5. 9 10
      src/test-fortran.f90

+ 2 - 1
include/sctl/mem_mgr.txx

@@ -122,6 +122,7 @@ inline MemoryManager::~MemoryManager() {
   {  // free buff
     SCTL_ASSERT(buff);
     ::free(buff - ((uint16_t*)buff)[-1]);
+    buff = nullptr;
   }
 }
 
@@ -289,7 +290,7 @@ inline void MemoryManager::free(Iterator<char> p) const {
         SCTL_ASSERT_MSG(p_[i] == init_mem_val, "memory corruption detected.");
       }
     }
-    {  // system_malloc.erase(p_)
+    if (buff != nullptr) {  // system_malloc.erase(p_)
       omp_set_lock(&omp_lock);
       SCTL_ASSERT_MSG(system_malloc.erase(p_) == 1, "double free or corruption.");
       omp_unset_lock(&omp_lock);

+ 61 - 68
include/sctl/tree.hpp

@@ -176,7 +176,7 @@ struct VTUData {
   };
 };
 
-template <class Real, Integer DIM> class Tree {
+template <Integer DIM> class Tree {
   public:
 
     struct NodeAttr {
@@ -194,12 +194,11 @@ template <class Real, Integer DIM> class Tree {
       return DIM;
     }
 
-    Tree(const Comm& comm_ = Comm::Self()) {
-      comm = comm_;
+    Tree(const Comm& comm_ = Comm::Self()) : comm(comm_) {
       Integer rank = comm.Rank();
       Integer np = comm.Size();
 
-      Vector<Real> coord;
+      Vector<double> coord;
       { // Set coord
         Long N0 = 1;
         while (sctl::pow<DIM,Long>(N0) < np) N0++;
@@ -210,7 +209,7 @@ template <class Real, Integer DIM> class Tree {
         for (Long i = start; i < end; i++) {
           Long  idx = i;
           for (Integer k = 0; k < DIM; k++) {
-            coord[(i-start)*DIM+k] = (idx % N0) / (Real)N0;
+            coord[(i-start)*DIM+k] = (idx % N0) / (double)N0;
             idx /= N0;
           }
         }
@@ -242,7 +241,7 @@ template <class Real, Integer DIM> class Tree {
       return comm;
     }
 
-    void UpdateRefinement(const Vector<Real>& coord, Long M = 1, bool balance21 = 0, bool periodic = 0) {
+    template <class Real> void UpdateRefinement(const Vector<Real>& coord, Long M = 1, bool balance21 = 0, bool periodic = 0) {
       Integer np = comm.Size();
       Integer rank = comm.Rank();
 
@@ -503,7 +502,7 @@ template <class Real, Integer DIM> class Tree {
           child_cnt[depth] = 0;
           if (depth) {
             Long p = ancestors[depth-1];
-            Integer& c = child_cnt[depth-1];
+            Long& c = child_cnt[depth-1];
             node_lst[i].parent = p;
             node_lst[p].child[c] = i;
             node_lst[p].p2n = c;
@@ -552,14 +551,14 @@ template <class Real, Integer DIM> class Tree {
         }
 
         Vector<Long> cnt_tmp;
-        Vector<Real> data_tmp;
+        Vector<char> data_tmp;
         for (const auto& pair : node_data) {
           const std::string& data_name = pair.first;
 
           Long dof;
-          Iterator<Vector<Real>> data_;
+          Iterator<Vector<char>> data_;
           Iterator<Vector<Long>> cnt_;
-          GetData(data_, cnt_, data_name);
+          GetData_(data_, cnt_, data_name);
           { // Set dof
             StaticArray<Long,2> Nl, Ng;
             Nl[0] = data_->Dim();
@@ -595,7 +594,7 @@ template <class Real, Integer DIM> class Tree {
       }
     }
 
-    void AddData(const std::string& name, const Vector<Real>& data, const Vector<Long>& cnt) {
+    template <class ValueType> void AddData(const std::string& name, const Vector<ValueType>& data, const Vector<Long>& cnt) {
       Long dof;
       { // Check dof
         StaticArray<Long,2> Nl, Ng;
@@ -609,36 +608,29 @@ template <class Real, Integer DIM> class Tree {
       if (dof) SCTL_ASSERT(cnt.Dim() == node_mid.Dim());
 
       SCTL_ASSERT(node_data.find(name) == node_data.end());
-      node_data[name] = data;
-      node_cnt [name] =  cnt;
+      node_data[name].ReInit(data.Dim()*sizeof(ValueType), (Iterator<char>)data.begin(), true);
+      node_cnt [name] = cnt;
     }
 
-    void GetData(Iterator<Vector<Real>>& data, ConstIterator<Vector<Long>>& cnt, const std::string& name) {
-      auto data_ = node_data.find(name);
-      const auto cnt_ = node_cnt.find(name);
-      SCTL_ASSERT(data_ != node_data.end());
-      SCTL_ASSERT( cnt_ != node_cnt .end());
-      data = Ptr2Itr<Vector<Real>>(&data_->second,1);
-      cnt  = Ptr2ConstItr<Vector<Long>>(& cnt_->second,1);
-    }
-    void GetData(ConstIterator<Vector<Real>>& data, ConstIterator<Vector<Long>>& cnt, const std::string& name) const {
+    template <class ValueType> void GetData(Vector<ValueType>& data, Vector<Long>& cnt, const std::string& name) const {
       const auto data_ = node_data.find(name);
       const auto cnt_ = node_cnt.find(name);
       SCTL_ASSERT(data_ != node_data.end());
       SCTL_ASSERT( cnt_ != node_cnt .end());
-      data = Ptr2ConstItr<Vector<Real>>(&data_->second,1);
-      cnt  = Ptr2ConstItr<Vector<Long>>(& cnt_->second,1);
+      data.ReInit(data_->second.Dim()/sizeof(ValueType), (Iterator<ValueType>)data_->second.begin(), false);
+      SCTL_ASSERT(data.Dim()*(Long)sizeof(ValueType) == data_->second.Dim());
+      cnt .ReInit( cnt_->second.Dim(), (Iterator<Long>)cnt_->second.begin(), false);
     }
 
-    void ReduceBroadcast(const std::string& name) {
+    template <class ValueType> void ReduceBroadcast(const std::string& name) {
       Integer np = comm.Size();
       Integer rank = comm.Rank();
 
       Vector<Long> dsp;
-      Iterator<Vector<Real>> data_;
+      Iterator<Vector<char>> data_;
       Iterator<Vector<Long>> cnt_;
-      GetData(data_, cnt_, name);
-      Vector<Real>& data = *data_;
+      GetData_(data_, cnt_, name);
+      Vector<ValueType> data(data_->Dim()/sizeof(ValueType), (Iterator<ValueType>)data_->begin(), false);
       Vector<Long>& cnt = *cnt_;
       scan(dsp, cnt);
 
@@ -693,7 +685,7 @@ template <class Real, Integer DIM> class Tree {
           scan(recv_data_dsp, recv_data_cnt);
         }
 
-        Vector<Real> send_buff, recv_buff;
+        Vector<ValueType> send_buff, recv_buff;
         Vector<Long> send_buff_cnt(np), send_buff_dsp(np);
         Vector<Long> recv_buff_cnt(np), recv_buff_dsp(np);
         { // Set send_buff, send_buff_cnt, send_buff_dsp, recv_buff, recv_buff_cnt, recv_buff_dsp
@@ -782,7 +774,7 @@ template <class Real, Integer DIM> class Tree {
           scan(recv_data_dsp, recv_data_cnt);
         }
 
-        Vector<Real> send_buff, recv_buff;
+        Vector<ValueType> send_buff, recv_buff;
         Vector<Long> send_buff_cnt(np), send_buff_dsp(np);
         Vector<Long> recv_buff_cnt(np), recv_buff_dsp(np);
         { // Set send_buff, send_buff_cnt, send_buff_dsp, recv_buff, recv_buff_cnt, recv_buff_dsp
@@ -837,9 +829,10 @@ template <class Real, Integer DIM> class Tree {
           Long N1 = (end_idx ? dsp[end_idx-1] + cnt[end_idx-1] : 0) * dof;
           Long Ns = (Nsplit ? recv_data_dsp[Nsplit-1] + recv_data_cnt[Nsplit-1] : 0) * dof;
           if (N0 != Ns || recv_buff.Dim() != N0+data.Dim()-N1) { // resize data and preserve non-ghost data
-            Vector<Real> data_new(recv_buff.Dim() + N1-N0);
-            memcopy(data_new.begin() + Ns, data.begin() + N0, N1-N0);
-            data.Swap(data_new);
+            Vector<char> data_new((recv_buff.Dim() + N1-N0) * sizeof(ValueType));
+            memcopy(data_new.begin() + Ns * sizeof(ValueType), data_->begin() + N0 * sizeof(ValueType), (N1-N0) * sizeof(ValueType));
+            data_->Swap(data_new);
+            data.ReInit(data_->Dim()/sizeof(ValueType), (Iterator<ValueType>)data_->begin(), false);
           }
 
           memcopy(cnt.begin(), recv_data_cnt.begin(), start_idx);
@@ -903,12 +896,12 @@ template <class Real, Integer DIM> class Tree {
 
   protected:
 
-    void GetData(Iterator<Vector<Real>>& data, Iterator<Vector<Long>>& cnt, const std::string& name) {
+    void GetData_(Iterator<Vector<char>>& data, Iterator<Vector<Long>>& cnt, const std::string& name) {
       auto data_ = node_data.find(name);
       const auto cnt_ = node_cnt.find(name);
       SCTL_ASSERT(data_ != node_data.end());
       SCTL_ASSERT( cnt_ != node_cnt .end());
-      data = Ptr2Itr<Vector<Real>>(&data_->second,1);
+      data = Ptr2Itr<Vector<char>>(&data_->second,1);
       cnt  = Ptr2Itr<Vector<Long>>(& cnt_->second,1);
     }
 
@@ -931,7 +924,7 @@ template <class Real, Integer DIM> class Tree {
     Vector<NodeAttr> node_attr;
     Vector<NodeLists> node_lst;
 
-    std::map<std::string, Vector<Real>> node_data;
+    std::map<std::string, Vector<char>> node_data;
     std::map<std::string, Vector<Long>> node_cnt;
 
     Vector<Morton<DIM>> user_mid;
@@ -940,7 +933,7 @@ template <class Real, Integer DIM> class Tree {
     Comm comm;
 };
 
-template <class Real, Integer DIM, class BaseTree = Tree<Real,DIM>> class PtTree : public BaseTree {
+template <class Real, Integer DIM, class BaseTree = Tree<DIM>> class PtTree : public BaseTree {
   public:
 
     PtTree(const Comm& comm = Comm::Self()) : BaseTree(comm) {}
@@ -948,8 +941,8 @@ template <class Real, Integer DIM, class BaseTree = Tree<Real,DIM>> class PtTree
     ~PtTree() {
       #ifdef SCTL_MEMDEBUG
       for (auto& pair : data_pt_name) {
-        ConstIterator<Vector<Real>> data;
-        ConstIterator<Vector<Long>> cnt;
+        Vector<Real> data;
+        Vector<Long> cnt;
         this->GetData(data, cnt, pair.second);
         SCTL_ASSERT(scatter_idx.find(pair.second) != scatter_idx.end());
       }
@@ -992,9 +985,9 @@ template <class Real, Integer DIM, class BaseTree = Tree<Real,DIM>> class PtTree
           if (pair.second == pt_name) {
             const auto& data_name = pair.first;
 
-            Iterator<Vector<Real>> data;
+            Iterator<Vector<char>> data;
             Iterator<Vector<Long>> cnt;
-            this->GetData(data, cnt, data_name);
+            this->GetData_(data, cnt, data_name);
 
             { // Update data
               Long dof = 0;
@@ -1014,7 +1007,7 @@ template <class Real, Integer DIM, class BaseTree = Tree<Real,DIM>> class PtTree
               offset *= dof;
               count *= dof;
 
-              Vector<Real> data_(count, data->begin() + offset);
+              Vector<char> data_(count, data->begin() + offset);
               comm.PartitionN(data_, pt_mid_.Dim());
               data->Swap(data_);
             }
@@ -1046,9 +1039,9 @@ template <class Real, Integer DIM, class BaseTree = Tree<Real,DIM>> class PtTree
       AddParticleData(name, name, coord);
 
       { // Set node_cnt
-        Iterator<Vector<Real>> data_;
+        Iterator<Vector<char>> data_;
         Iterator<Vector<Long>> cnt_;
-        this->GetData(data_,cnt_,name);
+        this->GetData_(data_,cnt_,name);
         cnt_[0].ReInit(node_mid.Dim());
         for (Long i = 0; i < node_mid.Dim(); i++) {
           Long start = std::lower_bound(pt_mid_.begin(), pt_mid_.end(), node_mid[i]) - pt_mid_.begin();
@@ -1065,19 +1058,19 @@ template <class Real, Integer DIM, class BaseTree = Tree<Real,DIM>> class PtTree
       SCTL_ASSERT(data_pt_name.find(data_name) == data_pt_name.end());
       data_pt_name[data_name] = particle_name;
 
-      Iterator<Vector<Real>> data_;
+      Iterator<Vector<char>> data_;
       Iterator<Vector<Long>> cnt_;
       this->AddData(data_name, Vector<Real>(), Vector<Long>());
-      this->GetData(data_,cnt_,data_name);
+      this->GetData_(data_,cnt_,data_name);
       { // Set data_[0]
-        data_[0] = data;
+        data_[0].ReInit(data.Dim()*sizeof(Real), (Iterator<char>)data.begin(), true);
         this->GetComm().ScatterForward(data_[0], scatter_idx[particle_name]);
       }
       if (data_name != particle_name) { // Set cnt_[0]
-        Iterator<Vector<Real>> pt_coord;
-        Iterator<Vector<Long>> pt_cnt;
+        Vector<Real> pt_coord;
+        Vector<Long> pt_cnt;
         this->GetData(pt_coord, pt_cnt, particle_name);
-        cnt_[0] = pt_cnt[0];
+        cnt_[0] = pt_cnt;
       }
     }
 
@@ -1094,14 +1087,14 @@ template <class Real, Integer DIM, class BaseTree = Tree<Real,DIM>> class PtTree
 
       Long dof;
       Vector<Long> dsp;
-      ConstIterator<Vector<Long>> cnt_;
-      ConstIterator<Vector<Real>> data_;
+      Vector<Long> cnt_;
+      Vector<Real> data_;
       this->GetData(data_, cnt_, data_name);
-      SCTL_ASSERT(cnt_->Dim() == node_mid.Dim());
-      BaseTree::scan(dsp, cnt_[0]);
+      SCTL_ASSERT(cnt_.Dim() == node_mid.Dim());
+      BaseTree::scan(dsp, cnt_);
       { // Set dof
         Long Nn = node_mid.Dim();
-        StaticArray<Long,2> Ng, Nl = {data_->Dim(), dsp[Nn-1]+cnt_[0][Nn-1]};
+        StaticArray<Long,2> Ng, Nl = {data_.Dim(), dsp[Nn-1]+cnt_[Nn-1]};
         comm.Allreduce((ConstIterator<Long>)Nl, (Iterator<Long>)Ng, 2, Comm::CommOp::SUM);
         dof = Ng[0] / std::max<Long>(Ng[1],1);
       }
@@ -1111,8 +1104,8 @@ template <class Real, Integer DIM, class BaseTree = Tree<Real,DIM>> class PtTree
         Long N0 = std::lower_bound(node_mid.begin(), node_mid.end(), mins[rank]) - node_mid.begin();
         Long N1 = std::lower_bound(node_mid.begin(), node_mid.end(), (rank+1==np ? Morton<DIM>().Next() : mins[rank+1])) - node_mid.begin();
         Long start = dsp[N0] * dof;
-        Long end = (N1<dsp.Dim() ? dsp[N1] : dsp[N1-1]+cnt_[0][N1-1]) * dof;
-        data.ReInit(end-start, (Iterator<Real>)data_->begin()+start, true);
+        Long end = (N1<dsp.Dim() ? dsp[N1] : dsp[N1-1]+cnt_[N1-1]) * dof;
+        data.ReInit(end-start, data_.begin()+start, true);
         comm.ScatterReverse(data, scatter_idx_, Nlocal_ * dof);
       }
     }
@@ -1148,20 +1141,20 @@ template <class Real, Integer DIM, class BaseTree = Tree<Real,DIM>> class PtTree
         SCTL_ASSERT(data_pt_name.find(data_name) != data_pt_name.end());
         std::string particle_name = data_pt_name.find(data_name)->second;
 
-        ConstIterator<Vector<Real>> pt_coord = NullIterator<Vector<Real>>();
-        ConstIterator<Vector<Real>> pt_value = NullIterator<Vector<Real>>();
-        ConstIterator<Vector<Long>> pt_cnt = NullIterator<Vector<Long>>();
+        Vector<Real> pt_coord;
+        Vector<Real> pt_value;
+        Vector<Long> pt_cnt;
         Vector<Long> pt_dsp;
         Long value_dof = 0;
         { // Set pt_coord, pt_cnt, pt_dsp
           this->GetData(pt_coord, pt_cnt, particle_name);
-          Tree<Real,DIM>::scan(pt_dsp, pt_cnt[0]);
+          Tree<DIM>::scan(pt_dsp, pt_cnt);
         }
         if (particle_name != data_name) { // Set pt_value, value_dof
-          ConstIterator<Vector<Long>> pt_cnt = NullIterator<Vector<Long>>();
+          Vector<Long> pt_cnt;
           this->GetData(pt_value, pt_cnt, data_name);
-          Long Npt = omp_par::reduce(pt_cnt->begin(), pt_cnt->Dim());
-          value_dof = pt_value->Dim() / std::max<Long>(Npt,1);
+          Long Npt = omp_par::reduce(pt_cnt.begin(), pt_cnt.Dim());
+          value_dof = pt_value.Dim() / std::max<Long>(Npt,1);
         }
 
         Vector<VTKReal> &coord = vtu_data.coord;
@@ -1177,14 +1170,14 @@ template <class Real, Integer DIM, class BaseTree = Tree<Real,DIM>> class PtTree
         value.SetZero();
 
         SCTL_ASSERT(node_mid.Dim() == node_attr.Dim());
-        SCTL_ASSERT(node_mid.Dim() == pt_cnt->Dim());
+        SCTL_ASSERT(node_mid.Dim() == pt_cnt.Dim());
         for (Long i = 0; i < node_mid.Dim(); i++) {
           if (!show_ghost && node_attr[i].Ghost) continue;
           if (!node_attr[i].Leaf) continue;
 
-          for (Long j = 0; j < pt_cnt[0][i]; j++) {
-            ConstIterator<Real> pt_coord_ = pt_coord->begin() + (pt_dsp[i] + j) * DIM;
-            ConstIterator<Real> pt_value_ = (value_dof ? pt_value->begin() + (pt_dsp[i] + j) * value_dof : NullIterator<Real>());
+          for (Long j = 0; j < pt_cnt[i]; j++) {
+            ConstIterator<Real> pt_coord_ = pt_coord.begin() + (pt_dsp[i] + j) * DIM;
+            ConstIterator<Real> pt_value_ = (value_dof ? pt_value.begin() + (pt_dsp[i] + j) * value_dof : NullIterator<Real>());
 
             for (Integer k = 0; k < DIM; k++) coord.PushBack((VTKReal)pt_coord_[k]);
             for (Integer k = DIM; k < 3; k++) coord.PushBack(0);

+ 7 - 7
src/libtree.cpp

@@ -37,7 +37,7 @@ template <class Real, Integer DIM> class ParticleTreeWrapper {
         child_cnt[depth] = 0;
         if (depth) {
           Long p = ancestors[depth-1];
-          Integer& c = child_cnt[depth-1];
+          Long& c = child_cnt[depth-1];
           parent_lst_[i] = p+1; // +1 because indexing starts from 1 in fortran
           child_lst_[p*MAX_CHILD+c] = i+1; // +1 because indexing starts from 1 in fortran
           c++;
@@ -66,16 +66,16 @@ template <class Real, Integer DIM> class ParticleTreeWrapper {
     }
 
     void GetData(Real*& data, int32_t& Ndata, const int32_t*& cnt, int32_t& Ncnt, const char* name) {
-      Iterator<Vector<Real>> data_;
-      ConstIterator<Vector<Long>> cnt_;
+      Vector<Real> data_;
+      Vector<Long> cnt_;
       tree_.GetData(data_, cnt_, name);
-      data = &data_[0][0];
-      Ndata = data_->Dim();
-      Ncnt = cnt_->Dim();
+      data = &data_[0];
+      Ndata = data_.Dim();
+      Ncnt = cnt_.Dim();
 
       Vector<int32_t>& data_cnt = data_cnt_[name];
       data_cnt.ReInit(Ncnt);
-      for (Long i = 0; i < Ncnt; i++) data_cnt[i] = cnt_[0][i];
+      for (Long i = 0; i < Ncnt; i++) data_cnt[i] = cnt_[i];
       cnt = &data_cnt[0];
     }
 

+ 6 - 6
src/test-cpp.cpp

@@ -41,6 +41,8 @@ int main(int argc, char** argv) {
     t.AddParticles("src_pts", coord);
     sctl::Profile::Toc();
 
+    t.ReduceBroadcast<Real>("src_pts");
+
     sctl::Profile::Tic("UpdateRefine",&comm);
     t.UpdateRefinement(coord, 100, true, false);
     sctl::Profile::Toc();
@@ -92,13 +94,11 @@ int main(int argc, char** argv) {
         }
         t.AddData("mult", node_data, node_cnt);
       }
-      t.ReduceBroadcast("mult");
+      t.ReduceBroadcast<Real>("mult");
       { // verify
-        sctl::ConstIterator<sctl::Vector<Real>> node_data_;
-        sctl::ConstIterator<sctl::Vector<sctl::Long>> node_cnt_;
-        t.GetData(node_data_, node_cnt_, "mult");
-        const sctl::Vector<Real>& node_data = *node_data_;
-        const sctl::Vector<sctl::Long>& node_cnt = *node_cnt_;
+        sctl::Vector<Real> node_data;
+        sctl::Vector<sctl::Long> node_cnt;
+        t.GetData(node_data, node_cnt, "mult");
         sctl::Vector<sctl::Long> node_dsp;
         scan(node_dsp, node_cnt);
 

+ 9 - 10
src/test-fortran.f90

@@ -6,7 +6,6 @@ program main
   type(c_ptr) :: tree_ctx
   real*8 :: pt_coord(100000*3)
   integer*4 :: Npt, mpi_rank, mpi_size, i, ierr
-  real*4 :: rand
 
   integer*4 :: Ndata, Ncnt
   real*8, pointer :: pt_data(:)
@@ -53,31 +52,31 @@ program main
 
 
   ! Add particles to the tree and associate some data with each particle
-  call AddParticles('src_pts', pt_coord, Npt, tree_ctx)
-  call AddParticleData('src_coord', 'src_pts', pt_coord, Npt*3, tree_ctx)
+  call AddParticles('src_pts'//C_NULL_CHAR, pt_coord, Npt, tree_ctx)
+  call AddParticleData('src_coord'//C_NULL_CHAR, 'src_pts'//C_NULL_CHAR, pt_coord, Npt*3, tree_ctx)
 
   ! Read the data and verify it is the same as original array
-  call GetParticleData(pt_data_, Ndata, 'src_coord', tree_ctx)
+  call GetParticleData(pt_data_, Ndata, 'src_coord'//C_NULL_CHAR, tree_ctx)
   call c_f_pointer(pt_data_, pt_data, [Ndata])
   call compare_arrays(pt_coord, Npt*3, pt_data, Ndata)
 
 
 
   ! Read node-wise particle data
-  call GetData(pt_data_, Ndata, pt_cnt_, Ncnt, 'src_coord', tree_ctx)
+  call GetData(pt_data_, Ndata, pt_cnt_, Ncnt, 'src_coord'//C_NULL_CHAR, tree_ctx)
   call c_f_pointer(pt_data_, pt_data, [Ndata])
   call c_f_pointer(pt_cnt_, pt_cnt, [Ncnt])
 
 
 
   ! Write tree node structure to VTK file
-  call WriteTreeVTK('tree', .false., tree_ctx)
+  call WriteTreeVTK('tree'//C_NULL_CHAR, .false., tree_ctx)
 
   ! Write points and associated data to VTK file
-  call WriteParticleVTK('pts', 'src_coord', .false., tree_ctx)
+  call WriteParticleVTK('pts'//C_NULL_CHAR, 'src_coord'//C_NULL_CHAR, .false., tree_ctx)
 
   ! Delete the particles (or the associated data)
-  call DeleteParticleData('src_pts', tree_ctx)
+  call DeleteParticleData('src_pts'//C_NULL_CHAR, tree_ctx)
 
 
 
@@ -87,8 +86,8 @@ program main
   do i = 1, Nnodes
     mult_cnt(i) = 1
   enddo
-  call AddData('mult_coeff', mult_coeff, Nnodes*10, mult_cnt, Nnodes, tree_ctx)
-  call DeleteData('mult_coeff', tree_ctx)
+  call AddData('mult_coeff'//C_NULL_CHAR, mult_coeff, Nnodes*10, mult_cnt, Nnodes, tree_ctx)
+  call DeleteData('mult_coeff'//C_NULL_CHAR, tree_ctx)