Sacado Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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_Core.hpp>
30 #include <Kokkos_Pair.hpp>
31 #include <Kokkos_UnorderedMap.hpp>
32 #include <Kokkos_StaticCrsGraph.hpp>
33 
34 #include <Kokkos_Timer.hpp>
35 
36 #include <BoxElemFixture.hpp>
37 #include <HexElement.hpp>
38 
39 #include "Sacado.hpp"
40 
41 //----------------------------------------------------------------------------
42 //----------------------------------------------------------------------------
43 
44 namespace Kokkos {
45 namespace Example {
46 namespace FENL {
47 
48 template< typename ValueType , class Space >
49 struct CrsMatrix {
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 ;
52 #else
53  typedef Kokkos::StaticCrsGraph< unsigned , Space , void , void , unsigned > StaticCrsGraphType ;
54 #endif
55  typedef View< ValueType * , Space > coeff_type ;
56 
59 
60  CrsMatrix() : graph(), coeff() {}
61 
62  CrsMatrix( const StaticCrsGraphType & arg_graph )
63  : graph( arg_graph )
64  , coeff( "crs_matrix_coeff" , arg_graph.entries.extent(0) )
65  {}
66 };
67 
68 template< class ElemNodeIdView , class CrsGraphType , unsigned ElemNode >
70 public:
71 
72  typedef typename ElemNodeIdView::execution_space execution_space ;
73  typedef pair<unsigned,unsigned> key_type ;
74 
75  typedef Kokkos::UnorderedMap< key_type, void , execution_space > SetType ;
76  typedef typename CrsGraphType::row_map_type::non_const_type RowMapType ;
77  typedef Kokkos::View< unsigned , execution_space > UnsignedValue ;
78 
79  // Static dimensions of 0 generate compiler warnings or errors.
80  typedef Kokkos::View< unsigned*[ElemNode][ElemNode] , execution_space >
82 
83  struct TagFillNodeSet {};
84  struct TagScanNodeCount {};
88 
89 private:
90 
96 
97  const unsigned node_count ;
98  const ElemNodeIdView elem_node_id ;
104 
105 public:
106 
107  CrsGraphType graph ;
109 
110  struct Times
111  {
112  double ratio;
118  };
119 
120  NodeNodeGraph( const ElemNodeIdView & arg_elem_node_id ,
121  const unsigned arg_node_count,
122  Times & results
123  )
124  : node_count(arg_node_count)
125  , elem_node_id( arg_elem_node_id )
126  , row_total( "row_total" )
127  , row_count(Kokkos::ViewAllocateWithoutInitializing("row_count") , node_count ) // will deep_copy to 0 inside loop
128  , row_map( "graph_row_map" , node_count + 1 )
129  , node_node_set()
130  , phase( FILL_NODE_SET )
131  , graph()
132  , elem_graph()
133  {
134  //--------------------------------
135  // Guess at capacity required for the map:
136 
137  Kokkos::Timer wall_clock ;
138 
139  wall_clock.reset();
140  phase = FILL_NODE_SET ;
141 
142  // upper bound on the capacity
143  size_t set_capacity = (28ull * node_count) / 2;
144  unsigned failed_insert_count = 0 ;
145 
146  do {
147  // Zero the row count to restart the fill
148  Kokkos::deep_copy( row_count , 0u );
149 
150  node_node_set = SetType( ( set_capacity += failed_insert_count ) );
151 
152  // May be larger that requested:
153  set_capacity = node_node_set.capacity();
154 
155  Kokkos::parallel_reduce( Kokkos::RangePolicy<execution_space,TagFillNodeSet>(0,elem_node_id.extent(0))
156  , *this
157  , failed_insert_count );
158 
159  } while ( failed_insert_count );
160 
161  execution_space().fence();
162  results.ratio = (double)node_node_set.size() / (double)node_node_set.capacity();
163  results.fill_node_set = wall_clock.seconds();
164  //--------------------------------
165 
166  wall_clock.reset();
168 
169  // Exclusive scan of row_count into row_map
170  // including the final total in the 'node_count + 1' position.
171  // Zero the 'row_count' values.
172  Kokkos::parallel_scan( node_count , *this );
173 
174  // Zero the row count for the fill:
175  Kokkos::deep_copy( row_count , 0u );
176 
177  unsigned graph_entry_count = 0 ;
178 
179  Kokkos::deep_copy( graph_entry_count , row_total );
180 
181  // Assign graph's row_map and allocate graph's entries
182  graph.row_map = row_map ;
183  graph.entries = typename CrsGraphType::entries_type( "graph_entries" , graph_entry_count );
184 
185  //--------------------------------
186  // Fill graph's entries from the (node,node) set.
187 
188  execution_space().fence();
189  results.scan_node_count = wall_clock.seconds();
190 
191  wall_clock.reset();
193  Kokkos::parallel_for( node_node_set.capacity() , *this );
194 
195  execution_space().fence();
196  results.fill_graph_entries = wall_clock.seconds();
197 
198  //--------------------------------
199  // Done with the temporary sets and arrays
200  wall_clock.reset();
202 
204  row_count = RowMapType();
205  row_map = RowMapType();
206  node_node_set.clear();
207 
208  //--------------------------------
209 
210  Kokkos::parallel_for( node_count , *this );
211 
212  execution_space().fence();
213  results.sort_graph_entries = wall_clock.seconds();
214 
215  //--------------------------------
216  // Element-to-graph mapping:
217  wall_clock.reset();
219  elem_graph = ElemGraphType("elem_graph", elem_node_id.extent(0) );
220  Kokkos::parallel_for( elem_node_id.extent(0) , *this );
221 
222  execution_space().fence();
223  results.fill_element_graph = wall_clock.seconds();
224  }
225 
226  //------------------------------------
227  // parallel_for: create map and count row length
228 
229  KOKKOS_INLINE_FUNCTION
230  void operator()( const TagFillNodeSet & , unsigned ielem , unsigned & count ) const
231  {
232  // Loop over element's (row_local_node,col_local_node) pairs:
233  for ( unsigned row_local_node = 0 ; row_local_node < elem_node_id.extent(1) ; ++row_local_node ) {
234 
235  const unsigned row_node = elem_node_id( ielem , row_local_node );
236 
237  for ( unsigned col_local_node = row_local_node ; col_local_node < elem_node_id.extent(1) ; ++col_local_node ) {
238 
239  const unsigned col_node = elem_node_id( ielem , col_local_node );
240 
241  // If either node is locally owned then insert the pair into the unordered map:
242 
243  if ( row_node < row_count.extent(0) || col_node < row_count.extent(0) ) {
244 
245  const key_type key = (row_node < col_node) ? make_pair( row_node, col_node ) : make_pair( col_node, row_node ) ;
246 
247  const typename SetType::insert_result result = node_node_set.insert( key );
248 
249  // A successfull insert: the first time this pair was added
250  if ( result.success() ) {
251 
252  // If row node is owned then increment count
253  if ( row_node < row_count.extent(0) ) { Kokkos::atomic_inc( & row_count( row_node ) ); }
254 
255  // If column node is owned and not equal to row node then increment count
256  if ( col_node < row_count.extent(0) && col_node != row_node ) { Kokkos::atomic_inc( & row_count( col_node ) ); }
257  }
258  else if ( result.failed() ) {
259  ++count ;
260  }
261  }
262  }
263  }
264  }
265 
266  KOKKOS_INLINE_FUNCTION
267  void fill_graph_entries( const unsigned iset ) const
268  {
269  typedef typename std::remove_reference< decltype( row_count(0) ) >::type atomic_incr_type;
270 
271  if ( node_node_set.valid_at(iset) ) {
272  // Add each entry to the graph entries.
273 
274  const key_type key = node_node_set.key_at(iset) ;
275  const unsigned row_node = key.first ;
276  const unsigned col_node = key.second ;
277 
278  if ( row_node < row_count.extent(0) ) {
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 ;
281  }
282 
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 ;
286  }
287  }
288  }
289 
290  KOKKOS_INLINE_FUNCTION
291  void sort_graph_entries( const unsigned irow ) const
292  {
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);
297  unsigned j = i ;
298  for ( ; row_beg < j && col < graph.entries(j-1) ; --j ) {
299  graph.entries(j) = graph.entries(j-1);
300  }
301  graph.entries(j) = col ;
302  }
303  }
304 
305  KOKKOS_INLINE_FUNCTION
306  void fill_elem_graph_map( const unsigned ielem ) const
307  {
308  for ( unsigned row_local_node = 0 ; row_local_node < elem_node_id.extent(1) ; ++row_local_node ) {
309 
310  const unsigned row_node = elem_node_id( ielem , row_local_node );
311 
312  for ( unsigned col_local_node = 0 ; col_local_node < elem_node_id.extent(1) ; ++col_local_node ) {
313 
314  const unsigned col_node = elem_node_id( ielem , col_local_node );
315 
316  unsigned entry = ~0u ;
317 
318  if ( row_node + 1 < graph.row_map.extent(0) ) {
319 
320  const unsigned entry_end = graph.row_map( row_node + 1 );
321 
322  entry = graph.row_map( row_node );
323 
324  for ( ; entry < entry_end && graph.entries(entry) != col_node ; ++entry );
325 
326  if ( entry == entry_end ) entry = ~0u ;
327  }
328 
329  elem_graph( ielem , row_local_node , col_local_node ) = entry ;
330  }
331  }
332  }
333 
334  KOKKOS_INLINE_FUNCTION
335  void operator()( const unsigned iwork ) const
336  {
337 /*
338  if ( phase == FILL_NODE_SET ) {
339  operator()( TagFillNodeSet() , iwork );
340  }
341  else */
342  if ( phase == FILL_GRAPH_ENTRIES ) {
343  fill_graph_entries( iwork );
344  }
345  else if ( phase == SORT_GRAPH_ENTRIES ) {
346  sort_graph_entries( iwork );
347  }
348  else if ( phase == FILL_ELEMENT_GRAPH ) {
349  fill_elem_graph_map( iwork );
350  }
351  }
352 
353  //------------------------------------
354  // parallel_scan: row offsets
355 
356  typedef unsigned value_type ;
357 
358  KOKKOS_INLINE_FUNCTION
359  void operator()( const unsigned irow , unsigned & update , const bool final ) const
360  {
361  // exclusive scan
362  if ( final ) { row_map( irow ) = update ; }
363 
364  update += row_count( irow );
365 
366  if ( final ) {
367  if ( irow + 1 == row_count.extent(0) ) {
368  row_map( irow + 1 ) = update ;
369  row_total() = update ;
370  }
371  }
372  }
373 
374  // For the reduce phase:
375  KOKKOS_INLINE_FUNCTION
376  void init( const TagFillNodeSet & , unsigned & update ) const { update = 0 ; }
377 
378  KOKKOS_INLINE_FUNCTION
379  void join( const TagFillNodeSet &
380  , unsigned & update
381  , const unsigned & input ) const { update += input ; }
382 
383  // For the scan phase::
384  KOKKOS_INLINE_FUNCTION
385  void init( unsigned & update ) const { update = 0 ; }
386 
387  KOKKOS_INLINE_FUNCTION
388  void join( unsigned & update
389  , const unsigned & input ) const { update += input ; }
390 
391  //------------------------------------
392 };
393 
394 } /* namespace FENL */
395 } /* namespace Example */
396 } /* namespace Kokkos */
397 
398 //----------------------------------------------------------------------------
399 //----------------------------------------------------------------------------
400 
401 namespace Kokkos {
402 namespace Example {
403 namespace FENL {
404 
405 template< class ExecutionSpace , BoxElemPart::ElemOrder Order ,
406  class CoordinateMap , typename ScalarType >
408 {
409 public:
410 
413 
414  //------------------------------------
415 
416  typedef ExecutionSpace execution_space ;
417  typedef ScalarType scalar_type ;
418 
421  typedef Kokkos::View< scalar_type* , Kokkos::LayoutLeft, execution_space > vector_type ;
422 
423  //------------------------------------
424 
426  static const unsigned TensorDim = SpatialDim * SpatialDim ;
430 
431  //------------------------------------
432 
435  typedef Kokkos::View< scalar_type*[FunctionCount][FunctionCount] , execution_space > elem_matrices_type ;
436  typedef Kokkos::View< scalar_type*[FunctionCount] , execution_space > elem_vectors_type ;
437 
439 
440  //------------------------------------
441 
442 
443  //------------------------------------
444  // Computational data:
445 
455 
457  : elem_data()
459  , node_coords( rhs.node_coords )
460  , elem_graph( rhs.elem_graph )
463  , solution( rhs.solution )
464  , residual( rhs.residual )
465  , jacobian( rhs.jacobian )
466  {}
467 
468  ElementComputationBase( const mesh_type & arg_mesh ,
469  const vector_type & arg_solution ,
470  const elem_graph_type & arg_elem_graph ,
471  const sparse_matrix_type & arg_jacobian ,
472  const vector_type & arg_residual )
473  : elem_data()
474  , elem_node_ids( arg_mesh.elem_node() )
475  , node_coords( arg_mesh.node_coord() )
476  , elem_graph( arg_elem_graph )
477  , elem_jacobians()
478  , elem_residuals()
479  , solution( arg_solution )
480  , residual( arg_residual )
481  , jacobian( arg_jacobian )
482  {}
483 
484  //------------------------------------
485 
486  KOKKOS_INLINE_FUNCTION
488  const double grad[][ FunctionCount ] , // Gradient of bases master element
489  const double x[] ,
490  const double y[] ,
491  const double z[] ,
492  double dpsidx[] ,
493  double dpsidy[] ,
494  double dpsidz[] ) const
495  {
496  enum { j11 = 0 , j12 = 1 , j13 = 2 ,
497  j21 = 3 , j22 = 4 , j23 = 5 ,
498  j31 = 6 , j32 = 7 , j33 = 8 };
499 
500  // Jacobian accumulation:
501 
502  double J[ TensorDim ] = { 0, 0, 0, 0, 0, 0, 0, 0, 0 };
503 
504  for( unsigned i = 0; i < FunctionCount ; ++i ) {
505  const double x1 = x[i] ;
506  const double x2 = y[i] ;
507  const double x3 = z[i] ;
508 
509  const double g1 = grad[0][i] ;
510  const double g2 = grad[1][i] ;
511  const double g3 = grad[2][i] ;
512 
513  J[j11] += g1 * x1 ;
514  J[j12] += g1 * x2 ;
515  J[j13] += g1 * x3 ;
516 
517  J[j21] += g2 * x1 ;
518  J[j22] += g2 * x2 ;
519  J[j23] += g2 * x3 ;
520 
521  J[j31] += g3 * x1 ;
522  J[j32] += g3 * x2 ;
523  J[j33] += g3 * x3 ;
524  }
525 
526  // Inverse jacobian:
527 
528  double invJ[ TensorDim ] = {
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] ) ,
532 
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] ) ,
536 
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] ) };
540 
541  const double detJ = J[j11] * invJ[j11] +
542  J[j21] * invJ[j12] +
543  J[j31] * invJ[j13] ;
544 
545  const double detJinv = 1.0 / detJ ;
546 
547  for ( unsigned i = 0 ; i < TensorDim ; ++i ) { invJ[i] *= detJinv ; }
548 
549  // Transform gradients:
550 
551  for( unsigned i = 0; i < FunctionCount ; ++i ) {
552  const double g0 = grad[0][i];
553  const double g1 = grad[1][i];
554  const double g2 = grad[2][i];
555 
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];
559  }
560 
561  return detJ ;
562  }
563 
564 };
565 
571 };
572 
573 template< class FiniteElementMeshType ,
574  class SparseMatrixType ,
576  >
578 
579 template< class ExecutionSpace , BoxElemPart::ElemOrder Order ,
580  class CoordinateMap , typename ScalarType >
582  < Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
583  CrsMatrix< ScalarType , ExecutionSpace > ,
584  Analytic > :
585  public ElementComputationBase<ExecutionSpace, Order, CoordinateMap,
586  ScalarType> {
587 public:
588 
589  typedef ElementComputationBase<ExecutionSpace, Order, CoordinateMap,
590  ScalarType> base_type;
591 
594 
595  static const unsigned FunctionCount = base_type::FunctionCount;
596  static const unsigned IntegrationCount = base_type::IntegrationCount;
597  static const unsigned ElemNodeCount = base_type::ElemNodeCount;
598 
600 
602  const typename base_type::mesh_type & arg_mesh ,
603  const typename base_type::vector_type & arg_solution ,
604  const typename base_type::elem_graph_type & arg_elem_graph ,
605  const typename base_type::sparse_matrix_type & arg_jacobian ,
606  const typename base_type::vector_type & arg_residual ) :
607  base_type(arg_mesh, arg_solution, arg_elem_graph,
608  arg_jacobian, arg_residual) {}
609 
610  //------------------------------------
611 
612  void apply() const
613  {
614  const size_t nelem = this->elem_node_ids.extent(0);
615  parallel_for( nelem , *this );
616  }
617 
618  KOKKOS_INLINE_FUNCTION
619  void gatherSolution(const unsigned ielem,
620  scalar_type val[],
621  unsigned node_index[],
622  double x[], double y[], double z[],
623  scalar_type res[],
624  scalar_type mat[][FunctionCount]) const
625  {
626  for ( unsigned i = 0 ; i < ElemNodeCount ; ++i ) {
627  const unsigned ni = this->elem_node_ids( ielem , i );
628 
629  node_index[i] = ni ;
630 
631  x[i] = this->node_coords( ni , 0 );
632  y[i] = this->node_coords( ni , 1 );
633  z[i] = this->node_coords( ni , 2 );
634 
635  val[i] = this->solution( ni ) ;
636  res[i] = 0 ;
637 
638  for( unsigned j = 0; j < FunctionCount ; j++){
639  mat[i][j] = 0 ;
640  }
641  }
642  }
643 
644  KOKKOS_INLINE_FUNCTION
645  void scatterResidual(const unsigned ielem,
646  const unsigned node_index[],
647  const scalar_type res[],
648  const scalar_type mat[][FunctionCount]) const
649  {
650  for( unsigned i = 0 ; i < FunctionCount ; i++ ) {
651  const unsigned row = node_index[i] ;
652  if ( row < this->residual.extent(0) ) {
653  Kokkos::atomic_add( & this->residual( row ) , res[i] );
654 
655  for( unsigned j = 0 ; j < FunctionCount ; j++ ) {
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] );
659  }
660  }
661  }
662  }
663  }
664 
665  KOKKOS_INLINE_FUNCTION
667  const scalar_type dof_values[] ,
668  const double x[],
669  const double y[],
670  const double z[],
671  scalar_type elem_res[] ,
672  scalar_type elem_mat[][FunctionCount] ) const
673  {
674  double coeff_k = 3.456;
675  double coeff_src = 1.234;
676  double advection[] = { 1.1, 1.2, 1.3 };
677  double dpsidx[ FunctionCount ] ;
678  double dpsidy[ FunctionCount ] ;
679  double dpsidz[ FunctionCount ] ;
680  for ( unsigned i = 0 ; i < IntegrationCount ; ++i ) {
681 
682  const double integ_weight = this->elem_data.weights[i];
683  const double* bases_vals = this->elem_data.values[i];
684  const double detJ =
685  this->transform_gradients( this->elem_data.gradients[i] ,
686  x , y , z ,
687  dpsidx , dpsidy , dpsidz );
688  const double detJ_weight = detJ * integ_weight;
689  const double detJ_weight_coeff_k = detJ_weight * coeff_k;
690 
691  scalar_type value_at_pt = 0 ;
692  scalar_type gradx_at_pt = 0 ;
693  scalar_type grady_at_pt = 0 ;
694  scalar_type gradz_at_pt = 0 ;
695  for ( unsigned m = 0 ; m < FunctionCount ; m++ ) {
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] ;
700  }
701 
702  const scalar_type source_term =
703  coeff_src * value_at_pt * value_at_pt ;
704  const scalar_type source_deriv =
705  2.0 * coeff_src * value_at_pt ;
706 
707  const scalar_type advection_x = advection[0];
708  const scalar_type advection_y = advection[1];
709  const scalar_type advection_z = advection[2];
710 
711  const scalar_type advection_term =
712  advection_x*gradx_at_pt +
713  advection_y*grady_at_pt +
714  advection_z*gradz_at_pt ;
715 
716  for ( unsigned m = 0; m < FunctionCount; ++m) {
717  scalar_type * const mat = elem_mat[m] ;
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] ;
722 
723  elem_res[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 ) ;
728 
729  for( unsigned n = 0; n < FunctionCount; n++) {
730  const double dpsidx_n = dpsidx[n] ;
731  const double dpsidy_n = dpsidy[n] ;
732  const double dpsidz_n = dpsidz[n] ;
733  mat[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] ) ;
741  }
742  }
743  }
744  }
745 
746  KOKKOS_INLINE_FUNCTION
747  void operator()( const unsigned ielem ) const
748  {
749  double x[ FunctionCount ] ;
750  double y[ FunctionCount ] ;
751  double z[ FunctionCount ] ;
752  unsigned node_index[ ElemNodeCount ];
753 
755  scalar_type elem_res[ FunctionCount ] ;
756  scalar_type elem_mat[ FunctionCount ][ FunctionCount ] ;
757 
758  // Gather nodal coordinates and solution vector:
759  gatherSolution(ielem, val, node_index, x, y, z, elem_res, elem_mat);
760 
761  // Compute nodal element residual vector and Jacobian matrix
762  computeElementResidualJacobian( val, x, y, z, elem_res , elem_mat );
763 
764  // Scatter nodal element residual and Jacobian in global vector and matrix:
765  scatterResidual( ielem, node_index, elem_res, elem_mat );
766  }
767 }; /* ElementComputation */
768 
769 template< class ExecutionSpace , BoxElemPart::ElemOrder Order ,
770  class CoordinateMap , typename ScalarType >
772  < Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
773  CrsMatrix< ScalarType , ExecutionSpace > ,
774  FadElement > : public ElementComputationBase<ExecutionSpace, Order, CoordinateMap,
775  ScalarType> {
776 public:
777 
778  typedef ElementComputationBase<ExecutionSpace, Order, CoordinateMap,
779  ScalarType> base_type;
780 
783 
784  static const unsigned FunctionCount = base_type::FunctionCount;
785  static const unsigned IntegrationCount = base_type::IntegrationCount;
786  static const unsigned ElemNodeCount = base_type::ElemNodeCount;
787 
789 
791 
793  const typename base_type::mesh_type & arg_mesh ,
794  const typename base_type::vector_type & arg_solution ,
795  const typename base_type::elem_graph_type & arg_elem_graph ,
796  const typename base_type::sparse_matrix_type & arg_jacobian ,
797  const typename base_type::vector_type & arg_residual ) :
798  base_type(arg_mesh, arg_solution, arg_elem_graph,
799  arg_jacobian, arg_residual) {}
800 
801  //------------------------------------
802 
803  void apply() const
804  {
805  const size_t nelem = this->elem_node_ids.extent(0);
806  parallel_for( nelem , *this );
807  }
808 
809  KOKKOS_INLINE_FUNCTION
810  void gatherSolution(const unsigned ielem,
812  unsigned node_index[],
813  double x[], double y[], double z[],
814  fad_scalar_type res[]) const
815  {
816  for ( unsigned i = 0 ; i < ElemNodeCount ; ++i ) {
817  const unsigned ni = this->elem_node_ids( ielem , i );
818 
819  node_index[i] = ni ;
820 
821  x[i] = this->node_coords( ni , 0 );
822  y[i] = this->node_coords( ni , 1 );
823  z[i] = this->node_coords( ni , 2 );
824 
825  val[i].val() = this->solution( ni );
826  val[i].diff( i, FunctionCount );
827  }
828  }
829 
830  KOKKOS_INLINE_FUNCTION
831  void scatterResidual(const unsigned ielem,
832  const unsigned node_index[],
833  fad_scalar_type res[]) const
834  {
835  for( unsigned i = 0 ; i < FunctionCount ; i++ ) {
836  const unsigned row = node_index[i] ;
837  if ( row < this->residual.extent(0) ) {
838  Kokkos::atomic_add( & this->residual( row ) , res[i].val() );
839 
840  for( unsigned j = 0 ; j < FunctionCount ; j++ ) {
841  const unsigned entry = this->elem_graph( ielem , i , j );
842  if ( entry != ~0u ) {
843  Kokkos::atomic_add( & this->jacobian.coeff( entry ) ,
844  res[i].fastAccessDx(j) );
845  }
846  }
847  }
848  }
849  }
850 
851  KOKKOS_INLINE_FUNCTION
852  void computeElementResidual(const fad_scalar_type dof_values[] ,
853  const double x[],
854  const double y[],
855  const double z[],
856  fad_scalar_type elem_res[] ) const
857  {
858  double coeff_k = 3.456;
859  double coeff_src = 1.234;
860  double advection[] = { 1.1, 1.2, 1.3 };
861  double dpsidx[ FunctionCount ] ;
862  double dpsidy[ FunctionCount ] ;
863  double dpsidz[ FunctionCount ] ;
864  for ( unsigned i = 0 ; i < IntegrationCount ; ++i ) {
865 
866  const double integ_weight = this->elem_data.weights[i];
867  const double* bases_vals = this->elem_data.values[i];
868  const double detJ =
869  this->transform_gradients( this->elem_data.gradients[i] ,
870  x , y , z ,
871  dpsidx , dpsidy , dpsidz );
872  const double detJ_weight = detJ * integ_weight;
873  const double detJ_weight_coeff_k = detJ_weight * coeff_k;
874 
875  fad_scalar_type value_at_pt = 0 ;
876  fad_scalar_type gradx_at_pt = 0 ;
877  fad_scalar_type grady_at_pt = 0 ;
878  fad_scalar_type gradz_at_pt = 0 ;
879  for ( unsigned m = 0 ; m < FunctionCount ; m++ ) {
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] ;
884  }
885 
886  const fad_scalar_type source_term =
887  coeff_src * value_at_pt * value_at_pt ;
888 
889  const fad_scalar_type advection_term =
890  advection[0]*gradx_at_pt +
891  advection[1]*grady_at_pt +
892  advection[2]*gradz_at_pt;
893 
894  for ( unsigned m = 0; m < FunctionCount; ++m) {
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] ;
899 
900  elem_res[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 ) ;
905  }
906  }
907  }
908 
909  KOKKOS_INLINE_FUNCTION
910  void operator()( const unsigned ielem ) const
911  {
912  double x[ FunctionCount ] ;
913  double y[ FunctionCount ] ;
914  double z[ FunctionCount ] ;
915  unsigned node_index[ ElemNodeCount ];
916 
918  fad_scalar_type elem_res[ FunctionCount ] ; // this zeros elem_res
919 
920  // Gather nodal coordinates and solution vector:
921  gatherSolution( ielem, val, node_index, x, y, z, elem_res );
922 
923  // Compute nodal element residual vector:
924  computeElementResidual( val, x, y, z, elem_res );
925 
926  // Scatter nodal element residual in global vector:
927  scatterResidual( ielem, node_index, elem_res );
928  }
929 }; /* ElementComputation */
930 
931 template< class ExecutionSpace , BoxElemPart::ElemOrder Order ,
932  class CoordinateMap , typename ScalarType >
934  < Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
935  CrsMatrix< ScalarType , ExecutionSpace > ,
937  public ElementComputation< Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
938  CrsMatrix< ScalarType , ExecutionSpace > ,
939  FadElement > {
940 public:
941 
945 
948 
949  static const unsigned FunctionCount = base_type::FunctionCount;
950  static const unsigned IntegrationCount = base_type::IntegrationCount;
951  static const unsigned ElemNodeCount = base_type::ElemNodeCount;
952 
954 
956 
958  const typename base_type::mesh_type & arg_mesh ,
959  const typename base_type::vector_type & arg_solution ,
960  const typename base_type::elem_graph_type & arg_elem_graph ,
961  const typename base_type::sparse_matrix_type & arg_jacobian ,
962  const typename base_type::vector_type & arg_residual ) :
963  base_type(arg_mesh, arg_solution, arg_elem_graph,
964  arg_jacobian, arg_residual) {}
965 
966  //------------------------------------
967 
968  void apply() const
969  {
970  const size_t nelem = this->elem_node_ids.extent(0);
971  parallel_for( nelem , *this );
972  }
973 
974  KOKKOS_INLINE_FUNCTION
975  void gatherSolution(const unsigned ielem,
976  scalar_type val[],
977  unsigned node_index[],
978  double x[], double y[], double z[],
979  fad_scalar_type res[]) const
980  {
981  for ( unsigned i = 0 ; i < ElemNodeCount ; ++i ) {
982  const unsigned ni = this->elem_node_ids( ielem , i );
983 
984  node_index[i] = ni ;
985 
986  x[i] = this->node_coords( ni , 0 );
987  y[i] = this->node_coords( ni , 1 );
988  z[i] = this->node_coords( ni , 2 );
989 
990  val[i] = this->solution( ni );
991  }
992  }
993 
994  KOKKOS_INLINE_FUNCTION
995  void computeElementResidual(const scalar_type dof_values[] ,
996  const double x[],
997  const double y[],
998  const double z[],
999  fad_scalar_type elem_res[] ) const
1000  {
1001  double coeff_k = 3.456;
1002  double coeff_src = 1.234;
1003  double advection[] = { 1.1, 1.2, 1.3 };
1004  double dpsidx[ FunctionCount ] ;
1005  double dpsidy[ FunctionCount ] ;
1006  double dpsidz[ FunctionCount ] ;
1007  for ( unsigned i = 0 ; i < IntegrationCount ; ++i ) {
1008 
1009  const double integ_weight = this->elem_data.weights[i];
1010  const double* bases_vals = this->elem_data.values[i];
1011  const double detJ =
1012  this->transform_gradients( this->elem_data.gradients[i] ,
1013  x , y , z ,
1014  dpsidx , dpsidy , dpsidz );
1015  const double detJ_weight = detJ * integ_weight;
1016  const double detJ_weight_coeff_k = detJ_weight * coeff_k;
1017 
1022  for ( unsigned m = 0 ; m < FunctionCount ; m++ ) {
1023  value_at_pt.val() += dof_values[m] * bases_vals[m] ;
1024  value_at_pt.fastAccessDx(m) = bases_vals[m] ;
1025 
1026  gradx_at_pt.val() += dof_values[m] * dpsidx[m] ;
1027  gradx_at_pt.fastAccessDx(m) = dpsidx[m] ;
1028 
1029  grady_at_pt.val() += dof_values[m] * dpsidy[m] ;
1030  grady_at_pt.fastAccessDx(m) = dpsidy[m] ;
1031 
1032  gradz_at_pt.val() += dof_values[m] * dpsidz[m] ;
1033  gradz_at_pt.fastAccessDx(m) = dpsidz[m] ;
1034  }
1035 
1036  const fad_scalar_type source_term =
1037  coeff_src * value_at_pt * value_at_pt ;
1038 
1039  const fad_scalar_type advection_term =
1040  advection[0]*gradx_at_pt +
1041  advection[1]*grady_at_pt +
1042  advection[2]*gradz_at_pt;
1043 
1044  for ( unsigned m = 0; m < FunctionCount; ++m) {
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] ;
1049 
1050  elem_res[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 ) ;
1055  }
1056  }
1057  }
1058 
1059  KOKKOS_INLINE_FUNCTION
1060  void operator()( const unsigned ielem ) const
1061  {
1062  double x[ FunctionCount ] ;
1063  double y[ FunctionCount ] ;
1064  double z[ FunctionCount ] ;
1065  unsigned node_index[ ElemNodeCount ];
1066 
1068  fad_scalar_type elem_res[ FunctionCount ] ;
1069 
1070  // Gather nodal coordinates and solution vector:
1071  gatherSolution( ielem, val, node_index, x, y, z, elem_res );
1072 
1073  // Compute nodal element residual vector:
1074  computeElementResidual( val, x, y, z, elem_res );
1075 
1076  // Scatter nodal element residual in global vector:
1077  this->scatterResidual( ielem, node_index, elem_res );
1078  }
1079 }; /* ElementComputation */
1080 
1081 template< class ExecutionSpace , BoxElemPart::ElemOrder Order ,
1082  class CoordinateMap , typename ScalarType >
1084  < Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
1085  CrsMatrix< ScalarType , ExecutionSpace > ,
1086  FadQuadPoint > :
1087  public ElementComputation< Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
1088  CrsMatrix< ScalarType , ExecutionSpace > ,
1089  Analytic > {
1090 public:
1091 
1095 
1098 
1099  static const unsigned FunctionCount = base_type::FunctionCount;
1100  static const unsigned IntegrationCount = base_type::IntegrationCount;
1101  static const unsigned ElemNodeCount = base_type::ElemNodeCount;
1102 
1104 
1106 
1108  const typename base_type::mesh_type & arg_mesh ,
1109  const typename base_type::vector_type & arg_solution ,
1110  const typename base_type::elem_graph_type & arg_elem_graph ,
1111  const typename base_type::sparse_matrix_type & arg_jacobian ,
1112  const typename base_type::vector_type & arg_residual ) :
1113  base_type(arg_mesh, arg_solution, arg_elem_graph,
1114  arg_jacobian, arg_residual) {}
1115 
1116  //------------------------------------
1117 
1118  void apply() const
1119  {
1120  const size_t nelem = this->elem_node_ids.extent(0);
1121  parallel_for( nelem , *this );
1122  }
1123 
1124  KOKKOS_INLINE_FUNCTION
1126  const scalar_type dof_values[] ,
1127  const double x[],
1128  const double y[],
1129  const double z[],
1130  scalar_type elem_res[] ,
1131  scalar_type elem_mat[][FunctionCount] ) const
1132  {
1133  double coeff_k = 3.456;
1134  double coeff_src = 1.234;
1135  double advection[] = { 1.1, 1.2, 1.3 };
1136  double dpsidx[ FunctionCount ] ;
1137  double dpsidy[ FunctionCount ] ;
1138  double dpsidz[ FunctionCount ] ;
1139 
1140  fad_scalar_type value_at_pt(4, 0, 0.0) ;
1141  fad_scalar_type gradx_at_pt(4, 1, 0.0) ;
1142  fad_scalar_type grady_at_pt(4, 2, 0.0) ;
1143  fad_scalar_type gradz_at_pt(4, 3, 0.0) ;
1144  for ( unsigned i = 0 ; i < IntegrationCount ; ++i ) {
1145 
1146  const double integ_weight = this->elem_data.weights[i];
1147  const double* bases_vals = this->elem_data.values[i];
1148  const double detJ =
1149  this->transform_gradients( this->elem_data.gradients[i] ,
1150  x , y , z ,
1151  dpsidx , dpsidy , dpsidz );
1152  const double detJ_weight = detJ * integ_weight;
1153  const double detJ_weight_coeff_k = detJ_weight * coeff_k;
1154 
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 ;
1159  for ( unsigned m = 0 ; m < FunctionCount ; m++ ) {
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] ;
1164  }
1165 
1166  const fad_scalar_type source_term =
1167  coeff_src * value_at_pt * value_at_pt ;
1168 
1169  const fad_scalar_type advection_term =
1170  advection[0]*gradx_at_pt +
1171  advection[1]*grady_at_pt +
1172  advection[2]*gradz_at_pt;
1173 
1174  for ( unsigned m = 0; m < FunctionCount; ++m) {
1175  const double bases_val_m = bases_vals[m] * detJ_weight ;
1176  fad_scalar_type res =
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 ) ;
1181 
1182  elem_res[m] += res.val();
1183 
1184  scalar_type * const mat = elem_mat[m] ;
1185  for( unsigned n = 0; n < FunctionCount; n++) {
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];
1190  }
1191  }
1192  }
1193  }
1194 
1195  KOKKOS_INLINE_FUNCTION
1196  void operator()( const unsigned ielem ) const
1197  {
1198  double x[ FunctionCount ] ;
1199  double y[ FunctionCount ] ;
1200  double z[ FunctionCount ] ;
1201  unsigned node_index[ ElemNodeCount ];
1202 
1204  scalar_type elem_res[ FunctionCount ] ;
1205  scalar_type elem_mat[ FunctionCount ][ FunctionCount ] ;
1206 
1207  // Gather nodal coordinates and solution vector:
1208  this->gatherSolution( ielem, val, node_index, x, y, z, elem_res, elem_mat );
1209 
1210  // Compute nodal element residual vector and Jacobian matrix:
1211  computeElementResidualJacobian( val, x, y, z, elem_res, elem_mat );
1212 
1213  // Scatter nodal element residual and Jacobian in global vector and matrix:
1214  this->scatterResidual( ielem, node_index, elem_res, elem_mat );
1215  }
1216 }; /* ElementComputation */
1217 
1218 } /* namespace FENL */
1219 } /* namespace Example */
1220 } /* namespace Kokkos */
1221 
1222 //----------------------------------------------------------------------------
1223 
1224 #endif /* #ifndef KOKKOS_EXAMPLE_FENLFUNCTORS_HPP */
double weights[integration_count]
Definition: HexElement.hpp:192
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
int * count
KOKKOS_INLINE_FUNCTION void fill_graph_entries(const unsigned iset) const
Kokkos::View< scalar_type *[FunctionCount], execution_space > elem_vectors_type
static const unsigned element_node_count
Definition: HexElement.hpp:188
CrsMatrix< ScalarType, ExecutionSpace > sparse_matrix_type
static const unsigned function_count
Definition: HexElement.hpp:190
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
double gradients[integration_count][spatial_dimension][function_count]
Definition: HexElement.hpp:194
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
Kokkos::View< scalar_type *, Kokkos::LayoutLeft, execution_space > vector_type
NodeNodeGraph< elem_node_type, sparse_graph_type, ElemNodeCount >::ElemGraphType elem_graph_type
expr val()
ElemNodeIdView::execution_space execution_space
KOKKOS_INLINE_FUNCTION void operator()(const unsigned irow, unsigned &update, const bool final) const
KOKKOS_INLINE_FUNCTION void init(unsigned &update) const
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
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 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]
Definition: HexElement.hpp:193
Do not initialize the derivative array.
KOKKOS_INLINE_FUNCTION void operator()(const unsigned iwork) const
static const unsigned integration_count
Definition: HexElement.hpp:189
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
Kokkos::Example::BoxElemFixture< ExecutionSpace, Order, CoordinateMap > mesh_type
NodeNodeGraph(const ElemNodeIdView &arg_elem_node_id, const unsigned arg_node_count, Times &results)
Uncopyable z
KOKKOS_INLINE_FUNCTION void init(const TagFillNodeSet &, unsigned &update) const
#define Method
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)
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 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 sort_graph_entries(const unsigned irow) const
Kokkos::View< scalar_type *[FunctionCount][FunctionCount], execution_space > elem_matrices_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
int n
Generate a distributed unstructured finite element mesh from a partitioned NX*NY*NZ box of elements...
Kokkos::UnorderedMap< key_type, void, execution_space > SetType
pair< unsigned, unsigned > key_type
const double y
static const unsigned spatial_dimension
Definition: HexElement.hpp:187