Tpetra parallel linear algebra  Version of the Day
TpetraExt_MatrixMatrix_OpenMP.hpp
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Tpetra: Templated Linear Algebra Services Package
5 // Copyright (2008) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 // @HEADER
41 
42 #ifndef TPETRA_MATRIXMATRIX_OPENMP_DEF_HPP
43 #define TPETRA_MATRIXMATRIX_OPENMP_DEF_HPP
44 
45 #ifdef HAVE_TPETRA_INST_OPENMP
46 namespace Tpetra {
47 namespace MMdetails {
48 
49 /*********************************************************************************************************/
50 // MMM KernelWrappers for Partial Specialization to OpenMP
51 template<class Scalar,
52  class LocalOrdinal,
53  class GlobalOrdinal, class LocalOrdinalViewType>
54 struct KernelWrappers<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode,LocalOrdinalViewType> {
55  static inline void mult_A_B_newmatrix_kernel_wrapper(CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& Aview,
56  CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& Bview,
57  const LocalOrdinalViewType & Acol2Brow,
58  const LocalOrdinalViewType & Acol2Irow,
59  const LocalOrdinalViewType & Bcol2Ccol,
60  const LocalOrdinalViewType & Icol2Ccol,
61  CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& C,
62  Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode> > Cimport,
63  const std::string& label = std::string(),
64  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
65 
66  static inline void mult_A_B_reuse_kernel_wrapper(CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& Aview,
67  CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& Bview,
68  const LocalOrdinalViewType & Acol2Brow,
69  const LocalOrdinalViewType & Acol2Irow,
70  const LocalOrdinalViewType & Bcol2Ccol,
71  const LocalOrdinalViewType & Icol2Ccol,
72  CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& C,
73  Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode> > Cimport,
74  const std::string& label = std::string(),
75  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
76 };
77 
78 
79 // Jacobi KernelWrappers for Partial Specialization to OpenMP
80 template<class Scalar,
81  class LocalOrdinal,
82  class GlobalOrdinal, class LocalOrdinalViewType>
83 struct KernelWrappers2<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode,LocalOrdinalViewType> {
84  static inline void jacobi_A_B_newmatrix_kernel_wrapper(Scalar omega,
85  const Vector<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode> & Dinv,
86  CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& Aview,
87  CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& Bview,
88  const LocalOrdinalViewType & Acol2Brow,
89  const LocalOrdinalViewType & Acol2Irow,
90  const LocalOrdinalViewType & Bcol2Ccol,
91  const LocalOrdinalViewType & Icol2Ccol,
92  CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& C,
93  Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode> > Cimport,
94  const std::string& label = std::string(),
95  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
96 
97  static inline void jacobi_A_B_reuse_kernel_wrapper(Scalar omega,
98  const Vector<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode> & Dinv,
99  CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& Aview,
100  CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& Bview,
101  const LocalOrdinalViewType & Acol2Brow,
102  const LocalOrdinalViewType & Acol2Irow,
103  const LocalOrdinalViewType & Bcol2Ccol,
104  const LocalOrdinalViewType & Icol2Ccol,
105  CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& C,
106  Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode> > Cimport,
107  const std::string& label = std::string(),
108  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
109 
110 };
111 
112 
113 /*********************************************************************************************************/
114 template<class Scalar,
115  class LocalOrdinal,
116  class GlobalOrdinal,
117  class LocalOrdinalViewType>
118 void KernelWrappers<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode,LocalOrdinalViewType>::mult_A_B_newmatrix_kernel_wrapper(CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& Aview,
119  CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& Bview,
120  const LocalOrdinalViewType & Acol2Brow,
121  const LocalOrdinalViewType & Acol2Irow,
122  const LocalOrdinalViewType & Bcol2Ccol,
123  const LocalOrdinalViewType & Icol2Ccol,
124  CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& C,
125  Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode> > Cimport,
126  const std::string& label,
127  const Teuchos::RCP<Teuchos::ParameterList>& params) {
128 
129 #ifdef HAVE_TPETRA_MMM_TIMINGS
130  std::string prefix_mmm = std::string("TpetraExt ") + label + std::string(": ");
131  using Teuchos::TimeMonitor;
132  Teuchos::RCP<TimeMonitor> MM;
133 #endif
134 
135  // Node-specific code
136  std::string nodename("OpenMP");
137 
138  // Lots and lots of typedefs
139  using Teuchos::RCP;
141  typedef typename KCRS::device_type device_t;
142  typedef typename KCRS::StaticCrsGraphType graph_t;
143  typedef typename graph_t::row_map_type::non_const_type lno_view_t;
144  typedef typename graph_t::entries_type::non_const_type lno_nnz_view_t;
145  typedef typename KCRS::values_type::non_const_type scalar_view_t;
146 
147  // Options
148  int team_work_size = 16; // Defaults to 16 as per Deveci 12/7/16 - csiefer
149  std::string myalg("SPGEMM_KK_MEMORY");
150 
151 
152  if(!params.is_null()) {
153  if(params->isParameter("openmp: algorithm"))
154  myalg = params->get("openmp: algorithm",myalg);
155  if(params->isParameter("openmp: team work size"))
156  team_work_size = params->get("openmp: team work size",team_work_size);
157  }
158 
159  if(myalg == "LTG") {
160  // Use the LTG kernel if requested
161  ::Tpetra::MatrixMatrix::ExtraKernels::mult_A_B_newmatrix_LowThreadGustavsonKernel(Aview,Bview,Acol2Brow,Acol2Irow,Bcol2Ccol,Icol2Ccol,C,Cimport,label,params);
162  }
163  else {
164  // Use the Kokkos-Kernels OpenMP Kernel
165 #ifdef HAVE_TPETRA_MMM_TIMINGS
166  MM = rcp(new TimeMonitor (*TimeMonitor::getNewTimer(prefix_mmm + std::string("MMM Newmatrix OpenMPWrapper"))));
167 #endif
168  // KokkosKernelsHandle
169  typedef KokkosKernels::Experimental::KokkosKernelsHandle<
170  typename lno_view_t::const_value_type,typename lno_nnz_view_t::const_value_type, typename scalar_view_t::const_value_type,
171  typename device_t::execution_space, typename device_t::memory_space,typename device_t::memory_space > KernelHandle;
172 
173  // Grab the Kokkos::SparseCrsMatrices
174  const KCRS & Ak = Aview.origMatrix->getLocalMatrix();
175  // const KCRS & Bk = Bview.origMatrix->getLocalMatrix();
176 
177  // Get the algorithm mode
178  std::string alg = nodename+std::string(" algorithm");
179  // printf("DEBUG: Using kernel: %s\n",myalg.c_str());
180  if(!params.is_null() && params->isParameter(alg)) myalg = params->get(alg,myalg);
181  KokkosSparse::SPGEMMAlgorithm alg_enum = KokkosSparse::StringToSPGEMMAlgorithm(myalg);
182 
183  // Merge the B and Bimport matrices
184  const KCRS Bmerged = Tpetra::MMdetails::merge_matrices(Aview,Bview,Acol2Brow,Acol2Irow,Bcol2Ccol,Icol2Ccol,C.getColMap()->getNodeNumElements());
185 
186 #ifdef HAVE_TPETRA_MMM_TIMINGS
187  MM = rcp(new TimeMonitor (*TimeMonitor::getNewTimer(prefix_mmm + std::string("MMM Newmatrix OpenMPCore"))));
188 #endif
189 
190  // Do the multiply on whatever we've got
191  typename KernelHandle::nnz_lno_t AnumRows = Ak.numRows();
192  // typename KernelHandle::nnz_lno_t BnumRows = Bmerged->numRows();
193  // typename KernelHandle::nnz_lno_t BnumCols = Bmerged->numCols();
194  typename KernelHandle::nnz_lno_t BnumRows = Bmerged.numRows();
195  typename KernelHandle::nnz_lno_t BnumCols = Bmerged.numCols();
196 
197 
198  lno_view_t row_mapC ("non_const_lnow_row", AnumRows + 1);
199  lno_nnz_view_t entriesC;
200  scalar_view_t valuesC;
201  KernelHandle kh;
202  kh.create_spgemm_handle(alg_enum);
203  kh.set_team_work_size(team_work_size);
204  // KokkosSparse::Experimental::spgemm_symbolic(&kh,AnumRows,BnumRows,BnumCols,Ak.graph.row_map,Ak.graph.entries,false,Bmerged->graph.row_map,Bmerged->graph.entries,false,row_mapC);
205  KokkosSparse::Experimental::spgemm_symbolic(&kh,AnumRows,BnumRows,BnumCols,Ak.graph.row_map,Ak.graph.entries,false,Bmerged.graph.row_map,Bmerged.graph.entries,false,row_mapC);
206 
207  size_t c_nnz_size = kh.get_spgemm_handle()->get_c_nnz();
208  if (c_nnz_size){
209  entriesC = lno_nnz_view_t (Kokkos::ViewAllocateWithoutInitializing("entriesC"), c_nnz_size);
210  valuesC = scalar_view_t (Kokkos::ViewAllocateWithoutInitializing("valuesC"), c_nnz_size);
211  }
212  // KokkosSparse::Experimental::spgemm_numeric(&kh,AnumRows,BnumRows,BnumCols,Ak.graph.row_map,Ak.graph.entries,Ak.values,false,Bmerged->graph.row_map,Bmerged->graph.entries,Bmerged->values,false,row_mapC,entriesC,valuesC);
213  KokkosSparse::Experimental::spgemm_numeric(&kh,AnumRows,BnumRows,BnumCols,Ak.graph.row_map,Ak.graph.entries,Ak.values,false,Bmerged.graph.row_map,Bmerged.graph.entries,Bmerged.values,false,row_mapC,entriesC,valuesC);
214  kh.destroy_spgemm_handle();
215 
216 #ifdef HAVE_TPETRA_MMM_TIMINGS
217  MM = rcp(new TimeMonitor (*TimeMonitor::getNewTimer(prefix_mmm + std::string("MMM Newmatrix OpenMPSort"))));
218 #endif
219  // Sort & set values
220  if (params.is_null() || params->get("sort entries",true))
221  Import_Util::sortCrsEntries(row_mapC, entriesC, valuesC);
222  C.setAllValues(row_mapC,entriesC,valuesC);
223 
224  }// end OMP KokkosKernels loop
225 
226 #ifdef HAVE_TPETRA_MMM_TIMINGS
227  MM = rcp(new TimeMonitor (*TimeMonitor::getNewTimer(prefix_mmm + std::string("MMM Newmatrix OpenMPESFC"))));
228 #endif
229 
230  // Final Fillcomplete
231  RCP<Teuchos::ParameterList> labelList = rcp(new Teuchos::ParameterList);
232  labelList->set("Timer Label",label);
233  if(!params.is_null()) labelList->set("compute global constants",params->get("compute global constants",true));
234  RCP<const Export<LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode> > dummyExport;
235  C.expertStaticFillComplete(Bview.origMatrix->getDomainMap(), Aview.origMatrix->getRangeMap(), Cimport,dummyExport,labelList);
236 
237 #if 0
238  {
239  Teuchos::ArrayRCP< const size_t > Crowptr;
240  Teuchos::ArrayRCP< const LocalOrdinal > Ccolind;
241  Teuchos::ArrayRCP< const Scalar > Cvalues;
242  C.getAllValues(Crowptr,Ccolind,Cvalues);
243 
244  // DEBUG
245  int MyPID = C->getComm()->getRank();
246  printf("[%d] Crowptr = ",MyPID);
247  for(size_t i=0; i<(size_t) Crowptr.size(); i++) {
248  printf("%3d ",(int)Crowptr.getConst()[i]);
249  }
250  printf("\n");
251  printf("[%d] Ccolind = ",MyPID);
252  for(size_t i=0; i<(size_t)Ccolind.size(); i++) {
253  printf("%3d ",(int)Ccolind.getConst()[i]);
254  }
255  printf("\n");
256  fflush(stdout);
257  // END DEBUG
258  }
259 #endif
260 }
261 
262 
263 /*********************************************************************************************************/
264 template<class Scalar,
265  class LocalOrdinal,
266  class GlobalOrdinal,
267  class LocalOrdinalViewType>
268 void KernelWrappers<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode,LocalOrdinalViewType>::mult_A_B_reuse_kernel_wrapper(CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& Aview,
269  CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& Bview,
270  const LocalOrdinalViewType & Acol2Brow,
271  const LocalOrdinalViewType & Acol2Irow,
272  const LocalOrdinalViewType & Bcol2Ccol,
273  const LocalOrdinalViewType & Icol2Ccol,
274  CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& C,
275  Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode> > Cimport,
276  const std::string& label,
277  const Teuchos::RCP<Teuchos::ParameterList>& params) {
278 #ifdef HAVE_TPETRA_MMM_TIMINGS
279  std::string prefix_mmm = std::string("TpetraExt ") + label + std::string(": ");
280  using Teuchos::TimeMonitor;
281  Teuchos::RCP<TimeMonitor> MM;
282 #endif
283 
284  // Lots and lots of typedefs
285  using Teuchos::RCP;
286 
287  // Options
288  int team_work_size = 16; // Defaults to 16 as per Deveci 12/7/16 - csiefer
289  std::string myalg("LTG");
290  if(!params.is_null()) {
291  if(params->isParameter("openmp: algorithm"))
292  myalg = params->get("openmp: algorithm",myalg);
293  if(params->isParameter("openmp: team work size"))
294  team_work_size = params->get("openmp: team work size",team_work_size);
295  }
296 
297  if(myalg == "LTG") {
298  // Use the LTG kernel if requested
299  ::Tpetra::MatrixMatrix::ExtraKernels::mult_A_B_reuse_LowThreadGustavsonKernel(Aview,Bview,Acol2Brow,Acol2Irow,Bcol2Ccol,Icol2Ccol,C,Cimport,label,params);
300  }
301  else {
302  throw std::runtime_error("Tpetra::MatrixMatrix::MMM reuse unknown kernel");
303  }
304 
305 #ifdef HAVE_TPETRA_MMM_TIMINGS
306  MM = rcp(new TimeMonitor (*TimeMonitor::getNewTimer(prefix_mmm + std::string("MMM Reuse OpenMPESFC"))));
307 #endif
308  C.fillComplete(C.getDomainMap(), C.getRangeMap());
309 }
310 
311 
312 /*********************************************************************************************************/
313 template<class Scalar,
314  class LocalOrdinal,
315  class GlobalOrdinal,
316  class LocalOrdinalViewType>
317 void KernelWrappers2<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode,LocalOrdinalViewType>::jacobi_A_B_newmatrix_kernel_wrapper(Scalar omega,
318  const Vector<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode> & Dinv,
319  CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& Aview,
320  CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& Bview,
321  const LocalOrdinalViewType & Acol2Brow,
322  const LocalOrdinalViewType & Acol2Irow,
323  const LocalOrdinalViewType & Bcol2Ccol,
324  const LocalOrdinalViewType & Icol2Ccol,
325  CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& C,
326  Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode> > Cimport,
327  const std::string& label,
328  const Teuchos::RCP<Teuchos::ParameterList>& params) {
329 
330 #ifdef HAVE_TPETRA_MMM_TIMINGS
331  std::string prefix_mmm = std::string("TpetraExt ") + label + std::string(": ");
332  using Teuchos::TimeMonitor;
333  Teuchos::RCP<TimeMonitor> MM;
334 #endif
335 
336  // Node-specific code
337  using Teuchos::RCP;
338 
339  // Options
340  int team_work_size = 16; // Defaults to 16 as per Deveci 12/7/16 - csiefer
341  std::string myalg("LTG");
342  if(!params.is_null()) {
343  if(params->isParameter("openmp: jacobi algorithm"))
344  myalg = params->get("openmp: jacobi algorithm",myalg);
345  if(params->isParameter("openmp: team work size"))
346  team_work_size = params->get("openmp: team work size",team_work_size);
347  }
348 
349  if(myalg == "LTG") {
350  // Use the LTG kernel if requested
351  ::Tpetra::MatrixMatrix::ExtraKernels::jacobi_A_B_newmatrix_LowThreadGustavsonKernel(omega,Dinv,Aview,Bview,Acol2Brow,Acol2Irow,Bcol2Ccol,Icol2Ccol,C,Cimport,label,params);
352  }
353  else if(myalg == "MSAK") {
354  ::Tpetra::MatrixMatrix::ExtraKernels::jacobi_A_B_newmatrix_MultiplyScaleAddKernel(omega,Dinv,Aview,Bview,Acol2Brow,Acol2Irow,Bcol2Ccol,Icol2Ccol,C,Cimport,label,params);
355  }
356  else {
357  throw std::runtime_error("Tpetra::MatrixMatrix::Jacobi newmatrix unknown kernel");
358  }
359 
360 #ifdef HAVE_TPETRA_MMM_TIMINGS
361  MM = rcp(new TimeMonitor (*TimeMonitor::getNewTimer(prefix_mmm + std::string("Jacobi Newmatrix OpenMPESFC"))));
362 #endif
363 
364  // Final Fillcomplete
365  RCP<Teuchos::ParameterList> labelList = rcp(new Teuchos::ParameterList);
366  labelList->set("Timer Label",label);
367  if(!params.is_null()) labelList->set("compute global constants",params->get("compute global constants",true));
368 
369  // NOTE: MSAK already fillCompletes, so we have to check here
370  if(!C.isFillComplete()) {
371  RCP<const Export<LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode> > dummyExport;
372  C.expertStaticFillComplete(Bview.origMatrix->getDomainMap(), Aview.origMatrix->getRangeMap(), Cimport,dummyExport,labelList);
373  }
374 
375 }
376 
377 
378 
379 /*********************************************************************************************************/
380 template<class Scalar,
381  class LocalOrdinal,
382  class GlobalOrdinal,
383  class LocalOrdinalViewType>
384 void KernelWrappers2<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode,LocalOrdinalViewType>::jacobi_A_B_reuse_kernel_wrapper(Scalar omega,
385  const Vector<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode> & Dinv,
386  CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& Aview,
387  CrsMatrixStruct<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& Bview,
388  const LocalOrdinalViewType & Acol2Brow,
389  const LocalOrdinalViewType & Acol2Irow,
390  const LocalOrdinalViewType & Bcol2Ccol,
391  const LocalOrdinalViewType & Icol2Ccol,
392  CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosOpenMPWrapperNode>& C,
393  Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosOpenMPWrapperNode> > Cimport,
394  const std::string& label,
395  const Teuchos::RCP<Teuchos::ParameterList>& params) {
396 
397 #ifdef HAVE_TPETRA_MMM_TIMINGS
398  std::string prefix_mmm = std::string("TpetraExt ") + label + std::string(": ");
399  using Teuchos::TimeMonitor;
400  Teuchos::RCP<TimeMonitor> MM;
401 #endif
402 
403  // Lots and lots of typedefs
404  using Teuchos::RCP;
405 
406  // Options
407  int team_work_size = 16; // Defaults to 16 as per Deveci 12/7/16 - csiefer
408  std::string myalg("LTG");
409  if(!params.is_null()) {
410  if(params->isParameter("openmp: jacobi algorithm"))
411  myalg = params->get("openmp: jacobi algorithm",myalg);
412  if(params->isParameter("openmp: team work size"))
413  team_work_size = params->get("openmp: team work size",team_work_size);
414  }
415 
416  if(myalg == "LTG") {
417  // Use the LTG kernel if requested
418  ::Tpetra::MatrixMatrix::ExtraKernels::jacobi_A_B_reuse_LowThreadGustavsonKernel(omega,Dinv,Aview,Bview,Acol2Brow,Acol2Irow,Bcol2Ccol,Icol2Ccol,C,Cimport,label,params);
419  }
420  else {
421  throw std::runtime_error("Tpetra::MatrixMatrix::Jacobi reuse unknown kernel");
422  }
423 
424 #ifdef HAVE_TPETRA_MMM_TIMINGS
425  MM = rcp(new TimeMonitor (*TimeMonitor::getNewTimer(prefix_mmm + std::string("Jacobi Reuse OpenMPESFC"))));
426 #endif
427  C.fillComplete(C.getDomainMap(), C.getRangeMap());
428 
429 }
430 
431 
432 
433 }//MMdetails
434 }//Tpetra
435 
436 #endif//OpenMP
437 
438 #endif
Tpetra::Classes::CrsMatrix::local_matrix_type
KokkosSparse::CrsMatrix< impl_scalar_type, LocalOrdinal, execution_space, void, typename local_graph_type::size_type > local_matrix_type
The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix on each MPI pro...
Definition: Tpetra_CrsMatrix_decl.hpp:483
Tpetra
Namespace Tpetra contains the class and methods constituting the Tpetra library.