Ignore:
Timestamp:
Apr 24, 2012, 2:26:14 PM (14 years ago)
Author:
Julian Iseringhausen <isering@…>
Children:
b51c3b
Parents:
e3dbbf
Message:

Fix energy calculation.

git-svn-id: https://svn.version.fz-juelich.de/scafacos/trunk@1729 5161e1c8-67bf-11de-9fd5-51895aff932f

File:
1 edited

Legend:

Unmodified
Added
Removed
  • src/comm/comm_mpi.cpp

    re3dbbf r716da7  
    232232  MPI_Alltoall(&send_sizes.front(), 1, MPI_INT, &recv_sizes.front(), 1, MPI_INT, comm_global);
    233233
     234  assert(RequestsPending() == 0);
     235
    234236  /*
    235237   * Send particles
     
    281283
    282284#ifdef VMG_ONE_SIDED
    283   // if (!win_created) {
    284   //   vmg_float* p = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_POTENTIAL_ARRAY");
    285   //   const vmg_int& num_particles_local = MG::GetFactory().GetObjectStorageVal<vmg_int>("PARTICLE_NUM_LOCAL");
    286   //   MPI_Win_create(p, num_particles_local*sizeof(vmg_float), sizeof(vmg_float), info, comm_global, &win);
    287   //   win_created = true;
    288   // }
    289 
    290   // MPI_Win_fence(MPI_MODE_NOPRECEDE, win);
    291 
    292   // for (iter=particles.begin(); iter!=particles.end(); ++iter)
    293   //   MPI_Put(&iter->Pot(), 1, MPI_DOUBLE, iter->Rank(), iter->Index(), 1, MPI_DOUBLE, win);
    294 
    295   // MPI_Win_fence(MPI_MODE_NOSTORE | MPI_MODE_NOSUCCEED, win);
     285  if (!win_created) {
     286    vmg_float* p = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_POTENTIAL_ARRAY");
     287    const vmg_int& num_particles_local = MG::GetFactory().GetObjectStorageVal<vmg_int>("PARTICLE_NUM_LOCAL");
     288    MPI_Win_create(p, num_particles_local*sizeof(vmg_float), sizeof(vmg_float), info, comm_global, &win);
     289    win_created = true;
     290  }
     291
     292  MPI_Win_fence(MPI_MODE_NOPRECEDE, win);
     293
     294  for (iter=particles.begin(); iter!=particles.end(); ++iter)
     295    MPI_Put(&iter->Pot(), 1, MPI_DOUBLE, iter->Rank(), iter->Index(), 1, MPI_DOUBLE, win);
     296
     297  MPI_Win_fence(MPI_MODE_NOSTORE | MPI_MODE_NOSUCCEED, win);
    296298#else
    297299  int rank, size;
     
    299301  MPI_Comm_size(comm_global, &size);
    300302
    301   std::vector< std::vector<vmg_float> > send_buffer_pot(size);
    302   std::vector< std::vector<vmg_float> > recv_buffer_pot(size);
     303  std::vector< std::vector<vmg_float> > send_buffer_float(size);
     304  std::vector< std::vector<vmg_float> > recv_buffer_float(size);
    303305  std::vector< std::vector<vmg_int> > send_buffer_index(size);
    304306  std::vector< std::vector<vmg_int> > recv_buffer_index(size);
     
    306308  vmg_int* size_receive = MG::GetFactory().GetObjectStorageArray<vmg_int>("PARTICLE_RECEIVER_ARRAY");
    307309  vmg_float* p = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_POTENTIAL_ARRAY");
     310  vmg_float* f = MG::GetFactory().GetObjectStorageArray<vmg_float>("PARTICLE_FIELD_ARRAY");
    308311
    309312  // Build send buffer
    310313  for (iter=particles.begin(); iter!=particles.end(); ++iter) {
    311     send_buffer_pot[iter->Rank()].push_back(iter->Pot());
     314    send_buffer_float[iter->Rank()].push_back(iter->Pot());
     315    send_buffer_float[iter->Rank()].push_back(iter->Field()[0]);
     316    send_buffer_float[iter->Rank()].push_back(iter->Field()[1]);
     317    send_buffer_float[iter->Rank()].push_back(iter->Field()[2]);
    312318    send_buffer_index[iter->Rank()].push_back(iter->Index());
    313319  }
     
    315321  // Send potentials
    316322  for (int i=0; i<size; ++i) {
    317     if (!send_buffer_pot[i].empty()) {
    318       MPI_Isend(&send_buffer_pot[i].front(), send_buffer_pot[i].size(), MPI_DOUBLE, i, 699+rank, comm_global, &Request());
     323    if (!send_buffer_float[i].empty()) {
     324      MPI_Isend(&send_buffer_float[i].front(), send_buffer_float[i].size(), MPI_DOUBLE, i, 699+rank, comm_global, &Request());
    319325      MPI_Isend(&send_buffer_index[i].front(), send_buffer_index[i].size(), MPI_INT, i, 32111+rank, comm_global, &Request());
    320326    }
     
    324330  for (int i=0; i<size; ++i) {
    325331    if (size_receive[i] > 0) {
    326       recv_buffer_pot[i].resize(size_receive[i]);
     332      recv_buffer_float[i].resize(4*size_receive[i]);
    327333      recv_buffer_index[i].resize(size_receive[i]);
    328       MPI_Irecv(&recv_buffer_pot[i].front(), size_receive[i], MPI_DOUBLE, i, 699+i, comm_global, &Request());
     334      MPI_Irecv(&recv_buffer_float[i].front(), 4*size_receive[i], MPI_DOUBLE, i, 699+i, comm_global, &Request());
    329335      MPI_Irecv(&recv_buffer_index[i].front(), size_receive[i], MPI_INT, i, 32111+i, comm_global, &Request());
    330336    }
     
    335341  // Add potential values
    336342  for (int i=0; i<size; ++i)
    337     for (unsigned int j=0; j<size_receive[i]; ++j)
    338       p[recv_buffer_index[i][j]] = recv_buffer_pot[i][j];
     343    for (unsigned int j=0; j<size_receive[i]; ++j) {
     344      p[recv_buffer_index[i][j]] = recv_buffer_float[i][4*j];
     345      std::memcpy(&f[recv_buffer_index[i][j]], &recv_buffer_float[i][4*j+1], 3*sizeof(vmg_float));
     346    }
    339347#endif
    340348
    341349}
    342350
    343 void CommMPI::CommLCListToGhosts(const Grid& grid, Particle::LinkedCellList& lc)
     351void CommMPI::CommLCListToGhosts(Particle::LinkedCellList& lc)
    344352{
    345353  VMG::MPI::DatatypesLocal types(lc, comm_global);
    346   std::vector< std::vector<vmg_float> > send_buffer(types.NB().size());
    347   std::vector< std::vector<vmg_float> > recv_buffer(types.Halo().size());
    348   std::vector<vmg_int> send_size(types.NB().size());
     354  std::vector<int> send_size(types.NB().size());
    349355  vmg_int recv_size;
    350356  std::list<Particle::Particle*>::iterator iter;
    351357  Index ind;
     358  Vector offset;
     359
     360  const Vector halo_length = lc.Local().HaloSize1() * lc.Extent().MeshWidth();
    352361
    353362  lc.ClearHalo();
    354 
    355   for (int i=0; i<3; ++i) {
    356     for (Grid::iterator iter=lc.Iterators().Halo1()[i].Begin(); iter!=lc.Iterators().Halo1()[i].End(); ++iter)
    357       assert(lc(*iter).size() == 0);
    358     for (Grid::iterator iter=lc.Iterators().Halo2()[i].Begin(); iter!=lc.Iterators().Halo2()[i].End(); ++iter)
    359       assert(lc(*iter).size() == 0);
    360   }
    361363
    362364  for (unsigned int i=0; i<types.NB().size(); ++i)
    363365    if (types.NB()[i].Feasible()) {
     366
     367      for (int j=0; j<3; ++j)
     368        if ((types.Offset()[i][j] < 0 && lc.Global().LocalBegin()[j] == 0) ||
     369            (types.Offset()[i][j] > 0 && lc.Global().LocalEnd()[j] == lc.Global().GlobalSize()[j]))
     370          offset[j] = -1 * types.Offset()[i][j] * lc.Extent().Size()[j];
     371        else
     372          offset[j] = 0;
     373
    364374      for (ind.X() = types.NB()[i].Starts().X(); ind.X() < types.NB()[i].Starts().X()+types.NB()[i].Subsizes().X(); ++ind.X())
    365375        for (ind.Y() = types.NB()[i].Starts().Y(); ind.Y() < types.NB()[i].Starts().Y()+types.NB()[i].Subsizes().Y(); ++ind.Y())
    366376          for (ind.Z() = types.NB()[i].Starts().Z(); ind.Z() < types.NB()[i].Starts().Z()+types.NB()[i].Subsizes().Z(); ++ind.Z())
    367377            for (iter=lc(ind).begin(); iter!=lc(ind).end(); ++iter) {
    368               if (grid.Global().LocalBegin()[i] == 0)
    369                 for (int j=0; j<3; ++j)
    370                   send_buffer[i].push_back((*iter)->Pos()[j] + (i==j ? grid.Extent().Size()[j] : 0.0));
    371               else
    372                 for (int j=0; j<3; ++j)
    373                   send_buffer[1].push_back((*iter)->Pos()[j]);
    374 
    375               send_buffer[1].push_back((*iter)->Charge());
     378              for (int j=0; j<3; ++j)
     379                types.NB()[i].Buffer().push_back((*iter)->Pos()[j] + offset[j]);
     380              types.NB()[i].Buffer().push_back((*iter)->Charge());
     381
     382              assert(lc.Extent().Begin().IsComponentwiseLessOrEqual((*iter)->Pos()));
     383              assert(lc.Extent().End().IsComponentwiseGreaterOrEqual((*iter)->Pos()));
     384              assert(lc.Extent().Begin().IsComponentwiseLessOrEqual((*iter)->Pos() + offset + halo_length));
     385              assert(lc.Extent().End().IsComponentwiseGreaterOrEqual((*iter)->Pos() + offset - halo_length));
    376386            }
    377387
     388      send_size[i] = types.NB()[i].Buffer().size();
    378389      MPI_Isend(&send_size[i], 1, MPI_INT, types.NB()[i].Rank(), 2048+types.NB()[i].TagSend(), comm_global, &Request());
    379390
    380391      if (send_size[i] > 0)
    381         MPI_Isend(&send_buffer[i].front(), send_size[i], MPI_DOUBLE,
     392        MPI_Isend(&types.NB()[i].Buffer().front(), send_size[i], MPI_DOUBLE,
    382393                  types.NB()[i].Rank(), 4096+types.NB()[i].TagSend(),
    383394                  comm_global, &Request());
     
    388399      MPI_Recv(&recv_size, 1, MPI_INT, types.Halo()[i].Rank(), 2048+types.Halo()[i].TagReceive(), comm_global, MPI_STATUS_IGNORE);
    389400      if (recv_size > 0) {
    390         recv_buffer[i].resize(recv_size);
    391         MPI_Irecv(&recv_buffer[i].front(), recv_size, MPI_DOUBLE,
     401        types.Halo()[i].Buffer().resize(recv_size);
     402        MPI_Irecv(&types.Halo()[i].Buffer().front(), recv_size, MPI_DOUBLE,
    392403                  types.Halo()[i].Rank(), 4096+types.Halo()[i].TagReceive(),
    393404                  comm_global, &Request());
     
    397408  WaitAll();
    398409
    399   for (unsigned int i=0; i<recv_buffer.size(); ++i)
    400     for (unsigned int j=0; j<recv_buffer[i].size(); j+=4)
    401       lc.AddParticleToHalo(&recv_buffer[i][j], recv_buffer[i][j+3]);
     410  for (unsigned int i=0; i<types.Halo().size(); ++i)
     411    for (unsigned int j=0; j<types.Halo()[i].Buffer().size(); j+=4)
     412      lc.AddParticleToHalo(&types.Halo()[i].Buffer()[j], types.Halo()[i].Buffer()[j+3]);
    402413}
    403414
     
    712723        << "    BOUNDARY_END_2:       " << grid.Local().BoundaryEnd2() << std::endl
    713724        << "    BOUNDARY_SIZE_2:      " << grid.Local().BoundarySize2() << std::endl
    714         << "    FINER_BEGIN:          " << grid.Local().FinerBeginFoo() << std::endl
    715         << "    FINER_END:            " << grid.Local().FinerEndFoo() << std::endl
    716         << "    FINER_SIZE:           " << grid.Local().FinerSizeFoo() << std::endl;
     725        << "    FINER_BEGIN:          " << grid.Local().FinerBegin() << std::endl
     726        << "    FINER_END:            " << grid.Local().FinerEnd() << std::endl
     727        << "    FINER_SIZE:           " << grid.Local().FinerSize() << std::endl;
    717728
    718729    if (rank == size-1)
Note: See TracChangeset for help on using the changeset viewer.