void EinsplineSetBuilder::ReadBands_ESHDF(int spin, EinsplineSetExtended<double>* orbitalSet)
{
  update_token(__FILE__,__LINE__,"ReadBands_ESHDF:double");
  ReportEngine PRE("EinsplineSetBuilder","ReadBands_ESHDF(EinsplineSetExtended<double>*");
  vector<AtomicOrbital<double> > realOrbs(AtomicOrbitals.size());
  for (int iat=0; iat<realOrbs.size(); iat++)
  {
    AtomicOrbital<complex<double> > &corb (AtomicOrbitals[iat]);
    realOrbs[iat].set_pos  (corb.Pos);
    realOrbs[iat].set_lmax (corb.lMax);
    realOrbs[iat].set_cutoff (corb.CutoffRadius);
    realOrbs[iat].set_spline (corb.SplineRadius, corb.SplinePoints);
    realOrbs[iat].set_polynomial (corb.PolyRadius, corb.PolyOrder);
    realOrbs[iat].Lattice = corb.Lattice;
  }
  bool root = myComm->rank()==0;
  // bcast other stuff
  myComm->bcast (NumDistinctOrbitals);
  myComm->bcast (NumValenceOrbs);
  myComm->bcast (NumCoreOrbs);
  int N = NumDistinctOrbitals;
  orbitalSet->kPoints.resize(N);
  orbitalSet->MakeTwoCopies.resize(N);
  orbitalSet->StorageValueVector.resize(N);
  orbitalSet->BlendValueVector.resize(N);
  orbitalSet->StorageLaplVector.resize(N);
  orbitalSet->BlendLaplVector.resize(N);
  orbitalSet->StorageGradVector.resize(N);
  orbitalSet->BlendGradVector.resize(N);
  orbitalSet->StorageHessVector.resize(N);
  orbitalSet->StorageGradHessVector.resize(N);
  orbitalSet->phase.resize(N);
  orbitalSet->eikr.resize(N);
  orbitalSet->NumValenceOrbs = NumValenceOrbs;
  orbitalSet->NumCoreOrbs    = NumCoreOrbs;
  orbitalSet->FirstOrderSplines.resize(IonPos.size());
  // Read in k-points
  int numOrbs = orbitalSet->getOrbitalSetSize();
  int num = 0;
  vector<BandInfo>& SortBands(*FullBands[spin]);
  if (root)
  {
    for (int iorb=0; iorb<N; iorb++)
    {
      int ti = SortBands[iorb].TwistIndex;
      PosType twist  = TwistAngles[ti];
      orbitalSet->kPoints[iorb] = orbitalSet->PrimLattice.k_cart(twist);
      orbitalSet->MakeTwoCopies[iorb] =
        (num < (numOrbs-1)) && SortBands[iorb].MakeTwoCopies;
      num += orbitalSet->MakeTwoCopies[iorb] ? 2 : 1;
    }
    PosType twist0 = TwistAngles[SortBands[0].TwistIndex];
    for (int i=0; i<OHMMS_DIM; i++)
      if (std::fabs(std::fabs(twist0[i]) - 0.5) < 1.0e-8)
        orbitalSet->HalfG[i] = 1;
      else
        orbitalSet->HalfG[i] = 0;
    EinsplineSetBuilder::RotateBands_ESHDF(spin, orbitalSet);
  }
  myComm->bcast(orbitalSet->kPoints);
  myComm->bcast(orbitalSet->MakeTwoCopies);
  myComm->bcast(orbitalSet->HalfG);
  // First, check to see if we have already read this in
  H5OrbSet set(H5FileName, spin, N);
  bool havePsir=!ReadGvectors_ESHDF();
  app_log() << "MeshSize = (" << MeshSize[0] << ", "
            << MeshSize[1] << ", " << MeshSize[2] << ")\n";
  //int nx, ny, nz, bi, ti;
  int nx, ny, nz;
  nx=MeshSize[0];
  ny=MeshSize[1];
  nz=MeshSize[2];
  Ugrid x_grid, y_grid, z_grid;
  BCtype_d xBC, yBC, zBC;
  if (orbitalSet->HalfG[0])
  {
    xBC.lCode = ANTIPERIODIC;
    xBC.rCode = ANTIPERIODIC;
  }
  else
  {
    xBC.lCode = PERIODIC;
    xBC.rCode = PERIODIC;
  }
  if (orbitalSet->HalfG[1])
  {
    yBC.lCode = ANTIPERIODIC;
    yBC.rCode = ANTIPERIODIC;
  }
  else
  {
    yBC.lCode = PERIODIC;
    yBC.rCode = PERIODIC;
  }
  if (orbitalSet->HalfG[2])
  {
    zBC.lCode = ANTIPERIODIC;
    zBC.rCode = ANTIPERIODIC;
  }
  else
  {
    zBC.lCode = PERIODIC;
    zBC.rCode = PERIODIC;
  }
  x_grid.start = 0.0;
  x_grid.end = 1.0;
  x_grid.num = nx;
  y_grid.start = 0.0;
  y_grid.end = 1.0;
  y_grid.num = ny;
  z_grid.start = 0.0;
  z_grid.end = 1.0;
  z_grid.num = nz;
  // Create the multiUBspline object
  orbitalSet->MultiSpline =
    create_multi_UBspline_3d_d (x_grid, y_grid, z_grid, xBC, yBC, zBC, NumValenceOrbs);
  if (HaveOrbDerivs)
  {
    orbitalSet->FirstOrderSplines.resize(IonPos.size());
    for (int ion=0; ion<IonPos.size(); ion++)
      for (int dir=0; dir<OHMMS_DIM; dir++)
        orbitalSet->FirstOrderSplines[ion][dir] =
          create_multi_UBspline_3d_d (x_grid, y_grid, z_grid, xBC, yBC, zBC, NumValenceOrbs);
  }
  //////////////////////////////////////
  // Create the MuffinTin APW splines //
  //////////////////////////////////////
  orbitalSet->MuffinTins.resize(NumMuffinTins);
  for (int tin=0; tin<NumMuffinTins; tin++)
  {
    orbitalSet->MuffinTins[tin].Atom = tin;
    orbitalSet->MuffinTins[tin].set_center (MT_centers[tin]);
    orbitalSet->MuffinTins[tin].set_lattice(Lattice);
    orbitalSet->MuffinTins[tin].init_APW
    (MT_APW_rgrids[tin], MT_APW_lmax[tin],
     NumValenceOrbs);
  }
  for (int iat=0; iat<realOrbs.size(); iat++)
  {
    realOrbs[iat].set_num_bands(NumValenceOrbs);
    realOrbs[iat].allocate();
  }
  int isComplex;
  if (root)
  {
    HDFAttribIO<int> h_isComplex(isComplex);
    h_isComplex.read(H5FileID, "/electrons/psi_r_is_complex");
  }
  myComm->bcast(isComplex);
  bool isCore = bcastSortBands(spin,N,root);
  if(isCore)
  {
    APP_ABORT("Core states not supported by ES-HDF yet.");
  }
  //this is common
  Array<double,3> splineData(nx,ny,nz);
  if(havePsir)
  {
    if(isComplex)
    {
      app_log() << "   Reading complex psi_r and convert to real" << endl;
      Array<complex<double>,3> rawData;
      for(int iorb=0,ival=0; iorb<N; ++iorb, ++ival)
      {
        int ti=SortBands[iorb].TwistIndex;
        if(root)
        {
          ostringstream path;
          path << "/electrons/kpoint_" << SortBands[iorb].TwistIndex
               << "/spin_" << spin << "/state_" << SortBands[iorb].BandIndex << "/psi_r";
          HDFAttribIO<Array<complex<double>,3> >  h_splineData(rawData);
          h_splineData.read(H5FileID, path.str().c_str());
        }
        myComm->bcast(rawData);
        //multiply twist factor and project on the real
        fix_phase_c2r(rawData,splineData,TwistAngles[ti]);
        set_multi_UBspline_3d_d (orbitalSet->MultiSpline, ival, splineData.data());
      }
    }
    else
    {
      app_log() << "   Reading real psi_r" << endl;
      for(int iorb=0,ival=0; iorb<N; ++iorb, ++ival)
      {
        if(root)
        {
          ostringstream path;
          path << "/electrons/kpoint_" << SortBands[iorb].TwistIndex
               << "/spin_" << spin << "/state_" << SortBands[iorb].BandIndex << "/psi_r";
          HDFAttribIO<Array<double,3> >  h_splineData(splineData);
          h_splineData.read(H5FileID, path.str().c_str());
        }
        myComm->bcast(splineData);
        set_multi_UBspline_3d_d (orbitalSet->MultiSpline, ival, splineData.data());
      }
    }
  }
  else
  {
    Array<ComplexType,3> FFTbox;
    FFTbox.resize(MeshSize[0], MeshSize[1], MeshSize[2]);
    fftw_plan FFTplan = fftw_plan_dft_3d
                        (MeshSize[0], MeshSize[1], MeshSize[2],
                         reinterpret_cast<fftw_complex*>(FFTbox.data()),
                         reinterpret_cast<fftw_complex*>(FFTbox.data()),
                         +1, FFTW_ESTIMATE);
    for(int iorb=0,ival=0; iorb<N; ++iorb, ++ival)
    {
      Vector<complex<double> > cG;
      int ncg=0;
      int ti=SortBands[iorb].TwistIndex;
      if(root)
      {
        ostringstream path;
        path << "/electrons/kpoint_" << SortBands[iorb].TwistIndex
             << "/spin_" << spin << "/state_" << SortBands[iorb].BandIndex << "/psi_g";
        HDFAttribIO<Vector<complex<double> > >  h_cG(cG);
        h_cG.read (H5FileID, path.str().c_str());
        ncg=cG.size();
      }
      myComm->bcast(ncg);
      if(ncg != Gvecs[0].size())
      {
        APP_ABORT("Failed : ncg != Gvecs[0].size()");
      }
      if(!root)
        cG.resize(ncg);
      myComm->bcast(cG);
      unpack4fftw(cG,Gvecs[0],MeshSize,FFTbox);
      fftw_execute (FFTplan);
      fix_phase_rotate_c2r(FFTbox,splineData,TwistAngles[ti]);
      set_multi_UBspline_3d_d (orbitalSet->MultiSpline, ival, splineData.data());
    }
    fftw_destroy_plan(FFTplan);
  }
  for(int iorb=0,ival=0; iorb<N; ++iorb, ++ival)
  {
    // Read atomic orbital information
    for (int iat=0; iat<realOrbs.size(); iat++)
    {
      app_log() << "Reading orbital " << iat << " for band " << ival << endl;
      AtomicOrbital<double> &orb = realOrbs[iat];
      //AtomicOrbital<complex<double> > &orb = realOrbs[iat];
      Array<complex<double>,2> radial_spline(orb.SplinePoints,orb.Numlm),
            poly_coefs(orb.PolyOrder+1,orb.Numlm);
      int ti   = SortBands[iorb].TwistIndex;
      if (root)
      {
        int bi   = SortBands[iorb].BandIndex;
        ostringstream path;
        path << "/electrons/kpoint_" << ti << "/spin_" << spin << "/state_" << bi << "/";
        ostringstream spline_path, poly_path;
        spline_path << path.str() << "radial_spline_" << iat;
        poly_path   << path.str() << "poly_coefs_"    << iat;
        HDFAttribIO<Array<complex<double>,2> > h_radial_spline(radial_spline);
        HDFAttribIO<Array<complex<double>,2> > h_poly_coefs(poly_coefs);
        h_radial_spline.read(H5FileID, spline_path.str().c_str());
        h_poly_coefs.read   (H5FileID, poly_path.str().c_str());
      }
      myComm->bcast(radial_spline);
      myComm->bcast(poly_coefs);
      realOrbs[iat].set_band (ival, radial_spline, poly_coefs, TwistAngles[ti]);
    }
  }
  for(int iorb=0,ival=0; iorb<N; ++iorb, ++ival)
  {
    // Now read muffin tin data
    for (int tin=0; tin<NumMuffinTins; tin++)
    {
      // app_log() << "Reading data for muffin tin " << tin << endl;
      PosType twist, k;
      int lmax = MT_APW_lmax[tin];
      int numYlm = (lmax+1)*(lmax+1);
      Array<complex<double>,2>
      u_lm_r(numYlm, MT_APW_num_radial_points[tin]);
      Array<complex<double>,1> du_lm_dr (numYlm);
      int ti   = SortBands[iorb].TwistIndex;
      if (root)
      {
        int bi   = SortBands[iorb].BandIndex;
        twist = TwistAngles[ti];
        k = orbitalSet->PrimLattice.k_cart(twist);
        string uName  = MuffinTinPath (ti, bi,tin) + "u_lm_r";
        string duName = MuffinTinPath (ti, bi,tin) + "du_lm_dr";
        HDFAttribIO<Array<complex<double>,2> > h_u_lm_r(u_lm_r);
        HDFAttribIO<Array<complex<double>,1> > h_du_lm_dr(du_lm_dr);
        h_u_lm_r.read(H5FileID, uName.c_str());
        h_du_lm_dr.read(H5FileID, duName.c_str());
      }
      myComm->bcast(u_lm_r);
      myComm->bcast(du_lm_dr);
      myComm->bcast(k);
      double Z = (double)IonTypes(tin);
      OrbitalSet->MuffinTins[tin].set_APW (ival, k, u_lm_r, du_lm_dr, Z);
    }
  }
  //FIX HaveOrbDerivs after debugging
//	// Now read orbital derivatives if we have them
//	if (HaveOrbDerivs) {
//	  for (int ion=0; ion<IonPos.size(); ion++)
//	    for (int dim=0; dim<OHMMS_DIM; dim++) {
//	      if (root) {
//		int ti   = SortBands[iorb].TwistIndex;
//		int bi   = SortBands[iorb].BandIndex;
//
//		app_log() << "Reading orbital derivative for ion " << ion
//			  << " dim " << dim << " spin " << spin << " band "
//			  << bi << " kpoint " << ti << endl;
//		ostringstream path;
//		path << "/electrons/kpoint_" << ti << "/spin_" << spin << "/state_" << bi << "/"
//		     << "dpsi_" << ion << "_" << dim << "_r";
//		string psirName = path.str();
//		if (isComplex) {
//		  HDFAttribIO<Array<complex<double>,3> > h_rawData(rawData);
//		  h_rawData.read(H5FileID, psirName.c_str());
//		  if ((rawData.size(0) != nx) ||
//		      (rawData.size(1) != ny) ||
//		      (rawData.size(2) != nz)) {
//		    fprintf (stderr, "Error in EinsplineSetBuilder::ReadBands.\n");
//		    fprintf (stderr, "Extended orbitals should all have the same dimensions\n");
//		    abort();
//		  }
//#pragma omp parallel for
//		  for (int ix=0; ix<nx; ix++) {
//		  PosType ru;
//		    ru[0] = (RealType)ix / (RealType)nx;
//		    for (int iy=0; iy<ny; iy++) {
//		      ru[1] = (RealType)iy / (RealType)ny;
//		      for (int iz=0; iz<nz; iz++) {
//			ru[2] = (RealType)iz / (RealType)nz;
//			double phi = -2.0*M_PI*dot (ru, TwistAngles[ti]);
//			double s, c;
//			sincos(phi, &s, &c);
//			complex<double> phase(c,s);
//			complex<double> z = phase*rawData(ix,iy,iz);
//			splineData(ix,iy,iz) = z.real();
//		      }
//		    }
//		  }
//		}
//		else {
//		  HDFAttribIO<Array<double,3> >  h_splineData(splineData);
//		  h_splineData.read(H5FileID, psirName.c_str());
//		  if ((splineData.size(0) != nx) ||
//		      (splineData.size(1) != ny) ||
//		      (splineData.size(2) != nz)) {
//		    fprintf (stderr, "Error in EinsplineSetBuilder::ReadBands.\n");
//		    fprintf (stderr, "Extended orbitals should all have the same dimensions\n");
//		    abort();
//		  }
//		}
//	      }
//	      myComm->bcast(splineData);
//	      set_multi_UBspline_3d_d
//		(orbitalSet->FirstOrderSplines[ion][dim], ival, splineData.data());
//	    }
//	}
//
//
//
  orbitalSet->AtomicOrbitals = realOrbs;
  for (int i=0; i<orbitalSet->AtomicOrbitals.size(); i++)
    orbitalSet->AtomicOrbitals[i].registerTimers();
  //ExtendedMap_d[set] = orbitalSet->MultiSpline;
}
SPOSetBase*
EinsplineSetBuilder::createSPOSet(xmlNodePtr cur)
{
  //use 2 bohr as the default when truncated orbitals are used based on the extend of the ions
  BufferLayer=2.0;
  OhmmsAttributeSet attribs;
  int numOrbs = 0;
  qafm=0;
  int sortBands(1);
  string sourceName;
  string spo_prec("double");
  string truncate("no");
#if defined(QMC_CUDA)
  string useGPU="yes";
#else
  string useGPU="no";
#endif
  attribs.add (H5FileName, "href");
  attribs.add (TileFactor, "tile");
  attribs.add (sortBands,  "sort");
  attribs.add (qafm,  "afmshift");
  attribs.add (TileMatrix, "tilematrix");
  attribs.add (TwistNum,   "twistnum");
  attribs.add (givenTwist,   "twist");
  attribs.add (sourceName, "source");
  attribs.add (MeshFactor, "meshfactor");
  attribs.add (useGPU,     "gpu");
  attribs.add (spo_prec,   "precision");
  attribs.add (truncate,   "truncate");
  attribs.add (BufferLayer, "buffer");
  attribs.put (XMLRoot);
  attribs.add (numOrbs,    "size");
  attribs.add (numOrbs,    "norbs");
  attribs.put (cur);
  ///////////////////////////////////////////////
  // Read occupation information from XML file //
  ///////////////////////////////////////////////
  cur = cur->children;
  int spinSet = -1;
  vector<int> Occ_Old(0,0);
  Occ.resize(0,0);
  bool NewOcc(false);
  while (cur != NULL)
  {
    string cname((const char*)(cur->name));
    if(cname == "occupation")
    {
      string occ_mode("ground");
      occ_format="energy";
      particle_hole_pairs=0;
      OhmmsAttributeSet oAttrib;
      oAttrib.add(occ_mode,"mode");
      oAttrib.add(spinSet,"spindataset");
      oAttrib.add(occ_format,"format");
      oAttrib.add(particle_hole_pairs,"pairs");
      oAttrib.put(cur);
      if(occ_mode == "excited")
      {
        putContent(Occ,cur);
      }
      else
        if(occ_mode != "ground")
        {
          app_error() << "Only ground state occupation currently supported "
                      << "in EinsplineSetBuilder.\n";
          APP_ABORT("EinsplineSetBuilder::createSPOSet");
        }
    }
    cur = cur->next;
  }
  if (Occ != Occ_Old)
  {
    NewOcc=true;
    Occ_Old = Occ;
  }
  else
    NewOcc=false;
#if defined(QMC_CUDA)
  app_log() << "\t  QMC_CUDA=1 Overwriting the precision of the einspline storage on the host.\n";
  spo_prec="double"; //overwrite
#endif
  H5OrbSet aset(H5FileName, spinSet, numOrbs);
  std::map<H5OrbSet,SPOSetBase*,H5OrbSet>::iterator iter;
  iter = SPOSetMap.find (aset);
  if ((iter != SPOSetMap.end() ) && (!NewOcc) && (qafm==0))
  {
    qafm=0;
    app_log() << "SPOSet parameters match in EinsplineSetBuilder:  "
              << "cloning EinsplineSet object.\n";
    return iter->second->makeClone();
  }
  // The tiling can be set by a simple vector, (e.g. 2x2x2), or by a
  // full 3x3 matrix of integers.  If the tilematrix was not set in
  // the input file...
  bool matrixNotSet = true;
  for (int i=0; i<3; i++)
    for (int j=0; j<3; j++)
      matrixNotSet = matrixNotSet && (TileMatrix(i,j) == 0);
  // then set the matrix to what may have been specified in the
  // tiling vector
  if (matrixNotSet)
    for (int i=0; i<3; i++)
      for (int j=0; j<3; j++)
        TileMatrix(i,j) = (i==j) ? TileFactor(i) : 0;
  if (myComm->rank() == 0)
    fprintf (stderr, " [ %2d %2d %2d\n   %2d %2d %2d\n   %2d %2d %2d ]\n",
             TileMatrix(0,0), TileMatrix(0,1), TileMatrix(0,2),
             TileMatrix(1,0), TileMatrix(1,1), TileMatrix(1,2),
             TileMatrix(2,0), TileMatrix(2,1), TileMatrix(2,2));
  if (numOrbs == 0)
  {
    app_error() << "You must specify the number of orbitals in the input file.\n";
    APP_ABORT("EinsplineSetBuilder::createSPOSet");
  }
  else
    app_log() << "  Reading " << numOrbs << " orbitals from HDF5 file.\n";
  Timer mytimer;
  mytimer.restart();
  /////////////////////////////////////////////////////////////////
  // Read the basic orbital information, without reading all the //
  // orbitals themselves.                                        //
  /////////////////////////////////////////////////////////////////
  if (myComm->rank() == 0)
    if (!ReadOrbitalInfo())
    {
      app_error() << "Error reading orbital info from HDF5 file.  Aborting.\n";
      APP_ABORT("EinsplineSetBuilder::createSPOSet");
    }
  app_log() <<  "TIMER  EinsplineSetBuilder::ReadOrbitalInfo " << mytimer.elapsed() << endl;
  myComm->barrier();
  mytimer.restart();
  BroadcastOrbitalInfo();
  app_log() <<  "TIMER  EinsplineSetBuilder::BroadcastOrbitalInfo " << mytimer.elapsed() << endl;
  app_log().flush();
  ///////////////////////////////////////////////////////////////////
  // Now, analyze the k-point mesh to figure out the what k-points //
  // are needed                                                    //
  ///////////////////////////////////////////////////////////////////
  PrimCell.set(Lattice);
  SuperCell.set(SuperLattice);
  for (int iat=0; iat<AtomicOrbitals.size(); iat++)
    AtomicOrbitals[iat].Lattice = Lattice;
  // Copy supercell into the ParticleSets
//     app_log() << "Overwriting XML lattice with that from the ESHDF file.\n";
//     PtclPoolType::iterator piter;
//     for(piter = ParticleSets.begin(); piter != ParticleSets.end(); piter++)
//       piter->second->Lattice.copy(SuperCell);
  AnalyzeTwists2();
  //////////////////////////////////
  // Create the OrbitalSet object
  //////////////////////////////////
  if (HaveLocalizedOrbs)
    OrbitalSet = new EinsplineSetLocal;
#ifdef QMC_CUDA
  else
    if (AtomicOrbitals.size() > 0)
    {
      if (UseRealOrbitals)
        OrbitalSet = new EinsplineSetHybrid<double>;
      else
        OrbitalSet = new EinsplineSetHybrid<complex<double> >;
    }
#endif
    else
    {
      if (UseRealOrbitals)
        OrbitalSet = new EinsplineSetExtended<double>;
      else
        OrbitalSet = new EinsplineSetExtended<complex<double> >;
    }
  //set the internal parameters
  setTiling(OrbitalSet,numOrbs);
  if (HaveLocalizedOrbs)
  {
    EinsplineSetLocal *restrict orbitalSet =
      dynamic_cast<EinsplineSetLocal*>(OrbitalSet);
    #pragma omp critical(read_einspline_orbs)
    {
      if ((spinSet == LastSpinSet) && (numOrbs <= NumOrbitalsRead) && (!NewOcc) )
        CopyBands(numOrbs);
      else
      {
        // Now, figure out occupation for the bands and read them
        OccupyBands(spinSet, sortBands);
        ReadBands (spinSet, orbitalSet);
      }
    }
    // Now, store what we have read
    LastOrbitalSet = OrbitalSet;
    LastSpinSet = spinSet;
    NumOrbitalsRead = numOrbs;
  }
  else // Otherwise, use EinsplineSetExtended
  {
    mytimer.restart();
    bool use_single= (spo_prec == "single" || spo_prec == "float");
    if (UseRealOrbitals)
    {
      OccupyBands(spinSet, sortBands);
      //check if a matching BsplineReaderBase exists
      BsplineReaderBase* spline_reader=0;
      //if(TargetPtcl.Lattice.SuperCellEnum != SUPERCELL_BULK && truncate=="yes")
      if(truncate=="yes")
      {
        if(use_single)
        {
          if(TargetPtcl.Lattice.SuperCellEnum == SUPERCELL_OPEN)
            spline_reader= new SplineMixedAdoptorReader<SplineOpenAdoptor<float,double,3> >(this);
          else
            if(TargetPtcl.Lattice.SuperCellEnum == SUPERCELL_SLAB)
              spline_reader= new SplineMixedAdoptorReader<SplineMixedAdoptor<float,double,3> >(this);
            else
              spline_reader= new SplineAdoptorReader<SplineR2RAdoptor<float,double,3> >(this);
        }
        else
        {
          if(TargetPtcl.Lattice.SuperCellEnum == SUPERCELL_OPEN)
            spline_reader= new SplineMixedAdoptorReader<SplineOpenAdoptor<double,double,3> >(this);
          else
            if(TargetPtcl.Lattice.SuperCellEnum == SUPERCELL_SLAB)
              spline_reader= new SplineMixedAdoptorReader<SplineMixedAdoptor<double,double,3> >(this);
            else
              spline_reader= new SplineAdoptorReader<SplineR2RAdoptor<double,double,3> >(this);
        }
      }
      else
      {
        if(use_single)
          spline_reader= new SplineAdoptorReader<SplineR2RAdoptor<float,double,3> >(this);
      }
      if(spline_reader)
      {
        HasCoreOrbs=bcastSortBands(NumDistinctOrbitals,myComm->rank()==0);
        SPOSetBase* bspline_zd=spline_reader->create_spline_set(spinSet,OrbitalSet);
        delete spline_reader;
        if(bspline_zd)
          SPOSetMap[aset] = bspline_zd;
        return bspline_zd;
      }
      else
      {
        app_log() << ">>>> Creating EinsplineSetExtended<double> <<<< " << endl;
        EinsplineSetExtended<double> *restrict orbitalSet =
          dynamic_cast<EinsplineSetExtended<double>* > (OrbitalSet);
        if (Format == ESHDF)
          ReadBands_ESHDF(spinSet,orbitalSet);
        else
          ReadBands(spinSet, orbitalSet);
      }
    }
    else
    {
      OccupyBands(spinSet, sortBands);
      BsplineReaderBase* spline_reader=0;
      if(truncate == "yes")
      {
        app_log() << "  Truncated orbitals with multiple kpoints are not supported yet!" << endl;
      }
      if(use_single)
      {
#if defined(QMC_COMPLEX)
        spline_reader= new SplineAdoptorReader<SplineC2CPackedAdoptor<float,double,3> >(this);
#else
        spline_reader= new SplineAdoptorReader<SplineC2RPackedAdoptor<float,double,3> >(this);
#endif
      }
      if(spline_reader)
      {
        RotateBands_ESHDF(spinSet, dynamic_cast<EinsplineSetExtended<complex<double> >*>(OrbitalSet));
        HasCoreOrbs=bcastSortBands(NumDistinctOrbitals,myComm->rank()==0);
        SPOSetBase* bspline_zd=spline_reader->create_spline_set(spinSet,OrbitalSet);
        delete spline_reader;
        if(bspline_zd)
          SPOSetMap[aset] = bspline_zd;
        return bspline_zd;
      }
      else
      {
        EinsplineSetExtended<complex<double> > *restrict orbitalSet =
          dynamic_cast<EinsplineSetExtended<complex<double> >*>(OrbitalSet);
        if (Format == ESHDF)
          ReadBands_ESHDF(spinSet,orbitalSet);
        else
          ReadBands(spinSet, orbitalSet);
      }
    }
    app_log() <<  "TIMER  EinsplineSetBuilder::ReadBands " << mytimer.elapsed() << endl;
  }
#ifndef QMC_COMPLEX
  if (myComm->rank()==0 && OrbitalSet->MuffinTins.size() > 0)
  {
    FILE *fout  = fopen ("TestMuffins.dat", "w");
    Vector<double> phi(numOrbs), lapl(numOrbs);
    Vector<PosType> grad(numOrbs);
    ParticleSet P;
    P.R.resize(6);
    for (int i=0; i<P.R.size(); i++)
      P.R[i] = PosType (0.0, 0.0, 0.0);
    PosType N = 0.25*PrimCell.a(0) + 0.25*PrimCell.a(1) + 0.25*PrimCell.a(2);
    for (double x=-1.0; x<=1.0; x+=0.0000500113412)
    {
      // for (double x=-0.003; x<=0.003; x+=0.0000011329343481381) {
      P.R[0] = x * (PrimCell.a(0) + 0.914*PrimCell.a(1) +
                    0.781413*PrimCell.a(2));
      double r = std::sqrt(dot(P.R[0], P.R[0]));
      double rN = std::sqrt(dot(P.R[0]-N, P.R[0]-N));
      OrbitalSet->evaluate(P, 0, phi, grad, lapl);
      // OrbitalSet->evaluate(P, 0, phi);
      fprintf (fout, "%1.12e ", r*x/std::fabs(x));
      for (int j=0; j<numOrbs; j++)
      {
        double gmag = std::sqrt(dot(grad[j],grad[j]));
        fprintf (fout, "%16.12e ",
                 /*phi[j]*phi[j]**/(-5.0/r  -0.5*lapl[j]/phi[j]));
        // double E = -5.0/r -0.5*lapl[j]/phi[j];
        fprintf (fout, "%16.12e ", phi[j]);
        fprintf (fout, "%16.12e ", gmag);
      }
      fprintf (fout, "\n");
    }
    fclose(fout);
  }
#endif
  SPOSetMap[aset] = OrbitalSet;
  if (sourceName.size() && (ParticleSets.find(sourceName) == ParticleSets.end()))
  {
    app_log() << "  EinsplineSetBuilder creates a ParticleSet " << sourceName << endl;
    ParticleSet* ions=new ParticleSet;
    ions->Lattice=TargetPtcl.Lattice;
    ESHDFIonsParser ap(*ions,H5FileID,myComm);
    ap.put(XMLRoot);
    ap.expand(TileMatrix);
    ions->setName(sourceName);
    ParticleSets[sourceName]=ions;
    //overwrite the lattice and assign random
    if(TargetPtcl.Lattice.SuperCellEnum)
    {
      TargetPtcl.Lattice=ions->Lattice;
      makeUniformRandom(TargetPtcl.R);
      TargetPtcl.R.setUnit(PosUnit::LatticeUnit);
      TargetPtcl.convert2Cart(TargetPtcl.R);
      TargetPtcl.createSK();
    }
  }
#ifdef QMC_CUDA
  if (useGPU == "yes" || useGPU == "1")
  {
    app_log() << "Initializing GPU data structures.\n";
    OrbitalSet->initGPU();
  }
#endif
  return OrbitalSet;
}
void EinsplineSetBuilder::ReadBands_ESHDF(int spin, EinsplineSetExtended<complex<double > >* orbitalSet)
{
  update_token(__FILE__,__LINE__,"ReadBands_ESHDF:complex");
  ReportEngine PRE("EinsplineSetBuilder","ReadBands_ESHDF(EinsplineSetExtended<complex<double > >*");
  Timer c_prep, c_unpack,c_fft, c_phase, c_spline, c_newphase, c_h5, c_init;
  double t_prep=0.0, t_unpack=0.0, t_fft=0.0, t_phase=0.0, t_spline=0.0, t_newphase=0.0, t_h5=0.0, t_init=0.0;
  c_prep.restart();
  bool root = myComm->rank()==0;
  vector<BandInfo>& SortBands(*FullBands[spin]);
  // bcast other stuff
  myComm->bcast (NumDistinctOrbitals);
  myComm->bcast (NumValenceOrbs);
  myComm->bcast (NumCoreOrbs);
  int N = NumDistinctOrbitals;
  orbitalSet->kPoints.resize(N);
  orbitalSet->MakeTwoCopies.resize(N);
  orbitalSet->StorageValueVector.resize(N);
  orbitalSet->BlendValueVector.resize(N);
  orbitalSet->StorageLaplVector.resize(N);
  orbitalSet->BlendLaplVector.resize(N);
  orbitalSet->StorageGradVector.resize(N);
  orbitalSet->BlendGradVector.resize(N);
  orbitalSet->StorageHessVector.resize(N);
  orbitalSet->StorageGradHessVector.resize(N);
  orbitalSet->phase.resize(N);
  orbitalSet->eikr.resize(N);
  orbitalSet->NumValenceOrbs = NumValenceOrbs;
  orbitalSet->NumCoreOrbs    = NumCoreOrbs;
  // Read in k-points
  int numOrbs = orbitalSet->getOrbitalSetSize();
  int num = 0;
  if (root)
  {
    for (int iorb=0; iorb<N; iorb++)
    {
      int ti = SortBands[iorb].TwistIndex;
      PosType twist  = TwistAngles[ti];
      orbitalSet->kPoints[iorb] = orbitalSet->PrimLattice.k_cart(twist);
      orbitalSet->MakeTwoCopies[iorb] =
        (num < (numOrbs-1)) && SortBands[iorb].MakeTwoCopies;
      num += orbitalSet->MakeTwoCopies[iorb] ? 2 : 1;
    }
  }
  myComm->bcast(orbitalSet->kPoints);
  myComm->bcast(orbitalSet->MakeTwoCopies);
  // First, check to see if we have already read this in
  H5OrbSet set(H5FileName, spin, N);
  ///check mesh or ready for FFT grid
  bool havePsig=ReadGvectors_ESHDF();
  app_log() << "MeshSize = (" << MeshSize[0] << ", " << MeshSize[1] << ", " << MeshSize[2] << ")\n";
  int nx, ny, nz, bi, ti;
  nx=MeshSize[0];
  ny=MeshSize[1];
  nz=MeshSize[2];
  Ugrid x_grid, y_grid, z_grid;
  BCtype_z xBC, yBC, zBC;
  xBC.lCode = PERIODIC;
  xBC.rCode = PERIODIC;
  yBC.lCode = PERIODIC;
  yBC.rCode = PERIODIC;
  zBC.lCode = PERIODIC;
  zBC.rCode = PERIODIC;
  x_grid.start = 0.0;
  x_grid.end = 1.0;
  x_grid.num = nx;
  y_grid.start = 0.0;
  y_grid.end = 1.0;
  y_grid.num = ny;
  z_grid.start = 0.0;
  z_grid.end = 1.0;
  z_grid.num = nz;
  // Create the multiUBspline object
  orbitalSet->MultiSpline =
    create_multi_UBspline_3d_z (x_grid, y_grid, z_grid, xBC, yBC, zBC, NumValenceOrbs);
  //////////////////////////////////////
  // Create the MuffinTin APW splines //
  //////////////////////////////////////
  orbitalSet->MuffinTins.resize(NumMuffinTins);
  for (int tin=0; tin<NumMuffinTins; tin++)
  {
    orbitalSet->MuffinTins[tin].Atom = tin;
    orbitalSet->MuffinTins[tin].set_center (MT_centers[tin]);
    orbitalSet->MuffinTins[tin].set_lattice(Lattice);
    orbitalSet->MuffinTins[tin].init_APW
    (MT_APW_rgrids[tin], MT_APW_lmax[tin],
     NumValenceOrbs);
  }
  for (int iat=0; iat<AtomicOrbitals.size(); iat++)
  {
    AtomicOrbitals[iat].set_num_bands(NumValenceOrbs);
    AtomicOrbitals[iat].allocate();
  }
  int isComplex=1;
  if (root)
  {
    HDFAttribIO<int> h_isComplex(isComplex);
    h_isComplex.read(H5FileID, "/electrons/psi_r_is_complex");
  }
  myComm->bcast(isComplex);
  if (!isComplex)
  {
    APP_ABORT("Expected complex orbitals in ES-HDF file, but found real ones.");
  }
  EinsplineSetBuilder::RotateBands_ESHDF(spin, orbitalSet);
  bool isCore = bcastSortBands(spin,N,root);
  if(isCore)
  {
    APP_ABORT("Core states not supported by ES-HDF yet.");
  }
  t_prep += c_prep.elapsed();
  /** For valence orbitals,
   * - extended orbitals either in G or in R
   * - localized orbitals
   */
  //this can potentially break
  Array<ComplexType,3> splineData(nx,ny,nz);
  if(havePsig)//perform FFT using FFTW
  {
    c_init.restart();
    Array<ComplexType,3> FFTbox;
    FFTbox.resize(MeshSize[0], MeshSize[1], MeshSize[2]);
    fftw_plan FFTplan = fftw_plan_dft_3d
                        (MeshSize[0], MeshSize[1], MeshSize[2],
                         reinterpret_cast<fftw_complex*>(FFTbox.data()),
                         reinterpret_cast<fftw_complex*>(FFTbox.data()),
                         +1, FFTW_ESTIMATE);
    Vector<complex<double> > cG(MaxNumGvecs);
    //this will be parallelized with OpenMP
    for(int iorb=0,ival=0; iorb<N; ++iorb, ++ival)
    {
      //Vector<complex<double> > cG;
      int ncg=0;
      int ti=SortBands[iorb].TwistIndex;
      c_h5.restart();
      if(root)
      {
        ostringstream path;
        path << "/electrons/kpoint_" << SortBands[iorb].TwistIndex
             << "/spin_" << spin << "/state_" << SortBands[iorb].BandIndex << "/psi_g";
        HDFAttribIO<Vector<complex<double> > >  h_cG(cG);
        h_cG.read (H5FileID, path.str().c_str());
        ncg=cG.size();
      }
      myComm->bcast(ncg);
      if(ncg != Gvecs[0].size())
      {
        APP_ABORT("Failed : ncg != Gvecs[0].size()");
      }
      if(!root)
        cG.resize(ncg);
      myComm->bcast(cG);
      t_h5 += c_h5.elapsed();
      c_unpack.restart();
      unpack4fftw(cG,Gvecs[0],MeshSize,FFTbox);
      t_unpack+= c_unpack.elapsed();
      c_fft.restart();
      fftw_execute (FFTplan);
      t_fft+= c_fft.elapsed();
      c_phase.restart();
      fix_phase_rotate_c2c(FFTbox,splineData,TwistAngles[ti]);
      t_phase+= c_phase.elapsed();
      c_spline.restart();
      set_multi_UBspline_3d_z(orbitalSet->MultiSpline, ival, splineData.data());
      t_spline+= c_spline.elapsed();
    }
    fftw_destroy_plan(FFTplan);
    t_init+=c_init.elapsed();
  }
  else
  {
    //this will be parallelized with OpenMP
    for(int iorb=0,ival=0; iorb<N; ++iorb, ++ival)
    {
      //check dimension
      if(root)
      {
        ostringstream path;
        path << "/electrons/kpoint_" << SortBands[iorb].TwistIndex
             << "/spin_" << spin << "/state_" << SortBands[iorb].BandIndex << "/psi_r";
        HDFAttribIO<Array<complex<double>,3> >  h_splineData(splineData);
        h_splineData.read(H5FileID, path.str().c_str());
      }
      myComm->bcast(splineData);
      set_multi_UBspline_3d_z(orbitalSet->MultiSpline, ival, splineData.data());
    }
    //return true;
  }
  app_log() << "    READBANDS::PREP   = " << t_prep << endl;
  app_log() << "    READBANDS::H5     = " << t_h5 << endl;
  app_log() << "    READBANDS::UNPACK = " << t_unpack << endl;
  app_log() << "    READBANDS::FFT    = " << t_fft << endl;
  app_log() << "    READBANDS::PHASE  = " << t_phase << endl;
  app_log() << "    READBANDS::SPLINE = " << t_spline << endl;
  app_log() << "    READBANDS::SUM    = " << t_init << endl;
  //now localized orbitals
  for(int iorb=0,ival=0; iorb<N; ++iorb, ++ival)
  {
    PosType twist=TwistAngles[SortBands[iorb].TwistIndex];
    // Read atomic orbital information
    for (int iat=0; iat<AtomicOrbitals.size(); iat++)
    {
      app_log() << "Reading orbital " << iat << " for band " << ival << endl;
      AtomicOrbital<complex<double> > &orb = AtomicOrbitals[iat];
      Array<complex<double>,2> radial_spline(orb.SplinePoints,orb.Numlm),
            poly_coefs(orb.PolyOrder+1,orb.Numlm);
      if (root)
      {
        int ti   = SortBands[iorb].TwistIndex;
        int bi   = SortBands[iorb].BandIndex;
        ostringstream path;
        path << "/electrons/kpoint_" << ti << "/spin_" << spin << "/state_" << bi << "/";
        AtomicOrbital<complex<double> > &orb = AtomicOrbitals[iat];
        ostringstream spline_path, poly_path;
        spline_path << path.str() << "radial_spline_" << iat;
        poly_path   << path.str() << "poly_coefs_"    << iat;
        HDFAttribIO<Array<complex<double>,2> > h_radial_spline(radial_spline);
        HDFAttribIO<Array<complex<double>,2> > h_poly_coefs(poly_coefs);
        h_radial_spline.read(H5FileID, spline_path.str().c_str());
        h_poly_coefs.read   (H5FileID, poly_path.str().c_str());
        // cerr << "radial_spline.size = (" << radial_spline.size(0)
        // 	 << ", " << radial_spline.size(1) << ")\n";
        // cerr << "poly_coefs.size = (" << poly_coefs.size(0)
        // 	 << ", " << poly_coefs.size(1) << ")\n";
      }
      myComm->bcast(radial_spline);
      myComm->bcast(poly_coefs);
      AtomicOrbitals[iat].set_band (ival, radial_spline, poly_coefs, twist);
    }
    // Now read muffin tin data
    for (int tin=0; tin<NumMuffinTins; tin++)
    {
      // app_log() << "Reading data for muffin tin " << tin << endl;
      PosType twist, k;
      int lmax = MT_APW_lmax[tin];
      int numYlm = (lmax+1)*(lmax+1);
      Array<complex<double>,2>
      u_lm_r(numYlm, MT_APW_num_radial_points[tin]);
      Array<complex<double>,1> du_lm_dr (numYlm);
      if (root)
      {
        int ti   = SortBands[iorb].TwistIndex;
        int bi   = SortBands[iorb].BandIndex;
        twist = TwistAngles[ti];
        k = orbitalSet->PrimLattice.k_cart(twist);
        string uName  = MuffinTinPath (ti, bi,tin) + "u_lm_r";
        string duName = MuffinTinPath (ti, bi,tin) + "du_lm_dr";
        HDFAttribIO<Array<complex<double>,2> > h_u_lm_r(u_lm_r);
        HDFAttribIO<Array<complex<double>,1> > h_du_lm_dr(du_lm_dr);
        h_u_lm_r.read(H5FileID, uName.c_str());
        h_du_lm_dr.read(H5FileID, duName.c_str());
      }
      myComm->bcast(u_lm_r);
      myComm->bcast(du_lm_dr);
      myComm->bcast(k);
      double Z = (double)IonTypes(tin);
      OrbitalSet->MuffinTins[tin].set_APW (ival, k, u_lm_r, du_lm_dr, Z);
    }
  }
  orbitalSet->AtomicOrbitals = AtomicOrbitals;
  for (int i=0; i<orbitalSet->AtomicOrbitals.size(); i++)
    orbitalSet->AtomicOrbitals[i].registerTimers();
  //ExtendedMap_z[set] = orbitalSet->MultiSpline;
}