IPCC  1.0
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
CKNLanczosLaunching Class Reference

This class includes functions for matrix debugging. More...

#include "KNLanczosLaunching.h"

Collaboration diagram for CKNLanczosLaunching:
Collaboration graph

Classes

struct  LANCZOS_PARAM
 Structure for Lanczos parameters. More...
 

Public Types

typedef struct
CKNLanczosLaunching::LANCZOS_PARAM
LPLANCZOS_PARAM
 

Public Member Functions

 CKNLanczosLaunching ()
 
 ~CKNLanczosLaunching ()
 

Static Public Member Functions

static
CKNLanczosMethod::LPEIGENVALUE_RESULT 
LaunchingLanczos (CKNMatrixOperation::CKNCSR *pA, bool bShowMsg, LPLANCZOS_PARAM lpParam, bool bMPI)
 Lanczos API entry function. More...
 

Detailed Description

This class includes functions for matrix debugging.

Date
10/July/2014
Author
Kyu Nam Cho(mysto.nosp@m.us@k.nosp@m.orea..nosp@m.ac.k.nosp@m.r), Hoon Ryu(elec1.nosp@m.020@.nosp@m.gmail.nosp@m..com)

Definition at line 19 of file KNLanczosLaunching.h.

Member Typedef Documentation

Constructor & Destructor Documentation

CKNLanczosLaunching::CKNLanczosLaunching ( )

Definition at line 26 of file KNLanczosLaunching.cpp.

27 {
28 }
CKNLanczosLaunching::~CKNLanczosLaunching ( )

Definition at line 30 of file KNLanczosLaunching.cpp.

31 {
32 }

Member Function Documentation

CKNLanczosMethod::LPEIGENVALUE_RESULT CKNLanczosLaunching::LaunchingLanczos ( CKNMatrixOperation::CKNCSR pA,
bool  bShowMsg,
LPLANCZOS_PARAM  lpParam,
bool  bMPI 
)
static

Lanczos API entry function.

Parameters
pACSR matrix for eigenvalue solving
bShowMsgFlag for display message
lpParamSeveral parameters for Lanczos method setting
bMPIApplying MPI environment or not

< Remove output directory first

< Launching lanczos method

Definition at line 40 of file KNLanczosLaunching.cpp.

References CKNMatrixOperation::AllocateLocalCSR(), CKNLanczosLaunching::LANCZOS_PARAM::bCalculateEigenVectors, CKNLanczosLaunching::LANCZOS_PARAM::bCalculateWaveFunction, CKNLanczosLaunching::LANCZOS_PARAM::bDoSelectiveReorthogonalization, CKNLanczosLaunching::LANCZOS_PARAM::bSaveResultToFile, CKNMatrixOperation::BuildLocalCSR(), CALCULATION_SUCCESS, CKNLanczosMethod::DoLanczosMethod(), CKNLanczosLaunching::LANCZOS_PARAM::fConvergeceCriteria, CKNLanczosLaunching::LANCZOS_PARAM::fevMax, CKNLanczosLaunching::LANCZOS_PARAM::fevMin, CKNMatrixOperation::FreeLocalCSR(), CKNMPIManager::GetCurrentRank(), CKNMatrixOperation::CKNCSR::GetNoneZeroCount(), CKNMPIManager::GetTotalNodeCount(), CKNMPIManager::InitCommunicationBufferMetric(), CKNMPIManager::InitLevel(), CKNTimeMeasurement::InitTimer(), CKNMPIManager::IsMultiLevelMPI(), CKNMPIManager::IsRootRank(), CKNLanczosLaunching::LANCZOS_PARAM::load_in_MIC, CKNMatrixOperation::CKNCSR::m_vectColumn, CKNMatrixOperation::CKNCSR::m_vectRow, CKNMatrixOperation::CKNCSR::m_vectValueImaginaryBuffer, CKNMatrixOperation::CKNCSR::m_vectValueRealBuffer, CKNLanczosMethod::MergeDegeneratedEigenvalues(), CKNLanczosLaunching::LANCZOS_PARAM::nCheckEigenvalueInterval, CKNLanczosLaunching::LANCZOS_PARAM::nFindingDegeneratedEVCount, CKNLanczosLaunching::LANCZOS_PARAM::nFindingEigenValueCount, CKNLanczosLaunching::LANCZOS_PARAM::nLanczosIterationCount, CKNLanczosLaunching::LANCZOS_PARAM::nMPILevel, phi_tid, CKNLanczosMethod::RecalcuWaveFunction(), CKNLanczosMethod::SaveLanczosResult(), CKNMPIManager::SetMPIEnviroment(), CKNIPCCUtility::SetShow(), SHOW_SIMPLE_MSG, CKNLanczosMethod::ShowLanczosResult(), CKNIPCCUtility::ShowMsg(), CKNTimeMeasurement::TotalMeasurementEnd(), CKNTimeMeasurement::TotalMeasurementStart(), and CKNMatrixOperation::UpdateLocalCSR().

41 {
42  CKNMatrixOperation::CKNCSR *pHamiltonian = pA;
43  CKNMatrixOperation::CKNCSR *plocalH = NULL;
44  CKNMatrixOperation::CKNCSR *plocalHm1 = NULL;
45  CKNMatrixOperation::CKNCSR *plocalHp1 = NULL;
46  unsigned int nRtnCode = CALCULATION_SUCCESS;
47  unsigned int i, nRepeatCount = 1;
49  CKNLanczosMethod lanczos;
50  char szMsg[1024];
51  double *pKValue[3] = {NULL, NULL, NULL};
52 
53  if (bShowMsg)
55 
56  int rank;
57  int world_size;
58 
59  if (bMPI)
60  {
61  MPI_Comm_size(MPI_COMM_WORLD, &world_size);
62  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
63 
64  CKNMPIManager::SetMPIEnviroment(rank, world_size);
65  if(!CKNMPIManager::InitLevel(lpParam->nMPILevel, lpParam->nFindingDegeneratedEVCount))
66  {
67  SHOW_SIMPLE_MSG("Adjusting MPI level and node count first!\n");
68  return lpResult;
69  }
71  }
72 
73 
75  {
76 #ifdef _WIN32
77  system("rmdir /s /q result");
78 #else //WIN32
79  system("rm -rf result");
80 #endif //WIN32
81 
82 #pragma omp parallel
83 {
84  if(omp_get_num_threads() > 1)
85  {
86  if(omp_get_thread_num() == 0)
87  {
88  sprintf(szMsg, "-Run with MPI (%d cores) and OpenMP (%d threads)\n", CKNMPIManager::GetTotalNodeCount(), omp_get_num_threads());
90  }
91  }
92  else
93  {
94  sprintf(szMsg, "-Run with MPI only (%d cores)\n", CKNMPIManager::GetTotalNodeCount());
96  }
97 }
98  }
99 
101  CKNMatrixOperation::AllocateLocalCSR(&plocalH, &plocalHm1, &plocalHp1);
102 
103 
104 #ifndef DISABLE_MPI_ROUTINE
106 // FIXME test code
107 #ifndef _WIN32
108  char host_name[256];
109  gethostname(host_name, 256);
110  //printf("[%s] MPI rank: %d | phi tid: %d\n", host_name, CKNMPIManager::GetCurrentRank(), phi_tid);
111 #endif //_WIN32
112 #else
113  phi_tid = 0;
114 #endif
115 
116  double *local_real = NULL;
117  double *local_imaginary = NULL;
118  unsigned int *local_row = NULL;
119  unsigned int *local_col = NULL;
120  unsigned int local_size = 0;
121  unsigned int local_row_size = 0;
122  unsigned int local_col_size = 0;
123 
124  double *left_real = NULL;
125  double *left_imaginary = NULL;
126  unsigned int *left_row = NULL;
127  unsigned int *left_col = NULL;
128  unsigned int left_size = 0;
129  unsigned int left_row_size = 0;
130  unsigned int left_col_size = 0;
131 
132  double *right_real = NULL;
133  double *right_imaginary = NULL;
134  unsigned int *right_row = NULL;
135  unsigned int *right_col = NULL;
136  unsigned int right_size = 0;
137  unsigned int right_row_size = 0;
138  unsigned int right_col_size = 0;
140  {
141  CKNMatrixOperation::BuildLocalCSR(pHamiltonian, plocalH, plocalHm1, plocalHp1);
142  SHOW_SIMPLE_MSG("-Local matrices built\n");
143 
144  local_real = plocalH->m_vectValueRealBuffer.data();
145  local_imaginary = plocalH->m_vectValueImaginaryBuffer.data();
146  local_row = plocalH->m_vectRow.data();
147  local_col = plocalH->m_vectColumn.data();
148  local_size = plocalH->GetNoneZeroCount();
149  local_row_size = plocalH->m_vectRow.size();
150  local_col_size = plocalH->m_vectColumn.size();
151 #pragma offload_transfer target(mic:phi_tid) \
152  nocopy(local_real[0:local_size] : ALLOC) \
153  nocopy(local_imaginary[0:local_size] : ALLOC) \
154  nocopy(local_row[0:local_row_size] : ALLOC) \
155  nocopy(local_col[0:local_col_size] : ALLOC)
156 
157  left_real = plocalHm1->m_vectValueRealBuffer.data();
158  left_imaginary = plocalHm1->m_vectValueImaginaryBuffer.data();
159  left_row = plocalHm1->m_vectRow.data();
160  left_col = plocalHm1->m_vectColumn.data();
161  left_size = plocalHm1->GetNoneZeroCount();
162  left_row_size = plocalHm1->m_vectRow.size();
163  left_col_size = plocalHm1->m_vectColumn.size();
164 #pragma offload_transfer target(mic:phi_tid) \
165  nocopy(left_real[0:left_size] : ALLOC) \
166  nocopy(left_imaginary[0:left_size] : ALLOC) \
167  nocopy(left_row[0:left_row_size] : ALLOC) \
168  nocopy(left_col[0:left_col_size] : ALLOC)
169 
170 // FIXME jinpil: error occurred without the following statement
171 // error message: "offload error: unexpected variable type"
172 // any statements with side effects can also avoid the error, e.g. MPI_Barrier(MPI_COMM_WORLD)
173  SHOW_SIMPLE_MSG("-Xeon Phi: mem alloc\n");
174 // the following statement cannot avoid the error
175 // compiler bug?
176 // int dummy = 0;
177 // -----
178 
179  right_real = plocalHp1->m_vectValueRealBuffer.data();
180  right_imaginary = plocalHp1->m_vectValueImaginaryBuffer.data();
181  right_row = plocalHp1->m_vectRow.data();
182  right_col = plocalHp1->m_vectColumn.data();
183  right_size = plocalHp1->GetNoneZeroCount();
184  right_row_size = plocalHp1->m_vectRow.size();
185  right_col_size = plocalHp1->m_vectColumn.size();
186 #pragma offload_transfer target(mic:phi_tid) \
187  nocopy(right_real[0:right_size] : ALLOC) \
188  nocopy(right_imaginary[0:right_size] : ALLOC) \
189  nocopy(right_row[0:right_row_size] : ALLOC) \
190  nocopy(right_col[0:right_col_size] : ALLOC)
191  }
192  else
193  {
194  local_real = pHamiltonian->m_vectValueRealBuffer.data();
195  local_imaginary = pHamiltonian->m_vectValueImaginaryBuffer.data();
196  local_row = pHamiltonian->m_vectRow.data();
197  local_col = pHamiltonian->m_vectColumn.data();
198  local_size = pHamiltonian->GetNoneZeroCount();
199  local_row_size = pHamiltonian->m_vectRow.size();
200  local_col_size = pHamiltonian->m_vectColumn.size();
201 #pragma offload_transfer target(mic:phi_tid) \
202  nocopy(local_real[0:local_size] : ALLOC) \
203  nocopy(local_imaginary[0:local_size] : ALLOC) \
204  nocopy(local_row[0:local_row_size] : ALLOC) \
205  nocopy(local_col[0:local_col_size] : ALLOC)
206  }
207 
209  {
210  CKNMatrixOperation::UpdateLocalCSR(pHamiltonian, plocalH, plocalHm1, plocalHp1);
211  SHOW_SIMPLE_MSG("-Local matrices updated\n");
212 
213 #pragma offload_transfer target(mic:phi_tid) \
214  in(local_real[0:local_size] : REUSE) \
215  in(local_imaginary[0:local_size] : REUSE) \
216  in(local_row[0:local_row_size] : REUSE) \
217  in(local_col[0:local_col_size] : REUSE)
218 
219 #pragma offload_transfer target(mic:phi_tid) \
220  in(left_real[0:left_size] : REUSE) \
221  in(left_imaginary[0:left_size] : REUSE) \
222  in(left_row[0:left_row_size] : REUSE) \
223  in(left_col[0:left_col_size] : REUSE)
224 
225 #pragma offload_transfer target(mic:phi_tid) \
226  in(right_real[0:right_size] : REUSE) \
227  in(right_imaginary[0:right_size] : REUSE) \
228  in(right_row[0:right_row_size] : REUSE) \
229  in(right_col[0:right_col_size] : REUSE)
230  }
231  else
232  {
233 #pragma offload_transfer target(mic:phi_tid) \
234  in(local_real[0:local_size] : REUSE) \
235  in(local_imaginary[0:local_size] : REUSE) \
236  in(local_row[0:local_row_size] : REUSE) \
237  in(local_col[0:local_col_size] : REUSE)
238  }
239 
243  lpResult = lanczos.DoLanczosMethod(pHamiltonian,
244  lpParam->nLanczosIterationCount,
245  lpParam->nCheckEigenvalueInterval,
246  lpParam->nFindingEigenValueCount,
247  lpParam->fevMin,
248  lpParam->fevMax,
249  lpParam->fConvergeceCriteria,
250  lpParam->bDoSelectiveReorthogonalization,
251  lpParam->bCalculateEigenVectors,
252  lpParam->bCalculateWaveFunction,
253  lpParam->load_in_MIC,
254  plocalH, // ADDED
255  plocalHm1, // ADDED
256  plocalHp1);// ADDED
257 
259  {
260  CKNLanczosMethod::MergeDegeneratedEigenvalues(lpResult, lpParam->nFindingDegeneratedEVCount, pHamiltonian, plocalH, plocalHm1, plocalHp1);
261  if( lpParam->bCalculateWaveFunction )
263  }
264 
265  if( lpParam->bSaveResultToFile )
266  {
268  CKNIPCCUtility::ShowMsg("\n-Save result into file\n");
269 
270  CKNLanczosMethod::SaveLanczosResult(lpResult, lpParam->bCalculateEigenVectors, lpParam->bCalculateWaveFunction, NULL, 0);
271  }
272 
274 
275  CKNLanczosMethod::ShowLanczosResult(lpResult, lpParam->bCalculateEigenVectors, lpParam->bCalculateWaveFunction, NULL, 0);
276 
278  SHOW_SIMPLE_MSG("\n-Finished Lanczos eigenvalue solver.\n\n");
279 
281  {
282  CKNMatrixOperation::FreeLocalCSR(plocalH, plocalHm1, plocalHp1);
283 
284 #pragma offload_transfer target(mic:phi_tid) \
285  nocopy(local_real : FREE) \
286  nocopy(local_imaginary : FREE) \
287  nocopy(local_row : FREE) \
288  nocopy(local_col : FREE)
289 
290 #pragma offload_transfer target(mic:phi_tid) \
291  nocopy(left_real : FREE) \
292  nocopy(left_imaginary : FREE) \
293  nocopy(left_row : FREE) \
294  nocopy(left_col : FREE)
295 
296 #pragma offload_transfer target(mic:phi_tid) \
297  nocopy(right_real : FREE) \
298  nocopy(right_imaginary : FREE) \
299  nocopy(right_row : FREE) \
300  nocopy(right_col : FREE)
301  }
302  else
303  {
304 #pragma offload_transfer target(mic:phi_tid) \
305  nocopy(local_real : FREE) \
306  nocopy(local_imaginary : FREE) \
307  nocopy(local_row : FREE) \
308  nocopy(local_col : FREE)
309  }
310 
311  return lpResult;
312 }
static void InitCommunicationBufferMetric()
Initializing MPI Communication buffer for MVMul.
static void ShowMsg(char *pszBuffer)
Show message.
static void SetMPIEnviroment(int nRank, int nTotalNode)
Set MPI Enviroment.
static void MergeDegeneratedEigenvalues(CKNLanczosMethod::LPEIGENVALUE_RESULT lpResult, unsigned int nFindingDegeneratedEVCount, CKNMatrixOperation::CKNCSR *pA, CKNMatrixOperation::CKNCSR *pLocalBlock, CKNMatrixOperation::CKNCSR *pLeftBlock, CKNMatrixOperation::CKNCSR *pRightBlock)
Merging eigenvalue into mater group.
static bool InitLevel(int nMPILevel, int nFindingDegeneratedEVCount)
Init MPI Level, most low level is for multi node cacluation for Lanczos.
unsigned int GetNoneZeroCount()
Getting numbers of none zero elements.
double_vector_t m_vectValueRealBuffer
A member variable for saving none zero elements.
Data and operation representation of CSR(Compressed Sparse Row)
This class for doing Lanczos method.
double_vector_t m_vectValueImaginaryBuffer
A member variable for saving none zero elements.
static int GetTotalNodeCount()
Definition: KNMPIManager.h:44
#define CALCULATION_SUCCESS
Return code at main loop, every steps completed.
Definition: CKNGlobal.h:70
#define SHOW_SIMPLE_MSG(message)
Definition: CKNGlobal.h:41
uint_vector_t m_vectColumn
A member variable for saving column information.
static void RecalcuWaveFunction(CKNLanczosMethod::LPEIGENVALUE_RESULT lpResult)
Recalculating wavefunction after merging degenerated eigenvalues.
static void SaveLanczosResult(CKNLanczosMethod::LPEIGENVALUE_RESULT lpResult, bool bCalcuEigenvalue, bool bWaveFunction, double *pKValue, int nRepeatCount)
Saving Lanczos computation result into file.
static void SetShow(bool bShow)
Definition: KNIPCCUtility.h:30
Structure for engienvalue computing.
static void UpdateLocalCSR(CKNMatrixOperation::CKNCSR *source, CKNMatrixOperation::CKNCSR *mine, CKNMatrixOperation::CKNCSR *left, CKNMatrixOperation::CKNCSR *right)
LPEIGENVALUE_RESULT DoLanczosMethod(CKNMatrixOperation::CKNCSR *pAMatrix, unsigned int nIterationCount, unsigned int nEigenValueCheckInterval, unsigned int nEigenValueCount, double fEigenvalueMin, double fEignevalueMax, double fConvergenceTolerance, bool bReorthogonalization, bool bCalcuEigVector, bool bWaveFunction, double load_in_MIC, CKNMatrixOperation::CKNCSR *pmylocalblock=NULL, CKNMatrixOperation::CKNCSR *leftlocalblock=NULL, CKNMatrixOperation::CKNCSR *rightlocalblock=NULL)
Doing lanczos method.
static int GetCurrentRank()
Definition: KNMPIManager.h:42
static void TotalMeasurementEnd()
Measurement end for total taken time.
static void BuildLocalCSR(CKNMatrixOperation::CKNCSR *source, CKNMatrixOperation::CKNCSR *mine, CKNMatrixOperation::CKNCSR *left, CKNMatrixOperation::CKNCSR *right)
static void TotalMeasurementStart()
Measurement start for total taken time.
static void AllocateLocalCSR(CKNMatrixOperation::CKNCSR **mine, CKNMatrixOperation::CKNCSR **left, CKNMatrixOperation::CKNCSR **right)
static bool IsRootRank()
Get Total node count.
uint_vector_t m_vectRow
A member variable for saving row information.
int phi_tid
static bool IsMultiLevelMPI()
Get MPI_Comm.
Definition: KNMPIManager.h:68
static void InitTimer()
Init time related variable.
static void ShowLanczosResult(CKNLanczosMethod::LPEIGENVALUE_RESULT lpResult, bool bCalculateEigenVectors, bool bCalculateWaveFunction, double *pKValue, int nRepeatCount)
Save calculating result into file.
static void FreeLocalCSR(CKNMatrixOperation::CKNCSR *mine, CKNMatrixOperation::CKNCSR *left, CKNMatrixOperation::CKNCSR *right)

Here is the call graph for this function:


The documentation for this class was generated from the following files: