44 #ifndef KOKKOS_EXAMPLE_FENLFUNCTORS_HPP
45 #define KOKKOS_EXAMPLE_FENLFUNCTORS_HPP
56 #include <Kokkos_Pair.hpp>
57 #include <Kokkos_UnorderedMap.hpp>
59 #include <impl/Kokkos_Timer.hpp>
80 Dim3(
const size_t x_,
const size_t y_ = 1,
const size_t z_ = 1) :
81 x(x_),
y(y_),
z(z_) {}
89 const size_t threads_per_block_x_ = 0,
90 const size_t threads_per_block_y_ = 0,
91 const size_t threads_per_block_z_ = 1) :
92 block_dim(threads_per_block_x_,threads_per_block_y_,threads_per_block_z_),
98 template<
typename ValueType ,
class Space >
100 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE // Don't remove this until Kokkos has removed the deprecated code path probably around September 2018
101 typedef Kokkos::StaticCrsGraph< unsigned , Space , void , unsigned >
StaticCrsGraphType ;
114 ,
values(
"crs_matrix_values" , arg_graph.entries.extent(0) )
120 template <
typename ViewType,
typename Enabled =
void>
130 KOKKOS_INLINE_FUNCTION
132 const unsigned local_rank)
136 #if defined( KOKKOS_ENABLE_CUDA )
138 template <
typename ViewType>
139 struct LocalViewTraits<
141 typename std::enable_if< std::is_same<typename ViewType::execution_space,
142 Kokkos::Cuda>::value &&
143 Kokkos::is_view_mp_vector<ViewType>::value
150 KOKKOS_INLINE_FUNCTION
152 const unsigned local_rank)
154 return Kokkos::partition<1>(v, local_rank);
161 template <
typename ScalarType>
171 template <
typename StorageType>
176 static const unsigned VectorSize = StorageType::static_size;
177 #if defined( KOKKOS_ENABLE_CUDA )
178 enum { is_cuda = Kokkos::Impl::is_same< execution_space, Kokkos::Cuda >::value };
180 enum { is_cuda =
false };
193 template<
class ElemNodeIdView ,
class CrsGraphType ,
unsigned ElemNode >
200 typedef Kokkos::UnorderedMap< key_type, void , execution_space >
SetType ;
201 typedef typename CrsGraphType::row_map_type::non_const_type
RowMapType ;
205 typedef Kokkos::View< unsigned*[ElemNode][ElemNode] , execution_space >
240 const unsigned arg_node_count,
256 Kokkos::Impl::Timer wall_clock ;
262 size_t set_capacity = (((28ull *
node_count) / 2ull)*4ull)/3ull;
278 execution_space::fence();
294 unsigned graph_entry_count = 0 ;
300 graph.entries =
typename CrsGraphType::entries_type(
"graph_entries" , graph_entry_count );
305 execution_space::fence();
312 execution_space::fence();
329 execution_space::fence();
339 execution_space::fence();
346 KOKKOS_INLINE_FUNCTION
350 for (
unsigned row_local_node = 0 ; row_local_node <
elem_node_id.extent(1) ; ++row_local_node ) {
352 const unsigned row_node =
elem_node_id( ielem , row_local_node );
354 for (
unsigned col_local_node = row_local_node ; col_local_node <
elem_node_id.extent(1) ; ++col_local_node ) {
356 const unsigned col_node =
elem_node_id( ielem , col_local_node );
362 const key_type key = (row_node < col_node) ? make_pair( row_node, col_node ) : make_pair( col_node, row_node ) ;
364 const typename SetType::insert_result result =
node_node_set.insert( key );
366 if ( result.success() ) {
375 KOKKOS_INLINE_FUNCTION
380 const unsigned row_node = key.first ;
381 const unsigned col_node = key.second ;
385 graph.entries( offset ) = col_node ;
388 if ( col_node <
row_count.extent(0) && col_node != row_node ) {
390 graph.entries( offset ) = row_node ;
395 KOKKOS_INLINE_FUNCTION
398 typedef typename CrsGraphType::size_type size_type;
399 const size_type row_beg =
graph.row_map( irow );
400 const size_type row_end =
graph.row_map( irow + 1 );
401 for ( size_type i = row_beg + 1 ; i < row_end ; ++i ) {
402 const typename CrsGraphType::data_type col =
graph.entries(i);
404 for ( ; row_beg < j && col <
graph.entries(j-1) ; --
j ) {
407 graph.entries(j) = col ;
411 KOKKOS_INLINE_FUNCTION
414 typedef typename CrsGraphType::data_type entry_type;
415 for (
unsigned row_local_node = 0 ; row_local_node <
elem_node_id.extent(1) ; ++row_local_node ) {
417 const unsigned row_node =
elem_node_id( ielem , row_local_node );
419 for (
unsigned col_local_node = 0 ; col_local_node <
elem_node_id.extent(1) ; ++col_local_node ) {
421 const unsigned col_node =
elem_node_id( ielem , col_local_node );
423 entry_type entry = 0 ;
425 if ( row_node + 1 <
graph.row_map.extent(0) ) {
427 const entry_type entry_end =
static_cast<entry_type
> (
graph.row_map( row_node + 1 ));
429 entry =
graph.row_map( row_node );
431 for ( ; entry < entry_end &&
graph.entries(entry) !=
static_cast<entry_type
> (col_node) ; ++entry );
433 if ( entry == entry_end ) entry = ~0u ;
436 elem_graph( ielem , row_local_node , col_local_node ) = entry ;
441 KOKKOS_INLINE_FUNCTION
463 KOKKOS_INLINE_FUNCTION
479 KOKKOS_INLINE_FUNCTION
482 KOKKOS_INLINE_FUNCTION
483 void join(
volatile unsigned &
update ,
const volatile unsigned & input )
const { update += input ; }
504 KOKKOS_INLINE_FUNCTION
505 double operator()(
double pt[],
unsigned ensemble_rank)
const
516 template <
typename Scalar,
typename MeshScalar,
typename Device >
521 template <
typename T1,
typename T2 = MeshScalar,
typename T3 = Device>
528 typedef typename RandomVariableView::size_type
size_type;
545 const MeshScalar mean ,
546 const MeshScalar variance ,
547 const MeshScalar correlation_length ,
556 solverParams.
set(
"Number of KL Terms",
int(num_rv));
557 solverParams.
set(
"Mean", mean);
558 solverParams.
set(
"Standard Deviation",
std::sqrt(variance));
561 correlation_lengths(ndim, correlation_length);
562 solverParams.
set(
"Domain Upper Bounds", domain_upper);
563 solverParams.
set(
"Domain Lower Bounds", domain_lower);
564 solverParams.
set(
"Correlation Lengths", correlation_lengths);
577 KOKKOS_INLINE_FUNCTION
580 KOKKOS_INLINE_FUNCTION
583 KOKKOS_INLINE_FUNCTION
597 class CoordinateMap ,
typename ScalarType >
612 typedef Kokkos::View< scalar_type* , Kokkos::LayoutLeft, execution_space >
vector_type ;
626 typedef Kokkos::View< scalar_type*[FunctionCount][FunctionCount] , execution_space >
elem_matrices_type ;
677 KOKKOS_INLINE_FUNCTION
685 double dpsidz[] )
const
687 enum { j11 = 0 , j12 = 1 , j13 = 2 ,
688 j21 = 3 , j22 = 4 , j23 = 5 ,
689 j31 = 6 , j32 = 7 , j33 = 8 };
693 double J[
TensorDim ] = { 0, 0, 0, 0, 0, 0, 0, 0, 0 };
696 const double x1 = x[i] ;
697 const double x2 = y[i] ;
698 const double x3 = z[i] ;
700 const double g1 = grad[0][i] ;
701 const double g2 = grad[1][i] ;
702 const double g3 = grad[2][i] ;
720 static_cast<double>( J[j22] * J[j33] - J[j23] * J[j32] ) ,
721 static_cast<double>( J[j13] * J[j32] - J[j12] * J[j33] ) ,
722 static_cast<double>( J[j12] * J[j23] - J[j13] * J[j22] ) ,
724 static_cast<double>( J[j23] * J[j31] - J[j21] * J[j33] ) ,
725 static_cast<double>( J[j11] * J[j33] - J[j13] * J[j31] ) ,
726 static_cast<double>( J[j13] * J[j21] - J[j11] * J[j23] ) ,
728 static_cast<double>( J[j21] * J[j32] - J[j22] * J[j31] ) ,
729 static_cast<double>( J[j12] * J[j31] - J[j11] * J[j32] ) ,
730 static_cast<double>( J[j11] * J[j22] - J[j12] * J[j21] ) };
732 const double detJ = J[j11] * invJ[j11] +
736 const double detJinv = 1.0 / detJ ;
738 for (
unsigned i = 0 ; i <
TensorDim ; ++i ) { invJ[i] *= detJinv ; }
743 const double g0 = grad[0][i];
744 const double g1 = grad[1][i];
745 const double g2 = grad[2][i];
747 dpsidx[i] = g0 * invJ[j11] + g1 * invJ[j12] + g2 * invJ[j13];
748 dpsidy[i] = g0 * invJ[j21] + g1 * invJ[j22] + g2 * invJ[j23];
749 dpsidz[i] = g0 * invJ[j31] + g1 * invJ[j32] + g2 * invJ[j33];
764 template<
class FiniteElementMeshType ,
765 class SparseMatrixType ,
767 class CoeffFunctionType = ElementComputationConstantCoefficient
772 class CoordinateMap ,
typename ScalarType ,
class CoeffFunctionType >
775 CrsMatrix< ScalarType , ExecutionSpace > ,
796 coeff_function( rhs.coeff_function ) ,
797 dev_config( rhs.dev_config ) {}
801 const CoeffFunctionType & arg_coeff_function ,
807 base_type(arg_mesh, arg_solution, arg_elem_graph, arg_jacobian,
809 coeff_function( arg_coeff_function ) ,
810 dev_config( arg_dev_config ) {}
817 parallel_for( nelem , *
this );
820 KOKKOS_INLINE_FUNCTION
823 unsigned node_index[],
824 double x[],
double y[],
double z[],
846 KOKKOS_INLINE_FUNCTION
848 const unsigned node_index[],
853 const unsigned row = node_index[i] ;
854 if ( row < this->
residual.extent(0) ) {
858 const unsigned entry = this->
elem_graph( ielem , i ,
j );
859 if ( entry != ~0u ) {
867 KOKKOS_INLINE_FUNCTION
877 double coeff_src = 1.234;
878 double advection[] = { 1.1, 1.2, 1.3 };
882 double pt[] = {0.0, 0.0, 0.0};
887 if ( ! coeff_function.is_constant ) {
894 coeff_k = coeff_function(pt, 0);
901 dpsidx , dpsidy , dpsidz );
902 const double detJ_weight = detJ * integ_weight;
903 const scalar_type detJ_weight_coeff_k = detJ_weight * coeff_k;
910 value_at_pt += dof_values[m] * bases_vals[m] ;
911 gradx_at_pt += dof_values[m] * dpsidx[m] ;
912 grady_at_pt += dof_values[m] * dpsidy[m] ;
913 gradz_at_pt += dof_values[m] * dpsidz[m] ;
917 coeff_src * value_at_pt * value_at_pt ;
919 2.0 * coeff_src * value_at_pt ;
926 advection_x*gradx_at_pt +
927 advection_y*grady_at_pt +
928 advection_z*gradz_at_pt ;
932 const double bases_val_m = bases_vals[m] * detJ_weight ;
933 const double dpsidx_m = dpsidx[m] ;
934 const double dpsidy_m = dpsidy[m] ;
935 const double dpsidz_m = dpsidz[m] ;
938 detJ_weight_coeff_k * ( dpsidx_m * gradx_at_pt +
939 dpsidy_m * grady_at_pt +
940 dpsidz_m * gradz_at_pt ) +
941 bases_val_m * ( advection_term + source_term ) ;
944 const double dpsidx_n = dpsidx[
n] ;
945 const double dpsidy_n = dpsidy[
n] ;
946 const double dpsidz_n = dpsidz[
n] ;
948 detJ_weight_coeff_k * ( dpsidx_m * dpsidx_n +
949 dpsidy_m * dpsidy_n +
950 dpsidz_m * dpsidz_n ) +
951 bases_val_m * ( advection_x * dpsidx_n +
952 advection_y * dpsidy_n +
953 advection_z * dpsidz_n +
954 source_deriv * bases_vals[
n] ) ;
960 KOKKOS_INLINE_FUNCTION
973 gatherSolution(ielem, val, node_index, x, y, z, elem_res, elem_mat);
976 computeElementResidualJacobian( val, x, y, z, elem_res , elem_mat );
979 scatterResidual( ielem, node_index, elem_res, elem_mat );
984 class CoordinateMap ,
typename ScalarType ,
class CoeffFunctionType >
987 CrsMatrix< ScalarType , ExecutionSpace > ,
1011 coeff_function( rhs.coeff_function ) ,
1012 dev_config( rhs.dev_config ) {}
1016 const CoeffFunctionType & arg_coeff_function ,
1022 base_type(arg_mesh, arg_solution, arg_elem_graph,
1023 arg_jacobian, arg_residual),
1024 coeff_function( arg_coeff_function ) ,
1025 dev_config( arg_dev_config ) {}
1032 parallel_for( nelem , *
this );
1035 KOKKOS_INLINE_FUNCTION
1038 unsigned node_index[],
1039 double x[],
double y[],
double z[],
1045 node_index[i] = ni ;
1051 val[i].val() = this->
solution( ni );
1056 KOKKOS_INLINE_FUNCTION
1058 const unsigned node_index[],
1062 const unsigned row = node_index[i] ;
1063 if ( row < this->
residual.extent(0) ) {
1067 const unsigned entry = this->
elem_graph( ielem , i ,
j );
1068 if ( entry != ~0u ) {
1070 res[i].fastAccessDx(
j) );
1077 template <
typename local_scalar_type>
1078 KOKKOS_INLINE_FUNCTION
1083 local_scalar_type elem_res[] )
const
1086 double coeff_src = 1.234;
1087 double advection[] = { 1.1, 1.2, 1.3 };
1091 double pt[] = {0.0, 0.0, 0.0};
1096 if ( ! coeff_function.is_constant ) {
1103 coeff_k = coeff_function(pt, 0);
1110 dpsidx , dpsidy , dpsidz );
1111 const double detJ_weight = detJ * integ_weight;
1112 const scalar_type detJ_weight_coeff_k = detJ_weight * coeff_k;
1114 local_scalar_type value_at_pt = 0 ;
1115 local_scalar_type gradx_at_pt = 0 ;
1116 local_scalar_type grady_at_pt = 0 ;
1117 local_scalar_type gradz_at_pt = 0 ;
1119 value_at_pt += dof_values[m] * bases_vals[m] ;
1120 gradx_at_pt += dof_values[m] * dpsidx[m] ;
1121 grady_at_pt += dof_values[m] * dpsidy[m] ;
1122 gradz_at_pt += dof_values[m] * dpsidz[m] ;
1125 const local_scalar_type source_term =
1126 coeff_src * value_at_pt * value_at_pt ;
1128 const local_scalar_type advection_term =
1129 advection[0]*gradx_at_pt +
1130 advection[1]*grady_at_pt +
1131 advection[2]*gradz_at_pt;
1134 const double bases_val_m = bases_vals[m] * detJ_weight ;
1135 const double dpsidx_m = dpsidx[m] ;
1136 const double dpsidy_m = dpsidy[m] ;
1137 const double dpsidz_m = dpsidz[m] ;
1140 detJ_weight_coeff_k * ( dpsidx_m * gradx_at_pt +
1141 dpsidy_m * grady_at_pt +
1142 dpsidz_m * gradz_at_pt ) +
1143 bases_val_m * ( advection_term + source_term ) ;
1148 KOKKOS_INLINE_FUNCTION
1160 gatherSolution( ielem, val, node_index, x, y, z, elem_res );
1163 computeElementResidual( val, x, y, z, elem_res );
1166 scatterResidual( ielem, node_index, elem_res );
1171 class CoordinateMap ,
typename ScalarType ,
class CoeffFunctionType >
1174 CrsMatrix< ScalarType , ExecutionSpace > ,
1176 public ElementComputation< Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
1177 CrsMatrix< ScalarType , ExecutionSpace > ,
1178 FadElement, CoeffFunctionType > {
1199 const CoeffFunctionType & arg_coeff_function ,
1205 base_type(arg_mesh, arg_coeff_function, arg_solution, arg_elem_graph,
1206 arg_jacobian, arg_residual, arg_dev_config) {}
1213 parallel_for( nelem , *
this );
1216 KOKKOS_INLINE_FUNCTION
1219 unsigned node_index[],
1220 double x[],
double y[],
double z[],
1226 node_index[i] = ni ;
1236 template <
typename local_scalar_type>
1237 KOKKOS_INLINE_FUNCTION
1242 local_scalar_type elem_res[] )
const
1245 double coeff_src = 1.234;
1246 double advection[] = { 1.1, 1.2, 1.3 };
1250 double pt[] = {0.0, 0.0, 0.0};
1255 if ( ! this->coeff_function.is_constant ) {
1262 coeff_k = this->coeff_function(pt, 0);
1269 dpsidx , dpsidy , dpsidz );
1270 const double detJ_weight = detJ * integ_weight;
1271 const scalar_type detJ_weight_coeff_k = detJ_weight * coeff_k;
1278 value_at_pt.val() += dof_values[m] * bases_vals[m] ;
1279 value_at_pt.fastAccessDx(m) = bases_vals[m] ;
1281 gradx_at_pt.val() += dof_values[m] * dpsidx[m] ;
1282 gradx_at_pt.fastAccessDx(m) = dpsidx[m] ;
1284 grady_at_pt.val() += dof_values[m] * dpsidy[m] ;
1285 grady_at_pt.fastAccessDx(m) = dpsidy[m] ;
1287 gradz_at_pt.val() += dof_values[m] * dpsidz[m] ;
1288 gradz_at_pt.fastAccessDx(m) = dpsidz[m] ;
1291 const local_scalar_type source_term =
1292 coeff_src * value_at_pt * value_at_pt ;
1294 const local_scalar_type advection_term =
1295 advection[0]*gradx_at_pt +
1296 advection[1]*grady_at_pt +
1297 advection[2]*gradz_at_pt;
1300 const double bases_val_m = bases_vals[m] * detJ_weight ;
1301 const double dpsidx_m = dpsidx[m] ;
1302 const double dpsidy_m = dpsidy[m] ;
1303 const double dpsidz_m = dpsidz[m] ;
1306 detJ_weight_coeff_k * ( dpsidx_m * gradx_at_pt +
1307 dpsidy_m * grady_at_pt +
1308 dpsidz_m * gradz_at_pt ) +
1309 bases_val_m * ( advection_term + source_term ) ;
1314 KOKKOS_INLINE_FUNCTION
1326 gatherSolution( ielem, val, node_index, x, y, z, elem_res );
1329 computeElementResidual( val, x, y, z, elem_res );
1332 this->scatterResidual( ielem, node_index, elem_res );
1337 class CoordinateMap ,
typename ScalarType ,
class CoeffFunctionType >
1340 CrsMatrix< ScalarType , ExecutionSpace > ,
1342 public ElementComputation< Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
1343 CrsMatrix< ScalarType , ExecutionSpace > ,
1344 Analytic , CoeffFunctionType > {
1365 const CoeffFunctionType & arg_coeff_function ,
1371 base_type(arg_mesh, arg_coeff_function, arg_solution, arg_elem_graph,
1372 arg_jacobian, arg_residual, arg_dev_config) {}
1379 parallel_for( nelem , *
this );
1382 KOKKOS_INLINE_FUNCTION
1392 double coeff_src = 1.234;
1393 double advection[] = { 1.1, 1.2, 1.3 };
1397 double pt[] = {0.0, 0.0, 0.0};
1407 if ( ! this->coeff_function.is_constant ) {
1414 coeff_k = this->coeff_function(pt, 0);
1421 dpsidx , dpsidy , dpsidz );
1422 const double detJ_weight = detJ * integ_weight;
1423 const scalar_type detJ_weight_coeff_k = detJ_weight * coeff_k;
1425 value_at_pt.val() = 0.0 ;
1426 gradx_at_pt.val() = 0.0 ;
1427 grady_at_pt.val() = 0.0 ;
1428 gradz_at_pt.val() = 0.0 ;
1430 value_at_pt.val() += dof_values[m] * bases_vals[m] ;
1431 gradx_at_pt.val() += dof_values[m] * dpsidx[m] ;
1432 grady_at_pt.val() += dof_values[m] * dpsidy[m] ;
1433 gradz_at_pt.val() += dof_values[m] * dpsidz[m] ;
1437 coeff_src * value_at_pt * value_at_pt ;
1440 advection[0]*gradx_at_pt +
1441 advection[1]*grady_at_pt +
1442 advection[2]*gradz_at_pt;
1445 const double bases_val_m = bases_vals[m] * detJ_weight ;
1447 detJ_weight_coeff_k * ( dpsidx[m] * gradx_at_pt +
1448 dpsidy[m] * grady_at_pt +
1449 dpsidz[m] * gradz_at_pt ) +
1450 bases_val_m * ( advection_term + source_term ) ;
1452 elem_res[m] += res.val();
1456 mat[
n] += res.fastAccessDx(0) * bases_vals[
n] +
1457 res.fastAccessDx(1) * dpsidx[
n] +
1458 res.fastAccessDx(2) * dpsidy[
n] +
1459 res.fastAccessDx(3) * dpsidz[
n];
1465 KOKKOS_INLINE_FUNCTION
1478 this->gatherSolution( ielem, val, node_index, x, y, z, elem_res, elem_mat );
1481 computeElementResidualJacobian( val, x, y, z, elem_res, elem_mat );
1484 this->scatterResidual( ielem, node_index, elem_res, elem_mat );
1489 template<
class FiniteElementMeshType ,
class SparseMatrixType
1490 ,
class CoeffFunctionType = ElementComputationConstantCoefficient
1492 class ElementComputation ;
1496 typename ScalarType ,
class CoeffFunctionType >
1497 class ElementComputation
1500 , CoeffFunctionType >
1515 typedef Kokkos::View< scalar_type* , Kokkos::LayoutLeft, execution_space >
vector_type ;
1519 typedef LocalViewTraits< vector_type > local_vector_view_traits;
1520 typedef LocalViewTraits< matrix_values_type> local_matrix_view_traits;
1521 typedef typename local_vector_view_traits::local_view_type local_vector_type;
1522 typedef typename local_matrix_view_traits::local_view_type local_matrix_type;
1523 typedef typename local_vector_view_traits::local_value_type local_scalar_type;
1524 static const bool use_team = local_vector_view_traits::use_team;
1536 typedef Kokkos::View< scalar_type*[FunctionCount][FunctionCount] , execution_space >
elem_matrices_type ;
1537 typedef Kokkos::View< scalar_type*[FunctionCount] , execution_space >
elem_vectors_type ;
1539 typedef LocalViewTraits< elem_matrices_type > local_elem_matrices_traits;
1540 typedef LocalViewTraits< elem_vectors_type > local_elem_vectors_traits;
1541 typedef typename local_elem_matrices_traits::local_view_type local_elem_matrices_type;
1542 typedef typename local_elem_vectors_traits::local_view_type local_elem_vectors_type;
1544 typedef typename NodeNodeGraph< elem_node_type , sparse_graph_type , ElemNodeCount >::ElemGraphType
elem_graph_type ;
1561 const CoeffFunctionType coeff_function ;
1564 ElementComputation(
const ElementComputation & rhs )
1574 , coeff_function( rhs.coeff_function )
1575 , dev_config( rhs.dev_config )
1580 ElementComputation(
const mesh_type & arg_mesh ,
1581 const CoeffFunctionType & arg_coeff_function ,
1596 , coeff_function( arg_coeff_function )
1597 , dev_config( arg_dev_config )
1606 const size_t team_size = dev_config.block_dim.x * dev_config.block_dim.y;
1607 const size_t league_size =
1608 (nelem + dev_config.block_dim.y-1) / dev_config.block_dim.y;
1609 Kokkos::TeamPolicy< execution_space > config( league_size, team_size );
1610 parallel_for( config , *
this );
1613 parallel_for( nelem , *
this );
1619 static const unsigned FLOPS_transform_gradients =
1624 KOKKOS_INLINE_FUNCTION
1632 double dpsidz[] )
const
1634 enum { j11 = 0 , j12 = 1 , j13 = 2 ,
1635 j21 = 3 , j22 = 4 , j23 = 5 ,
1636 j31 = 6 , j32 = 7 , j33 = 8 };
1640 double J[
TensorDim ] = { 0, 0, 0, 0, 0, 0, 0, 0, 0 };
1643 const double x1 = x[i] ;
1644 const double x2 = y[i] ;
1645 const double x3 = z[i] ;
1647 const double g1 = grad[0][i] ;
1648 const double g2 = grad[1][i] ;
1649 const double g3 = grad[2][i] ;
1667 static_cast<double>( J[j22] * J[j33] - J[j23] * J[j32] ) ,
1668 static_cast<double>( J[j13] * J[j32] - J[j12] * J[j33] ) ,
1669 static_cast<double>( J[j12] * J[j23] - J[j13] * J[j22] ) ,
1671 static_cast<double>( J[j23] * J[j31] - J[j21] * J[j33] ) ,
1672 static_cast<double>( J[j11] * J[j33] - J[j13] * J[j31] ) ,
1673 static_cast<double>( J[j13] * J[j21] - J[j11] * J[j23] ) ,
1675 static_cast<double>( J[j21] * J[j32] - J[j22] * J[j31] ) ,
1676 static_cast<double>( J[j12] * J[j31] - J[j11] * J[j32] ) ,
1677 static_cast<double>( J[j11] * J[j22] - J[j12] * J[j21] ) };
1679 const double detJ = J[j11] * invJ[j11] +
1680 J[j21] * invJ[j12] +
1681 J[j31] * invJ[j13] ;
1683 const double detJinv = 1.0 / detJ ;
1685 for (
unsigned i = 0 ; i <
TensorDim ; ++i ) { invJ[i] *= detJinv ; }
1690 const double g0 = grad[0][i];
1691 const double g1 = grad[1][i];
1692 const double g2 = grad[2][i];
1694 dpsidx[i] = g0 * invJ[j11] + g1 * invJ[j12] + g2 * invJ[j13];
1695 dpsidy[i] = g0 * invJ[j21] + g1 * invJ[j22] + g2 * invJ[j23];
1696 dpsidz[i] = g0 * invJ[j31] + g1 * invJ[j32] + g2 * invJ[j33];
1702 KOKKOS_INLINE_FUNCTION
1703 void contributeResidualJacobian(
1704 const local_scalar_type dof_values[] ,
1705 const double dpsidx[] ,
1706 const double dpsidy[] ,
1707 const double dpsidz[] ,
1709 const local_scalar_type coeff_k ,
1710 const double integ_weight ,
1711 const double bases_vals[] ,
1712 local_scalar_type elem_res[] ,
1713 local_scalar_type elem_mat[][ FunctionCount ] )
const
1715 local_scalar_type value_at_pt = 0 ;
1716 local_scalar_type gradx_at_pt = 0 ;
1717 local_scalar_type grady_at_pt = 0 ;
1718 local_scalar_type gradz_at_pt = 0 ;
1721 value_at_pt += dof_values[m] * bases_vals[m] ;
1722 gradx_at_pt += dof_values[m] * dpsidx[m] ;
1723 grady_at_pt += dof_values[m] * dpsidy[m] ;
1724 gradz_at_pt += dof_values[m] * dpsidz[m] ;
1727 const local_scalar_type k_detJ_weight = coeff_k * detJ * integ_weight ;
1728 const local_scalar_type res_val = value_at_pt * value_at_pt * detJ * integ_weight ;
1729 const local_scalar_type mat_val = 2.0 * value_at_pt * detJ * integ_weight ;
1735 local_scalar_type *
const mat = elem_mat[m] ;
1736 const double bases_val_m = bases_vals[m];
1737 const double dpsidx_m = dpsidx[m] ;
1738 const double dpsidy_m = dpsidy[m] ;
1739 const double dpsidz_m = dpsidz[m] ;
1741 elem_res[m] += k_detJ_weight * ( dpsidx_m * gradx_at_pt +
1742 dpsidy_m * grady_at_pt +
1743 dpsidz_m * gradz_at_pt ) +
1744 res_val * bases_val_m ;
1748 mat[
n] += k_detJ_weight * ( dpsidx_m * dpsidx[
n] +
1749 dpsidy_m * dpsidy[
n] +
1750 dpsidz_m * dpsidz[
n] ) +
1751 mat_val * bases_val_m * bases_vals[n];
1756 KOKKOS_INLINE_FUNCTION
1757 void operator()(
const typename TeamPolicy< execution_space >::member_type & dev )
const
1760 const unsigned num_ensemble_threads = dev_config.block_dim.x ;
1761 const unsigned num_element_threads = dev_config.block_dim.y ;
1762 const unsigned element_rank = dev.team_rank() / num_ensemble_threads ;
1763 const unsigned ensemble_rank = dev.team_rank() % num_ensemble_threads ;
1765 const unsigned ielem =
1766 dev.league_rank() * num_element_threads + element_rank;
1771 (*this)( ielem, ensemble_rank );
1775 KOKKOS_INLINE_FUNCTION
1776 void operator()(
const unsigned ielem ,
1777 const unsigned ensemble_rank = 0 )
const
1779 local_vector_type local_solution =
1780 local_vector_view_traits::create_local_view(
solution,
1782 local_vector_type local_residual =
1783 local_vector_view_traits::create_local_view(
residual,
1785 local_matrix_type local_jacobian_values =
1800 node_index[i] = ni ;
1806 val[i] = local_solution( ni );
1816 elem_mat[i][
j] = 0 ;
1826 local_scalar_type coeff_k = 0 ;
1829 double pt[] = {0.0, 0.0, 0.0};
1833 if ( ! coeff_function.is_constant ) {
1842 coeff_k = coeff_function(pt, ensemble_rank);
1847 dpsidx , dpsidy , dpsidz );
1849 contributeResidualJacobian( val , dpsidx , dpsidy , dpsidz ,
1853 elem_vec , elem_mat );
1857 const unsigned row = node_index[i] ;
1859 atomic_add( & local_residual( row ) , elem_vec[i] );
1862 const unsigned entry =
elem_graph( ielem , i ,
j );
1863 if ( entry != ~0u ) {
1864 atomic_add( & local_jacobian_values( entry ) , elem_mat[i][
j] );
1875 template<
class FixtureType ,
class SparseMatrixType >
1879 typename ScalarType >
1904 static const bool use_team = local_vector_view_traits::use_team;
1929 const unsigned arg_bc_plane ,
1937 , bc_lower_value( arg_bc_lower_value )
1938 , bc_upper_value( arg_bc_upper_value )
1941 , bc_plane( arg_bc_plane )
1942 , node_count( arg_mesh.node_count_owned() )
1944 , dev_config( arg_dev_config )
1946 parallel_for( node_count , *
this );
1953 const size_t team_size = dev_config.block_dim.x * dev_config.block_dim.y;
1954 const size_t league_size =
1955 (node_count + dev_config.block_dim.y-1) / dev_config.block_dim.y;
1956 Kokkos::TeamPolicy< execution_space > config( league_size, team_size );
1957 parallel_for( config , *
this );
1960 parallel_for( node_count , *
this );
1965 KOKKOS_INLINE_FUNCTION
1966 void operator()(
const typename TeamPolicy< execution_space >::member_type & dev )
const
1969 const unsigned num_ensemble_threads = dev_config.block_dim.x ;
1970 const unsigned num_node_threads = dev_config.block_dim.y ;
1971 const unsigned node_rank = dev.team_rank() / num_ensemble_threads ;
1972 const unsigned ensemble_rank = dev.team_rank() % num_ensemble_threads ;
1974 const unsigned inode = dev.league_rank() * num_node_threads + node_rank;
1976 if (inode >= node_count)
1979 (*this)( inode, ensemble_rank );
1983 KOKKOS_INLINE_FUNCTION
1985 const unsigned ensemble_rank = 0)
const
1988 local_vector_view_traits::create_local_view(
residual,
2003 const bool bc_lower = c <= bc_lower_limit ;
2004 const bool bc_upper = bc_upper_limit <= c ;
2007 solution(inode) = bc_lower ? bc_lower_value : (
2008 bc_upper ? bc_upper_value : 0 );
2011 if ( bc_lower || bc_upper ) {
2013 local_residual(inode) = 0 ;
2018 for(
unsigned i = iBeg ; i < iEnd ; ++i ) {
2027 for(
unsigned i = iBeg ; i < iEnd ; ++i ) {
2031 if ( ( cc <= bc_lower_limit ) || ( bc_upper_limit <= cc ) ) {
2032 local_jacobian_values(i) = 0 ;
2047 #if defined( __CUDACC__ )
double weights[integration_count]
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)
LocalViewTraits< vector_type > local_vector_view_traits
ExponentialKLCoefficient(const ExponentialKLCoefficient &rhs)
const scalar_coord_type bc_upper_limit
static const unsigned SpatialDim
Kokkos::View< scalar_type *, execution_space > vector_type
const Kokkos::Example::FENL::DeviceConfig dev_config
KOKKOS_INLINE_FUNCTION void setRandomVariables(const RandomVariableView &rv)
LocalViewTraits< matrix_values_type > local_matrix_view_traits
Sacado::Fad::SLFad< scalar_type, FunctionCount > fad_scalar_type
const unsigned node_count
KOKKOS_INLINE_FUNCTION void computeElementResidualJacobian(const scalar_type dof_values[], const double x[], const double y[], const double z[], scalar_type elem_res[], scalar_type elem_mat[][FunctionCount]) const
ElementComputation(const ElementComputation &rhs)
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.
const elem_vectors_type elem_residuals
KOKKOS_INLINE_FUNCTION void fill_graph_entries(const unsigned iset) const
const vector_type solution
static const bool use_team
Kokkos::View< scalar_type *[FunctionCount], execution_space > elem_vectors_type
static const unsigned element_node_count
CrsMatrix< ScalarType, ExecutionSpace > sparse_matrix_type
const elem_matrices_type elem_jacobians
view_type::value_type local_value_type
static const unsigned function_count
static void eval(Kokkos::Example::FENL::DeviceConfig &dev_config_elem, Kokkos::Example::FENL::DeviceConfig &dev_config_bc)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
ElementComputation(const typename base_type::mesh_type &arg_mesh, const CoeffFunctionType &arg_coeff_function, const typename base_type::vector_type &arg_solution, const typename base_type::elem_graph_type &arg_elem_graph, const typename base_type::sparse_matrix_type &arg_jacobian, const typename base_type::vector_type &arg_residual, const Kokkos::Example::FENL::DeviceConfig arg_dev_config)
base_type::execution_space execution_space
ElementComputationBase< ExecutionSpace, Order, CoordinateMap, ScalarType > base_type
Kokkos::View< const double *[SpaceDim], Device > node_coord_type
Kokkos::DefaultExecutionSpace execution_space
double gradients[integration_count][spatial_dimension][function_count]
Kokkos::View< const unsigned *[ElemNode], Device > elem_node_type
static const unsigned ElemNodeCount
Sacado::Fad::SLFad< scalar_type, 4 > fad_scalar_type
KOKKOS_INLINE_FUNCTION RandomVariableView getRandomVariables() const
static void eval(Kokkos::Example::FENL::DeviceConfig &dev_config_elem, Kokkos::Example::FENL::DeviceConfig &dev_config_bc)
const Kokkos::Example::FENL::DeviceConfig dev_config
const CoeffFunctionType coeff_function
ExecutionSpace execution_space
ElementComputation(const ElementComputation &rhs)
KOKKOS_INLINE_FUNCTION void scatterResidual(const unsigned ielem, const unsigned node_index[], fad_scalar_type res[]) const
Kokkos::View< scalar_type *, Kokkos::LayoutLeft, execution_space > vector_type
KOKKOS_INLINE_FUNCTION void atomic_add(volatile Sacado::UQ::PCE< Storage > *const dest, const Sacado::UQ::PCE< Storage > &src)
const vector_type residual
const bc_scalar_type bc_lower_value
ElementComputation(const typename base_type::mesh_type &arg_mesh, const CoeffFunctionType &arg_coeff_function, const typename base_type::vector_type &arg_solution, const typename base_type::elem_graph_type &arg_elem_graph, const typename base_type::sparse_matrix_type &arg_jacobian, const typename base_type::vector_type &arg_residual, const Kokkos::Example::FENL::DeviceConfig arg_dev_config)
Kokkos::View< Scalar *, Kokkos::LayoutLeft, Device > RandomVariableView
const elem_node_type elem_node_ids
const vector_type residual
const unsigned node_count
KOKKOS_INLINE_FUNCTION void computeElementResidualJacobian(const scalar_type dof_values[], const double x[], const double y[], const double z[], scalar_type elem_res[], scalar_type elem_mat[][FunctionCount]) const
const bc_scalar_type bc_upper_value
base_type::scalar_type scalar_type
base_type::execution_space execution_space
KOKKOS_INLINE_FUNCTION local_scalar_type operator()(const MeshScalar point[], const size_type ensemble_rank) const
KOKKOS_INLINE_FUNCTION void operator()(const unsigned ielem) const
NodeNodeGraph< elem_node_type, sparse_graph_type, ElemNodeCount >::ElemGraphType elem_graph_type
const vector_type solution
const MeshScalar m_variance
local_rv_view_traits::local_view_type local_rv_view_type
ElementComputation< Kokkos::Example::BoxElemFixture< ExecutionSpace, Order, CoordinateMap >, CrsMatrix< ScalarType, ExecutionSpace >, Analytic, CoeffFunctionType > base_type
mesh_type::node_coord_type node_coord_type
ElemNodeIdView::execution_space execution_space
node_coord_type::value_type scalar_coord_type
static const unsigned TensorDim
KOKKOS_INLINE_FUNCTION void operator()(const unsigned irow, unsigned &update, const bool final) const
ElementComputationConstantCoefficient(const ElementComputationConstantCoefficient &rhs)
RandomVariableView::size_type size_type
KOKKOS_INLINE_FUNCTION void init(unsigned &update) const
ExecutionSpace execution_space
double sort_graph_entries
size_t num_threads_per_block
double fill_graph_entries
KOKKOS_INLINE_FUNCTION void operator()(const unsigned ielem) const
ElementComputationBase(const ElementComputationBase &rhs)
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)
Kokkos::Example::HexElement_Data< mesh_type::ElemNode > element_data_type
CrsGraphType::row_map_type::non_const_type RowMapType
base_type::scalar_type scalar_type
const element_data_type elem_data
const MeshScalar m_corr_len
KOKKOS_INLINE_FUNCTION void join(volatile unsigned &update, const volatile unsigned &input) const
ElementComputationBase(const mesh_type &arg_mesh, const vector_type &arg_solution, const elem_graph_type &arg_elem_graph, const sparse_matrix_type &arg_jacobian, const vector_type &arg_residual)
KOKKOS_INLINE_FUNCTION void operator()(const typename TeamPolicy< execution_space >::member_type &dev) const
ElementComputation(const ElementComputation &rhs)
base_type::execution_space execution_space
CrsMatrix(const StaticCrsGraphType &arg_graph)
const sparse_matrix_type jacobian
Kokkos::StaticCrsGraph< unsigned, Space, void, void, unsigned > StaticCrsGraphType
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
double values[integration_count][function_count]
mesh_type::elem_node_type elem_node_type
Kokkos::Example::FENL::CrsMatrix< ScalarType, ExecutionSpace > sparse_matrix_type
void deep_copy(const Stokhos::CrsMatrix< ValueType, DstDevice, Layout > &dst, const Stokhos::CrsMatrix< ValueType, SrcDevice, Layout > &src)
KOKKOS_INLINE_FUNCTION void operator()(const unsigned inode, const unsigned ensemble_rank=0) const
KOKKOS_INLINE_FUNCTION void operator()(const unsigned iwork) const
base_type::execution_space execution_space
KOKKOS_INLINE_FUNCTION double operator()(double pt[], unsigned ensemble_rank) const
base_type::scalar_type scalar_type
KOKKOS_INLINE_FUNCTION void gatherSolution(const unsigned ielem, scalar_type val[], unsigned node_index[], double x[], double y[], double z[], fad_scalar_type res[]) const
const CoeffFunctionType coeff_function
const scalar_coord_type bc_lower_limit
const ElemNodeIdView elem_node_id
static const unsigned integration_count
base_type::scalar_type scalar_type
Kokkos::Example::BoxElemFixture< ExecutionSpace, Order, CoordinateMap > mesh_type
mesh_type::node_coord_type node_coord_type
const elem_graph_type elem_graph
sparse_matrix_type::StaticCrsGraphType sparse_graph_type
const Kokkos::Example::FENL::DeviceConfig dev_config
Kokkos::Example::BoxElemFixture< ExecutionSpace, Order, CoordinateMap > mesh_type
const view_type & local_view_type
NodeNodeGraph(const ElemNodeIdView &arg_elem_node_id, const unsigned arg_node_count, Times &results)
const sparse_matrix_type jacobian
const node_coord_type node_coords
static const unsigned FunctionCount
KOKKOS_INLINE_FUNCTION void gatherSolution(const unsigned ielem, scalar_type val[], unsigned node_index[], double x[], double y[], double z[], scalar_type res[], scalar_type mat[][FunctionCount]) const
ElementComputation< Kokkos::Example::BoxElemFixture< ExecutionSpace, Order, CoordinateMap >, CrsMatrix< ScalarType, ExecutionSpace >, FadElement, CoeffFunctionType > base_type
sparse_matrix_type::StaticCrsGraphType sparse_graph_type
Kokkos::View< unsigned, execution_space > UnsignedValue
const node_coord_type node_coords
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)
ExponentialKLCoefficient< T1, T2, T3 > type
KOKKOS_INLINE_FUNCTION void operator()(const unsigned ielem) const
ElementComputation(const typename base_type::mesh_type &arg_mesh, const CoeffFunctionType &arg_coeff_function, const typename base_type::vector_type &arg_solution, const typename base_type::elem_graph_type &arg_elem_graph, const typename base_type::sparse_matrix_type &arg_jacobian, const typename base_type::vector_type &arg_residual, const Kokkos::Example::FENL::DeviceConfig arg_dev_config)
KOKKOS_INLINE_FUNCTION void gatherSolution(const unsigned ielem, fad_scalar_type val[], unsigned node_index[], double x[], double y[], double z[], fad_scalar_type res[]) const
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)
ElementComputationBase< ExecutionSpace, Order, CoordinateMap, ScalarType > base_type
KOKKOS_INLINE_FUNCTION PCE< Storage > sin(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION void scatterResidual(const unsigned ielem, const unsigned node_index[], const scalar_type res[], const scalar_type mat[][FunctionCount]) const
ElementComputation(const typename base_type::mesh_type &arg_mesh, const CoeffFunctionType &arg_coeff_function, const typename base_type::vector_type &arg_solution, const typename base_type::elem_graph_type &arg_elem_graph, const typename base_type::sparse_matrix_type &arg_jacobian, const typename base_type::vector_type &arg_residual, const Kokkos::Example::FENL::DeviceConfig arg_dev_config)
View< ValueType *, Space > values_type
StorageType::execution_space execution_space
local_matrix_view_traits::local_view_type local_matrix_type
KOKKOS_INLINE_FUNCTION void sort_graph_entries(const unsigned irow) const
double fill_element_graph
Kokkos::View< scalar_type *[FunctionCount][FunctionCount], execution_space > elem_matrices_type
const unsigned VectorSize
KOKKOS_INLINE_FUNCTION void fill_set(const unsigned ielem) const
ElementComputationConstantCoefficient(const double val)
Kokkos::View< unsigned *[ElemNode][ElemNode], execution_space > ElemGraphType
Sacado::Fad::SLFad< scalar_type, FunctionCount > fad_scalar_type
Generate a distributed unstructured finite element mesh from a partitioned NX*NY*NZ box of elements...
Stokhos::KL::ExponentialRandomField< MeshScalar, Device > rf_type
sparse_matrix_type::values_type matrix_values_type
ExponentialKLCoefficient(const MeshScalar mean, const MeshScalar variance, const MeshScalar correlation_length, const size_type num_rv)
KOKKOS_INLINE_FUNCTION void operator()(const unsigned ielem) const
local_vector_view_traits::local_view_type local_vector_type
KOKKOS_INLINE_FUNCTION void computeElementResidual(const local_scalar_type dof_values[], const double x[], const double y[], const double z[], local_scalar_type elem_res[]) const
ElementComputation(const ElementComputation &rhs)
void update(const ValueType &alpha, VectorType &x, const ValueType &beta, const VectorType &y)
Kokkos::UnorderedMap< key_type, void, execution_space > SetType
pair< unsigned, unsigned > key_type
KOKKOS_INLINE_FUNCTION void computeElementResidual(const scalar_type dof_values[], const double x[], const double y[], const double z[], local_scalar_type elem_res[]) const
static const unsigned IntegrationCount
static const unsigned spatial_dimension