Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MPAssembly/fenl_functors.hpp
Go to the documentation of this file.
1 //@HEADER
2 // ************************************************************************
3 //
4 // Kokkos v. 4.0
5 // Copyright (2022) National Technology & Engineering
6 // Solutions of Sandia, LLC (NTESS).
7 //
8 // Under the terms of Contract DE-NA0003525 with NTESS,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Part of Kokkos, under the Apache License v2.0 with LLVM Exceptions.
12 // See https://kokkos.org/LICENSE for license information.
13 // SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
14 //
15 //@HEADER
16 
17 #ifndef KOKKOS_EXAMPLE_FENLFUNCTORS_HPP
18 #define KOKKOS_EXAMPLE_FENLFUNCTORS_HPP
19 
20 #include <stdio.h>
21 
22 #include <iostream>
23 #include <fstream>
24 #include <iomanip>
25 #include <cstdlib>
26 #include <cmath>
27 #include <limits>
28 
29 #include <Kokkos_Pair.hpp>
30 #include <Kokkos_UnorderedMap.hpp>
31 
32 #include <Kokkos_Timer.hpp>
33 
34 #include <BoxElemFixture.hpp>
35 #include <HexElement.hpp>
36 
38 
39 //----------------------------------------------------------------------------
40 //----------------------------------------------------------------------------
41 
42 namespace Kokkos {
43 namespace Example {
44 namespace FENL {
45 
46 struct DeviceConfig {
47  struct Dim3 {
48  size_t x, y, z;
49  Dim3(const size_t x_, const size_t y_ = 1, const size_t z_ = 1) :
50  x(x_), y(y_), z(z_) {}
51  };
52 
53  Dim3 block_dim;
54  size_t num_blocks;
55  size_t num_threads_per_block;
56 
57  DeviceConfig(const size_t num_blocks_ = 0,
58  const size_t threads_per_block_x_ = 0,
59  const size_t threads_per_block_y_ = 0,
60  const size_t threads_per_block_z_ = 1) :
61  block_dim(threads_per_block_x_,threads_per_block_y_,threads_per_block_z_),
62  num_blocks(num_blocks_),
64  {}
65 };
66 
67 template< typename ValueType , class Space >
68 struct CrsMatrix {
69 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE // Don't remove this until Kokkos has removed the deprecated code path probably around September 2018
70  typedef Kokkos::StaticCrsGraph< unsigned , Space , void , unsigned > StaticCrsGraphType ;
71 #else
72  typedef Kokkos::StaticCrsGraph< unsigned , Space , void , void , unsigned > StaticCrsGraphType ;
73 #endif
74  typedef View< ValueType * , Space > values_type ;
75 
78 
79  CrsMatrix() : graph(), values() {}
80 
81  CrsMatrix( const StaticCrsGraphType & arg_graph )
82  : graph( arg_graph )
83  , values( "crs_matrix_values" , arg_graph.entries.extent(0) )
84  {}
85 };
86 
87 // Traits class for creating strided local views for embedded ensemble-based,
88 // specialized for ensemble UQ scalar type
89 template <typename ViewType, typename Enabled = void>
90 struct LocalViewTraits {
91  typedef ViewType view_type;
92  // typedef Kokkos::View<typename view_type::data_type,
93  // typename view_type::array_layout,
94  // typename view_type::execution_space,
95  // Kokkos::MemoryUnmanaged> local_view_type;
96  typedef const view_type& local_view_type;
98  static const bool use_team = false;
99  KOKKOS_INLINE_FUNCTION
101  const unsigned local_rank)
102  { return v; }
103 };
104 
105 #if defined( KOKKOS_ENABLE_CUDA )
106 
107 template <typename ViewType>
108 struct LocalViewTraits<
109  ViewType,
110  typename std::enable_if< std::is_same<typename ViewType::execution_space,
111  Kokkos::Cuda>::value &&
112  Kokkos::is_view_mp_vector<ViewType>::value
113  >::type > {
114  typedef ViewType view_type;
117  static const bool use_team = true;
118 
119  KOKKOS_INLINE_FUNCTION
121  const unsigned local_rank)
122  {
123  return Kokkos::partition<1>(v, local_rank);
124  }
125 };
126 
127 #endif /* #if defined( KOKKOS_ENABLE_CUDA ) */
128 
129 // Compute DeviceConfig struct's based on scalar type
130 template <typename ScalarType>
131 struct CreateDeviceConfigs {
132  static void eval( Kokkos::Example::FENL::DeviceConfig& dev_config_elem,
133  Kokkos::Example::FENL::DeviceConfig& dev_config_bc ) {
134  dev_config_elem = Kokkos::Example::FENL::DeviceConfig( 0 , 1 , 1 );
135  dev_config_bc = Kokkos::Example::FENL::DeviceConfig( 0 , 1 , 1 );
136  }
137 };
138 
139 // Compute DeviceConfig struct's based on scalar type
140 template <typename StorageType>
141 struct CreateDeviceConfigs< Sacado::MP::Vector<StorageType> > {
143  static void eval( Kokkos::Example::FENL::DeviceConfig& dev_config_elem,
144  Kokkos::Example::FENL::DeviceConfig& dev_config_bc ) {
145  static const unsigned VectorSize = StorageType::static_size;
146 #if defined( KOKKOS_ENABLE_CUDA )
147  enum { is_cuda = std::is_same< execution_space, Kokkos::Cuda >::value };
148 #else
149  enum { is_cuda = false };
150 #endif /* #if defined( KOKKOS_ENABLE_CUDA ) */
151  if ( is_cuda ) {
152  dev_config_elem = Kokkos::Example::FENL::DeviceConfig( 0 , VectorSize , 64/VectorSize );
153  dev_config_bc = Kokkos::Example::FENL::DeviceConfig( 0 , VectorSize , 256/VectorSize );
154  }
155  else {
156  dev_config_elem = Kokkos::Example::FENL::DeviceConfig( 0 , 1 , 1 );
157  dev_config_bc = Kokkos::Example::FENL::DeviceConfig( 0 , 1 , 1 );
158  }
159  }
160 };
161 
162 template< class ElemNodeIdView , class CrsGraphType , unsigned ElemNode >
163 class NodeNodeGraph {
164 public:
165 
167  typedef pair<unsigned,unsigned> key_type ;
168 
169  typedef Kokkos::UnorderedMap< key_type, void , execution_space > SetType ;
170  typedef typename CrsGraphType::row_map_type::non_const_type RowMapType ;
171  typedef Kokkos::View< unsigned , execution_space > UnsignedValue ;
172 
173  // Static dimensions of 0 generate compiler warnings or errors.
174  typedef Kokkos::View< unsigned*[ElemNode][ElemNode] , execution_space >
176 
177 private:
178 
184 
185  const unsigned node_count ;
186  const ElemNodeIdView elem_node_id ;
191  PhaseType phase ;
192 
193 public:
194 
195  CrsGraphType graph ;
197 
198  struct Times
199  {
200  double ratio;
201  double fill_node_set;
202  double scan_node_count;
203  double fill_graph_entries;
204  double sort_graph_entries;
205  double fill_element_graph;
206  };
207 
208  NodeNodeGraph( const ElemNodeIdView & arg_elem_node_id ,
209  const unsigned arg_node_count,
210  Times & results
211  )
212  : node_count(arg_node_count)
213  , elem_node_id( arg_elem_node_id )
214  , row_total( "row_total" )
215  , row_count(Kokkos::ViewAllocateWithoutInitializing("row_count") , node_count ) // will deep_copy to 0 inside loop
216  , row_map( "graph_row_map" , node_count + 1 )
217  , node_node_set()
218  , phase( FILL_NODE_SET )
219  , graph()
220  , elem_graph()
221  {
222  //--------------------------------
223  // Guess at capacity required for the map:
224 
225  Kokkos::Timer wall_clock ;
226 
227  wall_clock.reset();
228  phase = FILL_NODE_SET ;
229 
230  // upper bound on the capacity
231  size_t set_capacity = (((28ull * node_count) / 2ull)*4ull)/3ull;
232 
233 
234  // Increase capacity until the (node,node) map is successfully filled.
235  {
236  // Zero the row count to restart the fill
238 
239  node_node_set = SetType( set_capacity );
240 
241  // May be larger that requested:
242  set_capacity = node_node_set.capacity();
243 
244  Kokkos::parallel_for( elem_node_id.extent(0) , *this );
245  }
246 
247  execution_space().fence();
248  results.ratio = (double)node_node_set.size() / (double)node_node_set.capacity();
249  results.fill_node_set = wall_clock.seconds();
250  //--------------------------------
251 
252  wall_clock.reset();
254 
255  // Exclusive scan of row_count into row_map
256  // including the final total in the 'node_count + 1' position.
257  // Zero the 'row_count' values.
258  Kokkos::parallel_scan( node_count , *this );
259 
260  // Zero the row count for the fill:
262 
263  unsigned graph_entry_count = 0 ;
264 
265  Kokkos::deep_copy( graph_entry_count , row_total );
266 
267  // Assign graph's row_map and allocate graph's entries
268  graph.row_map = row_map ;
269  graph.entries = typename CrsGraphType::entries_type( "graph_entries" , graph_entry_count );
270 
271  //--------------------------------
272  // Fill graph's entries from the (node,node) set.
273 
274  execution_space().fence();
275  results.scan_node_count = wall_clock.seconds();
276 
277  wall_clock.reset();
279  Kokkos::parallel_for( node_node_set.capacity() , *this );
280 
281  execution_space().fence();
282  results.fill_graph_entries = wall_clock.seconds();
283 
284  //--------------------------------
285  // Done with the temporary sets and arrays
286  wall_clock.reset();
288 
290  row_count = RowMapType();
291  row_map = RowMapType();
292  node_node_set.clear();
293 
294  //--------------------------------
295 
296  Kokkos::parallel_for( node_count , *this );
297 
298  execution_space().fence();
299  results.sort_graph_entries = wall_clock.seconds();
300 
301  //--------------------------------
302  // Element-to-graph mapping:
303  wall_clock.reset();
305  elem_graph = ElemGraphType("elem_graph", elem_node_id.extent(0) );
306  Kokkos::parallel_for( elem_node_id.extent(0) , *this );
307 
308  execution_space().fence();
309  results.fill_element_graph = wall_clock.seconds();
310  }
311 
312  //------------------------------------
313  // parallel_for: create map and count row length
314 
315  KOKKOS_INLINE_FUNCTION
316  void fill_set( const unsigned ielem ) const
317  {
318  // Loop over element's (row_local_node,col_local_node) pairs:
319  for ( unsigned row_local_node = 0 ; row_local_node < elem_node_id.extent(1) ; ++row_local_node ) {
320 
321  const unsigned row_node = elem_node_id( ielem , row_local_node );
322 
323  for ( unsigned col_local_node = row_local_node ; col_local_node < elem_node_id.extent(1) ; ++col_local_node ) {
324 
325  const unsigned col_node = elem_node_id( ielem , col_local_node );
326 
327  // If either node is locally owned then insert the pair into the unordered map:
328 
329  if ( row_node < row_count.extent(0) || col_node < row_count.extent(0) ) {
330 
331  const key_type key = (row_node < col_node) ? make_pair( row_node, col_node ) : make_pair( col_node, row_node ) ;
332 
333  const typename SetType::insert_result result = node_node_set.insert( key );
334 
335  if ( result.success() ) {
336  if ( row_node < row_count.extent(0) ) { atomic_fetch_add( & row_count( row_node ) , (typename RowMapType::value_type)1 ); }
337  if ( col_node < row_count.extent(0) && col_node != row_node ) { atomic_fetch_add( & row_count( col_node ) , (typename RowMapType::value_type)1 ); }
338  }
339  }
340  }
341  }
342  }
343 
344  KOKKOS_INLINE_FUNCTION
345  void fill_graph_entries( const unsigned iset ) const
346  {
347  if ( node_node_set.valid_at(iset) ) {
348  const key_type key = node_node_set.key_at(iset) ;
349  const unsigned row_node = key.first ;
350  const unsigned col_node = key.second ;
351 
352  if ( row_node < row_count.extent(0) ) {
353  const unsigned offset = graph.row_map( row_node ) + atomic_fetch_add( & row_count( row_node ) , (typename RowMapType::value_type)1 );
354  graph.entries( offset ) = col_node ;
355  }
356 
357  if ( col_node < row_count.extent(0) && col_node != row_node ) {
358  const unsigned offset = graph.row_map( col_node ) + atomic_fetch_add( & row_count( col_node ) , (typename RowMapType::value_type)1 );
359  graph.entries( offset ) = row_node ;
360  }
361  }
362  }
363 
364  KOKKOS_INLINE_FUNCTION
365  void sort_graph_entries( const unsigned irow ) const
366  {
367  typedef typename CrsGraphType::size_type size_type;
368  const size_type row_beg = graph.row_map( irow );
369  const size_type row_end = graph.row_map( irow + 1 );
370  for ( size_type i = row_beg + 1 ; i < row_end ; ++i ) {
371  const typename CrsGraphType::data_type col = graph.entries(i);
372  size_type j = i ;
373  for ( ; row_beg < j && col < graph.entries(j-1) ; --j ) {
374  graph.entries(j) = graph.entries(j-1);
375  }
376  graph.entries(j) = col ;
377  }
378  }
379 
380  KOKKOS_INLINE_FUNCTION
381  void fill_elem_graph_map( const unsigned ielem ) const
382  {
383  typedef typename CrsGraphType::data_type entry_type;
384  for ( unsigned row_local_node = 0 ; row_local_node < elem_node_id.extent(1) ; ++row_local_node ) {
385 
386  const unsigned row_node = elem_node_id( ielem , row_local_node );
387 
388  for ( unsigned col_local_node = 0 ; col_local_node < elem_node_id.extent(1) ; ++col_local_node ) {
389 
390  const unsigned col_node = elem_node_id( ielem , col_local_node );
391 
392  entry_type entry = 0 ;
393 
394  if ( row_node + 1 < graph.row_map.extent(0) ) {
395 
396  const entry_type entry_end = static_cast<entry_type> (graph.row_map( row_node + 1 ));
397 
398  entry = graph.row_map( row_node );
399 
400  for ( ; entry < entry_end && graph.entries(entry) != static_cast<entry_type> (col_node) ; ++entry );
401 
402  if ( entry == entry_end ) entry = ~0u ;
403  }
404 
405  elem_graph( ielem , row_local_node , col_local_node ) = entry ;
406  }
407  }
408  }
409 
410  KOKKOS_INLINE_FUNCTION
411  void operator()( const unsigned iwork ) const
412  {
413  if ( phase == FILL_NODE_SET ) {
414  fill_set( iwork );
415  }
416  else if ( phase == FILL_GRAPH_ENTRIES ) {
417  fill_graph_entries( iwork );
418  }
419  else if ( phase == SORT_GRAPH_ENTRIES ) {
420  sort_graph_entries( iwork );
421  }
422  else if ( phase == FILL_ELEMENT_GRAPH ) {
423  fill_elem_graph_map( iwork );
424  }
425  }
426 
427  //------------------------------------
428  // parallel_scan: row offsets
429 
430  typedef unsigned value_type ;
431 
432  KOKKOS_INLINE_FUNCTION
433  void operator()( const unsigned irow , unsigned & update , const bool final ) const
434  {
435  // exclusive scan
436  if ( final ) { row_map( irow ) = update ; }
437 
438  update += row_count( irow );
439 
440  if ( final ) {
441  if ( irow + 1 == row_count.extent(0) ) {
442  row_map( irow + 1 ) = update ;
443  row_total() = update ;
444  }
445  }
446  }
447 
448  KOKKOS_INLINE_FUNCTION
449  void init( unsigned & update ) const { update = 0 ; }
450 
451  KOKKOS_INLINE_FUNCTION
452  void join( unsigned & update , const unsigned & input ) const { update += input ; }
453 
454  //------------------------------------
455 };
456 
457 } /* namespace FENL */
458 } /* namespace Example */
459 } /* namespace Kokkos */
460 
461 //----------------------------------------------------------------------------
462 //----------------------------------------------------------------------------
463 
464 namespace Kokkos {
465 namespace Example {
466 namespace FENL {
467 
468 struct ElementComputationConstantCoefficient {
469  enum { is_constant = false };
470 
471  const double coeff_k ;
472 
473  KOKKOS_INLINE_FUNCTION
474  double operator()( double pt[], unsigned ensemble_rank) const
475  { return coeff_k * std::sin(pt[0]) * std::sin(pt[1]) * std::sin(pt[2]); }
476 
478  : coeff_k( val ) {}
479 
481  : coeff_k( rhs.coeff_k ) {}
482 };
483 
484 // Exponential KL from Stokhos
485 template < typename Scalar, typename MeshScalar, typename Device >
486 class ExponentialKLCoefficient {
487 public:
488 
489  // Turn into a meta-function class usable with Sacado::mpl::apply
490  template <typename T1, typename T2 = MeshScalar, typename T3 = Device>
491  struct apply {
493  };
494 
495  enum { is_constant = false };
496  typedef Kokkos::View<Scalar*, Kokkos::LayoutLeft, Device> RandomVariableView;
497  typedef typename RandomVariableView::size_type size_type;
498 
503 
504  rf_type m_rf; // Exponential random field
505  const MeshScalar m_mean; // Mean of random field
506  const MeshScalar m_variance; // Variance of random field
507  const MeshScalar m_corr_len; // Correlation length of random field
508  const size_type m_num_rv; // Number of random variables
509  RandomVariableView m_rv; // KL random variables
510 
511 public:
512 
514  const MeshScalar mean ,
515  const MeshScalar variance ,
516  const MeshScalar correlation_length ,
517  const size_type num_rv ) :
518  m_mean( mean ),
519  m_variance( variance ),
520  m_corr_len( correlation_length ),
521  m_num_rv( num_rv ),
522  m_rv( "KL Random Variables", m_num_rv )
523  {
524  Teuchos::ParameterList solverParams;
525  solverParams.set("Number of KL Terms", int(num_rv));
526  solverParams.set("Mean", mean);
527  solverParams.set("Standard Deviation", std::sqrt(variance));
528  int ndim = 3;
529  Teuchos::Array<double> domain_upper(ndim, 1.0), domain_lower(ndim, 0.0),
530  correlation_lengths(ndim, correlation_length);
531  solverParams.set("Domain Upper Bounds", domain_upper);
532  solverParams.set("Domain Lower Bounds", domain_lower);
533  solverParams.set("Correlation Lengths", correlation_lengths);
534 
535  m_rf = rf_type(solverParams);
536  }
537 
539  m_rf( rhs.m_rf ) ,
540  m_mean( rhs.m_mean ) ,
541  m_variance( rhs.m_variance ) ,
542  m_corr_len( rhs.m_corr_len ) ,
543  m_num_rv( rhs.m_num_rv ) ,
544  m_rv( rhs.m_rv ) {}
545 
546  KOKKOS_INLINE_FUNCTION
547  void setRandomVariables( const RandomVariableView& rv) { m_rv = rv; }
548 
549  KOKKOS_INLINE_FUNCTION
551 
552  KOKKOS_INLINE_FUNCTION
553  local_scalar_type operator() ( const MeshScalar point[],
554  const size_type ensemble_rank ) const
555  {
556  local_rv_view_type local_rv =
558 
559  local_scalar_type val = m_rf.evaluate(point, local_rv);
560 
561  return val;
562  }
563 };
564 
565 template< class FiniteElementMeshType , class SparseMatrixType
566  , class CoeffFunctionType = ElementComputationConstantCoefficient
567  >
568 class ElementComputation ;
569 
570 
571 template< class ExecutionSpace , BoxElemPart::ElemOrder Order , class CoordinateMap ,
572  typename ScalarType , class CoeffFunctionType >
574  < Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap >
575  , Kokkos::Example::FENL::CrsMatrix< ScalarType , ExecutionSpace >
576  , CoeffFunctionType >
577 {
578 public:
579 
582 
583  //------------------------------------
584 
585  typedef ExecutionSpace execution_space ;
586  typedef ScalarType scalar_type ;
587 
591  typedef Kokkos::View< scalar_type* , Kokkos::LayoutLeft, execution_space > vector_type ;
592 
593  //------------------------------------
594 
600  static const bool use_team = local_vector_view_traits::use_team;
601 
602  static const unsigned SpatialDim = element_data_type::spatial_dimension ;
603  static const unsigned TensorDim = SpatialDim * SpatialDim ;
604  static const unsigned ElemNodeCount = element_data_type::element_node_count ;
605  static const unsigned FunctionCount = element_data_type::function_count ;
606  static const unsigned IntegrationCount = element_data_type::integration_count ;
607 
608  //------------------------------------
609 
612  typedef Kokkos::View< scalar_type*[FunctionCount][FunctionCount] , execution_space > elem_matrices_type ;
613  typedef Kokkos::View< scalar_type*[FunctionCount] , execution_space > elem_vectors_type ;
614 
619 
621 
622  //------------------------------------
623 
624 
625  //------------------------------------
626  // Computational data:
627 
637  const CoeffFunctionType coeff_function ;
639 
641  : elem_data()
642  , elem_node_ids( rhs.elem_node_ids )
643  , node_coords( rhs.node_coords )
644  , elem_graph( rhs.elem_graph )
645  , elem_jacobians( rhs.elem_jacobians )
646  , elem_residuals( rhs.elem_residuals )
647  , solution( rhs.solution )
648  , residual( rhs.residual )
649  , jacobian( rhs.jacobian )
650  , coeff_function( rhs.coeff_function )
651  , dev_config( rhs.dev_config )
652  {}
653 
654  // If the element->sparse_matrix graph is provided then perform atomic updates
655  // Otherwise fill per-element contributions for subequent gather-add into a residual and jacobian.
656  ElementComputation( const mesh_type & arg_mesh ,
657  const CoeffFunctionType & arg_coeff_function ,
658  const vector_type & arg_solution ,
659  const elem_graph_type & arg_elem_graph ,
660  const sparse_matrix_type & arg_jacobian ,
661  const vector_type & arg_residual ,
662  const Kokkos::Example::FENL::DeviceConfig arg_dev_config )
663  : elem_data()
664  , elem_node_ids( arg_mesh.elem_node() )
665  , node_coords( arg_mesh.node_coord() )
666  , elem_graph( arg_elem_graph )
667  , elem_jacobians()
668  , elem_residuals()
669  , solution( arg_solution )
670  , residual( arg_residual )
671  , jacobian( arg_jacobian )
672  , coeff_function( arg_coeff_function )
673  , dev_config( arg_dev_config )
674  {}
675 
676  //------------------------------------
677 
678  void apply() const
679  {
680  const size_t nelem = elem_node_ids.extent(0);
681  if ( use_team ) {
682  const size_t team_size = dev_config.block_dim.x * dev_config.block_dim.y;
683  const size_t league_size =
684  (nelem + dev_config.block_dim.y-1) / dev_config.block_dim.y;
685  Kokkos::TeamPolicy< execution_space > config( league_size, team_size );
686  parallel_for( config , *this );
687  }
688  else {
689  parallel_for( nelem , *this );
690  }
691  }
692 
693  //------------------------------------
694 
695  static const unsigned FLOPS_transform_gradients =
696  /* Jacobian */ FunctionCount * TensorDim * 2 +
697  /* Inverse jacobian */ TensorDim * 6 + 6 +
698  /* Gradient transform */ FunctionCount * 15 ;
699 
700  KOKKOS_INLINE_FUNCTION
702  const double grad[][ FunctionCount ] , // Gradient of bases master element
703  const double x[] ,
704  const double y[] ,
705  const double z[] ,
706  double dpsidx[] ,
707  double dpsidy[] ,
708  double dpsidz[] ) const
709  {
710  enum { j11 = 0 , j12 = 1 , j13 = 2 ,
711  j21 = 3 , j22 = 4 , j23 = 5 ,
712  j31 = 6 , j32 = 7 , j33 = 8 };
713 
714  // Jacobian accumulation:
715 
716  double J[ TensorDim ] = { 0, 0, 0, 0, 0, 0, 0, 0, 0 };
717 
718  for( unsigned i = 0; i < FunctionCount ; ++i ) {
719  const double x1 = x[i] ;
720  const double x2 = y[i] ;
721  const double x3 = z[i] ;
722 
723  const double g1 = grad[0][i] ;
724  const double g2 = grad[1][i] ;
725  const double g3 = grad[2][i] ;
726 
727  J[j11] += g1 * x1 ;
728  J[j12] += g1 * x2 ;
729  J[j13] += g1 * x3 ;
730 
731  J[j21] += g2 * x1 ;
732  J[j22] += g2 * x2 ;
733  J[j23] += g2 * x3 ;
734 
735  J[j31] += g3 * x1 ;
736  J[j32] += g3 * x2 ;
737  J[j33] += g3 * x3 ;
738  }
739 
740  // Inverse jacobian:
741 
742  double invJ[ TensorDim ] = {
743  static_cast<double>( J[j22] * J[j33] - J[j23] * J[j32] ) ,
744  static_cast<double>( J[j13] * J[j32] - J[j12] * J[j33] ) ,
745  static_cast<double>( J[j12] * J[j23] - J[j13] * J[j22] ) ,
746 
747  static_cast<double>( J[j23] * J[j31] - J[j21] * J[j33] ) ,
748  static_cast<double>( J[j11] * J[j33] - J[j13] * J[j31] ) ,
749  static_cast<double>( J[j13] * J[j21] - J[j11] * J[j23] ) ,
750 
751  static_cast<double>( J[j21] * J[j32] - J[j22] * J[j31] ) ,
752  static_cast<double>( J[j12] * J[j31] - J[j11] * J[j32] ) ,
753  static_cast<double>( J[j11] * J[j22] - J[j12] * J[j21] ) };
754 
755  const double detJ = J[j11] * invJ[j11] +
756  J[j21] * invJ[j12] +
757  J[j31] * invJ[j13] ;
758 
759  const double detJinv = 1.0 / detJ ;
760 
761  for ( unsigned i = 0 ; i < TensorDim ; ++i ) { invJ[i] *= detJinv ; }
762 
763  // Transform gradients:
764 
765  for( unsigned i = 0; i < FunctionCount ; ++i ) {
766  const double g0 = grad[0][i];
767  const double g1 = grad[1][i];
768  const double g2 = grad[2][i];
769 
770  dpsidx[i] = g0 * invJ[j11] + g1 * invJ[j12] + g2 * invJ[j13];
771  dpsidy[i] = g0 * invJ[j21] + g1 * invJ[j22] + g2 * invJ[j23];
772  dpsidz[i] = g0 * invJ[j31] + g1 * invJ[j32] + g2 * invJ[j33];
773  }
774 
775  return detJ ;
776  }
777 
778  KOKKOS_INLINE_FUNCTION
780  const local_scalar_type dof_values[] ,
781  const double dpsidx[] ,
782  const double dpsidy[] ,
783  const double dpsidz[] ,
784  const double detJ ,
785  const local_scalar_type coeff_k ,
786  const double integ_weight ,
787  const double bases_vals[] ,
788  local_scalar_type elem_res[] ,
789  local_scalar_type elem_mat[][ FunctionCount ] ) const
790  {
791  local_scalar_type value_at_pt = 0 ;
792  local_scalar_type gradx_at_pt = 0 ;
793  local_scalar_type grady_at_pt = 0 ;
794  local_scalar_type gradz_at_pt = 0 ;
795 
796  for ( unsigned m = 0 ; m < FunctionCount ; m++ ) {
797  value_at_pt += dof_values[m] * bases_vals[m] ;
798  gradx_at_pt += dof_values[m] * dpsidx[m] ;
799  grady_at_pt += dof_values[m] * dpsidy[m] ;
800  gradz_at_pt += dof_values[m] * dpsidz[m] ;
801  }
802 
803  const local_scalar_type k_detJ_weight = coeff_k * detJ * integ_weight ;
804  const local_scalar_type res_val = value_at_pt * value_at_pt * detJ * integ_weight ;
805  const local_scalar_type mat_val = 2.0 * value_at_pt * detJ * integ_weight ;
806 
807  // $$ R_i = \int_{\Omega} \nabla \phi_i \cdot (k \nabla T) + \phi_i T^2 d \Omega $$
808  // $$ J_{i,j} = \frac{\partial R_i}{\partial T_j} = \int_{\Omega} k \nabla \phi_i \cdot \nabla \phi_j + 2 \phi_i \phi_j T d \Omega $$
809 
810  for ( unsigned m = 0; m < FunctionCount; ++m) {
811  local_scalar_type * const mat = elem_mat[m] ;
812  const double bases_val_m = bases_vals[m];
813  const double dpsidx_m = dpsidx[m] ;
814  const double dpsidy_m = dpsidy[m] ;
815  const double dpsidz_m = dpsidz[m] ;
816 
817  elem_res[m] += k_detJ_weight * ( dpsidx_m * gradx_at_pt +
818  dpsidy_m * grady_at_pt +
819  dpsidz_m * gradz_at_pt ) +
820  res_val * bases_val_m ;
821 
822  for( unsigned n = 0; n < FunctionCount; n++) {
823 
824  mat[n] += k_detJ_weight * ( dpsidx_m * dpsidx[n] +
825  dpsidy_m * dpsidy[n] +
826  dpsidz_m * dpsidz[n] ) +
827  mat_val * bases_val_m * bases_vals[n];
828  }
829  }
830  }
831 
832  KOKKOS_INLINE_FUNCTION
833  void operator()( const typename TeamPolicy< execution_space >::member_type & dev ) const
834  {
835 
836  const unsigned num_ensemble_threads = dev_config.block_dim.x ;
837  const unsigned num_element_threads = dev_config.block_dim.y ;
838  const unsigned element_rank = dev.team_rank() / num_ensemble_threads ;
839  const unsigned ensemble_rank = dev.team_rank() % num_ensemble_threads ;
840 
841  const unsigned ielem =
842  dev.league_rank() * num_element_threads + element_rank;
843 
844  if (ielem >= elem_node_ids.extent(0))
845  return;
846 
847  (*this)( ielem, ensemble_rank );
848 
849  }
850 
851  KOKKOS_INLINE_FUNCTION
852  void operator()( const unsigned ielem ,
853  const unsigned ensemble_rank = 0 ) const
854  {
855  local_vector_type local_solution =
856  local_vector_view_traits::create_local_view(solution,
857  ensemble_rank);
858  local_vector_type local_residual =
859  local_vector_view_traits::create_local_view(residual,
860  ensemble_rank);
861  local_matrix_type local_jacobian_values =
862  local_matrix_view_traits::create_local_view(jacobian.values,
863  ensemble_rank);
864 
865  // Gather nodal coordinates and solution vector:
866 
867  double x[ FunctionCount ] ;
868  double y[ FunctionCount ] ;
869  double z[ FunctionCount ] ;
870  local_scalar_type val[ FunctionCount ] ;
871  unsigned node_index[ ElemNodeCount ];
872 
873  for ( unsigned i = 0 ; i < ElemNodeCount ; ++i ) {
874  const unsigned ni = elem_node_ids( ielem , i );
875 
876  node_index[i] = ni ;
877 
878  x[i] = node_coords( ni , 0 );
879  y[i] = node_coords( ni , 1 );
880  z[i] = node_coords( ni , 2 );
881 
882  val[i] = local_solution( ni );
883  }
884 
885 
886  local_scalar_type elem_vec[ FunctionCount ] ;
887  local_scalar_type elem_mat[ FunctionCount ][ FunctionCount ] ;
888 
889  for( unsigned i = 0; i < FunctionCount ; i++ ) {
890  elem_vec[i] = 0 ;
891  for( unsigned j = 0; j < FunctionCount ; j++){
892  elem_mat[i][j] = 0 ;
893  }
894  }
895 
896 
897  for ( unsigned i = 0 ; i < IntegrationCount ; ++i ) {
898  double dpsidx[ FunctionCount ] ;
899  double dpsidy[ FunctionCount ] ;
900  double dpsidz[ FunctionCount ] ;
901 
902  local_scalar_type coeff_k = 0 ;
903 
904  {
905  double pt[] = {0.0, 0.0, 0.0};
906 
907  // If function is not constant
908  // then compute physical coordinates of integration point
909  if ( ! coeff_function.is_constant ) {
910  for ( unsigned j = 0 ; j < FunctionCount ; ++j ) {
911  pt[0] += x[j] * elem_data.values[i][j] ;
912  pt[1] += y[j] * elem_data.values[i][j] ;
913  pt[2] += z[j] * elem_data.values[i][j] ;
914  }
915  }
916 
917  // Need to fix this for local_scalar_type!!!!!!
918  coeff_k = coeff_function(pt, ensemble_rank);
919  }
920 
921  const double detJ =
922  transform_gradients( elem_data.gradients[i] , x , y , z ,
923  dpsidx , dpsidy , dpsidz );
924 
925  contributeResidualJacobian( val , dpsidx , dpsidy , dpsidz ,
926  detJ , coeff_k ,
927  elem_data.weights[i] ,
928  elem_data.values[i] ,
929  elem_vec , elem_mat );
930  }
931 
932  for( unsigned i = 0 ; i < FunctionCount ; i++ ) {
933  const unsigned row = node_index[i] ;
934  if ( row < residual.extent(0) ) {
935  atomic_add( & local_residual( row ) , elem_vec[i] );
936 
937  for( unsigned j = 0 ; j < FunctionCount ; j++ ) {
938  const unsigned entry = elem_graph( ielem , i , j );
939  if ( entry != ~0u ) {
940  atomic_add( & local_jacobian_values( entry ) , elem_mat[i][j] );
941  }
942  }
943  }
944  }
945  }
946 }; /* ElementComputation */
947 
948 //----------------------------------------------------------------------------
949 
950 template< class FixtureType , class SparseMatrixType >
951 class DirichletComputation ;
952 
953 template< class ExecutionSpace , BoxElemPart::ElemOrder Order , class CoordinateMap ,
954  typename ScalarType >
955 class DirichletComputation<
956  Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
957  Kokkos::Example::FENL::CrsMatrix< ScalarType , ExecutionSpace > >
958 {
959 public:
960 
964 
965  typedef ExecutionSpace execution_space ;
966  typedef ScalarType scalar_type ;
967 
971  typedef Kokkos::View< scalar_type* , execution_space > vector_type ;
972 
973  //------------------------------------
974 
979  static const bool use_team = local_vector_view_traits::use_team;
980 
981  typedef double bc_scalar_type ;
982 
983  //------------------------------------
984  // Computational data:
985 
986  const node_coord_type node_coords ;
987  const vector_type solution ;
988  const sparse_matrix_type jacobian ;
989  const vector_type residual ;
990  const bc_scalar_type bc_lower_value ;
991  const bc_scalar_type bc_upper_value ;
992  const scalar_coord_type bc_lower_limit ;
993  const scalar_coord_type bc_upper_limit ;
994  const unsigned bc_plane ;
995  const unsigned node_count ;
996  bool init ;
997  const Kokkos::Example::FENL::DeviceConfig dev_config ;
998 
999 
1000  DirichletComputation( const mesh_type & arg_mesh ,
1001  const vector_type & arg_solution ,
1002  const sparse_matrix_type & arg_jacobian ,
1003  const vector_type & arg_residual ,
1004  const unsigned arg_bc_plane ,
1005  const bc_scalar_type arg_bc_lower_value ,
1006  const bc_scalar_type arg_bc_upper_value ,
1007  const Kokkos::Example::FENL::DeviceConfig arg_dev_config )
1008  : node_coords( arg_mesh.node_coord() )
1009  , solution( arg_solution )
1010  , jacobian( arg_jacobian )
1011  , residual( arg_residual )
1012  , bc_lower_value( arg_bc_lower_value )
1013  , bc_upper_value( arg_bc_upper_value )
1014  , bc_lower_limit( std::numeric_limits<scalar_coord_type>::epsilon() )
1015  , bc_upper_limit( scalar_coord_type(1) - std::numeric_limits<scalar_coord_type>::epsilon() )
1016  , bc_plane( arg_bc_plane )
1017  , node_count( arg_mesh.node_count_owned() )
1018  , init( false )
1019  , dev_config( arg_dev_config )
1020  {
1021  parallel_for( node_count , *this );
1022  init = true ;
1023  }
1024 
1025  void apply() const
1026  {
1027  if ( use_team ) {
1028  const size_t team_size = dev_config.block_dim.x * dev_config.block_dim.y;
1029  const size_t league_size =
1030  (node_count + dev_config.block_dim.y-1) / dev_config.block_dim.y;
1031  Kokkos::TeamPolicy< execution_space > config( league_size, team_size );
1032  parallel_for( config , *this );
1033  }
1034  else
1035  parallel_for( node_count , *this );
1036  }
1037 
1038  //------------------------------------
1039 
1040  KOKKOS_INLINE_FUNCTION
1041  void operator()( const typename TeamPolicy< execution_space >::member_type & dev ) const
1042  {
1043 
1044  const unsigned num_ensemble_threads = dev_config.block_dim.x ;
1045  const unsigned num_node_threads = dev_config.block_dim.y ;
1046  const unsigned node_rank = dev.team_rank() / num_ensemble_threads ;
1047  const unsigned ensemble_rank = dev.team_rank() % num_ensemble_threads ;
1048 
1049  const unsigned inode = dev.league_rank() * num_node_threads + node_rank;
1050 
1051  if (inode >= node_count)
1052  return;
1053 
1054  (*this)( inode, ensemble_rank );
1055 
1056  }
1057 
1058  KOKKOS_INLINE_FUNCTION
1059  void operator()( const unsigned inode ,
1060  const unsigned ensemble_rank = 0) const
1061  {
1062  local_vector_type local_residual =
1063  local_vector_view_traits::create_local_view(residual,
1064  ensemble_rank);
1065  local_matrix_type local_jacobian_values =
1066  local_matrix_view_traits::create_local_view(jacobian.values,
1067  ensemble_rank);
1068 
1069  // Apply dirichlet boundary condition on the Solution and Residual vectors.
1070  // To maintain the symmetry of the original global stiffness matrix,
1071  // zero out the columns that correspond to boundary conditions, and
1072  // update the residual vector accordingly
1073 
1074  const unsigned iBeg = jacobian.graph.row_map[inode];
1075  const unsigned iEnd = jacobian.graph.row_map[inode+1];
1076 
1077  const scalar_coord_type c = node_coords(inode,bc_plane);
1078  const bool bc_lower = c <= bc_lower_limit ;
1079  const bool bc_upper = bc_upper_limit <= c ;
1080 
1081  if ( ! init ) {
1082  solution(inode) = bc_lower ? bc_lower_value : (
1083  bc_upper ? bc_upper_value : 0 );
1084  }
1085  else {
1086  if ( bc_lower || bc_upper ) {
1087 
1088  local_residual(inode) = 0 ;
1089 
1090  // zero each value on the row, and leave a one
1091  // on the diagonal
1092 
1093  for( unsigned i = iBeg ; i < iEnd ; ++i ) {
1094  local_jacobian_values(i) = int(inode) == int(jacobian.graph.entries(i)) ? 1 : 0 ;
1095  }
1096  }
1097  else {
1098 
1099  // Find any columns that are boundary conditions.
1100  // Clear them and adjust the residual vector
1101 
1102  for( unsigned i = iBeg ; i < iEnd ; ++i ) {
1103  const unsigned cnode = jacobian.graph.entries(i) ;
1104  const scalar_coord_type cc = node_coords(cnode,bc_plane);
1105 
1106  if ( ( cc <= bc_lower_limit ) || ( bc_upper_limit <= cc ) ) {
1107  local_jacobian_values(i) = 0 ;
1108  }
1109  }
1110  }
1111  }
1112  }
1113 };
1114 
1115 template< typename FixtureType , typename VectorType >
1117 {
1118 public:
1119 
1120  typedef FixtureType fixture_type ;
1121  typedef VectorType vector_type ;
1124 
1126  static const unsigned SpatialDim = element_data_type::spatial_dimension ;
1127  static const unsigned TensorDim = SpatialDim * SpatialDim ;
1128  static const unsigned ElemNodeCount = element_data_type::element_node_count ;
1129  static const unsigned IntegrationCount = element_data_type::integration_count ;
1130 
1131  //------------------------------------
1132  // Computational data:
1133 
1137 
1139  : elem_data()
1140  , fixture( rhs.fixture )
1141  , solution( rhs.solution )
1142  {}
1143 
1144  ResponseComputation( const fixture_type& arg_fixture ,
1145  const vector_type & arg_solution )
1146  : elem_data()
1147  , fixture( arg_fixture )
1148  , solution( arg_solution )
1149  {}
1150 
1151  //------------------------------------
1152 
1154  {
1155  value_type response = 0;
1156  //Kokkos::parallel_reduce( fixture.elem_count() , *this , response );
1157  Kokkos::parallel_reduce( solution.extent(0) , *this , response );
1158  return response;
1159  }
1160 
1161  //------------------------------------
1162 
1163  KOKKOS_INLINE_FUNCTION
1165  const double grad[][ ElemNodeCount ] , // Gradient of bases master element
1166  const double x[] ,
1167  const double y[] ,
1168  const double z[] ) const
1169  {
1170  enum { j11 = 0 , j12 = 1 , j13 = 2 ,
1171  j21 = 3 , j22 = 4 , j23 = 5 ,
1172  j31 = 6 , j32 = 7 , j33 = 8 };
1173 
1174  // Jacobian accumulation:
1175 
1176  double J[ TensorDim ] = { 0, 0, 0, 0, 0, 0, 0, 0, 0 };
1177 
1178  for( unsigned i = 0; i < ElemNodeCount ; ++i ) {
1179  const double x1 = x[i] ;
1180  const double x2 = y[i] ;
1181  const double x3 = z[i] ;
1182 
1183  const double g1 = grad[0][i] ;
1184  const double g2 = grad[1][i] ;
1185  const double g3 = grad[2][i] ;
1186 
1187  J[j11] += g1 * x1 ;
1188  J[j12] += g1 * x2 ;
1189  J[j13] += g1 * x3 ;
1190 
1191  J[j21] += g2 * x1 ;
1192  J[j22] += g2 * x2 ;
1193  J[j23] += g2 * x3 ;
1194 
1195  J[j31] += g3 * x1 ;
1196  J[j32] += g3 * x2 ;
1197  J[j33] += g3 * x3 ;
1198  }
1199 
1200  // Inverse jacobian:
1201 
1202  double invJ[ TensorDim ] = {
1203  static_cast<double>( J[j22] * J[j33] - J[j23] * J[j32] ) ,
1204  static_cast<double>( J[j13] * J[j32] - J[j12] * J[j33] ) ,
1205  static_cast<double>( J[j12] * J[j23] - J[j13] * J[j22] ) ,
1206 
1207  static_cast<double>( J[j23] * J[j31] - J[j21] * J[j33] ) ,
1208  static_cast<double>( J[j11] * J[j33] - J[j13] * J[j31] ) ,
1209  static_cast<double>( J[j13] * J[j21] - J[j11] * J[j23] ) ,
1210 
1211  static_cast<double>( J[j21] * J[j32] - J[j22] * J[j31] ) ,
1212  static_cast<double>( J[j12] * J[j31] - J[j11] * J[j32] ) ,
1213  static_cast<double>( J[j11] * J[j22] - J[j12] * J[j21] ) };
1214 
1215  const double detJ = J[j11] * invJ[j11] +
1216  J[j21] * invJ[j12] +
1217  J[j31] * invJ[j13] ;
1218 
1219  return detJ ;
1220  }
1221 
1222  KOKKOS_INLINE_FUNCTION
1224  const value_type dof_values[] ,
1225  const double detJ ,
1226  const double integ_weight ,
1227  const double bases_vals[] ) const
1228  {
1229  // $$ g_i = \int_{\Omega} T^2 d \Omega $$
1230 
1231  value_type value_at_pt = 0 ;
1232  for ( unsigned m = 0 ; m < ElemNodeCount ; m++ ) {
1233  value_at_pt += dof_values[m] * bases_vals[m] ;
1234  }
1235 
1236  value_type elem_response =
1237  value_at_pt * value_at_pt * detJ * integ_weight ;
1238 
1239  return elem_response;
1240  }
1241 
1242  /*
1243  KOKKOS_INLINE_FUNCTION
1244  void operator()( const unsigned ielem , value_type& response ) const
1245  {
1246  // Gather nodal coordinates and solution vector:
1247 
1248  double x[ ElemNodeCount ] ;
1249  double y[ ElemNodeCount ] ;
1250  double z[ ElemNodeCount ] ;
1251  value_type val[ ElemNodeCount ] ;
1252 
1253  for ( unsigned i = 0 ; i < ElemNodeCount ; ++i ) {
1254  const unsigned ni = fixture.elem_node( ielem , i );
1255 
1256  x[i] = fixture.node_coord( ni , 0 );
1257  y[i] = fixture.node_coord( ni , 1 );
1258  z[i] = fixture.node_coord( ni , 2 );
1259 
1260  val[i] = solution( ni );
1261  }
1262 
1263  for ( unsigned i = 0 ; i < IntegrationCount ; ++i ) {
1264 
1265  const double detJ = compute_detJ( elem_data.gradients[i] , x , y , z );
1266 
1267  response += contributeResponse( val , detJ , elem_data.weights[i] ,
1268  elem_data.values[i] );
1269  }
1270  }
1271  */
1272 
1273  KOKKOS_INLINE_FUNCTION
1274  void operator()( const unsigned i , value_type& response ) const
1275  {
1276  value_type u = solution(i);
1277  response += (u * u) / fixture.node_count_global();
1278  }
1279 
1280  KOKKOS_INLINE_FUNCTION
1281  void init( value_type & response ) const
1282  { response = 0 ; }
1283 
1284  KOKKOS_INLINE_FUNCTION
1285  void join( value_type & response ,
1286  const value_type & input ) const
1287  { response += input ; }
1288 
1289 }; /* ResponseComputation */
1290 
1291 } /* namespace FENL */
1292 } /* namespace Example */
1293 } /* namespace Kokkos */
1294 
1295 //----------------------------------------------------------------------------
1296 
1297 /* A Cuda-specific specialization for the element computation functor. */
1298 #if defined( __CUDACC__ )
1299 // #include <NonlinearElement_Cuda.hpp>
1300 #endif
1301 
1302 //----------------------------------------------------------------------------
1303 
1304 #endif /* #ifndef KOKKOS_EXAMPLE_FENLFUNCTORS_HPP */
static KOKKOS_INLINE_FUNCTION local_view_type create_local_view(const view_type &v, const unsigned local_rank)
KOKKOS_INLINE_FUNCTION PCE< Storage > sqrt(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION double compute_detJ(const double grad[][ElemNodeCount], const double x[], const double y[], const double z[]) const
KOKKOS_INLINE_FUNCTION void operator()(const unsigned i, value_type &response) const
ExponentialKLCoefficient(const ExponentialKLCoefficient &rhs)
KOKKOS_INLINE_FUNCTION void setRandomVariables(const RandomVariableView &rv)
LocalViewTraits< RandomVariableView > local_rv_view_traits
KOKKOS_INLINE_FUNCTION Teuchos::PromotionTraits< typename rv_type::value_type, value_type >::promote evaluate(const point_type &point, const rv_type &random_variables) const
Evaluate random field at a point.
KOKKOS_INLINE_FUNCTION void fill_graph_entries(const unsigned iset) const
KOKKOS_INLINE_FUNCTION void join(unsigned &update, const unsigned &input) const
static void eval(Kokkos::Example::FENL::DeviceConfig &dev_config_elem, Kokkos::Example::FENL::DeviceConfig &dev_config_bc)
KOKKOS_INLINE_FUNCTION value_type contributeResponse(const value_type dof_values[], const double detJ, const double integ_weight, const double bases_vals[]) const
Kokkos::View< const double *[SpaceDim], Device > node_coord_type
Kokkos::DefaultExecutionSpace execution_space
Kokkos::View< const unsigned *[ElemNode], Device > elem_node_type
KOKKOS_INLINE_FUNCTION RandomVariableView getRandomVariables() const
static void eval(Kokkos::Example::FENL::DeviceConfig &dev_config_elem, Kokkos::Example::FENL::DeviceConfig &dev_config_bc)
KOKKOS_INLINE_FUNCTION void atomic_add(volatile Sacado::UQ::PCE< Storage > *const dest, const Sacado::UQ::PCE< Storage > &src)
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Kokkos::View< Scalar *, Kokkos::LayoutLeft, Device > RandomVariableView
KOKKOS_INLINE_FUNCTION local_scalar_type operator()(const MeshScalar point[], const size_type ensemble_rank) const
local_rv_view_traits::local_view_type local_rv_view_type
ElemNodeIdView::execution_space execution_space
ResponseComputation(const fixture_type &arg_fixture, const vector_type &arg_solution)
KOKKOS_INLINE_FUNCTION void operator()(const unsigned irow, unsigned &update, const bool final) const
ElementComputationConstantCoefficient(const ElementComputationConstantCoefficient &rhs)
KOKKOS_INLINE_FUNCTION void init(unsigned &update) const
KOKKOS_INLINE_FUNCTION void fill_elem_graph_map(const unsigned ielem) const
Dim3(const size_t x_, const size_t y_=1, const size_t z_=1)
CrsGraphType::row_map_type::non_const_type RowMapType
KOKKOS_INLINE_FUNCTION void init(value_type &response) const
CrsMatrix(const StaticCrsGraphType &arg_graph)
Kokkos::StaticCrsGraph< unsigned, Space, void, void, unsigned > StaticCrsGraphType
void deep_copy(const Stokhos::CrsMatrix< ValueType, DstDevice, Layout > &dst, const Stokhos::CrsMatrix< ValueType, SrcDevice, Layout > &src)
KOKKOS_INLINE_FUNCTION void operator()(const unsigned iwork) const
KOKKOS_INLINE_FUNCTION void contributeResidualJacobian(const local_scalar_type dof_values[], const double dpsidx[], const double dpsidy[], const double dpsidz[], const double detJ, const local_scalar_type coeff_k, const double integ_weight, const double bases_vals[], local_scalar_type elem_res[], local_scalar_type elem_mat[][FunctionCount]) const
KOKKOS_INLINE_FUNCTION double operator()(double pt[], unsigned ensemble_rank) const
NodeNodeGraph(const ElemNodeIdView &arg_elem_node_id, const unsigned arg_node_count, Times &results)
expr val()
KOKKOS_INLINE_FUNCTION void join(value_type &response, const value_type &input) const
ElementComputation(const mesh_type &arg_mesh, const CoeffFunctionType &arg_coeff_function, const vector_type &arg_solution, const elem_graph_type &arg_elem_graph, const sparse_matrix_type &arg_jacobian, const vector_type &arg_residual, const Kokkos::Example::FENL::DeviceConfig arg_dev_config)
KOKKOS_INLINE_FUNCTION double transform_gradients(const double grad[][FunctionCount], const double x[], const double y[], const double z[], double dpsidx[], double dpsidy[], double dpsidz[]) const
Kokkos::View< unsigned, execution_space > UnsignedValue
DeviceConfig(const size_t num_blocks_=0, const size_t threads_per_block_x_=0, const size_t threads_per_block_y_=0, const size_t threads_per_block_z_=1)
Kokkos::Example::HexElement_Data< fixture_type::ElemNode > element_data_type
local_rv_view_traits::local_value_type local_scalar_type
DirichletComputation(const mesh_type &arg_mesh, const vector_type &arg_solution, const sparse_matrix_type &arg_jacobian, const vector_type &arg_residual, const unsigned arg_bc_plane, const bc_scalar_type arg_bc_lower_value, const bc_scalar_type arg_bc_upper_value, const Kokkos::Example::FENL::DeviceConfig arg_dev_config)
KOKKOS_INLINE_FUNCTION PCE< Storage > sin(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION void sort_graph_entries(const unsigned irow) const
KOKKOS_INLINE_FUNCTION void fill_set(const unsigned ielem) const
Kokkos::View< unsigned *[ElemNode][ElemNode], execution_space > ElemGraphType
int n
Generate a distributed unstructured finite element mesh from a partitioned NX*NY*NZ box of elements...
Stokhos::KL::ExponentialRandomField< MeshScalar, Device > rf_type
ExponentialKLCoefficient(const MeshScalar mean, const MeshScalar variance, const MeshScalar correlation_length, const size_type num_rv)
void update(const ValueType &alpha, VectorType &x, const ValueType &beta, const VectorType &y)
Kokkos::UnorderedMap< key_type, void, execution_space > SetType