44 #ifndef KOKKOS_EXAMPLE_FENLFUNCTORS_HPP
45 #define KOKKOS_EXAMPLE_FENLFUNCTORS_HPP
56 #include <Kokkos_Core.hpp>
57 #include <Kokkos_Pair.hpp>
58 #include <Kokkos_UnorderedMap.hpp>
59 #include <Kokkos_StaticCrsGraph.hpp>
61 #include <impl/Kokkos_Timer.hpp>
75 template<
typename ValueType ,
class Space >
77 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE // Don't remove this until Kokkos has removed the deprecated code path probably around September 2018
78 typedef Kokkos::StaticCrsGraph< unsigned , Space , void , unsigned >
StaticCrsGraphType ;
91 ,
coeff(
"crs_matrix_coeff" , arg_graph.entries.extent(0) )
95 template<
class ElemNodeIdView ,
class CrsGraphType ,
unsigned ElemNode >
102 typedef Kokkos::UnorderedMap< key_type, void , execution_space >
SetType ;
103 typedef typename CrsGraphType::row_map_type::non_const_type
RowMapType ;
107 typedef Kokkos::View< unsigned*[ElemNode][ElemNode] , execution_space >
148 const unsigned arg_node_count,
164 Kokkos::Impl::Timer wall_clock ;
170 size_t set_capacity = (28ull *
node_count) / 2;
171 unsigned failed_insert_count = 0 ;
182 Kokkos::parallel_reduce( Kokkos::RangePolicy<execution_space,TagFillNodeSet>(0,
elem_node_id.extent(0))
184 , failed_insert_count );
186 }
while ( failed_insert_count );
204 unsigned graph_entry_count = 0 ;
206 Kokkos::deep_copy( graph_entry_count ,
row_total );
210 graph.entries =
typename CrsGraphType::entries_type(
"graph_entries" , graph_entry_count );
260 for (
unsigned row_local_node = 0 ; row_local_node <
elem_node_id.extent(1) ; ++row_local_node ) {
262 const unsigned row_node =
elem_node_id( ielem , row_local_node );
264 for (
unsigned col_local_node = row_local_node ; col_local_node <
elem_node_id.extent(1) ; ++col_local_node ) {
266 const unsigned col_node =
elem_node_id( ielem , col_local_node );
272 const key_type key = (row_node < col_node) ? make_pair( row_node, col_node ) : make_pair( col_node, row_node ) ;
274 const typename SetType::insert_result result =
node_node_set.insert( key );
277 if ( result.success() ) {
280 if ( row_node <
row_count.extent(0) ) { atomic_increment( &
row_count( row_node ) ); }
283 if ( col_node <
row_count.extent(0) && col_node != row_node ) { atomic_increment( &
row_count( col_node ) ); }
285 else if ( result.failed() ) {
300 const unsigned row_node = key.first ;
301 const unsigned col_node = key.second ;
304 typedef typename std::remove_reference< decltype( row_count(0) ) >::type atomic_incr_type;
305 const unsigned offset =
graph.row_map( row_node ) + atomic_fetch_add( &
row_count( row_node ) , atomic_incr_type(1) );
306 graph.entries( offset ) = col_node ;
309 if ( col_node <
row_count.extent(0) && col_node != row_node ) {
310 typedef typename std::remove_reference< decltype( row_count(0) ) >::type atomic_incr_type;
311 const unsigned offset =
graph.row_map( col_node ) + atomic_fetch_add( &
row_count( col_node ) , atomic_incr_type(1) );
312 graph.entries( offset ) = row_node ;
320 const unsigned row_beg =
graph.row_map( irow );
321 const unsigned row_end =
graph.row_map( irow + 1 );
322 for (
unsigned i = row_beg + 1 ; i < row_end ; ++i ) {
323 const unsigned col =
graph.entries(i);
325 for ( ; row_beg < j && col <
graph.entries(j-1) ; --j ) {
328 graph.entries(j) = col ;
335 for (
unsigned row_local_node = 0 ; row_local_node <
elem_node_id.extent(1) ; ++row_local_node ) {
337 const unsigned row_node =
elem_node_id( ielem , row_local_node );
339 for (
unsigned col_local_node = 0 ; col_local_node <
elem_node_id.extent(1) ; ++col_local_node ) {
341 const unsigned col_node =
elem_node_id( ielem , col_local_node );
343 unsigned entry = ~0u ;
345 if ( row_node + 1 <
graph.row_map.extent(0) ) {
347 const unsigned entry_end =
graph.row_map( row_node + 1 );
349 entry =
graph.row_map( row_node );
351 for ( ; entry < entry_end &&
graph.entries(entry) != col_node ; ++entry );
353 if ( entry == entry_end ) entry = ~0u ;
356 elem_graph( ielem , row_local_node , col_local_node ) = entry ;
386 void operator()(
const unsigned irow ,
unsigned & update ,
const bool final )
const
389 if (
final ) {
row_map( irow ) = update ; }
407 ,
volatile unsigned & update
408 ,
volatile const unsigned & input )
const { update += input ; }
412 void init(
unsigned & update )
const { update = 0 ; }
415 void join(
volatile unsigned & update
416 ,
volatile const unsigned & input )
const { update += input ; }
433 class CoordinateMap ,
typename ScalarType >
448 typedef Kokkos::View< scalar_type* , Kokkos::LayoutLeft, execution_space >
vector_type ;
462 typedef Kokkos::View< scalar_type*[FunctionCount][FunctionCount] , execution_space >
elem_matrices_type ;
521 double dpsidz[] )
const
523 enum { j11 = 0 , j12 = 1 , j13 = 2 ,
524 j21 = 3 , j22 = 4 , j23 = 5 ,
525 j31 = 6 , j32 = 7 , j33 = 8 };
529 double J[
TensorDim ] = { 0, 0, 0, 0, 0, 0, 0, 0, 0 };
532 const double x1 = x[i] ;
533 const double x2 = y[i] ;
534 const double x3 = z[i] ;
536 const double g1 = grad[0][i] ;
537 const double g2 = grad[1][i] ;
538 const double g3 = grad[2][i] ;
556 static_cast<double>( J[j22] * J[j33] - J[j23] * J[j32] ) ,
557 static_cast<double>( J[j13] * J[j32] - J[j12] * J[j33] ) ,
558 static_cast<double>( J[j12] * J[j23] - J[j13] * J[j22] ) ,
560 static_cast<double>( J[j23] * J[j31] - J[j21] * J[j33] ) ,
561 static_cast<double>( J[j11] * J[j33] - J[j13] * J[j31] ) ,
562 static_cast<double>( J[j13] * J[j21] - J[j11] * J[j23] ) ,
564 static_cast<double>( J[j21] * J[j32] - J[j22] * J[j31] ) ,
565 static_cast<double>( J[j12] * J[j31] - J[j11] * J[j32] ) ,
566 static_cast<double>( J[j11] * J[j22] - J[j12] * J[j21] ) };
568 const double detJ = J[j11] * invJ[j11] +
572 const double detJinv = 1.0 / detJ ;
574 for (
unsigned i = 0 ; i <
TensorDim ; ++i ) { invJ[i] *= detJinv ; }
579 const double g0 = grad[0][i];
580 const double g1 = grad[1][i];
581 const double g2 = grad[2][i];
583 dpsidx[i] = g0 * invJ[j11] + g1 * invJ[j12] + g2 * invJ[j13];
584 dpsidy[i] = g0 * invJ[j21] + g1 * invJ[j22] + g2 * invJ[j23];
585 dpsidz[i] = g0 * invJ[j31] + g1 * invJ[j32] + g2 * invJ[j33];
600 template<
class FiniteElementMeshType ,
601 class SparseMatrixType ,
607 class CoordinateMap ,
typename ScalarType >
610 CrsMatrix< ScalarType , ExecutionSpace > ,
634 base_type(arg_mesh, arg_solution, arg_elem_graph,
635 arg_jacobian, arg_residual) {}
642 parallel_for( nelem , *
this );
648 unsigned node_index[],
649 double x[],
double y[],
double z[],
673 const unsigned node_index[],
678 const unsigned row = node_index[i] ;
679 if ( row < this->
residual.extent(0) ) {
680 atomic_add( & this->
residual( row ) , res[i] );
683 const unsigned entry = this->
elem_graph( ielem , i , j );
684 if ( entry != ~0u ) {
701 double coeff_k = 3.456;
702 double coeff_src = 1.234;
703 double advection[] = { 1.1, 1.2, 1.3 };
714 dpsidx , dpsidy , dpsidz );
715 const double detJ_weight = detJ * integ_weight;
716 const double detJ_weight_coeff_k = detJ_weight * coeff_k;
723 value_at_pt += dof_values[m] * bases_vals[m] ;
724 gradx_at_pt += dof_values[m] * dpsidx[m] ;
725 grady_at_pt += dof_values[m] * dpsidy[m] ;
726 gradz_at_pt += dof_values[m] * dpsidz[m] ;
730 coeff_src * value_at_pt * value_at_pt ;
732 2.0 * coeff_src * value_at_pt ;
739 advection_x*gradx_at_pt +
740 advection_y*grady_at_pt +
741 advection_z*gradz_at_pt ;
745 const double bases_val_m = bases_vals[m] * detJ_weight ;
746 const double dpsidx_m = dpsidx[m] ;
747 const double dpsidy_m = dpsidy[m] ;
748 const double dpsidz_m = dpsidz[m] ;
751 detJ_weight_coeff_k * ( dpsidx_m * gradx_at_pt +
752 dpsidy_m * grady_at_pt +
753 dpsidz_m * gradz_at_pt ) +
754 bases_val_m * ( advection_term + source_term ) ;
757 const double dpsidx_n = dpsidx[
n] ;
758 const double dpsidy_n = dpsidy[
n] ;
759 const double dpsidz_n = dpsidz[
n] ;
761 detJ_weight_coeff_k * ( dpsidx_m * dpsidx_n +
762 dpsidy_m * dpsidy_n +
763 dpsidz_m * dpsidz_n ) +
764 bases_val_m * ( advection_x * dpsidx_n +
765 advection_y * dpsidy_n +
766 advection_z * dpsidz_n +
767 source_deriv * bases_vals[
n] ) ;
786 gatherSolution(ielem, val, node_index, x, y, z, elem_res, elem_mat);
789 computeElementResidualJacobian( val, x, y, z, elem_res , elem_mat );
792 scatterResidual( ielem, node_index, elem_res, elem_mat );
797 class CoordinateMap ,
typename ScalarType >
800 CrsMatrix< ScalarType , ExecutionSpace > ,
825 base_type(arg_mesh, arg_solution, arg_elem_graph,
826 arg_jacobian, arg_residual) {}
833 parallel_for( nelem , *
this );
839 unsigned node_index[],
840 double x[],
double y[],
double z[],
852 val[i].val() = this->
solution( ni );
859 const unsigned node_index[],
863 const unsigned row = node_index[i] ;
864 if ( row < this->
residual.extent(0) ) {
865 atomic_add( & this->
residual( row ) , res[i].
val() );
868 const unsigned entry = this->
elem_graph( ielem , i , j );
869 if ( entry != ~0u ) {
871 res[i].fastAccessDx(j) );
885 double coeff_k = 3.456;
886 double coeff_src = 1.234;
887 double advection[] = { 1.1, 1.2, 1.3 };
898 dpsidx , dpsidy , dpsidz );
899 const double detJ_weight = detJ * integ_weight;
900 const double detJ_weight_coeff_k = detJ_weight * coeff_k;
907 value_at_pt += dof_values[m] * bases_vals[m] ;
908 gradx_at_pt += dof_values[m] * dpsidx[m] ;
909 grady_at_pt += dof_values[m] * dpsidy[m] ;
910 gradz_at_pt += dof_values[m] * dpsidz[m] ;
914 coeff_src * value_at_pt * value_at_pt ;
917 advection[0]*gradx_at_pt +
918 advection[1]*grady_at_pt +
919 advection[2]*gradz_at_pt;
922 const double bases_val_m = bases_vals[m] * detJ_weight ;
923 const double dpsidx_m = dpsidx[m] ;
924 const double dpsidy_m = dpsidy[m] ;
925 const double dpsidz_m = dpsidz[m] ;
928 detJ_weight_coeff_k * ( dpsidx_m * gradx_at_pt +
929 dpsidy_m * grady_at_pt +
930 dpsidz_m * gradz_at_pt ) +
931 bases_val_m * ( advection_term + source_term ) ;
948 gatherSolution( ielem, val, node_index, x, y, z, elem_res );
951 computeElementResidual( val, x, y, z, elem_res );
954 scatterResidual( ielem, node_index, elem_res );
959 class CoordinateMap ,
typename ScalarType >
962 CrsMatrix< ScalarType , ExecutionSpace > ,
964 public ElementComputation< Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
965 CrsMatrix< ScalarType , ExecutionSpace > ,
990 base_type(arg_mesh, arg_solution, arg_elem_graph,
991 arg_jacobian, arg_residual) {}
998 parallel_for( nelem , *
this );
1004 unsigned node_index[],
1005 double x[],
double y[],
double z[],
1011 node_index[i] = ni ;
1028 double coeff_k = 3.456;
1029 double coeff_src = 1.234;
1030 double advection[] = { 1.1, 1.2, 1.3 };
1041 dpsidx , dpsidy , dpsidz );
1042 const double detJ_weight = detJ * integ_weight;
1043 const double detJ_weight_coeff_k = detJ_weight * coeff_k;
1050 value_at_pt.val() += dof_values[m] * bases_vals[m] ;
1051 value_at_pt.fastAccessDx(m) = bases_vals[m] ;
1053 gradx_at_pt.val() += dof_values[m] * dpsidx[m] ;
1054 gradx_at_pt.fastAccessDx(m) = dpsidx[m] ;
1056 grady_at_pt.val() += dof_values[m] * dpsidy[m] ;
1057 grady_at_pt.fastAccessDx(m) = dpsidy[m] ;
1059 gradz_at_pt.val() += dof_values[m] * dpsidz[m] ;
1060 gradz_at_pt.fastAccessDx(m) = dpsidz[m] ;
1064 coeff_src * value_at_pt * value_at_pt ;
1067 advection[0]*gradx_at_pt +
1068 advection[1]*grady_at_pt +
1069 advection[2]*gradz_at_pt;
1072 const double bases_val_m = bases_vals[m] * detJ_weight ;
1073 const double dpsidx_m = dpsidx[m] ;
1074 const double dpsidy_m = dpsidy[m] ;
1075 const double dpsidz_m = dpsidz[m] ;
1078 detJ_weight_coeff_k * ( dpsidx_m * gradx_at_pt +
1079 dpsidy_m * grady_at_pt +
1080 dpsidz_m * gradz_at_pt ) +
1081 bases_val_m * ( advection_term + source_term ) ;
1098 gatherSolution( ielem, val, node_index, x, y, z, elem_res );
1101 computeElementResidual( val, x, y, z, elem_res );
1104 this->scatterResidual( ielem, node_index, elem_res );
1109 class CoordinateMap ,
typename ScalarType >
1112 CrsMatrix< ScalarType , ExecutionSpace > ,
1114 public ElementComputation< Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
1115 CrsMatrix< ScalarType , ExecutionSpace > ,
1140 base_type(arg_mesh, arg_solution, arg_elem_graph,
1141 arg_jacobian, arg_residual) {}
1148 parallel_for( nelem , *
this );
1160 double coeff_k = 3.456;
1161 double coeff_src = 1.234;
1162 double advection[] = { 1.1, 1.2, 1.3 };
1178 dpsidx , dpsidy , dpsidz );
1179 const double detJ_weight = detJ * integ_weight;
1180 const double detJ_weight_coeff_k = detJ_weight * coeff_k;
1182 value_at_pt.val() = 0.0 ;
1183 gradx_at_pt.val() = 0.0 ;
1184 grady_at_pt.val() = 0.0 ;
1185 gradz_at_pt.val() = 0.0 ;
1187 value_at_pt.val() += dof_values[m] * bases_vals[m] ;
1188 gradx_at_pt.val() += dof_values[m] * dpsidx[m] ;
1189 grady_at_pt.val() += dof_values[m] * dpsidy[m] ;
1190 gradz_at_pt.val() += dof_values[m] * dpsidz[m] ;
1194 coeff_src * value_at_pt * value_at_pt ;
1197 advection[0]*gradx_at_pt +
1198 advection[1]*grady_at_pt +
1199 advection[2]*gradz_at_pt;
1202 const double bases_val_m = bases_vals[m] * detJ_weight ;
1204 detJ_weight_coeff_k * ( dpsidx[m] * gradx_at_pt +
1205 dpsidy[m] * grady_at_pt +
1206 dpsidz[m] * gradz_at_pt ) +
1207 bases_val_m * ( advection_term + source_term ) ;
1209 elem_res[m] += res.val();
1213 mat[
n] += res.fastAccessDx(0) * bases_vals[
n] +
1214 res.fastAccessDx(1) * dpsidx[
n] +
1215 res.fastAccessDx(2) * dpsidy[
n] +
1216 res.fastAccessDx(3) * dpsidz[
n];
1235 this->gatherSolution( ielem, val, node_index, x, y, z, elem_res, elem_mat );
1238 computeElementResidualJacobian( val, x, y, z, elem_res, elem_mat );
1241 this->scatterResidual( ielem, node_index, elem_res, elem_mat );
double weights[integration_count]
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 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)
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
static const unsigned TensorDim
base_type::scalar_type scalar_type
base_type::execution_space execution_space
const elem_vectors_type elem_residuals
KOKKOS_INLINE_FUNCTION void fill_graph_entries(const unsigned iset) const
static const unsigned SpatialDim
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
static const unsigned function_count
ElementComputation< Kokkos::Example::BoxElemFixture< ExecutionSpace, Order, CoordinateMap >, CrsMatrix< ScalarType, ExecutionSpace >, Analytic > base_type
Kokkos::View< const double *[SpaceDim], Device > node_coord_type
KOKKOS_INLINE_FUNCTION void operator()(const unsigned ielem) const
double gradients[integration_count][spatial_dimension][function_count]
Kokkos::View< const unsigned *[ElemNode], Device > elem_node_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
ElementComputationBase< ExecutionSpace, Order, CoordinateMap, ScalarType > base_type
KOKKOS_INLINE_FUNCTION void join(volatile unsigned &update, volatile const unsigned &input) const
ExecutionSpace execution_space
Kokkos::View< scalar_type *, Kokkos::LayoutLeft, execution_space > vector_type
base_type::execution_space execution_space
const vector_type residual
Sacado::Fad::SFad< scalar_type, 4 > fad_scalar_type
static const unsigned IntegrationCount
const elem_node_type elem_node_ids
const unsigned node_count
NodeNodeGraph< elem_node_type, sparse_graph_type, ElemNodeCount >::ElemGraphType elem_graph_type
const vector_type solution
base_type::execution_space execution_space
#define KOKKOS_INLINE_FUNCTION
ElemNodeIdView::execution_space execution_space
KOKKOS_INLINE_FUNCTION void operator()(const unsigned irow, unsigned &update, const bool final) const
ElementComputationBase< ExecutionSpace, Order, CoordinateMap, ScalarType > base_type
KOKKOS_INLINE_FUNCTION void init(unsigned &update) const
double sort_graph_entries
double fill_graph_entries
ElementComputationBase(const ElementComputationBase &rhs)
KOKKOS_INLINE_FUNCTION void fill_elem_graph_map(const unsigned ielem) const
Kokkos::Example::HexElement_Data< mesh_type::ElemNode > element_data_type
CrsGraphType::row_map_type::non_const_type RowMapType
KOKKOS_INLINE_FUNCTION void computeElementResidual(const scalar_type dof_values[], const double x[], const double y[], const double z[], fad_scalar_type elem_res[]) const
const element_data_type elem_data
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)
ElementComputation(const ElementComputation &rhs)
KOKKOS_INLINE_FUNCTION void operator()(const TagFillNodeSet &, unsigned ielem, unsigned &count) const
CrsMatrix(const StaticCrsGraphType &arg_graph)
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
Do not initialize the derivative array.
KOKKOS_INLINE_FUNCTION void operator()(const unsigned ielem) const
KOKKOS_INLINE_FUNCTION void operator()(const unsigned iwork) const
KOKKOS_INLINE_FUNCTION void operator()(const unsigned ielem) const
static const unsigned ElemNodeCount
const ElemNodeIdView elem_node_id
static const unsigned integration_count
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
mesh_type::node_coord_type node_coord_type
const elem_graph_type elem_graph
ElementComputation(const ElementComputation &rhs)
Kokkos::Example::BoxElemFixture< ExecutionSpace, Order, CoordinateMap > mesh_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
base_type::execution_space execution_space
Sacado::Fad::SFad< scalar_type, FunctionCount > fad_scalar_type
KOKKOS_INLINE_FUNCTION void init(const TagFillNodeSet &, unsigned &update) const
View< ValueType *, Space > coeff_type
sparse_matrix_type::StaticCrsGraphType sparse_graph_type
Kokkos::View< unsigned, execution_space > UnsignedValue
ElementComputation(const typename base_type::mesh_type &arg_mesh, 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)
ElementComputation< Kokkos::Example::BoxElemFixture< ExecutionSpace, Order, CoordinateMap >, CrsMatrix< ScalarType, ExecutionSpace >, FadElement > base_type
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(const typename base_type::mesh_type &arg_mesh, 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)
ElementComputation(const ElementComputation &rhs)
KOKKOS_INLINE_FUNCTION void computeElementResidual(const fad_scalar_type dof_values[], const double x[], const double y[], const double z[], fad_scalar_type elem_res[]) const
ElementComputation(const ElementComputation &rhs)
ElementComputation(const typename base_type::mesh_type &arg_mesh, 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)
KOKKOS_INLINE_FUNCTION void operator()(const unsigned ielem) const
base_type::scalar_type scalar_type
KOKKOS_INLINE_FUNCTION void sort_graph_entries(const unsigned irow) const
base_type::scalar_type scalar_type
double fill_element_graph
Kokkos::View< scalar_type *[FunctionCount][FunctionCount], execution_space > elem_matrices_type
Sacado::Fad::SFad< scalar_type, FunctionCount > fad_scalar_type
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
KOKKOS_INLINE_FUNCTION void join(const TagFillNodeSet &, volatile unsigned &update, volatile const unsigned &input) const
Kokkos::View< unsigned *[ElemNode][ElemNode], execution_space > ElemGraphType
static const unsigned FunctionCount
Generate a distributed unstructured finite element mesh from a partitioned NX*NY*NZ box of elements...
KOKKOS_INLINE_FUNCTION void scatterResidual(const unsigned ielem, const unsigned node_index[], fad_scalar_type res[]) const
base_type::scalar_type scalar_type
Kokkos::UnorderedMap< key_type, void, execution_space > SetType
pair< unsigned, unsigned > key_type
static const unsigned spatial_dimension