Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_IlukGraph.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4//
5// Copyright 2009 NTESS and the Ifpack2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
15
16#ifndef IFPACK2_ILUK_GRAPH_HPP
17#define IFPACK2_ILUK_GRAPH_HPP
18
19#include <algorithm>
20#include <vector>
21
22#include "KokkosSparse_spiluk.hpp"
23
24#include <Ifpack2_ConfigDefs.hpp>
25#include <Teuchos_ParameterList.hpp>
26#include <Teuchos_CommHelpers.hpp>
27#include <Tpetra_CrsGraph.hpp>
28#include <Tpetra_Details_WrappedDualView.hpp>
29#include <Tpetra_Import.hpp>
30#include <Ifpack2_CreateOverlapGraph.hpp>
31#include <Ifpack2_Parameters.hpp>
32
33namespace Ifpack2 {
34
66template<class GraphType, class KKHandleType>
67class IlukGraph : public Teuchos::Describable {
68public:
69 typedef typename GraphType::local_ordinal_type local_ordinal_type;
70 typedef typename GraphType::global_ordinal_type global_ordinal_type;
71 typedef typename GraphType::node_type node_type;
72
74 typedef Tpetra::RowGraph<local_ordinal_type,
75 global_ordinal_type,
76 node_type> row_graph_type;
78 typedef Tpetra::CrsGraph<local_ordinal_type,
79 global_ordinal_type,
80 node_type> crs_graph_type;
81
82
83
84 typedef typename crs_graph_type::nonconst_global_inds_host_view_type nonconst_global_inds_host_view_type;
85 typedef typename crs_graph_type::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
86 typedef typename crs_graph_type::global_inds_host_view_type global_inds_host_view_type;
87 typedef typename crs_graph_type::local_inds_host_view_type local_inds_host_view_type;
88
100 IlukGraph (const Teuchos::RCP<const GraphType>& G,
101 const int levelFill,
102 const int levelOverlap,
103 const double overalloc = 2.);
104
106 virtual ~IlukGraph ();
107
113 void setParameters (const Teuchos::ParameterList& parameterlist);
114
127 void initialize(const Teuchos::RCP<KKHandleType>& KernelHandle);
128
130 int getLevelFill () const { return LevelFill_; }
131
133 int getLevelOverlap () const { return LevelOverlap_; }
134
136 Teuchos::RCP<const GraphType> getA_Graph () const {
137 return Graph_;
138 }
139
141 Teuchos::RCP<crs_graph_type> getL_Graph () const {
142 return L_Graph_;
143 }
144
146 Teuchos::RCP<crs_graph_type> getU_Graph () const {
147 return U_Graph_;
148 }
149
151 Teuchos::RCP<const crs_graph_type> getOverlapGraph () const {
152 return OverlapGraph_;
153 }
154
156 size_t getNumGlobalDiagonals() const { return NumGlobalDiagonals_; }
157
158private:
159 typedef typename GraphType::map_type map_type;
160
170
180
182 void constructOverlapGraph();
183
184 Teuchos::RCP<const GraphType> Graph_;
185 Teuchos::RCP<const crs_graph_type> OverlapGraph_;
186 int LevelFill_;
187 int LevelOverlap_;
188 const double Overalloc_;
189 Teuchos::RCP<crs_graph_type> L_Graph_;
190 Teuchos::RCP<crs_graph_type> U_Graph_;
191 size_t NumMyDiagonals_;
192 size_t NumGlobalDiagonals_;
193};
194
195
196template<class GraphType, class KKHandleType>
198IlukGraph (const Teuchos::RCP<const GraphType>& G,
199 const int levelFill,
200 const int levelOverlap,
201 const double overalloc)
202 : Graph_ (G),
203 LevelFill_ (levelFill),
204 LevelOverlap_ (levelOverlap),
205 Overalloc_ (overalloc),
206 NumMyDiagonals_ (0),
207 NumGlobalDiagonals_ (0)
208{
209 TEUCHOS_TEST_FOR_EXCEPTION(Overalloc_ <= 1., std::runtime_error,
210 "Ifpack2::IlukGraph: FATAL: overalloc must be greater than 1.")
211}
212
213
214template<class GraphType, class KKHandleType>
217
218
219template<class GraphType, class KKHandleType>
221setParameters (const Teuchos::ParameterList& parameterlist)
222{
223 getParameter (parameterlist, "iluk level-of-fill", LevelFill_);
224 getParameter (parameterlist, "iluk level-of-overlap", LevelOverlap_);
225}
226
227
228template<class GraphType, class KKHandleType>
229void IlukGraph<GraphType, KKHandleType>::constructOverlapGraph () {
230 // FIXME (mfh 22 Dec 2013) This won't do if we want
231 // RILUK::initialize() to do the right thing (that is,
232 // unconditionally recompute the "symbolic factorization").
233 if (OverlapGraph_ == Teuchos::null) {
234 OverlapGraph_ = createOverlapGraph<GraphType> (Graph_, LevelOverlap_);
235 }
236}
237
238
239template<class GraphType, class KKHandleType>
241{
242 using Teuchos::Array;
243 using Teuchos::ArrayView;
244 using Teuchos::RCP;
245 using Teuchos::rcp;
246 using Teuchos::REDUCE_SUM;
247 using Teuchos::reduceAll;
248
249 size_t NumIn, NumL, NumU;
250 bool DiagFound;
251
252 constructOverlapGraph();
253
254 // Get Maximum Row length
255 const int MaxNumIndices = OverlapGraph_->getLocalMaxNumRowEntries ();
256
257 // FIXME (mfh 23 Dec 2013) Use size_t or whatever
258 // getLocalNumElements() returns, instead of ptrdiff_t.
259 const int NumMyRows = OverlapGraph_->getRowMap ()->getLocalNumElements ();
260
261 using device_type = typename node_type::device_type;
262 using execution_space = typename device_type::execution_space;
263 using dual_view_type = Kokkos::DualView<size_t*,device_type>;
264 dual_view_type numEntPerRow_dv("numEntPerRow",NumMyRows);
265 Tpetra::Details::WrappedDualView<dual_view_type> numEntPerRow(numEntPerRow_dv);
266
267 const auto overalloc = Overalloc_;
268 const auto levelfill = LevelFill_;
269 {
270 // Scoping for the localOverlapGraph access
271 auto numEntPerRow_d = numEntPerRow.getDeviceView(Tpetra::Access::OverwriteAll);
272 auto localOverlapGraph = OverlapGraph_->getLocalGraphDevice();
273 Kokkos::parallel_for("CountOverlapGraphRowEntries",
274 Kokkos::RangePolicy<execution_space>(0, NumMyRows),
275 KOKKOS_LAMBDA(const int i)
276 {
277 // Heuristic to get the maximum number of entries per row.
278 int RowMaxNumIndices = localOverlapGraph.rowConst(i).length;
279 numEntPerRow_d(i) = (levelfill == 0) ? RowMaxNumIndices // No additional storage needed
280 : Kokkos::ceil(static_cast<double>(RowMaxNumIndices)
281 * Kokkos::pow(overalloc, levelfill));
282 });
283
284 };
285
286 bool insertError; // No error found yet while inserting entries
287 do {
288 insertError = false;
289 Teuchos::ArrayView<const size_t> a_numEntPerRow(numEntPerRow.getHostView(Tpetra::Access::ReadOnly).data(),NumMyRows);
290 L_Graph_ = rcp (new crs_graph_type (OverlapGraph_->getRowMap (),
291 OverlapGraph_->getRowMap (),
292 a_numEntPerRow));
293 U_Graph_ = rcp (new crs_graph_type (OverlapGraph_->getRowMap (),
294 OverlapGraph_->getRowMap (),
295 a_numEntPerRow));
296
297 Array<local_ordinal_type> L (MaxNumIndices);
298 Array<local_ordinal_type> U (MaxNumIndices);
299
300 // First we copy the user's graph into L and U, regardless of fill level
301
302 NumMyDiagonals_ = 0;
303
304 for (int i = 0; i< NumMyRows; ++i) {
305 local_inds_host_view_type my_indices;
306 OverlapGraph_->getLocalRowView (i, my_indices);
307
308 // Split into L and U (we don't assume that indices are ordered).
309
310 NumL = 0;
311 NumU = 0;
312 DiagFound = false;
313 NumIn = my_indices.size();
314
315 for (size_t j = 0; j < NumIn; ++j) {
316 const local_ordinal_type k = my_indices[j];
317
318 if (k<NumMyRows) { // Ignore column elements that are not in the square matrix
319
320 if (k==i) {
321 DiagFound = true;
322 }
323 else if (k < i) {
324 L[NumL] = k;
325 NumL++;
326 }
327 else {
328 U[NumU] = k;
329 NumU++;
330 }
331 }
332 }
333
334 // Check in things for this row of L and U
335
336 if (DiagFound) {
337 ++NumMyDiagonals_;
338 }
339 if (NumL) {
340 L_Graph_->insertLocalIndices (i, NumL, L.data());
341 }
342 if (NumU) {
343 U_Graph_->insertLocalIndices (i, NumU, U.data());
344 }
345 }
346
347 if (LevelFill_ > 0) {
348 // Complete Fill steps
349 RCP<const map_type> L_DomainMap = OverlapGraph_->getRowMap ();
350 RCP<const map_type> L_RangeMap = Graph_->getRangeMap ();
351 RCP<const map_type> U_DomainMap = Graph_->getDomainMap ();
352 RCP<const map_type> U_RangeMap = OverlapGraph_->getRowMap ();
353 RCP<Teuchos::ParameterList> params = Teuchos::parameterList ();
354 params->set ("Optimize Storage",false);
355 L_Graph_->fillComplete (L_DomainMap, L_RangeMap, params);
356 U_Graph_->fillComplete (U_DomainMap, U_RangeMap, params);
357 L_Graph_->resumeFill ();
358 U_Graph_->resumeFill ();
359
360 // At this point L_Graph and U_Graph are filled with the pattern of input graph,
361 // sorted and have redundant indices (if any) removed. Indices are zero based.
362 // LevelFill is greater than zero, so continue...
363
364 int MaxRC = NumMyRows;
365 std::vector<std::vector<int> > Levels(MaxRC);
366 std::vector<int> LinkList(MaxRC);
367 std::vector<int> CurrentLevel(MaxRC);
368 Array<local_ordinal_type> CurrentRow (MaxRC + 1);
369 std::vector<int> LevelsRowU(MaxRC);
370
371 try {
372 for (int i = 0; i < NumMyRows; ++i) {
373 int First, Next;
374
375 // copy column indices of row into workspace and sort them
376
377 size_t LenL = L_Graph_->getNumEntriesInLocalRow(i);
378 size_t LenU = U_Graph_->getNumEntriesInLocalRow(i);
379 size_t Len = LenL + LenU + 1;
380 CurrentRow.resize(Len);
381 nonconst_local_inds_host_view_type CurrentRow_view(CurrentRow.data(),CurrentRow.size());
382 L_Graph_->getLocalRowCopy(i, CurrentRow_view, LenL); // Get L Indices
383 CurrentRow[LenL] = i; // Put in Diagonal
384 if (LenU > 0) {
385 ArrayView<local_ordinal_type> URowView = CurrentRow.view (LenL+1,LenU);
386 nonconst_local_inds_host_view_type URowView_v(URowView.data(),URowView.size());
387
388 // Get U Indices
389 U_Graph_->getLocalRowCopy (i, URowView_v, LenU);
390 }
391
392 // Construct linked list for current row
393
394 for (size_t j=0; j<Len-1; j++) {
395 LinkList[CurrentRow[j]] = CurrentRow[j+1];
396 CurrentLevel[CurrentRow[j]] = 0;
397 }
398
399 LinkList[CurrentRow[Len-1]] = NumMyRows;
400 CurrentLevel[CurrentRow[Len-1]] = 0;
401
402 // Merge List with rows in U
403
404 First = CurrentRow[0];
405 Next = First;
406 while (Next < i) {
407 int PrevInList = Next;
408 int NextInList = LinkList[Next];
409 int RowU = Next;
410 // Get Indices for this row of U
411 local_inds_host_view_type IndicesU;
412 U_Graph_->getLocalRowView (RowU, IndicesU);
413 // FIXME (mfh 23 Dec 2013) size() returns ptrdiff_t, not int.
414 int LengthRowU = IndicesU.size ();
415
416 int ii;
417
418 // Scan RowU
419
420 for (ii = 0; ii < LengthRowU; /*nop*/) {
421 int CurInList = IndicesU[ii];
422 if (CurInList < NextInList) {
423 // new fill-in
424 int NewLevel = CurrentLevel[RowU] + Levels[RowU][ii+1] + 1;
425 if (NewLevel <= LevelFill_) {
426 LinkList[PrevInList] = CurInList;
427 LinkList[CurInList] = NextInList;
428 PrevInList = CurInList;
429 CurrentLevel[CurInList] = NewLevel;
430 }
431 ii++;
432 }
433 else if (CurInList == NextInList) {
434 PrevInList = NextInList;
435 NextInList = LinkList[PrevInList];
436 int NewLevel = CurrentLevel[RowU] + Levels[RowU][ii+1] + 1;
437 CurrentLevel[CurInList] = std::min (CurrentLevel[CurInList],
438 NewLevel);
439 ii++;
440 }
441 else { // (CurInList > NextInList)
442 PrevInList = NextInList;
443 NextInList = LinkList[PrevInList];
444 }
445 }
446 Next = LinkList[Next];
447 }
448
449 // Put pattern into L and U
450 CurrentRow.resize(0);
451
452 Next = First;
453
454 // Lower
455 while (Next < i) {
456 CurrentRow.push_back(Next);
457 Next = LinkList[Next];
458 }
459
460 // FIXME (mfh 23 Dec 2013) It's not clear to me that
461 // removeLocalIndices works like people expect it to work. In
462 // particular, it does not actually change the column Map.
463 L_Graph_->removeLocalIndices (i); // Delete current set of Indices
464 if (CurrentRow.size() > 0) {
465 L_Graph_->insertLocalIndices (i, CurrentRow.size(),CurrentRow.data());
466 }
467
468 // Diagonal
469
470 TEUCHOS_TEST_FOR_EXCEPTION(
471 Next != i, std::runtime_error,
472 "Ifpack2::IlukGraph::initialize: FATAL: U has zero diagonal")
473
474 LevelsRowU[0] = CurrentLevel[Next];
475 Next = LinkList[Next];
476
477 // Upper
478 CurrentRow.resize(0);
479 LenU = 0;
480
481 while (Next < NumMyRows) {
482 LevelsRowU[LenU+1] = CurrentLevel[Next];
483 CurrentRow.push_back (Next);
484 ++LenU;
485 Next = LinkList[Next];
486 }
487
488 // FIXME (mfh 23 Dec 2013) It's not clear to me that
489 // removeLocalIndices works like people expect it to work. In
490 // particular, it does not actually change the column Map.
491
492 U_Graph_->removeLocalIndices (i); // Delete current set of Indices
493 if (LenU > 0) {
494 U_Graph_->insertLocalIndices (i, CurrentRow.size(),CurrentRow.data());
495 }
496
497 // Allocate and fill Level info for this row
498 Levels[i] = std::vector<int> (LenU+1);
499 for (size_t jj=0; jj<LenU+1; jj++) {
500 Levels[i][jj] = LevelsRowU[jj];
501 }
502 }
503 }
504 catch (std::runtime_error &e) {
505 insertError = true;
506 auto numEntPerRow_d = numEntPerRow.getDeviceView(Tpetra::Access::OverwriteAll);
507 Kokkos::parallel_for("CountOverlapGraphRowEntries",
508 Kokkos::RangePolicy<execution_space>(0, NumMyRows),
509 KOKKOS_LAMBDA(const int i)
510 {
511 const auto numRowEnt = numEntPerRow_d(i);
512 numEntPerRow_d(i) = ceil(static_cast<double>((numRowEnt != 0 ? numRowEnt : 1)) * overalloc);
513 });
514 }
515 const int localInsertError = insertError ? 1 : 0;
516 int globalInsertError = 0;
517 reduceAll (* (OverlapGraph_->getRowMap ()->getComm ()), REDUCE_SUM, 1,
518 &localInsertError, &globalInsertError);
519 insertError = globalInsertError > 0;
520 }
521 } while (insertError); // do until all insertions complete successfully
522
523 // Complete Fill steps
524 RCP<const map_type> L_DomainMap = OverlapGraph_->getRowMap ();
525 RCP<const map_type> L_RangeMap = Graph_->getRangeMap ();
526 RCP<const map_type> U_DomainMap = Graph_->getDomainMap ();
527 RCP<const map_type> U_RangeMap = OverlapGraph_->getRowMap ();
528 L_Graph_->fillComplete (L_DomainMap, L_RangeMap);//DoOptimizeStorage is default here...
529 U_Graph_->fillComplete (U_DomainMap, U_RangeMap);//DoOptimizeStorage is default here...
530
531 reduceAll<int, size_t> (* (L_DomainMap->getComm ()), REDUCE_SUM, 1,
532 &NumMyDiagonals_, &NumGlobalDiagonals_);
533}
534
535
536template<class GraphType, class KKHandleType>
537void IlukGraph<GraphType, KKHandleType>::initialize(const Teuchos::RCP<KKHandleType>& KernelHandle)
538{
539 using Teuchos::Array;
540 using Teuchos::ArrayView;
541 using Teuchos::RCP;
542 using Teuchos::rcp;
543 using Teuchos::REDUCE_SUM;
544 using Teuchos::reduceAll;
545
546 typedef typename crs_graph_type::local_graph_device_type local_graph_device_type;
547 typedef typename local_graph_device_type::size_type size_type;
548 typedef typename local_graph_device_type::data_type data_type;
549 typedef typename local_graph_device_type::array_layout array_layout;
550 typedef typename local_graph_device_type::device_type device_type;
551
552 typedef typename Kokkos::View<size_type*, array_layout, device_type> lno_row_view_t;
553 typedef typename Kokkos::View<data_type*, array_layout, device_type> lno_nonzero_view_t;
554
555 constructOverlapGraph();
556
557 // FIXME (mfh 23 Dec 2013) Use size_t or whatever
558 // getLocalNumElements() returns, instead of ptrdiff_t.
559 const int NumMyRows = OverlapGraph_->getRowMap()->getLocalNumElements();
560 auto localOverlapGraph = OverlapGraph_->getLocalGraphDevice();
561
562 if (KernelHandle->get_spiluk_handle()->get_nrows() < static_cast<size_type>(NumMyRows)) {
563 KernelHandle->get_spiluk_handle()->reset_handle(NumMyRows,
564 KernelHandle->get_spiluk_handle()->get_nnzL(),
565 KernelHandle->get_spiluk_handle()->get_nnzU());
566 }
567
568 lno_row_view_t L_row_map("L_row_map", NumMyRows + 1);
569 lno_nonzero_view_t L_entries("L_entries", KernelHandle->get_spiluk_handle()->get_nnzL());
570 lno_row_view_t U_row_map("U_row_map", NumMyRows + 1);
571 lno_nonzero_view_t U_entries("U_entries", KernelHandle->get_spiluk_handle()->get_nnzU());
572
573 bool symbolicError;
574 do {
575 symbolicError = false;
576 try {
577 KokkosSparse::Experimental::spiluk_symbolic( KernelHandle.getRawPtr(), LevelFill_,
578 localOverlapGraph.row_map, localOverlapGraph.entries,
579 L_row_map, L_entries, U_row_map, U_entries );
580 }
581 catch (std::runtime_error &e) {
582 symbolicError = true;
583 data_type nnzL = static_cast<data_type>(Overalloc_)*L_entries.extent(0);
584 data_type nnzU = static_cast<data_type>(Overalloc_)*U_entries.extent(0);
585 KernelHandle->get_spiluk_handle()->reset_handle(NumMyRows, nnzL, nnzU);
586 Kokkos::resize(L_entries, KernelHandle->get_spiluk_handle()->get_nnzL());
587 Kokkos::resize(U_entries, KernelHandle->get_spiluk_handle()->get_nnzU());
588 }
589 const int localSymbolicError = symbolicError ? 1 : 0;
590 int globalSymbolicError = 0;
591 reduceAll (* (OverlapGraph_->getRowMap ()->getComm ()), REDUCE_SUM, 1,
592 &localSymbolicError, &globalSymbolicError);
593 symbolicError = globalSymbolicError > 0;
594 } while (symbolicError);
595
596 Kokkos::resize(L_entries, KernelHandle->get_spiluk_handle()->get_nnzL());
597 Kokkos::resize(U_entries, KernelHandle->get_spiluk_handle()->get_nnzU());
598
599 RCP<Teuchos::ParameterList> params = Teuchos::parameterList ();
600 params->set ("Optimize Storage",false);
601
602 L_Graph_ = rcp (new crs_graph_type (OverlapGraph_->getRowMap (),
603 OverlapGraph_->getRowMap (),
604 L_row_map, L_entries));
605 U_Graph_ = rcp (new crs_graph_type (OverlapGraph_->getRowMap (),
606 OverlapGraph_->getRowMap (),
607 U_row_map, U_entries));
608
609 RCP<const map_type> L_DomainMap = OverlapGraph_->getRowMap ();
610 RCP<const map_type> L_RangeMap = Graph_->getRangeMap ();
611 RCP<const map_type> U_DomainMap = Graph_->getDomainMap ();
612 RCP<const map_type> U_RangeMap = OverlapGraph_->getRowMap ();
613 L_Graph_->fillComplete (L_DomainMap, L_RangeMap, params);
614 U_Graph_->fillComplete (U_DomainMap, U_RangeMap, params);
615}
616
617} // namespace Ifpack2
618
619#endif /* IFPACK2_ILUK_GRAPH_HPP */
Construct a level filled graph for use in computing an ILU(k) incomplete factorization.
Definition Ifpack2_IlukGraph.hpp:67
size_t getNumGlobalDiagonals() const
Returns the global number of diagonals in the ILU(k) graph.
Definition Ifpack2_IlukGraph.hpp:156
Tpetra::CrsGraph< local_ordinal_type, global_ordinal_type, node_type > crs_graph_type
Tpetra::CrsGraph specialization used by this class.
Definition Ifpack2_IlukGraph.hpp:80
Teuchos::RCP< crs_graph_type > getU_Graph() const
Returns the graph of upper triangle of the ILU(k) graph as a Tpetra::CrsGraph.
Definition Ifpack2_IlukGraph.hpp:146
int getLevelFill() const
The level of fill used to construct this graph.
Definition Ifpack2_IlukGraph.hpp:130
int getLevelOverlap() const
The level of overlap used to construct this graph.
Definition Ifpack2_IlukGraph.hpp:133
Teuchos::RCP< const GraphType > getA_Graph() const
Returns the original graph given.
Definition Ifpack2_IlukGraph.hpp:136
void initialize()
Set up the graph structure of the L and U factors.
Definition Ifpack2_IlukGraph.hpp:240
Teuchos::RCP< const crs_graph_type > getOverlapGraph() const
Returns the the overlapped graph.
Definition Ifpack2_IlukGraph.hpp:151
void setParameters(const Teuchos::ParameterList &parameterlist)
Set parameters.
Definition Ifpack2_IlukGraph.hpp:221
Tpetra::RowGraph< local_ordinal_type, global_ordinal_type, node_type > row_graph_type
Tpetra::RowGraph specialization used by this class.
Definition Ifpack2_IlukGraph.hpp:76
Teuchos::RCP< crs_graph_type > getL_Graph() const
Returns the graph of lower triangle of the ILU(k) graph as a Tpetra::CrsGraph.
Definition Ifpack2_IlukGraph.hpp:141
IlukGraph(const Teuchos::RCP< const GraphType > &G, const int levelFill, const int levelOverlap, const double overalloc=2.)
Constructor.
Definition Ifpack2_IlukGraph.hpp:198
virtual ~IlukGraph()
IlukGraph Destructor.
Definition Ifpack2_IlukGraph.hpp:215
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41
void getParameter(const Teuchos::ParameterList &params, const std::string &name, T &value)
Set a value from a ParameterList if a parameter with the specified name exists.
Definition Ifpack2_Parameters.hpp:26
Teuchos::RCP< const GraphType > createOverlapGraph(const Teuchos::RCP< const GraphType > &inputGraph, const int overlapLevel)
Construct an overlapped graph for use with Ifpack2 preconditioners.
Definition Ifpack2_CreateOverlapGraph.hpp:40