Skip to content

Commit

Permalink
Implemented overloaded dealii's constraints.distribute function to re…
Browse files Browse the repository at this point in the history
…duce the overhead costs

associated with native dealii implementation
  • Loading branch information
phanimotamarri committed Feb 27, 2018
1 parent b68056b commit 9c2ca58
Show file tree
Hide file tree
Showing 7 changed files with 162 additions and 6 deletions.
Original file line number Diff line number Diff line change
Expand Up @@ -54,7 +54,7 @@ set LATTICE VECTORS FILE = latticeVectors.inp
#kPoint Quadrature Rule file for sampling the Brillouin zone (used only in periodic calculation)
#(Located in $DFT PATH/data/kPointList)
#
set kPOINT RULE FILE = fccBandkPointRule.inp
#set kPOINT RULE FILE = fccBandkPointRule.inp
#set kPOINT RULE FILE = GammaPoint.inp
set BZ SAMPLING POINTS ALONG X = 4
set BZ SAMPLING POINTS ALONG Y = 4
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -78,7 +78,7 @@ set EXCHANGE CORRELATION TYPE = 1
#

#Parameter indicating number of eigenfunctions to be solved
set NUMBER OF KOHN-SHAM WAVEFUNCTIONS = 10
set NUMBER OF KOHN-SHAM WAVEFUNCTIONS = 9
#Approximate estimate of lower bound of the wanted spectrum
set LOWER BOUND WANTED SPECTRUM = -10.0
#Chebyshev Filter polynomial degree
Expand Down
60 changes: 60 additions & 0 deletions include/dealiiOverloadedFunc.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,60 @@
// ---------------------------------------------------------------------
//
// Copyright (c) 2017 The Regents of the University of Michigan and DFT-FE authors.
//
// This file is part of the DFT-FE code.
//
// The DFT-FE code is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE at
// the top level of the DFT-FE distribution.
//
// ---------------------------------------------------------------------
//
// @author Phani Motamarri (2018)
//

#ifndef dealiiOverloadedFunc_H_
#define dealiiOverloadedFunc_H_

#include <vector>

#include "headers.h"

//
//Declare dftUtils functions
//
namespace dftUtils
{

struct constraintMatrixInfo
{
std::vector<unsigned int> rowIdsGlobal;
std::vector<unsigned int> rowIdsLocal;
std::vector<unsigned int> columnIdsLocal;
std::vector<double> columnValues;
std::vector<double> inhomogenities;
std::vector<unsigned int> rowSizes;
};

/**
* convert a given constraintMatrix to simple arrays (STL) for fast access
*/
void convertConstraintMatrixToSTLVector(dealii::parallel::distributed::Vector<double> &fieldVector,
dealii::ConstraintMatrix & constraintMatrixData,
dealii::IndexSet & locally_owned_dofs,
constraintMatrixInfo & constraintMatrixDataInVector);


/**
* overload dealii internal function distribute
*/
void distribute(constraintMatrixInfo &constraintMatrixDataInVector,
dealii::parallel::distributed::Vector<double> &fieldVector);


};

#endif
2 changes: 2 additions & 0 deletions setup.sh
Original file line number Diff line number Diff line change
Expand Up @@ -5,6 +5,8 @@ set -o pipefail
#Provide paths for external libraries and optimization flag (0 for Debug, 1 for Release)
dealiiPetscRealDir="/home/vikramg/DFT-FE-softwares/softwareCentos/dealiiDev/intel_18.0.1_petscReal_noavx"
dealiiPetscComplexDir="/home/vikramg/DFT-FE-softwares/softwareCentos/dealiiDev/intel_18.0.1_petscComplex_noavx"
#dealiiPetscRealDir="/home/vikramg/DFT-FE-softwares/softwareCentos/dealii8.5.1/intel_18.0.1_petscDouble_noavx"
#dealiiPetscComplexDir="/home/vikramg/DFT-FE-softwares/softwareCentos/dealii8.5.1/intel_18.0.1_petscComplex_noavx"
alglibDir="/nfs/mcfs_comp/home/rudraa/software/alglib/cpp/src"
libxcDir="/nfs/mcfs_comp/home/rudraa/software/libxc/libxc-2.2.0/installDir"
spglibDir="/home/vikramg/DFT-FE-softwares/softwareCentos/spglib"
Expand Down
2 changes: 1 addition & 1 deletion src/dft/dft.cc
Original file line number Diff line number Diff line change
Expand Up @@ -391,7 +391,7 @@ void dftClass<FEOrder>::init ()
//
dftUtils::convertConstraintMatrixToSTLVector(vChebyshev,
constraintsNoneEigen,
locally_owned_dofs,
locally_owned_dofsEigen,
constraintsNoneEigenDataInVector);


Expand Down
6 changes: 3 additions & 3 deletions src/eigen/eigen.cc
Original file line number Diff line number Diff line change
Expand Up @@ -875,9 +875,9 @@ void eigenClass<FEOrder>::HX(const std::vector<vectorType*> &src,
*(dftPtr->tempPSI2[i]) = *src[i];
dftPtr->tempPSI2[i]->scale(massVector); //MX
dftPtr->tempPSI2[i]->update_ghost_values();
dftPtr->constraintsNoneEigen.distribute(*(dftPtr->tempPSI2[i]));
//dftUtils::distribute(dftPtr->constraintsNoneEigenDataInVector,
// *(dftPtr->tempPSI2[i]));
//dftPtr->constraintsNoneEigen.distribute(*(dftPtr->tempPSI2[i]));
dftUtils::distribute(dftPtr->constraintsNoneEigenDataInVector,
*(dftPtr->tempPSI2[i]));
dftPtr->tempPSI2[i]->update_ghost_values();
*dst[i] = 0.0;
}
Expand Down
94 changes: 94 additions & 0 deletions utils/dealiiOverloadedFunc.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,94 @@
// ---------------------------------------------------------------------
//
// Copyright (c) 2017 The Regents of the University of Michigan and DFT-FE authors.
//
// This file is part of the DFT-FE code.
//
// The DFT-FE code is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE at
// the top level of the DFT-FE distribution.
//
// ---------------------------------------------------------------------
//
// @author Phani Motamarri (2018)
//
#include "../include/dealiiOverloadedFunc.h"




namespace dftUtils{

void convertConstraintMatrixToSTLVector(dealii::parallel::distributed::Vector<double> &fieldVector,
dealii::ConstraintMatrix & constraintMatrixData,
dealii::IndexSet & locally_owned_dofs,
constraintMatrixInfo & constraintMatrixDataInVector)
{
const std::shared_ptr< const dealii::Utilities::MPI::Partitioner > & partitioner = fieldVector.get_partitioner();
constraintMatrixDataInVector.rowIdsGlobal.clear();
constraintMatrixDataInVector.rowIdsLocal.clear();
constraintMatrixDataInVector.columnIdsLocal.clear();
constraintMatrixDataInVector.columnValues.clear();
constraintMatrixDataInVector.inhomogenities.clear();
constraintMatrixDataInVector.rowSizes.clear();

for(dealii::IndexSet::ElementIterator it = locally_owned_dofs.begin(); it != locally_owned_dofs.end();++it)
{
if(constraintMatrixData.is_constrained(*it))
{
(constraintMatrixDataInVector.rowIdsGlobal).push_back(*it);
(constraintMatrixDataInVector.rowIdsLocal).push_back(partitioner->global_to_local(*it));
}
}

for(unsigned int i = 0; i < (constraintMatrixDataInVector.rowIdsGlobal).size(); ++i)
{
const int lineDof = constraintMatrixDataInVector.rowIdsGlobal[i];
(constraintMatrixDataInVector.inhomogenities).push_back(constraintMatrixData.get_inhomogeneity(lineDof));
const std::vector<std::pair<dealii::types::global_dof_index, double > > * rowData=constraintMatrixData.get_constraint_entries(lineDof);
(constraintMatrixDataInVector.rowSizes).push_back(rowData->size());
for(unsigned int j = 0; j < rowData->size();++j)
{
(constraintMatrixDataInVector.columnIdsLocal).push_back(partitioner->global_to_local((*rowData)[j].first));
(constraintMatrixDataInVector.columnValues).push_back((*rowData)[j].second);
}
}

return;

}



void distribute(constraintMatrixInfo &constraintMatrixDataInVector,
dealii::parallel::distributed::Vector<double> &fieldVector)
{
unsigned int count = 0;
//std::cout<<"Size Local: "<< (constraintMatrixDataInVector.rowIdsLocal).size()<<std::endl;
//std::cout<<"Size Global: "<< (constraintMatrixDataInVector.rowIdsGlobal).size()<<std::endl;

for(unsigned int i = 0; i < (constraintMatrixDataInVector.rowIdsLocal).size();++i)
{
double new_value = constraintMatrixDataInVector.inhomogenities[i];

for(unsigned int j = 0; j < constraintMatrixDataInVector.rowSizes[i]; ++j)
{

new_value += fieldVector.local_element(constraintMatrixDataInVector.columnIdsLocal[count])*constraintMatrixDataInVector.columnValues[count];
count++;
}

fieldVector.local_element(constraintMatrixDataInVector.rowIdsLocal[i]) = new_value;

}

return;
}

}



0 comments on commit 9c2ca58

Please sign in to comment.