Commit a5d38c08 authored by Axel Kohlmeyer's avatar Axel Kohlmeyer
Browse files

prototype implementation for extensible imbalance scheme

(cherry picked from commit 362a26a3def0803546b986f951bf9868acc0db15)
parent eb273ab9
Loading
Loading
Loading
Loading
+64 −3
Original line number Diff line number Diff line
@@ -30,6 +30,11 @@
#include "group.h"
#include "timer.h"

#include "imbalance_group.h"
#include "imbalance_time.h"
#include "imbalance_neigh.h"
#include "imbalance_var.h"

using namespace LAMMPS_NS;

enum{XYZ,SHIFT,BISECTION};
@@ -55,6 +60,10 @@ Balance::Balance(LAMMPS *lmp) : Pointers(lmp)
  fp = NULL;
  firststep = 1;

  nimbalance = 0;
  imbalance = NULL;
  weight = NULL;

  ngroup = 0;
  group_id = NULL;
  group_weight = NULL;
@@ -86,11 +95,16 @@ Balance::~Balance()
  }

  delete rcb;
  for (int i; i < nimbalance; ++i)
    delete imbalance[i];
  delete [] imbalance;
  delete [] weight;

#if 1
  delete [] group_id;
  delete [] group_weight;
  delete [] clock_imbalance;

#endif
  if (fp) fclose(fp);
}

@@ -100,6 +114,8 @@ Balance::~Balance()

void Balance::command(int narg, char **arg)
{
  double start_time = MPI_Wtime();

  if (domain->box_exist == 0)
    error->all(FLERR,"Balance command before simulation box is defined");

@@ -199,10 +215,15 @@ void Balance::command(int narg, char **arg)
    } else break;
  }

  // optional keywords
  // process optional keywords

  outflag = 0;
  // get max number of imbalance weight flags/classes
  nimbalance = 0;
  for (int i=iarg; i < narg; ++i)
    if (strcmp(arg[iarg],"weight") == 0) ++nimbalance;
  imbalance = new Imbalance*[nimbalance];

  nimbalance = outflag = 0;
  while (iarg < narg) {
    if (strcmp(arg[iarg],"out") == 0) {
      if (iarg+2 > narg) error->all(FLERR,"Illegal balance command");
@@ -213,6 +234,30 @@ void Balance::command(int narg, char **arg)
        if (fp == NULL) error->one(FLERR,"Cannot open balance output file");
      }
      iarg += 2;
    } else if (strcmp(arg[iarg],"weight") == 0) {
      Imbalance *imb;
      int nopt = 0;
      if (strcmp(arg[iarg+1],"group") == 0) {
        imb = new ImbalanceGroup;
        nopt = imb->options(lmp,narg-iarg-1,arg+iarg+1);
        imbalance[nimbalance] = imb;
      } else if (strcmp(arg[iarg+1],"time") == 0) {
        imb = new ImbalanceTime;
        nopt = imb->options(lmp,narg-iarg-1,arg+iarg+1);
        imbalance[nimbalance] = imb;
      } else if (strcmp(arg[iarg+1],"neigh") == 0) {
        imb = new ImbalanceNeigh;
        nopt = imb->options(lmp,narg-iarg-1,arg+iarg+1);
        imbalance[nimbalance] = imb;
      } else if (strcmp(arg[iarg+1],"var") == 0) {
        imb = new ImbalanceVar;
        nopt = imb->options(lmp,narg-iarg-1,arg+iarg+1);
        imbalance[nimbalance] = imb;
      } else {
        error->all(FLERR,"Unknown balance weight method");
      }
      iarg += 2+nopt;
#if 1
    } else if (strcmp(arg[iarg],"clock") == 0) {
      if (iarg+2 > narg) error->all(FLERR,"Illegal balance command");
      double factor = force->numeric(FLERR,arg[iarg+1]);
@@ -223,6 +268,7 @@ void Balance::command(int narg, char **arg)
    } else if (strcmp(arg[iarg],"group") == 0) {
      group_setup(narg-iarg-1,arg+iarg+1);
      iarg += 2*ngroup + 2;
#endif
    } else error->all(FLERR,"Illegal balance command");
  }

@@ -275,6 +321,17 @@ void Balance::command(int narg, char **arg)
  comm->exchange();
  if (domain->triclinic) domain->lamda2x(atom->nlocal);

  // compute and apply imbalance weights
  if (nimbalance > 0) {
    int i;
    const int nlocal = atom->nlocal;
    weight = new double[nlocal];
    for (i = 0; i < nlocal; ++i)
      weight[i] = 1.0;
    for (i = 0; i < nimbalance; ++i)
      imbalance[i]->compute(lmp,weight);
  } else weight = NULL;

  // imbinit = initial imbalance

  int maxinit;
@@ -385,6 +442,8 @@ void Balance::command(int narg, char **arg)

  if (me == 0) {
    if (screen) {
      fprintf(screen,"  rebalancing time: %g seconds\n",
              MPI_Wtime()-start_time);
      fprintf(screen,"  iteration count = %d\n",niter);
      if (ngroup > 0) {
        fprintf(screen,"  group weights:");
@@ -641,12 +700,14 @@ int *Balance::bisection(int sortflag)
  // then invert() to create list of proc assignements for my atoms
  // Use specified weightings for each atom rather than atom count

#if 1
  double factor = 1.0;
  if (clock_imbalance) factor = clock_imbalance[me];

  double *weights = new double[nlocal];
  for (int i = 0; i < nlocal; i++)
    weights[i] = getcost(i)*factor;
#endif

  rcb->compute(dim,atom->nlocal,atom->x,weights,shrinklo,shrinkhi);
  rcb->invert(sortflag);
+6 −0
Original line number Diff line number Diff line
@@ -67,11 +67,17 @@ class Balance : protected Pointers {
  int *proccount;            // particle count per processor
  int *allproccount;

  int nimbalance;
  class Imbalance **imbalance; // list of imbalance compute classes
  double *weight;            // per (local) atom weight factor or NULL

#if 1
  int    ngroup;             // number of groups weights
  int    *group_id;          // group ids for weights
  double *group_weight;      // weights of groups

  double *clock_imbalance;   // computed wall clock imbalance, NULL if not available
#endif

  int outflag;               // for output of balance results to file
  FILE *fp;

src/imbalance.h

0 → 100644
+40 −0
Original line number Diff line number Diff line
/* -*- c++ -*- ----------------------------------------------------------
   LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
   http://lammps.sandia.gov, Sandia National Laboratories
   Steve Plimpton, sjplimp@sandia.gov

   Copyright (2003) Sandia Corporation.  Under the terms of Contract
   DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
   certain rights in this software.  This software is distributed under
   the GNU General Public License.

   See the README file in the top-level LAMMPS directory.
------------------------------------------------------------------------- */

#ifndef LMP_IMBALANCE_H
#define LMP_IMBALANCE_H

namespace LAMMPS_NS {
 class LAMMPS;

class Imbalance {
 public:
  Imbalance() {};
  virtual ~Imbalance() {};

  // disallow copy constructor and assignment operator
 private:
  Imbalance(const Imbalance &) {};
  Imbalance &operator=(const Imbalance &) {return *this;};

  // required member functions
 public:
  // parse options. return number of arguments consumed.
  virtual int options(LAMMPS *lmp, int narg, char **arg) = 0;
  // compute per-atom imbalance and apply to weight array
  virtual void compute(LAMMPS *lmp, double *weights) = 0;
};

}

#endif
+27 −0
Original line number Diff line number Diff line
/* ----------------------------------------------------------------------
   LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
   http://lammps.sandia.gov, Sandia National Laboratories
   Steve Plimpton, sjplimp@sandia.gov

   Copyright (2003) Sandia Corporation.  Under the terms of Contract
   DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
   certain rights in this software.  This software is distributed under
   the GNU General Public License.

   See the README file in the top-level LAMMPS directory.
------------------------------------------------------------------------- */


#include "pointers.h"
#include "imbalance_group.h"

using namespace LAMMPS_NS;

int ImbalanceGroup::options(LAMMPS *lmp, int narg, char **arg)
{
  return 0;
}
 
void ImbalanceGroup::compute(LAMMPS *lmp, double *weight)
{
}

src/imbalance_group.h

0 → 100644
+41 −0
Original line number Diff line number Diff line
/* -*- c++ -*- ----------------------------------------------------------
   LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
   http://lammps.sandia.gov, Sandia National Laboratories
   Steve Plimpton, sjplimp@sandia.gov

   Copyright (2003) Sandia Corporation.  Under the terms of Contract
   DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
   certain rights in this software.  This software is distributed under
   the GNU General Public License.

   See the README file in the top-level LAMMPS directory.
------------------------------------------------------------------------- */

#ifndef LMP_IMBALANCE_GROUP_H
#define LMP_IMBALANCE_GROUP_H

#include "imbalance.h"

namespace LAMMPS_NS {

class ImbalanceGroup : public Imbalance {
 public:
  ImbalanceGroup() : Imbalance() {};
  virtual ~ImbalanceGroup() {};

  // disallow copy constructor and assignment operator
 private:
  ImbalanceGroup(const ImbalanceGroup &) {};
  ImbalanceGroup &operator=(const ImbalanceGroup &) {return *this;};

  // required member functions
 public:
  // parse options. return number of arguments consumed.
  virtual int options(LAMMPS *lmp, int narg, char **arg);
  // compute per-atom imbalance and apply to weight array
  virtual void compute(LAMMPS *lmp, double *weight);
};

}

#endif
Loading