Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add OMP parallelization for cusp correction #1643

Merged
merged 10 commits into from
Jun 19, 2019
117 changes: 79 additions & 38 deletions src/QMCWaveFunctions/lcao/CuspCorrectionConstruction.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -27,6 +27,11 @@ void applyCuspCorrection(const Matrix<CuspCorrectionParameters>& info,
{
typedef QMCTraits::RealType RealType;

NewTimer* cuspApplyTimer =
TimerManager.createTimer("CuspCorrectionConstruction::applyCuspCorrection", timer_level_medium);

ScopedTimer cuspApplyTimerWrapper(cuspApplyTimer);

LCAOrbitalSet phi = LCAOrbitalSet(lcwc.myBasisSet);
phi.setOrbitalSetSize(lcwc.OrbitalSetSize);
phi.BasisSetSize = lcwc.BasisSetSize;
Expand Down Expand Up @@ -188,6 +193,13 @@ void generateCuspInfo(int orbital_set_size,
{
typedef QMCTraits::RealType RealType;

NewTimer* cuspCreateTimer =
TimerManager.createTimer("CuspCorrectionConstruction::createCuspParameters", timer_level_medium);
NewTimer* splitPhiEtaTimer = TimerManager.createTimer("CuspCorrectionConstruction::splitPhiEta", timer_level_fine);
NewTimer* computeTimer = TimerManager.createTimer("CuspCorrectionConstruction::computeCorrection", timer_level_fine);

ScopedTimer createCuspTimerWrapper(cuspCreateTimer);

LCAOrbitalSet phi = LCAOrbitalSet(lcwc.myBasisSet);
phi.setOrbitalSetSize(lcwc.OrbitalSetSize);
phi.BasisSetSize = lcwc.BasisSetSize;
Expand All @@ -211,55 +223,84 @@ void generateCuspInfo(int orbital_set_size,
int start_mo = offset[Comm.rank()];
int end_mo = offset[Comm.rank() + 1];
app_log() << " Number of molecular orbitals to compute correction on this rank: " << end_mo - start_mo << std::endl;
for (int mo_idx = start_mo; mo_idx < end_mo; mo_idx++)

#pragma omp parallel
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I'd like to see us move toward abstract or C++ concurrency constructions. The implicit capture of variables by the omp parallel pragma contributes to hazy scope and fears about race conditions. It makes it very convenient to write monster methods that build up a huge local namespace.
Is the potential confusion about thread scope worth it just to have what I think is initialization code be faster?

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Background info: The Cusp Correction construction/calculation is distressingly slow, so maintainable improvements are definitely welcome imo.

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Integrated the parallel scope with the for loop scope. This will result in additional allocations/deallocations with the objects, but it shouldn't affect the run time.

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Sorry I still have question with respect to the last change.

#pragma omp parallel
{
  ParticleSet localTargetPtcl(targetPtcl);
  #pragma omp for
  for(..)
  {
    do something on localTargetPtcl;
  }
}

This is a valid OpenMP code and this is a recommended optimization to avoid allocation and deallocation. I had a couple times corrected by OpenMP experts that my expectation of "#omp for" of concurrent execution of the loop iterations was wrong. "#omp for" is a workshare construct and indicating the loop can be distributed among the threads spanned by the parallel region. It didn't say that the loop iterations are independent and can be executed concurrently. Instead in OpenMP 5.0, the loop construct is introduced indicating that the loop iterations are independent and this definition aligns with the concept of the concurrent loop in C++ and Fortran.

Going back to the example, the thread scope is the parallel region and the localTargetPtcl is defined in the right scope.
It should have no difference with a code without OpenMP.

{
  ParticleSet localTargetPtcl;
  for(..) // The effect of "omp for" is changing the lower and upper bound of the loop, although it depends on the scheduling
  {
    do something on localTargetPtcl;
  }
}

The actual code did suffer from allcation/deallication and the imbalance between iterations make all the overhead amplified. I remember that when I profiled the cusp correction, the constructor and destructor took a lot of time. For this reason, I think the old way from @markdewing of doing parallel region is preferred.

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@PDoakORNL can we be a bit more flexible allowing this optimization to reduce overhead? it does introduce more data in the thread scope but scope and lifetime is quite clear.

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I did some profiling on my laptop, and there doesn't seem to be significant overhead in creating the LCAOrbitalSet copies.

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@ye-luo I would rather leave the code as is less dependent on openmp syntax. Is that acceptable?

I also think that we'd be better off refactoring egregious design issues when they are the cause of performance issues. Although according to Mark this isn't really the hotspot. @markdewing Can you tell what it is?

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

As @markdewing confirming that overhead is small. Putting object within the innermost scope is the cleanest way. If this really impact a workload, we can revisit this.

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Bulk of time (85%) is in DGEMM called from LCAOrbitalSet::evaluate, called from OneMolecularOrbital::phi_vgl, called from getCurrentLocalEnergy.
This is for a larger system, and I killed it part way through the cusp correction (so pretty much all the time in the run was spent doing cusp correction)

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

This is not the type of overhead I worried. I mean measuring the time of cusp construction with and without the optimization of copy.

{
app_log() << " Working on MO: " << mo_idx << std::endl;
ParticleSet localTargetPtcl(targetPtcl);
ParticleSet localSourcePtcl(sourcePtcl);

LCAOrbitalSet local_phi(phi);
local_phi.myBasisSet = phi.myBasisSet->makeClone();
local_phi.IsCloned = true;
local_phi.C = nullptr;
local_phi.setIdentity(false);

LCAOrbitalSet local_eta(eta);
local_eta.myBasisSet = eta.myBasisSet->makeClone();
local_eta.IsCloned = true;
local_eta.C = nullptr;
local_eta.setIdentity(false);

// Specify dynamic scheduling explicitly for load balancing. Each iteration should take enough
// time that scheduling overhead is not an issue.
#pragma omp for schedule(dynamic) collapse(2)
for (int center_idx = 0; center_idx < num_centers; center_idx++)
PDoakORNL marked this conversation as resolved.
Show resolved Hide resolved
{
*(eta.C) = *(lcwc.C);
*(phi.C) = *(lcwc.C);
splitPhiEta(center_idx, corrCenter, phi, eta);

bool corrO = false;
auto& cref(*(phi.C));
for (int ip = 0; ip < cref.cols(); ip++)
{
if (std::abs(cref(mo_idx, ip)) > 0)
{
corrO = true;
break;
}
}

if (corrO)
for (int mo_idx = start_mo; mo_idx < end_mo; mo_idx++)
{
OneMolecularOrbital etaMO(&targetPtcl, &sourcePtcl, &eta);
etaMO.changeOrbital(center_idx, mo_idx);
app_log() << " Working on MO: " << mo_idx << " Center: " << center_idx << std::endl;
PDoakORNL marked this conversation as resolved.
Show resolved Hide resolved

OneMolecularOrbital phiMO(&targetPtcl, &sourcePtcl, &phi);
phiMO.changeOrbital(center_idx, mo_idx);
splitPhiEtaTimer->start();

SpeciesSet& tspecies(sourcePtcl.getSpeciesSet());
int iz = tspecies.addAttribute("charge");
RealType Z = tspecies(iz, sourcePtcl.GroupID[center_idx]);
*(local_eta.C) = *(lcwc.C);
*(local_phi.C) = *(lcwc.C);
splitPhiEta(center_idx, corrCenter, local_phi, local_eta);

RealType Rc_max = 0.2;
RealType rc = 0.1;
splitPhiEtaTimer->stop();

RealType dx = rc * 1.2 / npts;
ValueVector_t pos(npts);
ValueVector_t ELideal(npts);
ValueVector_t ELcurr(npts);
for (int i = 0; i < npts; i++)
bool corrO = false;
auto& cref(*(local_phi.C));
for (int ip = 0; ip < cref.cols(); ip++)
{
pos[i] = (i + 1.0) * dx;
if (std::abs(cref(mo_idx, ip)) > 0)
{
corrO = true;
break;
}
}

RealType eta0 = etaMO.phi(0.0);
ValueVector_t ELorig(npts);
CuspCorrection cusp(info(center_idx, mo_idx));
minimizeForRc(cusp, phiMO, Z, rc, Rc_max, eta0, pos, ELcurr, ELideal);
info(center_idx, mo_idx) = cusp.cparam;
if (corrO)
{
OneMolecularOrbital etaMO(&localTargetPtcl, &localSourcePtcl, &local_eta);
etaMO.changeOrbital(center_idx, mo_idx);

OneMolecularOrbital phiMO(&localTargetPtcl, &localSourcePtcl, &local_phi);
phiMO.changeOrbital(center_idx, mo_idx);

SpeciesSet& tspecies(localSourcePtcl.getSpeciesSet());
int iz = tspecies.addAttribute("charge");
RealType Z = tspecies(iz, localSourcePtcl.GroupID[center_idx]);

RealType Rc_max = 0.2;
RealType rc = 0.1;

RealType dx = rc * 1.2 / npts;
ValueVector_t pos(npts);
ValueVector_t ELideal(npts);
ValueVector_t ELcurr(npts);
for (int i = 0; i < npts; i++)
{
pos[i] = (i + 1.0) * dx;
}

RealType eta0 = etaMO.phi(0.0);
ValueVector_t ELorig(npts);
CuspCorrection cusp(info(center_idx, mo_idx));
computeTimer->start();
minimizeForRc(cusp, phiMO, Z, rc, Rc_max, eta0, pos, ELcurr, ELideal);
computeTimer->stop();
info(center_idx, mo_idx) = cusp.cparam;
PDoakORNL marked this conversation as resolved.
Show resolved Hide resolved
}
}
}
}
Expand Down