17 #ifndef KOKKOS_EXAMPLE_FENLFUNCTORS_HPP
18 #define KOKKOS_EXAMPLE_FENLFUNCTORS_HPP
29 #include <Kokkos_Core.hpp>
30 #include <Kokkos_Pair.hpp>
31 #include <Kokkos_UnorderedMap.hpp>
32 #include <Kokkos_StaticCrsGraph.hpp>
34 #include <Kokkos_Timer.hpp>
48 template<
typename ValueType ,
class Space >
50 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE // Don't remove this until Kokkos has removed the deprecated code path probably around September 2018
51 typedef Kokkos::StaticCrsGraph< unsigned , Space , void , unsigned >
StaticCrsGraphType ;
64 ,
coeff(
"crs_matrix_coeff" , arg_graph.entries.extent(0) )
68 template<
class ElemNodeIdView ,
class CrsGraphType ,
unsigned ElemNode >
75 typedef Kokkos::UnorderedMap< key_type, void , execution_space >
SetType ;
76 typedef typename CrsGraphType::row_map_type::non_const_type
RowMapType ;
80 typedef Kokkos::View< unsigned*[ElemNode][ElemNode] , execution_space >
121 const unsigned arg_node_count,
137 Kokkos::Timer wall_clock ;
143 size_t set_capacity = (28ull *
node_count) / 2;
144 unsigned failed_insert_count = 0 ;
155 Kokkos::parallel_reduce( Kokkos::RangePolicy<execution_space,TagFillNodeSet>(0,
elem_node_id.extent(0))
157 , failed_insert_count );
159 }
while ( failed_insert_count );
177 unsigned graph_entry_count = 0 ;
179 Kokkos::deep_copy( graph_entry_count ,
row_total );
183 graph.entries =
typename CrsGraphType::entries_type(
"graph_entries" , graph_entry_count );
229 KOKKOS_INLINE_FUNCTION
233 for (
unsigned row_local_node = 0 ; row_local_node <
elem_node_id.extent(1) ; ++row_local_node ) {
235 const unsigned row_node =
elem_node_id( ielem , row_local_node );
237 for (
unsigned col_local_node = row_local_node ; col_local_node <
elem_node_id.extent(1) ; ++col_local_node ) {
239 const unsigned col_node =
elem_node_id( ielem , col_local_node );
245 const key_type key = (row_node < col_node) ? make_pair( row_node, col_node ) : make_pair( col_node, row_node ) ;
247 const typename SetType::insert_result result =
node_node_set.insert( key );
250 if ( result.success() ) {
253 if ( row_node <
row_count.extent(0) ) { Kokkos::atomic_inc( &
row_count( row_node ) ); }
256 if ( col_node <
row_count.extent(0) && col_node != row_node ) { Kokkos::atomic_inc( &
row_count( col_node ) ); }
258 else if ( result.failed() ) {
266 KOKKOS_INLINE_FUNCTION
269 typedef typename std::remove_reference< decltype( row_count(0) ) >::type atomic_incr_type;
275 const unsigned row_node = key.first ;
276 const unsigned col_node = key.second ;
279 const unsigned offset =
graph.row_map( row_node ) + Kokkos::atomic_fetch_add( &
row_count( row_node ) , atomic_incr_type(1) );
280 graph.entries( offset ) = col_node ;
283 if ( col_node <
row_count.extent(0) && col_node != row_node ) {
284 const unsigned offset =
graph.row_map( col_node ) + Kokkos::atomic_fetch_add( &
row_count( col_node ) , atomic_incr_type(1) );
285 graph.entries( offset ) = row_node ;
290 KOKKOS_INLINE_FUNCTION
293 const unsigned row_beg =
graph.row_map( irow );
294 const unsigned row_end =
graph.row_map( irow + 1 );
295 for (
unsigned i = row_beg + 1 ;
i < row_end ; ++
i ) {
296 const unsigned col =
graph.entries(
i);
298 for ( ; row_beg < j && col <
graph.entries(j-1) ; --j ) {
301 graph.entries(j) = col ;
305 KOKKOS_INLINE_FUNCTION
308 for (
unsigned row_local_node = 0 ; row_local_node <
elem_node_id.extent(1) ; ++row_local_node ) {
310 const unsigned row_node =
elem_node_id( ielem , row_local_node );
312 for (
unsigned col_local_node = 0 ; col_local_node <
elem_node_id.extent(1) ; ++col_local_node ) {
314 const unsigned col_node =
elem_node_id( ielem , col_local_node );
316 unsigned entry = ~0u ;
318 if ( row_node + 1 <
graph.row_map.extent(0) ) {
320 const unsigned entry_end =
graph.row_map( row_node + 1 );
322 entry =
graph.row_map( row_node );
324 for ( ; entry < entry_end &&
graph.entries(entry) != col_node ; ++entry );
326 if ( entry == entry_end ) entry = ~0u ;
329 elem_graph( ielem , row_local_node , col_local_node ) = entry ;
334 KOKKOS_INLINE_FUNCTION
358 KOKKOS_INLINE_FUNCTION
359 void operator()(
const unsigned irow ,
unsigned & update ,
const bool final )
const
362 if (
final ) {
row_map( irow ) = update ; }
375 KOKKOS_INLINE_FUNCTION
378 KOKKOS_INLINE_FUNCTION
381 ,
const unsigned & input )
const { update += input ; }
384 KOKKOS_INLINE_FUNCTION
385 void init(
unsigned & update )
const { update = 0 ; }
387 KOKKOS_INLINE_FUNCTION
389 ,
const unsigned & input )
const { update += input ; }
406 class CoordinateMap ,
typename ScalarType >
421 typedef Kokkos::View< scalar_type* , Kokkos::LayoutLeft, execution_space >
vector_type ;
435 typedef Kokkos::View< scalar_type*[FunctionCount][FunctionCount] , execution_space >
elem_matrices_type ;
486 KOKKOS_INLINE_FUNCTION
494 double dpsidz[] )
const
496 enum { j11 = 0 , j12 = 1 , j13 = 2 ,
497 j21 = 3 , j22 = 4 , j23 = 5 ,
498 j31 = 6 , j32 = 7 , j33 = 8 };
502 double J[
TensorDim ] = { 0, 0, 0, 0, 0, 0, 0, 0, 0 };
505 const double x1 = x[
i] ;
506 const double x2 = y[
i] ;
507 const double x3 = z[
i] ;
509 const double g1 = grad[0][
i] ;
510 const double g2 = grad[1][
i] ;
511 const double g3 = grad[2][
i] ;
529 static_cast<double>( J[j22] * J[j33] - J[j23] * J[j32] ) ,
530 static_cast<double>( J[j13] * J[j32] - J[j12] * J[j33] ) ,
531 static_cast<double>( J[j12] * J[j23] - J[j13] * J[j22] ) ,
533 static_cast<double>( J[j23] * J[j31] - J[j21] * J[j33] ) ,
534 static_cast<double>( J[j11] * J[j33] - J[j13] * J[j31] ) ,
535 static_cast<double>( J[j13] * J[j21] - J[j11] * J[j23] ) ,
537 static_cast<double>( J[j21] * J[j32] - J[j22] * J[j31] ) ,
538 static_cast<double>( J[j12] * J[j31] - J[j11] * J[j32] ) ,
539 static_cast<double>( J[j11] * J[j22] - J[j12] * J[j21] ) };
541 const double detJ = J[j11] * invJ[j11] +
545 const double detJinv = 1.0 / detJ ;
547 for (
unsigned i = 0 ;
i <
TensorDim ; ++
i ) { invJ[
i] *= detJinv ; }
552 const double g0 = grad[0][
i];
553 const double g1 = grad[1][
i];
554 const double g2 = grad[2][
i];
556 dpsidx[
i] = g0 * invJ[j11] + g1 * invJ[j12] + g2 * invJ[j13];
557 dpsidy[
i] = g0 * invJ[j21] + g1 * invJ[j22] + g2 * invJ[j23];
558 dpsidz[
i] = g0 * invJ[j31] + g1 * invJ[j32] + g2 * invJ[j33];
573 template<
class FiniteElementMeshType ,
574 class SparseMatrixType ,
580 class CoordinateMap ,
typename ScalarType >
583 CrsMatrix< ScalarType , ExecutionSpace > ,
607 base_type(arg_mesh, arg_solution, arg_elem_graph,
608 arg_jacobian, arg_residual) {}
615 parallel_for( nelem , *
this );
618 KOKKOS_INLINE_FUNCTION
621 unsigned node_index[],
622 double x[],
double y[],
double z[],
644 KOKKOS_INLINE_FUNCTION
646 const unsigned node_index[],
651 const unsigned row = node_index[
i] ;
652 if ( row < this->
residual.extent(0) ) {
653 Kokkos::atomic_add( & this->
residual( row ) , res[
i] );
656 const unsigned entry = this->
elem_graph( ielem , i , j );
657 if ( entry != ~0u ) {
658 Kokkos::atomic_add( & this->
jacobian.
coeff( entry ) , mat[
i][j] );
665 KOKKOS_INLINE_FUNCTION
674 double coeff_k = 3.456;
675 double coeff_src = 1.234;
676 double advection[] = { 1.1, 1.2, 1.3 };
687 dpsidx , dpsidy , dpsidz );
688 const double detJ_weight = detJ * integ_weight;
689 const double detJ_weight_coeff_k = detJ_weight * coeff_k;
696 value_at_pt += dof_values[m] * bases_vals[m] ;
697 gradx_at_pt += dof_values[m] * dpsidx[m] ;
698 grady_at_pt += dof_values[m] * dpsidy[m] ;
699 gradz_at_pt += dof_values[m] * dpsidz[m] ;
703 coeff_src * value_at_pt * value_at_pt ;
705 2.0 * coeff_src * value_at_pt ;
712 advection_x*gradx_at_pt +
713 advection_y*grady_at_pt +
714 advection_z*gradz_at_pt ;
718 const double bases_val_m = bases_vals[m] * detJ_weight ;
719 const double dpsidx_m = dpsidx[m] ;
720 const double dpsidy_m = dpsidy[m] ;
721 const double dpsidz_m = dpsidz[m] ;
724 detJ_weight_coeff_k * ( dpsidx_m * gradx_at_pt +
725 dpsidy_m * grady_at_pt +
726 dpsidz_m * gradz_at_pt ) +
727 bases_val_m * ( advection_term + source_term ) ;
730 const double dpsidx_n = dpsidx[
n] ;
731 const double dpsidy_n = dpsidy[
n] ;
732 const double dpsidz_n = dpsidz[
n] ;
734 detJ_weight_coeff_k * ( dpsidx_m * dpsidx_n +
735 dpsidy_m * dpsidy_n +
736 dpsidz_m * dpsidz_n ) +
737 bases_val_m * ( advection_x * dpsidx_n +
738 advection_y * dpsidy_n +
739 advection_z * dpsidz_n +
740 source_deriv * bases_vals[
n] ) ;
746 KOKKOS_INLINE_FUNCTION
759 gatherSolution(ielem, val, node_index, x, y, z, elem_res, elem_mat);
762 computeElementResidualJacobian( val, x, y, z, elem_res , elem_mat );
765 scatterResidual( ielem, node_index, elem_res, elem_mat );
770 class CoordinateMap ,
typename ScalarType >
773 CrsMatrix< ScalarType , ExecutionSpace > ,
798 base_type(arg_mesh, arg_solution, arg_elem_graph,
799 arg_jacobian, arg_residual) {}
806 parallel_for( nelem , *
this );
809 KOKKOS_INLINE_FUNCTION
812 unsigned node_index[],
813 double x[],
double y[],
double z[],
830 KOKKOS_INLINE_FUNCTION
832 const unsigned node_index[],
836 const unsigned row = node_index[
i] ;
837 if ( row < this->
residual.extent(0) ) {
838 Kokkos::atomic_add( & this->
residual( row ) , res[
i].
val() );
841 const unsigned entry = this->
elem_graph( ielem ,
i , j );
842 if ( entry != ~0u ) {
844 res[
i].fastAccessDx(j) );
851 KOKKOS_INLINE_FUNCTION
858 double coeff_k = 3.456;
859 double coeff_src = 1.234;
860 double advection[] = { 1.1, 1.2, 1.3 };
871 dpsidx , dpsidy , dpsidz );
872 const double detJ_weight = detJ * integ_weight;
873 const double detJ_weight_coeff_k = detJ_weight * coeff_k;
880 value_at_pt += dof_values[m] * bases_vals[m] ;
881 gradx_at_pt += dof_values[m] * dpsidx[m] ;
882 grady_at_pt += dof_values[m] * dpsidy[m] ;
883 gradz_at_pt += dof_values[m] * dpsidz[m] ;
887 coeff_src * value_at_pt * value_at_pt ;
890 advection[0]*gradx_at_pt +
891 advection[1]*grady_at_pt +
892 advection[2]*gradz_at_pt;
895 const double bases_val_m = bases_vals[m] * detJ_weight ;
896 const double dpsidx_m = dpsidx[m] ;
897 const double dpsidy_m = dpsidy[m] ;
898 const double dpsidz_m = dpsidz[m] ;
901 detJ_weight_coeff_k * ( dpsidx_m * gradx_at_pt +
902 dpsidy_m * grady_at_pt +
903 dpsidz_m * gradz_at_pt ) +
904 bases_val_m * ( advection_term + source_term ) ;
909 KOKKOS_INLINE_FUNCTION
921 gatherSolution( ielem, val, node_index, x, y, z, elem_res );
924 computeElementResidual( val, x, y, z, elem_res );
927 scatterResidual( ielem, node_index, elem_res );
932 class CoordinateMap ,
typename ScalarType >
935 CrsMatrix< ScalarType , ExecutionSpace > ,
937 public ElementComputation< Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
938 CrsMatrix< ScalarType , ExecutionSpace > ,
963 base_type(arg_mesh, arg_solution, arg_elem_graph,
964 arg_jacobian, arg_residual) {}
971 parallel_for( nelem , *
this );
974 KOKKOS_INLINE_FUNCTION
977 unsigned node_index[],
978 double x[],
double y[],
double z[],
994 KOKKOS_INLINE_FUNCTION
1001 double coeff_k = 3.456;
1002 double coeff_src = 1.234;
1003 double advection[] = { 1.1, 1.2, 1.3 };
1014 dpsidx , dpsidy , dpsidz );
1015 const double detJ_weight = detJ * integ_weight;
1016 const double detJ_weight_coeff_k = detJ_weight * coeff_k;
1023 value_at_pt.val() += dof_values[m] * bases_vals[m] ;
1024 value_at_pt.fastAccessDx(m) = bases_vals[m] ;
1026 gradx_at_pt.val() += dof_values[m] * dpsidx[m] ;
1027 gradx_at_pt.fastAccessDx(m) = dpsidx[m] ;
1029 grady_at_pt.val() += dof_values[m] * dpsidy[m] ;
1030 grady_at_pt.fastAccessDx(m) = dpsidy[m] ;
1032 gradz_at_pt.val() += dof_values[m] * dpsidz[m] ;
1033 gradz_at_pt.fastAccessDx(m) = dpsidz[m] ;
1037 coeff_src * value_at_pt * value_at_pt ;
1040 advection[0]*gradx_at_pt +
1041 advection[1]*grady_at_pt +
1042 advection[2]*gradz_at_pt;
1045 const double bases_val_m = bases_vals[m] * detJ_weight ;
1046 const double dpsidx_m = dpsidx[m] ;
1047 const double dpsidy_m = dpsidy[m] ;
1048 const double dpsidz_m = dpsidz[m] ;
1051 detJ_weight_coeff_k * ( dpsidx_m * gradx_at_pt +
1052 dpsidy_m * grady_at_pt +
1053 dpsidz_m * gradz_at_pt ) +
1054 bases_val_m * ( advection_term + source_term ) ;
1059 KOKKOS_INLINE_FUNCTION
1071 gatherSolution( ielem, val, node_index, x, y, z, elem_res );
1074 computeElementResidual( val, x, y, z, elem_res );
1077 this->scatterResidual( ielem, node_index, elem_res );
1082 class CoordinateMap ,
typename ScalarType >
1085 CrsMatrix< ScalarType , ExecutionSpace > ,
1087 public ElementComputation< Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
1088 CrsMatrix< ScalarType , ExecutionSpace > ,
1113 base_type(arg_mesh, arg_solution, arg_elem_graph,
1114 arg_jacobian, arg_residual) {}
1121 parallel_for( nelem , *
this );
1124 KOKKOS_INLINE_FUNCTION
1133 double coeff_k = 3.456;
1134 double coeff_src = 1.234;
1135 double advection[] = { 1.1, 1.2, 1.3 };
1151 dpsidx , dpsidy , dpsidz );
1152 const double detJ_weight = detJ * integ_weight;
1153 const double detJ_weight_coeff_k = detJ_weight * coeff_k;
1155 value_at_pt.val() = 0.0 ;
1156 gradx_at_pt.val() = 0.0 ;
1157 grady_at_pt.val() = 0.0 ;
1158 gradz_at_pt.val() = 0.0 ;
1160 value_at_pt.val() += dof_values[m] * bases_vals[m] ;
1161 gradx_at_pt.val() += dof_values[m] * dpsidx[m] ;
1162 grady_at_pt.val() += dof_values[m] * dpsidy[m] ;
1163 gradz_at_pt.val() += dof_values[m] * dpsidz[m] ;
1167 coeff_src * value_at_pt * value_at_pt ;
1170 advection[0]*gradx_at_pt +
1171 advection[1]*grady_at_pt +
1172 advection[2]*gradz_at_pt;
1175 const double bases_val_m = bases_vals[m] * detJ_weight ;
1177 detJ_weight_coeff_k * ( dpsidx[m] * gradx_at_pt +
1178 dpsidy[m] * grady_at_pt +
1179 dpsidz[m] * gradz_at_pt ) +
1180 bases_val_m * ( advection_term + source_term ) ;
1182 elem_res[m] += res.val();
1186 mat[
n] += res.fastAccessDx(0) * bases_vals[
n] +
1187 res.fastAccessDx(1) * dpsidx[
n] +
1188 res.fastAccessDx(2) * dpsidy[
n] +
1189 res.fastAccessDx(3) * dpsidz[
n];
1195 KOKKOS_INLINE_FUNCTION
1208 this->gatherSolution( ielem, val, node_index, x, y, z, elem_res, elem_mat );
1211 computeElementResidualJacobian( val, x, y, z, elem_res, elem_mat );
1214 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
KOKKOS_INLINE_FUNCTION void join(unsigned &update, const unsigned &input) const
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
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
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 join(const TagFillNodeSet &, unsigned &update, const unsigned &input) const
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::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