Exemplo n.º 1
0
Particle *move_indexed_particle(ParticleList *dl, ParticleList *sl, int i)
{
  int re = realloc_particlelist(dl, ++dl->n);
  Particle *dst = &dl->part[dl->n - 1];
  Particle *src = &sl->part[i];
  Particle *end = &sl->part[sl->n - 1];

  memcpy(dst, src, sizeof(Particle));
  if (re) {
    //fprintf(stderr, "%d: m_i_p: update destination list after realloc\n",this_node);
    update_local_particles(dl); }
  else {
    //fprintf(stderr, "%d: m_i_p: update loc_part entry for moved particle (id %d)\n",this_node,dst->p.identity);
    local_particles[dst->p.identity] = dst;
  }
  if ( src != end ) {
    //fprintf(stderr, "%d: m_i_p: copy end particle in source list (id %d)\n",this_node,end->p.identity);
    memcpy(src, end, sizeof(Particle));

  }
  if (realloc_particlelist(sl, --sl->n)) {
    //fprintf(stderr, "%d: m_i_p: update source list after realloc\n",this_node);
    update_local_particles(sl); }
  else if ( src != end ) {
    //fprintf(stderr, "%d: m_i_p: update loc_part entry for end particle (id %d)\n",this_node,src->p.identity);
    local_particles[src->p.identity] = src; }
  return dst;
}
Exemplo n.º 2
0
Particle *move_unindexed_particle(ParticleList *dl, ParticleList *sl, int i)
{
  Particle *dst, *src, *end;

  realloc_particlelist(dl, ++dl->n);
  dst = &dl->part[dl->n - 1];
  src = &sl->part[i];
  end = &sl->part[sl->n - 1];
  memcpy(dst, src, sizeof(Particle));
  if ( src != end )
    memcpy(src, end, sizeof(Particle));
  sl->n -= 1;
  realloc_particlelist(sl, sl->n);
  return dst;
}
Exemplo n.º 3
0
static void prepare_ghost_cell(Cell *cell, int size)
{
#ifdef GHOSTS_HAVE_BONDS
  // free all allocated information, will be resent
  {
    int np   = cell->n;
    Particle *part = cell->part;
    for (int p = 0; p < np; p++) {
      free_particle(part + p);
    }
  }          
#endif
  realloc_particlelist(cell, cell->n = size);
  // invalidate pointers etc
  {
    int np   = cell->n;
    Particle *part = cell->part;
    for (int p = 0; p < np; p++) {
      Particle *pt = &part[p];
      // no bonds or exclusions
      pt->bl.e = 0;
      pt->bl.n = 0;
      pt->bl.max = 0;
#ifdef EXCLUSIONS
      pt->el.e = 0;
      pt->el.n = 0;
      pt->el.max = 0;
#endif
#ifdef GHOST_FLAG
      //init ghost variable
      pt->l.ghost=1;
#endif
    }
  }
}
Exemplo n.º 4
0
void cells_re_init(int new_cs)
{
  CellPList tmp_local;
  Cell *tmp_cells;
  int tmp_n_cells,i;

  CELL_TRACE(fprintf(stderr, "%d: cells_re_init: convert type (%d->%d)\n", this_node, cell_structure.type, new_cs));

  invalidate_ghosts();

  /* 
     CELL_TRACE({
     int p;
     for (p = 0; p < n_total_particles; p++)
     if (local_particles[p])
     fprintf(stderr, "%d: cells_re_init: got particle %d\n", this_node, p);
     }
     );
  */

  topology_release(cell_structure.type);
  /* MOVE old local_cell list to temporary buffer */
  memcpy(&tmp_local,&local_cells,sizeof(CellPList));
  init_cellplist(&local_cells);

  /* MOVE old cells to temporary buffer */
  tmp_cells   = cells;
  tmp_n_cells = n_cells;
  cells   = NULL;
  n_cells = 0;

  topology_init(new_cs, &tmp_local);

  particle_invalidate_part_node();

  /* finally deallocate the old cells */
  realloc_cellplist(&tmp_local,0);
  for(i=0;i<tmp_n_cells;i++)
    realloc_particlelist(&tmp_cells[i],0);

  free(tmp_cells);
  CELL_TRACE(fprintf(stderr, "%d: old cells deallocated\n",this_node));

  /*
    CELL_TRACE({
    int p;
    for (p = 0; p < n_total_particles; p++)
    if (local_particles[p])
    fprintf(stderr, "%d: cells_re_init: now got particle %d\n", this_node, p);
    }
    );
  */

  /* to enforce initialization of the ghost cells */
  resort_particles = 1;

#ifdef ADDITIONAL_CHECKS
  check_cells_consistency();
#endif
}
Exemplo n.º 5
0
void local_place_particle(int part, double p[3], int _new)
{
  Cell *cell;
  double pp[3];
  int i[3], rl;
  Particle *pt;

  i[0] = 0;
  i[1] = 0;
  i[2] = 0;
  pp[0] = p[0];
  pp[1] = p[1];
  pp[2] = p[2];

#ifdef LEES_EDWARDS
  double vv[3]={0.,0.,0.};
  fold_position(pp, vv, i);
#else
  fold_position(pp, i);
#endif 
  
  if (_new) {
    /* allocate particle anew */
    cell = cell_structure.position_to_cell(pp);
    if (!cell) {
      fprintf(stderr, "%d: INTERNAL ERROR: particle %d at %f(%f) %f(%f) %f(%f) does not belong on this node\n",
	      this_node, part, p[0], pp[0], p[1], pp[1], p[2], pp[2]);
      errexit();
    }
    rl = realloc_particlelist(cell, ++cell->n);
    pt = &cell->part[cell->n - 1];
    init_particle(pt);

    pt->p.identity = part;
    if (rl)
      update_local_particles(cell);
    else
      local_particles[pt->p.identity] = pt;
  }
  else
    pt = local_particles[part];

  PART_TRACE(fprintf(stderr, "%d: local_place_particle: got particle id=%d @ %f %f %f\n",
		     this_node, part, p[0], p[1], p[2]));

#ifdef LEES_EDWARDS
  pt->m.v[0] += vv[0];  
  pt->m.v[1] += vv[1];  
  pt->m.v[2] += vv[2];  
#endif

  memcpy(pt->r.p, pp, 3*sizeof(double));
  memcpy(pt->l.i, i, 3*sizeof(int));
#ifdef BOND_CONSTRAINT
  memcpy(pt->r.p_old, pp, 3*sizeof(double));
#endif
}
Exemplo n.º 6
0
Particle *append_unindexed_particle(ParticleList *l, Particle *part)
{
  Particle *p;

  realloc_particlelist(l, ++l->n);
  p = &l->part[l->n - 1];

  memcpy(p, part, sizeof(Particle));
  return p;
}
Exemplo n.º 7
0
void cell_cell_transfer(GhostCommunication *gc, int data_parts)
{
    int pl, p, np, offset;
    Particle *part1, *part2, *pt1, *pt2;

    GHOST_TRACE(fprintf(stderr, "%d: local_transfer: type %d data_parts %d\n", this_node,gc->type,data_parts));

    /* transfer data */
    offset = gc->n_part_lists/2;
    for (pl = 0; pl < offset; pl++) {
        np   = gc->part_lists[pl]->n;
        if (data_parts == GHOSTTRANS_PARTNUM) {
            realloc_particlelist(gc->part_lists[pl + offset],
                                 gc->part_lists[pl + offset]->n = gc->part_lists[pl]->n);
#ifdef GHOST_FLAG
            {
                //init ghost variable
                int i;
                for (i=0; i<gc->part_lists[pl + offset]->n; i++) {
                    gc->part_lists[pl + offset]->part[i].l.ghost=1;
                }
            }
#endif
        }
        else {
            part1 = gc->part_lists[pl]->part;
            part2 = gc->part_lists[pl + offset]->part;
            for (p = 0; p < np; p++) {
                pt1 = &part1[p];
                pt2 = &part2[p];
                if (data_parts & GHOSTTRANS_PROPRTS)
                    memcpy(&pt2->p, &pt1->p, sizeof(ParticleProperties));
                if (data_parts & GHOSTTRANS_POSSHFTD) {
                    /* ok, this is not nice, but perhaps fast */
                    int i;
                    memcpy(&pt2->r, &pt1->r, sizeof(ParticlePosition));
                    for (i = 0; i < 3; i++)
                        pt2->r.p[i] += gc->shift[i];
                }
                else if (data_parts & GHOSTTRANS_POSITION)
                    memcpy(&pt2->r, &pt1->r, sizeof(ParticlePosition));
                if (data_parts & GHOSTTRANS_MOMENTUM)
                    memcpy(&pt2->m, &pt1->m, sizeof(ParticleMomentum));
                if (data_parts & GHOSTTRANS_FORCE)
                    add_force(&pt2->f, &pt1->f);
#ifdef LB
                if (data_parts & GHOSTTRANS_COUPLING)
                    memcpy(&pt2->lc, &pt1->lc, sizeof(ParticleLatticeCoupling));
#endif
            }
        }
    }
}
Exemplo n.º 8
0
Particle *append_indexed_particle(ParticleList *l, Particle *part)
{
  int re;
  Particle *p;

  re = realloc_particlelist(l, ++l->n);
  p  = &l->part[l->n - 1];

  memcpy(p, part, sizeof(Particle));

  if (re)
    update_local_particles(l);
  else
    local_particles[p->p.identity] = p;
  return p;
}
Exemplo n.º 9
0
void realloc_cells(int size)
{
  int i;
  CELL_TRACE(fprintf(stderr, "%d: realloc_cells %d\n", this_node, size));
  /* free all memory associated with cells to be deleted. */
  for(i=size; i<n_cells; i++) {
    realloc_particlelist(&cells[i],0);
  }
  /* resize the cell list */
  if(size != n_cells) {
    cells = (Cell *) realloc(cells, sizeof(Cell)*size);
  }
  /* initialize new cells */
  for(i=n_cells; i<size; i++) {
    init_particlelist(&cells[i]);
  }
  n_cells = size;
}  
Exemplo n.º 10
0
void send_particles(ParticleList *particles, int node)
{
  int pc;
  /* Dynamic data, bonds and exclusions */
  IntList local_dyn;

  PART_TRACE(fprintf(stderr, "%d: send_particles %d to %d\n", this_node, particles->n, node));

  MPI_Send(&particles->n, 1, MPI_INT, node, REQ_SNDRCV_PART, comm_cart);
  MPI_Send(particles->part, particles->n*sizeof(Particle),
	   MPI_BYTE, node, REQ_SNDRCV_PART, comm_cart);

  init_intlist(&local_dyn);
  for (pc = 0; pc < particles->n; pc++) {
    Particle *p = &particles->part[pc];
    int size =  local_dyn.n + p->bl.n;
#ifdef EXCLUSIONS
    size += p->el.n;
#endif
    realloc_intlist(&local_dyn, size);
    memcpy(local_dyn.e + local_dyn.n, p->bl.e, p->bl.n*sizeof(int));
    local_dyn.n += p->bl.n;
#ifdef EXCLUSIONS
    memcpy(local_dyn.e + local_dyn.n, p->el.e, p->el.n*sizeof(int));
    local_dyn.n += p->el.n;
#endif
  }

  PART_TRACE(fprintf(stderr, "%d: send_particles sending %d bond ints\n", this_node, local_dyn.n));
  if (local_dyn.n > 0) {
    MPI_Send(local_dyn.e, local_dyn.n*sizeof(int),
	     MPI_BYTE, node, REQ_SNDRCV_PART, comm_cart);
    realloc_intlist(&local_dyn, 0);
  }

  /* remove particles from this nodes local list and free data */
  for (pc = 0; pc < particles->n; pc++) {
    local_particles[particles->part[pc].p.identity] = NULL;
    free_particle(&particles->part[pc]);
  }

  realloc_particlelist(particles, particles->n = 0);
}
Exemplo n.º 11
0
void put_recv_buffer(GhostCommunication *gc, int data_parts)
{
  int pl, p, np;
  Particle *part, *pt;
  char *retrieve;

  /* put back data */
  retrieve = r_buffer;
  for (pl = 0; pl < gc->n_part_lists; pl++) {
    if (data_parts == GHOSTTRANS_PARTNUM) {
      GHOST_TRACE(fprintf(stderr, "%d: reallocating cell %p to size %d, assigned to node %d\n",
			  this_node, gc->part_lists[pl], *(int *)retrieve, gc->node));
      realloc_particlelist(gc->part_lists[pl], gc->part_lists[pl]->n = *(int *)retrieve);
#ifdef GHOST_FLAG
      //init ghost variable
      int i;
      for (i=0;i<gc->part_lists[pl]->n;i++){
           gc->part_lists[pl]->part[i].l.ghost=1;
      }
#endif
      retrieve += sizeof(int);
    }
    else {
      np   = gc->part_lists[pl]->n;
      part = gc->part_lists[pl]->part;
      for (p = 0; p < np; p++) {
	pt = &part[p];
	if (data_parts & GHOSTTRANS_PROPRTS) {
	  memcpy(&pt->p, retrieve, sizeof(ParticleProperties));
	  retrieve +=  sizeof(ParticleProperties);
	  /* GHOST_TRACE(fprintf(stderr, "%d: received ghost %d", this_node, pt->p.identity)); */
	  if (local_particles[pt->p.identity] == NULL) {
	    /* GHOST_TRACE(fprintf(stderr, ", using.\n")); */
	    local_particles[pt->p.identity] = pt;
	  }
	  /*
	    else {
	    GHOST_TRACE(fprintf(stderr, ", already here.\n"));
	    }
	  */
	}
	if (data_parts & GHOSTTRANS_POSITION) {
	  memcpy(&pt->r, retrieve, sizeof(ParticlePosition));
	  retrieve +=  sizeof(ParticlePosition);
	}
	if (data_parts & GHOSTTRANS_MOMENTUM) {
	  memcpy(&pt->m, retrieve, sizeof(ParticleMomentum));
	  retrieve +=  sizeof(ParticleMomentum);
	}
	if (data_parts & GHOSTTRANS_FORCE) {
	  memcpy(&pt->f, retrieve, sizeof(ParticleForce));
	  retrieve +=  sizeof(ParticleForce);
	}
#ifdef LB
	if (data_parts & GHOSTTRANS_COUPLING) {
	  memcpy(&pt->lc, retrieve, sizeof(ParticleLatticeCoupling));
	  retrieve +=  sizeof(ParticleLatticeCoupling);
	}
#endif
      }
    }
  }
#ifdef ADDITIONAL_CHECKS
  if (retrieve - r_buffer != n_r_buffer) {
    fprintf(stderr, "%d: recv buffer size %d differs from what I put in %d\n", this_node, n_r_buffer, retrieve - r_buffer);
    errexit();
  }
#endif
}
Exemplo n.º 12
0
void nsq_balance_particles()
{
  int i, n, surplus, s_node, tmp, lack, l_node, transfer;

  int pp = cells_get_n_particles();
  int *ppnode = malloc(n_nodes*sizeof(int));
  /* minimal difference between node shares */
  int minshare = n_total_particles/n_nodes;
  int maxshare = minshare + 1;

  CELL_TRACE(fprintf(stderr, "%d: nsq_balance_particles: load %d-%d\n", this_node, minshare, maxshare));

  MPI_Allgather(&pp, 1, MPI_INT, ppnode, 1, MPI_INT, MPI_COMM_WORLD);
  for (;;) {
    /* find node with most excessive particles */
    surplus = -1;
    s_node = -1;
    for (n = 0; n < n_nodes; n++) {
      tmp = ppnode[n] - minshare;
      CELL_TRACE(fprintf(stderr, "%d: nsq_balance_particles: node %d has %d\n", this_node, n, ppnode[n]));
      if (tmp > surplus) {
	surplus = tmp;
	s_node = n;
      }
    }
    CELL_TRACE(fprintf(stderr, "%d: nsq_balance_particles: excess %d on node %d\n", this_node, surplus, s_node));

    /* find node with most lacking particles */
    lack = -1;
    l_node = -1;
    for (n = 0; n < n_nodes; n++) {
      tmp = maxshare - ppnode[n];
      if (tmp > lack) {
	lack = tmp;
	l_node = n;
      }
    }
    CELL_TRACE(fprintf(stderr, "%d: nsq_balance_particles: lack %d on node %d\n", this_node, lack, l_node));

    /* should not happen: minshare or maxshare wrong or more likely,
       the algorithm */
    if (s_node == -1 || l_node == -1) {
      fprintf(stderr, "%d: Particle load balancing failed\n", this_node);
      break;
    }

    /* exit if all nodes load is withing min and max share */
    if (lack <= 1 && surplus <= 1)
      break;

    transfer = lack < surplus ? lack : surplus;

    if (s_node == this_node) {
      ParticleList send_buf;
      init_particlelist(&send_buf);
      realloc_particlelist(&send_buf, send_buf.n = transfer);
      for (i = 0; i < transfer; i++) {
	memcpy(&send_buf.part[i], &local->part[--local->n], sizeof(Particle));
      }
      realloc_particlelist(local, local->n);
      update_local_particles(local);

      send_particles(&send_buf, l_node);
#ifdef ADDITIONAL_CHECKS
      check_particle_consistency();
#endif
    }
    else if (l_node == this_node) {
      recv_particles(local, s_node);
#ifdef ADDITIONAL_CHECKS
      check_particle_consistency();
#endif
    }
    ppnode[s_node] -= transfer;
    ppnode[l_node] += transfer;
  }
  CELL_TRACE(fprintf(stderr, "%d: nsq_balance_particles: done\n", this_node));

  free(ppnode);
}
Exemplo n.º 13
0
void recv_particles(ParticleList *particles, int node)
{
  int transfer=0, read, pc;
  IntList local_dyn;

  PART_TRACE(fprintf(stderr, "%d: recv_particles from %d\n", this_node, node));

  MPI_Recv(&transfer, 1, MPI_INT, node,
	   REQ_SNDRCV_PART, comm_cart, MPI_STATUS_IGNORE);

  PART_TRACE(fprintf(stderr, "%d: recv_particles get %d\n", this_node, transfer));

  realloc_particlelist(particles, particles->n + transfer);
  MPI_Recv(&particles->part[particles->n], transfer*sizeof(Particle), MPI_BYTE, node,
	   REQ_SNDRCV_PART, comm_cart, MPI_STATUS_IGNORE);
  particles->n += transfer;

  init_intlist(&local_dyn);
  for (pc = particles->n - transfer; pc < particles->n; pc++) {
    Particle *p = &particles->part[pc];
    local_dyn.n += p->bl.n;
#ifdef EXCLUSIONS
    local_dyn.n += p->el.n;
#endif

    PART_TRACE(fprintf(stderr, "%d: recv_particles got particle %d\n", this_node, p->p.identity));
    if (local_particles[p->p.identity] != NULL) {
      fprintf(stderr, "%d: transmitted particle %d is already here...\n", this_node, p->p.identity);
      errexit();
    }
  }

  update_local_particles(particles);

  PART_TRACE(fprintf(stderr, "%d: recv_particles expecting %d bond ints\n", this_node, local_dyn.n));
  if (local_dyn.n > 0) {
    alloc_intlist(&local_dyn, local_dyn.n);
    MPI_Recv(local_dyn.e, local_dyn.n*sizeof(int), MPI_BYTE, node,
	     REQ_SNDRCV_PART, comm_cart, MPI_STATUS_IGNORE);
  }
  read = 0;
  for (pc = particles->n - transfer; pc < particles->n; pc++) {
    Particle *p = &particles->part[pc];
    if (p->bl.n > 0) {
      alloc_intlist(&p->bl, p->bl.n);
      memcpy(p->bl.e, &local_dyn.e[read], p->bl.n*sizeof(int));
      read += p->bl.n;
    }
    else
      p->bl.e = NULL;
#ifdef EXCLUSIONS
    if (p->el.n > 0) {
      alloc_intlist(&p->el, p->el.n);
      memcpy(p->el.e, &local_dyn.e[read], p->el.n*sizeof(int));
      read += p->el.n;
    }
    else
      p->el.e = NULL;
#endif
  }
  if (local_dyn.n > 0)
    realloc_intlist(&local_dyn, 0);
}