Sacado Package Browser (Single Doxygen Collection)  Version of the Day
fenl_functors.hpp
Go to the documentation of this file.
1 /*
2 //@HEADER
3 // ************************************************************************
4 //
5 // Kokkos: Manycore Performance-Portable Multidimensional Arrays
6 // Copyright (2012) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact H. Carter Edwards (hcedwar@sandia.gov)
39 //
40 // ************************************************************************
41 //@HEADER
42 */
43 
44 #ifndef KOKKOS_EXAMPLE_FENLFUNCTORS_HPP
45 #define KOKKOS_EXAMPLE_FENLFUNCTORS_HPP
46 
47 #include <stdio.h>
48 
49 #include <iostream>
50 #include <fstream>
51 #include <iomanip>
52 #include <cstdlib>
53 #include <cmath>
54 #include <limits>
55 
56 #include <Kokkos_Core.hpp>
57 #include <Kokkos_Pair.hpp>
58 #include <Kokkos_UnorderedMap.hpp>
59 #include <Kokkos_StaticCrsGraph.hpp>
60 
61 #include <impl/Kokkos_Timer.hpp>
62 
63 #include <BoxElemFixture.hpp>
64 #include <HexElement.hpp>
65 
66 #include "Sacado.hpp"
67 
68 //----------------------------------------------------------------------------
69 //----------------------------------------------------------------------------
70 
71 namespace Kokkos {
72 namespace Example {
73 namespace FENL {
74 
75 template< typename ValueType , class Space >
76 struct CrsMatrix {
77  typedef Kokkos::StaticCrsGraph< unsigned , Space , void , unsigned > StaticCrsGraphType ;
78  typedef View< ValueType * , Space > coeff_type ;
79 
82 
83  CrsMatrix() : graph(), coeff() {}
84 
85  CrsMatrix( const StaticCrsGraphType & arg_graph )
86  : graph( arg_graph )
87  , coeff( "crs_matrix_coeff" , arg_graph.entries.dimension_0() )
88  {}
89 };
90 
91 template< class ElemNodeIdView , class CrsGraphType , unsigned ElemNode >
93 public:
94 
95  typedef typename ElemNodeIdView::execution_space execution_space ;
96  typedef pair<unsigned,unsigned> key_type ;
97 
98  typedef Kokkos::UnorderedMap< key_type, void , execution_space > SetType ;
99  typedef typename CrsGraphType::row_map_type::non_const_type RowMapType ;
100  typedef Kokkos::View< unsigned , execution_space > UnsignedValue ;
101 
102  // Static dimensions of 0 generate compiler warnings or errors.
103  typedef Kokkos::View< unsigned*[ElemNode][ElemNode] , execution_space >
105 
106  struct TagFillNodeSet {};
107  struct TagScanNodeCount {};
111 
112 private:
113 
119 
120  const unsigned node_count ;
121  const ElemNodeIdView elem_node_id ;
127 
128 public:
129 
130  CrsGraphType graph ;
132 
133  struct Times
134  {
135  double ratio;
141  };
142 
143  NodeNodeGraph( const ElemNodeIdView & arg_elem_node_id ,
144  const unsigned arg_node_count,
145  Times & results
146  )
147  : node_count(arg_node_count)
148  , elem_node_id( arg_elem_node_id )
149  , row_total( "row_total" )
150  , row_count(Kokkos::ViewAllocateWithoutInitializing("row_count") , node_count ) // will deep_copy to 0 inside loop
151  , row_map( "graph_row_map" , node_count + 1 )
152  , node_node_set()
153  , phase( FILL_NODE_SET )
154  , graph()
155  , elem_graph()
156  {
157  //--------------------------------
158  // Guess at capacity required for the map:
159 
160  Kokkos::Impl::Timer wall_clock ;
161 
162  wall_clock.reset();
163  phase = FILL_NODE_SET ;
164 
165  // upper bound on the capacity
166  size_t set_capacity = (28ull * node_count) / 2;
167  unsigned failed_insert_count = 0 ;
168 
169  do {
170  // Zero the row count to restart the fill
171  Kokkos::deep_copy( row_count , 0u );
172 
173  node_node_set = SetType( ( set_capacity += failed_insert_count ) );
174 
175  // May be larger that requested:
176  set_capacity = node_node_set.capacity();
177 
178  Kokkos::parallel_reduce( Kokkos::RangePolicy<execution_space,TagFillNodeSet>(0,elem_node_id.dimension_0())
179  , *this
180  , failed_insert_count );
181 
182  } while ( failed_insert_count );
183 
184  execution_space::fence();
185  results.ratio = (double)node_node_set.size() / (double)node_node_set.capacity();
186  results.fill_node_set = wall_clock.seconds();
187  //--------------------------------
188 
189  wall_clock.reset();
191 
192  // Exclusive scan of row_count into row_map
193  // including the final total in the 'node_count + 1' position.
194  // Zero the 'row_count' values.
195  Kokkos::parallel_scan( node_count , *this );
196 
197  // Zero the row count for the fill:
198  Kokkos::deep_copy( row_count , 0u );
199 
200  unsigned graph_entry_count = 0 ;
201 
202  Kokkos::deep_copy( graph_entry_count , row_total );
203 
204  // Assign graph's row_map and allocate graph's entries
205  graph.row_map = row_map ;
206  graph.entries = typename CrsGraphType::entries_type( "graph_entries" , graph_entry_count );
207 
208  //--------------------------------
209  // Fill graph's entries from the (node,node) set.
210 
211  execution_space::fence();
212  results.scan_node_count = wall_clock.seconds();
213 
214  wall_clock.reset();
216  Kokkos::parallel_for( node_node_set.capacity() , *this );
217 
218  execution_space::fence();
219  results.fill_graph_entries = wall_clock.seconds();
220 
221  //--------------------------------
222  // Done with the temporary sets and arrays
223  wall_clock.reset();
225 
227  row_count = RowMapType();
228  row_map = RowMapType();
229  node_node_set.clear();
230 
231  //--------------------------------
232 
233  Kokkos::parallel_for( node_count , *this );
234 
235  execution_space::fence();
236  results.sort_graph_entries = wall_clock.seconds();
237 
238  //--------------------------------
239  // Element-to-graph mapping:
240  wall_clock.reset();
242  elem_graph = ElemGraphType("elem_graph", elem_node_id.dimension_0() );
243  Kokkos::parallel_for( elem_node_id.dimension_0() , *this );
244 
245  execution_space::fence();
246  results.fill_element_graph = wall_clock.seconds();
247  }
248 
249  //------------------------------------
250  // parallel_for: create map and count row length
251 
253  void operator()( const TagFillNodeSet & , unsigned ielem , unsigned & count ) const
254  {
255  // Loop over element's (row_local_node,col_local_node) pairs:
256  for ( unsigned row_local_node = 0 ; row_local_node < elem_node_id.dimension_1() ; ++row_local_node ) {
257 
258  const unsigned row_node = elem_node_id( ielem , row_local_node );
259 
260  for ( unsigned col_local_node = row_local_node ; col_local_node < elem_node_id.dimension_1() ; ++col_local_node ) {
261 
262  const unsigned col_node = elem_node_id( ielem , col_local_node );
263 
264  // If either node is locally owned then insert the pair into the unordered map:
265 
266  if ( row_node < row_count.dimension_0() || col_node < row_count.dimension_0() ) {
267 
268  const key_type key = (row_node < col_node) ? make_pair( row_node, col_node ) : make_pair( col_node, row_node ) ;
269 
270  const typename SetType::insert_result result = node_node_set.insert( key );
271 
272  // A successfull insert: the first time this pair was added
273  if ( result.success() ) {
274 
275  // If row node is owned then increment count
276  if ( row_node < row_count.dimension_0() ) { atomic_fetch_add( & row_count( row_node ) , 1 ); }
277 
278  // If column node is owned and not equal to row node then increment count
279  if ( col_node < row_count.dimension_0() && col_node != row_node ) { atomic_fetch_add( & row_count( col_node ) , 1 ); }
280  }
281  else if ( result.failed() ) {
282  ++count ;
283  }
284  }
285  }
286  }
287  }
288 
290  void fill_graph_entries( const unsigned iset ) const
291  {
292  if ( node_node_set.valid_at(iset) ) {
293  // Add each entry to the graph entries.
294 
295  const key_type key = node_node_set.key_at(iset) ;
296  const unsigned row_node = key.first ;
297  const unsigned col_node = key.second ;
298 
299  if ( row_node < row_count.dimension_0() ) {
300  const unsigned offset = graph.row_map( row_node ) + atomic_fetch_add( & row_count( row_node ) , 1 );
301  graph.entries( offset ) = col_node ;
302  }
303 
304  if ( col_node < row_count.dimension_0() && col_node != row_node ) {
305  const unsigned offset = graph.row_map( col_node ) + atomic_fetch_add( & row_count( col_node ) , 1 );
306  graph.entries( offset ) = row_node ;
307  }
308  }
309  }
310 
312  void sort_graph_entries( const unsigned irow ) const
313  {
314  const unsigned row_beg = graph.row_map( irow );
315  const unsigned row_end = graph.row_map( irow + 1 );
316  for ( unsigned i = row_beg + 1 ; i < row_end ; ++i ) {
317  const unsigned col = graph.entries(i);
318  unsigned j = i ;
319  for ( ; row_beg < j && col < graph.entries(j-1) ; --j ) {
320  graph.entries(j) = graph.entries(j-1);
321  }
322  graph.entries(j) = col ;
323  }
324  }
325 
327  void fill_elem_graph_map( const unsigned ielem ) const
328  {
329  for ( unsigned row_local_node = 0 ; row_local_node < elem_node_id.dimension_1() ; ++row_local_node ) {
330 
331  const unsigned row_node = elem_node_id( ielem , row_local_node );
332 
333  for ( unsigned col_local_node = 0 ; col_local_node < elem_node_id.dimension_1() ; ++col_local_node ) {
334 
335  const unsigned col_node = elem_node_id( ielem , col_local_node );
336 
337  unsigned entry = ~0u ;
338 
339  if ( row_node + 1 < graph.row_map.dimension_0() ) {
340 
341  const unsigned entry_end = graph.row_map( row_node + 1 );
342 
343  entry = graph.row_map( row_node );
344 
345  for ( ; entry < entry_end && graph.entries(entry) != col_node ; ++entry );
346 
347  if ( entry == entry_end ) entry = ~0u ;
348  }
349 
350  elem_graph( ielem , row_local_node , col_local_node ) = entry ;
351  }
352  }
353  }
354 
356  void operator()( const unsigned iwork ) const
357  {
358 /*
359  if ( phase == FILL_NODE_SET ) {
360  operator()( TagFillNodeSet() , iwork );
361  }
362  else */
363  if ( phase == FILL_GRAPH_ENTRIES ) {
364  fill_graph_entries( iwork );
365  }
366  else if ( phase == SORT_GRAPH_ENTRIES ) {
367  sort_graph_entries( iwork );
368  }
369  else if ( phase == FILL_ELEMENT_GRAPH ) {
370  fill_elem_graph_map( iwork );
371  }
372  }
373 
374  //------------------------------------
375  // parallel_scan: row offsets
376 
377  typedef unsigned value_type ;
378 
380  void operator()( const unsigned irow , unsigned & update , const bool final ) const
381  {
382  // exclusive scan
383  if ( final ) { row_map( irow ) = update ; }
384 
385  update += row_count( irow );
386 
387  if ( final ) {
388  if ( irow + 1 == row_count.dimension_0() ) {
389  row_map( irow + 1 ) = update ;
390  row_total() = update ;
391  }
392  }
393  }
394 
395  // For the reduce phase:
397  void init( const TagFillNodeSet & , unsigned & update ) const { update = 0 ; }
398 
400  void join( const TagFillNodeSet &
401  , volatile unsigned & update
402  , volatile const unsigned & input ) const { update += input ; }
403 
404  // For the scan phase::
406  void init( unsigned & update ) const { update = 0 ; }
407 
409  void join( volatile unsigned & update
410  , volatile const unsigned & input ) const { update += input ; }
411 
412  //------------------------------------
413 };
414 
415 } /* namespace FENL */
416 } /* namespace Example */
417 } /* namespace Kokkos */
418 
419 //----------------------------------------------------------------------------
420 //----------------------------------------------------------------------------
421 
422 namespace Kokkos {
423 namespace Example {
424 namespace FENL {
425 
426 template< class ExecutionSpace , BoxElemPart::ElemOrder Order ,
427  class CoordinateMap , typename ScalarType >
429 {
430 public:
431 
434 
435  //------------------------------------
436 
437  typedef ExecutionSpace execution_space ;
438  typedef ScalarType scalar_type ;
439 
442  typedef Kokkos::View< scalar_type* , Kokkos::LayoutLeft, execution_space > vector_type ;
443 
444  //------------------------------------
445 
447  static const unsigned TensorDim = SpatialDim * SpatialDim ;
451 
452  //------------------------------------
453 
456  typedef Kokkos::View< scalar_type*[FunctionCount][FunctionCount] , execution_space > elem_matrices_type ;
457  typedef Kokkos::View< scalar_type*[FunctionCount] , execution_space > elem_vectors_type ;
458 
460 
461  //------------------------------------
462 
463 
464  //------------------------------------
465  // Computational data:
466 
476 
478  : elem_data()
480  , node_coords( rhs.node_coords )
481  , elem_graph( rhs.elem_graph )
484  , solution( rhs.solution )
485  , residual( rhs.residual )
486  , jacobian( rhs.jacobian )
487  {}
488 
489  ElementComputationBase( const mesh_type & arg_mesh ,
490  const vector_type & arg_solution ,
491  const elem_graph_type & arg_elem_graph ,
492  const sparse_matrix_type & arg_jacobian ,
493  const vector_type & arg_residual )
494  : elem_data()
495  , elem_node_ids( arg_mesh.elem_node() )
496  , node_coords( arg_mesh.node_coord() )
497  , elem_graph( arg_elem_graph )
498  , elem_jacobians()
499  , elem_residuals()
500  , solution( arg_solution )
501  , residual( arg_residual )
502  , jacobian( arg_jacobian )
503  {}
504 
505  //------------------------------------
506 
509  const double grad[][ FunctionCount ] , // Gradient of bases master element
510  const double x[] ,
511  const double y[] ,
512  const double z[] ,
513  double dpsidx[] ,
514  double dpsidy[] ,
515  double dpsidz[] ) const
516  {
517  enum { j11 = 0 , j12 = 1 , j13 = 2 ,
518  j21 = 3 , j22 = 4 , j23 = 5 ,
519  j31 = 6 , j32 = 7 , j33 = 8 };
520 
521  // Jacobian accumulation:
522 
523  double J[ TensorDim ] = { 0, 0, 0, 0, 0, 0, 0, 0, 0 };
524 
525  for( unsigned i = 0; i < FunctionCount ; ++i ) {
526  const double x1 = x[i] ;
527  const double x2 = y[i] ;
528  const double x3 = z[i] ;
529 
530  const double g1 = grad[0][i] ;
531  const double g2 = grad[1][i] ;
532  const double g3 = grad[2][i] ;
533 
534  J[j11] += g1 * x1 ;
535  J[j12] += g1 * x2 ;
536  J[j13] += g1 * x3 ;
537 
538  J[j21] += g2 * x1 ;
539  J[j22] += g2 * x2 ;
540  J[j23] += g2 * x3 ;
541 
542  J[j31] += g3 * x1 ;
543  J[j32] += g3 * x2 ;
544  J[j33] += g3 * x3 ;
545  }
546 
547  // Inverse jacobian:
548 
549  double invJ[ TensorDim ] = {
550  static_cast<double>( J[j22] * J[j33] - J[j23] * J[j32] ) ,
551  static_cast<double>( J[j13] * J[j32] - J[j12] * J[j33] ) ,
552  static_cast<double>( J[j12] * J[j23] - J[j13] * J[j22] ) ,
553 
554  static_cast<double>( J[j23] * J[j31] - J[j21] * J[j33] ) ,
555  static_cast<double>( J[j11] * J[j33] - J[j13] * J[j31] ) ,
556  static_cast<double>( J[j13] * J[j21] - J[j11] * J[j23] ) ,
557 
558  static_cast<double>( J[j21] * J[j32] - J[j22] * J[j31] ) ,
559  static_cast<double>( J[j12] * J[j31] - J[j11] * J[j32] ) ,
560  static_cast<double>( J[j11] * J[j22] - J[j12] * J[j21] ) };
561 
562  const double detJ = J[j11] * invJ[j11] +
563  J[j21] * invJ[j12] +
564  J[j31] * invJ[j13] ;
565 
566  const double detJinv = 1.0 / detJ ;
567 
568  for ( unsigned i = 0 ; i < TensorDim ; ++i ) { invJ[i] *= detJinv ; }
569 
570  // Transform gradients:
571 
572  for( unsigned i = 0; i < FunctionCount ; ++i ) {
573  const double g0 = grad[0][i];
574  const double g1 = grad[1][i];
575  const double g2 = grad[2][i];
576 
577  dpsidx[i] = g0 * invJ[j11] + g1 * invJ[j12] + g2 * invJ[j13];
578  dpsidy[i] = g0 * invJ[j21] + g1 * invJ[j22] + g2 * invJ[j23];
579  dpsidz[i] = g0 * invJ[j31] + g1 * invJ[j32] + g2 * invJ[j33];
580  }
581 
582  return detJ ;
583  }
584 
585 };
586 
592 };
593 
594 template< class FiniteElementMeshType ,
595  class SparseMatrixType ,
596  AssemblyMethod Method
597  >
599 
600 template< class ExecutionSpace , BoxElemPart::ElemOrder Order ,
601  class CoordinateMap , typename ScalarType >
603  < Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
604  CrsMatrix< ScalarType , ExecutionSpace > ,
605  Analytic > :
606  public ElementComputationBase<ExecutionSpace, Order, CoordinateMap,
607  ScalarType> {
608 public:
609 
610  typedef ElementComputationBase<ExecutionSpace, Order, CoordinateMap,
611  ScalarType> base_type;
612 
615 
616  static const unsigned FunctionCount = base_type::FunctionCount;
617  static const unsigned IntegrationCount = base_type::IntegrationCount;
618  static const unsigned ElemNodeCount = base_type::ElemNodeCount;
619 
621 
623  const typename base_type::mesh_type & arg_mesh ,
624  const typename base_type::vector_type & arg_solution ,
625  const typename base_type::elem_graph_type & arg_elem_graph ,
626  const typename base_type::sparse_matrix_type & arg_jacobian ,
627  const typename base_type::vector_type & arg_residual ) :
628  base_type(arg_mesh, arg_solution, arg_elem_graph,
629  arg_jacobian, arg_residual) {}
630 
631  //------------------------------------
632 
633  void apply() const
634  {
635  const size_t nelem = this->elem_node_ids.dimension_0();
636  parallel_for( nelem , *this );
637  }
638 
640  void gatherSolution(const unsigned ielem,
641  scalar_type val[],
642  unsigned node_index[],
643  double x[], double y[], double z[],
644  scalar_type res[],
645  scalar_type mat[][FunctionCount]) const
646  {
647  for ( unsigned i = 0 ; i < ElemNodeCount ; ++i ) {
648  const unsigned ni = this->elem_node_ids( ielem , i );
649 
650  node_index[i] = ni ;
651 
652  x[i] = this->node_coords( ni , 0 );
653  y[i] = this->node_coords( ni , 1 );
654  z[i] = this->node_coords( ni , 2 );
655 
656  val[i] = this->solution( ni ) ;
657  res[i] = 0 ;
658 
659  for( unsigned j = 0; j < FunctionCount ; j++){
660  mat[i][j] = 0 ;
661  }
662  }
663  }
664 
666  void scatterResidual(const unsigned ielem,
667  const unsigned node_index[],
668  const scalar_type res[],
669  const scalar_type mat[][FunctionCount]) const
670  {
671  for( unsigned i = 0 ; i < FunctionCount ; i++ ) {
672  const unsigned row = node_index[i] ;
673  if ( row < this->residual.dimension_0() ) {
674  atomic_add( & this->residual( row ) , res[i] );
675 
676  for( unsigned j = 0 ; j < FunctionCount ; j++ ) {
677  const unsigned entry = this->elem_graph( ielem , i , j );
678  if ( entry != ~0u ) {
679  atomic_add( & this->jacobian.coeff( entry ) , mat[i][j] );
680  }
681  }
682  }
683  }
684  }
685 
688  const scalar_type dof_values[] ,
689  const double x[],
690  const double y[],
691  const double z[],
692  scalar_type elem_res[] ,
693  scalar_type elem_mat[][FunctionCount] ) const
694  {
695  double coeff_k = 3.456;
696  double coeff_src = 1.234;
697  double advection[] = { 1.1, 1.2, 1.3 };
698  double dpsidx[ FunctionCount ] ;
699  double dpsidy[ FunctionCount ] ;
700  double dpsidz[ FunctionCount ] ;
701  for ( unsigned i = 0 ; i < IntegrationCount ; ++i ) {
702 
703  const double integ_weight = this->elem_data.weights[i];
704  const double* bases_vals = this->elem_data.values[i];
705  const double detJ =
706  this->transform_gradients( this->elem_data.gradients[i] ,
707  x , y , z ,
708  dpsidx , dpsidy , dpsidz );
709  const double detJ_weight = detJ * integ_weight;
710  const double detJ_weight_coeff_k = detJ_weight * coeff_k;
711 
712  scalar_type value_at_pt = 0 ;
713  scalar_type gradx_at_pt = 0 ;
714  scalar_type grady_at_pt = 0 ;
715  scalar_type gradz_at_pt = 0 ;
716  for ( unsigned m = 0 ; m < FunctionCount ; m++ ) {
717  value_at_pt += dof_values[m] * bases_vals[m] ;
718  gradx_at_pt += dof_values[m] * dpsidx[m] ;
719  grady_at_pt += dof_values[m] * dpsidy[m] ;
720  gradz_at_pt += dof_values[m] * dpsidz[m] ;
721  }
722 
723  const scalar_type source_term =
724  coeff_src * value_at_pt * value_at_pt ;
725  const scalar_type source_deriv =
726  2.0 * coeff_src * value_at_pt ;
727 
728  const scalar_type advection_x = advection[0];
729  const scalar_type advection_y = advection[1];
730  const scalar_type advection_z = advection[2];
731 
732  const scalar_type advection_term =
733  advection_x*gradx_at_pt +
734  advection_y*grady_at_pt +
735  advection_z*gradz_at_pt ;
736 
737  for ( unsigned m = 0; m < FunctionCount; ++m) {
738  scalar_type * const mat = elem_mat[m] ;
739  const double bases_val_m = bases_vals[m] * detJ_weight ;
740  const double dpsidx_m = dpsidx[m] ;
741  const double dpsidy_m = dpsidy[m] ;
742  const double dpsidz_m = dpsidz[m] ;
743 
744  elem_res[m] +=
745  detJ_weight_coeff_k * ( dpsidx_m * gradx_at_pt +
746  dpsidy_m * grady_at_pt +
747  dpsidz_m * gradz_at_pt ) +
748  bases_val_m * ( advection_term + source_term ) ;
749 
750  for( unsigned n = 0; n < FunctionCount; n++) {
751  const double dpsidx_n = dpsidx[n] ;
752  const double dpsidy_n = dpsidy[n] ;
753  const double dpsidz_n = dpsidz[n] ;
754  mat[n] +=
755  detJ_weight_coeff_k * ( dpsidx_m * dpsidx_n +
756  dpsidy_m * dpsidy_n +
757  dpsidz_m * dpsidz_n ) +
758  bases_val_m * ( advection_x * dpsidx_n +
759  advection_y * dpsidy_n +
760  advection_z * dpsidz_n +
761  source_deriv * bases_vals[n] ) ;
762  }
763  }
764  }
765  }
766 
768  void operator()( const unsigned ielem ) const
769  {
770  double x[ FunctionCount ] ;
771  double y[ FunctionCount ] ;
772  double z[ FunctionCount ] ;
773  unsigned node_index[ ElemNodeCount ];
774 
775  scalar_type val[ FunctionCount ] ;
776  scalar_type elem_res[ FunctionCount ] ;
777  scalar_type elem_mat[ FunctionCount ][ FunctionCount ] ;
778 
779  // Gather nodal coordinates and solution vector:
780  gatherSolution(ielem, val, node_index, x, y, z, elem_res, elem_mat);
781 
782  // Compute nodal element residual vector and Jacobian matrix
783  computeElementResidualJacobian( val, x, y, z, elem_res , elem_mat );
784 
785  // Scatter nodal element residual and Jacobian in global vector and matrix:
786  scatterResidual( ielem, node_index, elem_res, elem_mat );
787  }
788 }; /* ElementComputation */
789 
790 template< class ExecutionSpace , BoxElemPart::ElemOrder Order ,
791  class CoordinateMap , typename ScalarType >
793  < Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
794  CrsMatrix< ScalarType , ExecutionSpace > ,
795  FadElement > : public ElementComputationBase<ExecutionSpace, Order, CoordinateMap,
796  ScalarType> {
797 public:
798 
799  typedef ElementComputationBase<ExecutionSpace, Order, CoordinateMap,
800  ScalarType> base_type;
801 
804 
805  static const unsigned FunctionCount = base_type::FunctionCount;
806  static const unsigned IntegrationCount = base_type::IntegrationCount;
807  static const unsigned ElemNodeCount = base_type::ElemNodeCount;
808 
810 
812 
814  const typename base_type::mesh_type & arg_mesh ,
815  const typename base_type::vector_type & arg_solution ,
816  const typename base_type::elem_graph_type & arg_elem_graph ,
817  const typename base_type::sparse_matrix_type & arg_jacobian ,
818  const typename base_type::vector_type & arg_residual ) :
819  base_type(arg_mesh, arg_solution, arg_elem_graph,
820  arg_jacobian, arg_residual) {}
821 
822  //------------------------------------
823 
824  void apply() const
825  {
826  const size_t nelem = this->elem_node_ids.dimension_0();
827  parallel_for( nelem , *this );
828  }
829 
831  void gatherSolution(const unsigned ielem,
833  unsigned node_index[],
834  double x[], double y[], double z[],
835  fad_scalar_type res[]) const
836  {
837  for ( unsigned i = 0 ; i < ElemNodeCount ; ++i ) {
838  const unsigned ni = this->elem_node_ids( ielem , i );
839 
840  node_index[i] = ni ;
841 
842  x[i] = this->node_coords( ni , 0 );
843  y[i] = this->node_coords( ni , 1 );
844  z[i] = this->node_coords( ni , 2 );
845 
846  val[i].val() = this->solution( ni );
847  val[i].diff( i, FunctionCount );
848  }
849  }
850 
852  void scatterResidual(const unsigned ielem,
853  const unsigned node_index[],
854  fad_scalar_type res[]) const
855  {
856  for( unsigned i = 0 ; i < FunctionCount ; i++ ) {
857  const unsigned row = node_index[i] ;
858  if ( row < this->residual.dimension_0() ) {
859  atomic_add( & this->residual( row ) , res[i].val() );
860 
861  for( unsigned j = 0 ; j < FunctionCount ; j++ ) {
862  const unsigned entry = this->elem_graph( ielem , i , j );
863  if ( entry != ~0u ) {
864  atomic_add( & this->jacobian.coeff( entry ) ,
865  res[i].fastAccessDx(j) );
866  }
867  }
868  }
869  }
870  }
871 
873  void computeElementResidual(const fad_scalar_type dof_values[] ,
874  const double x[],
875  const double y[],
876  const double z[],
877  fad_scalar_type elem_res[] ) const
878  {
879  double coeff_k = 3.456;
880  double coeff_src = 1.234;
881  double advection[] = { 1.1, 1.2, 1.3 };
882  double dpsidx[ FunctionCount ] ;
883  double dpsidy[ FunctionCount ] ;
884  double dpsidz[ FunctionCount ] ;
885  for ( unsigned i = 0 ; i < IntegrationCount ; ++i ) {
886 
887  const double integ_weight = this->elem_data.weights[i];
888  const double* bases_vals = this->elem_data.values[i];
889  const double detJ =
890  this->transform_gradients( this->elem_data.gradients[i] ,
891  x , y , z ,
892  dpsidx , dpsidy , dpsidz );
893  const double detJ_weight = detJ * integ_weight;
894  const double detJ_weight_coeff_k = detJ_weight * coeff_k;
895 
896  fad_scalar_type value_at_pt = 0 ;
897  fad_scalar_type gradx_at_pt = 0 ;
898  fad_scalar_type grady_at_pt = 0 ;
899  fad_scalar_type gradz_at_pt = 0 ;
900  for ( unsigned m = 0 ; m < FunctionCount ; m++ ) {
901  value_at_pt += dof_values[m] * bases_vals[m] ;
902  gradx_at_pt += dof_values[m] * dpsidx[m] ;
903  grady_at_pt += dof_values[m] * dpsidy[m] ;
904  gradz_at_pt += dof_values[m] * dpsidz[m] ;
905  }
906 
907  const fad_scalar_type source_term =
908  coeff_src * value_at_pt * value_at_pt ;
909 
910  const fad_scalar_type advection_term =
911  advection[0]*gradx_at_pt +
912  advection[1]*grady_at_pt +
913  advection[2]*gradz_at_pt;
914 
915  for ( unsigned m = 0; m < FunctionCount; ++m) {
916  const double bases_val_m = bases_vals[m] * detJ_weight ;
917  const double dpsidx_m = dpsidx[m] ;
918  const double dpsidy_m = dpsidy[m] ;
919  const double dpsidz_m = dpsidz[m] ;
920 
921  elem_res[m] +=
922  detJ_weight_coeff_k * ( dpsidx_m * gradx_at_pt +
923  dpsidy_m * grady_at_pt +
924  dpsidz_m * gradz_at_pt ) +
925  bases_val_m * ( advection_term + source_term ) ;
926  }
927  }
928  }
929 
931  void operator()( const unsigned ielem ) const
932  {
933  double x[ FunctionCount ] ;
934  double y[ FunctionCount ] ;
935  double z[ FunctionCount ] ;
936  unsigned node_index[ ElemNodeCount ];
937 
938  fad_scalar_type val[ FunctionCount ] ;
939  fad_scalar_type elem_res[ FunctionCount ] ; // this zeros elem_res
940 
941  // Gather nodal coordinates and solution vector:
942  gatherSolution( ielem, val, node_index, x, y, z, elem_res );
943 
944  // Compute nodal element residual vector:
945  computeElementResidual( val, x, y, z, elem_res );
946 
947  // Scatter nodal element residual in global vector:
948  scatterResidual( ielem, node_index, elem_res );
949  }
950 }; /* ElementComputation */
951 
952 template< class ExecutionSpace , BoxElemPart::ElemOrder Order ,
953  class CoordinateMap , typename ScalarType >
955  < Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
956  CrsMatrix< ScalarType , ExecutionSpace > ,
958  public ElementComputation< Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
959  CrsMatrix< ScalarType , ExecutionSpace > ,
960  FadElement > {
961 public:
962 
966 
969 
970  static const unsigned FunctionCount = base_type::FunctionCount;
971  static const unsigned IntegrationCount = base_type::IntegrationCount;
972  static const unsigned ElemNodeCount = base_type::ElemNodeCount;
973 
975 
977 
979  const typename base_type::mesh_type & arg_mesh ,
980  const typename base_type::vector_type & arg_solution ,
981  const typename base_type::elem_graph_type & arg_elem_graph ,
982  const typename base_type::sparse_matrix_type & arg_jacobian ,
983  const typename base_type::vector_type & arg_residual ) :
984  base_type(arg_mesh, arg_solution, arg_elem_graph,
985  arg_jacobian, arg_residual) {}
986 
987  //------------------------------------
988 
989  void apply() const
990  {
991  const size_t nelem = this->elem_node_ids.dimension_0();
992  parallel_for( nelem , *this );
993  }
994 
996  void gatherSolution(const unsigned ielem,
997  scalar_type val[],
998  unsigned node_index[],
999  double x[], double y[], double z[],
1000  fad_scalar_type res[]) const
1001  {
1002  for ( unsigned i = 0 ; i < ElemNodeCount ; ++i ) {
1003  const unsigned ni = this->elem_node_ids( ielem , i );
1004 
1005  node_index[i] = ni ;
1006 
1007  x[i] = this->node_coords( ni , 0 );
1008  y[i] = this->node_coords( ni , 1 );
1009  z[i] = this->node_coords( ni , 2 );
1010 
1011  val[i] = this->solution( ni );
1012  }
1013  }
1014 
1016  void computeElementResidual(const scalar_type dof_values[] ,
1017  const double x[],
1018  const double y[],
1019  const double z[],
1020  fad_scalar_type elem_res[] ) const
1021  {
1022  double coeff_k = 3.456;
1023  double coeff_src = 1.234;
1024  double advection[] = { 1.1, 1.2, 1.3 };
1025  double dpsidx[ FunctionCount ] ;
1026  double dpsidy[ FunctionCount ] ;
1027  double dpsidz[ FunctionCount ] ;
1028  for ( unsigned i = 0 ; i < IntegrationCount ; ++i ) {
1029 
1030  const double integ_weight = this->elem_data.weights[i];
1031  const double* bases_vals = this->elem_data.values[i];
1032  const double detJ =
1033  this->transform_gradients( this->elem_data.gradients[i] ,
1034  x , y , z ,
1035  dpsidx , dpsidy , dpsidz );
1036  const double detJ_weight = detJ * integ_weight;
1037  const double detJ_weight_coeff_k = detJ_weight * coeff_k;
1038 
1039  fad_scalar_type value_at_pt(FunctionCount, 0.0, Sacado::NoInitDerivArray) ;
1040  fad_scalar_type gradx_at_pt(FunctionCount, 0.0, Sacado::NoInitDerivArray) ;
1041  fad_scalar_type grady_at_pt(FunctionCount, 0.0, Sacado::NoInitDerivArray) ;
1042  fad_scalar_type gradz_at_pt(FunctionCount, 0.0, Sacado::NoInitDerivArray) ;
1043  for ( unsigned m = 0 ; m < FunctionCount ; m++ ) {
1044  value_at_pt.val() += dof_values[m] * bases_vals[m] ;
1045  value_at_pt.fastAccessDx(m) = bases_vals[m] ;
1046 
1047  gradx_at_pt.val() += dof_values[m] * dpsidx[m] ;
1048  gradx_at_pt.fastAccessDx(m) = dpsidx[m] ;
1049 
1050  grady_at_pt.val() += dof_values[m] * dpsidy[m] ;
1051  grady_at_pt.fastAccessDx(m) = dpsidy[m] ;
1052 
1053  gradz_at_pt.val() += dof_values[m] * dpsidz[m] ;
1054  gradz_at_pt.fastAccessDx(m) = dpsidz[m] ;
1055  }
1056 
1057  const fad_scalar_type source_term =
1058  coeff_src * value_at_pt * value_at_pt ;
1059 
1060  const fad_scalar_type advection_term =
1061  advection[0]*gradx_at_pt +
1062  advection[1]*grady_at_pt +
1063  advection[2]*gradz_at_pt;
1064 
1065  for ( unsigned m = 0; m < FunctionCount; ++m) {
1066  const double bases_val_m = bases_vals[m] * detJ_weight ;
1067  const double dpsidx_m = dpsidx[m] ;
1068  const double dpsidy_m = dpsidy[m] ;
1069  const double dpsidz_m = dpsidz[m] ;
1070 
1071  elem_res[m] +=
1072  detJ_weight_coeff_k * ( dpsidx_m * gradx_at_pt +
1073  dpsidy_m * grady_at_pt +
1074  dpsidz_m * gradz_at_pt ) +
1075  bases_val_m * ( advection_term + source_term ) ;
1076  }
1077  }
1078  }
1079 
1081  void operator()( const unsigned ielem ) const
1082  {
1083  double x[ FunctionCount ] ;
1084  double y[ FunctionCount ] ;
1085  double z[ FunctionCount ] ;
1086  unsigned node_index[ ElemNodeCount ];
1087 
1088  scalar_type val[ FunctionCount ] ;
1089  fad_scalar_type elem_res[ FunctionCount ] ;
1090 
1091  // Gather nodal coordinates and solution vector:
1092  gatherSolution( ielem, val, node_index, x, y, z, elem_res );
1093 
1094  // Compute nodal element residual vector:
1095  computeElementResidual( val, x, y, z, elem_res );
1096 
1097  // Scatter nodal element residual in global vector:
1098  this->scatterResidual( ielem, node_index, elem_res );
1099  }
1100 }; /* ElementComputation */
1101 
1102 template< class ExecutionSpace , BoxElemPart::ElemOrder Order ,
1103  class CoordinateMap , typename ScalarType >
1105  < Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
1106  CrsMatrix< ScalarType , ExecutionSpace > ,
1107  FadQuadPoint > :
1108  public ElementComputation< Kokkos::Example::BoxElemFixture< ExecutionSpace , Order , CoordinateMap > ,
1109  CrsMatrix< ScalarType , ExecutionSpace > ,
1110  Analytic > {
1111 public:
1112 
1116 
1119 
1120  static const unsigned FunctionCount = base_type::FunctionCount;
1121  static const unsigned IntegrationCount = base_type::IntegrationCount;
1122  static const unsigned ElemNodeCount = base_type::ElemNodeCount;
1123 
1125 
1127 
1129  const typename base_type::mesh_type & arg_mesh ,
1130  const typename base_type::vector_type & arg_solution ,
1131  const typename base_type::elem_graph_type & arg_elem_graph ,
1132  const typename base_type::sparse_matrix_type & arg_jacobian ,
1133  const typename base_type::vector_type & arg_residual ) :
1134  base_type(arg_mesh, arg_solution, arg_elem_graph,
1135  arg_jacobian, arg_residual) {}
1136 
1137  //------------------------------------
1138 
1139  void apply() const
1140  {
1141  const size_t nelem = this->elem_node_ids.dimension_0();
1142  parallel_for( nelem , *this );
1143  }
1144 
1147  const scalar_type dof_values[] ,
1148  const double x[],
1149  const double y[],
1150  const double z[],
1151  scalar_type elem_res[] ,
1152  scalar_type elem_mat[][FunctionCount] ) const
1153  {
1154  double coeff_k = 3.456;
1155  double coeff_src = 1.234;
1156  double advection[] = { 1.1, 1.2, 1.3 };
1157  double dpsidx[ FunctionCount ] ;
1158  double dpsidy[ FunctionCount ] ;
1159  double dpsidz[ FunctionCount ] ;
1160 
1161  fad_scalar_type value_at_pt(4, 0, 0.0) ;
1162  fad_scalar_type gradx_at_pt(4, 1, 0.0) ;
1163  fad_scalar_type grady_at_pt(4, 2, 0.0) ;
1164  fad_scalar_type gradz_at_pt(4, 3, 0.0) ;
1165  for ( unsigned i = 0 ; i < IntegrationCount ; ++i ) {
1166 
1167  const double integ_weight = this->elem_data.weights[i];
1168  const double* bases_vals = this->elem_data.values[i];
1169  const double detJ =
1170  this->transform_gradients( this->elem_data.gradients[i] ,
1171  x , y , z ,
1172  dpsidx , dpsidy , dpsidz );
1173  const double detJ_weight = detJ * integ_weight;
1174  const double detJ_weight_coeff_k = detJ_weight * coeff_k;
1175 
1176  value_at_pt.val() = 0.0 ;
1177  gradx_at_pt.val() = 0.0 ;
1178  grady_at_pt.val() = 0.0 ;
1179  gradz_at_pt.val() = 0.0 ;
1180  for ( unsigned m = 0 ; m < FunctionCount ; m++ ) {
1181  value_at_pt.val() += dof_values[m] * bases_vals[m] ;
1182  gradx_at_pt.val() += dof_values[m] * dpsidx[m] ;
1183  grady_at_pt.val() += dof_values[m] * dpsidy[m] ;
1184  gradz_at_pt.val() += dof_values[m] * dpsidz[m] ;
1185  }
1186 
1187  const fad_scalar_type source_term =
1188  coeff_src * value_at_pt * value_at_pt ;
1189 
1190  const fad_scalar_type advection_term =
1191  advection[0]*gradx_at_pt +
1192  advection[1]*grady_at_pt +
1193  advection[2]*gradz_at_pt;
1194 
1195  for ( unsigned m = 0; m < FunctionCount; ++m) {
1196  const double bases_val_m = bases_vals[m] * detJ_weight ;
1197  fad_scalar_type res =
1198  detJ_weight_coeff_k * ( dpsidx[m] * gradx_at_pt +
1199  dpsidy[m] * grady_at_pt +
1200  dpsidz[m] * gradz_at_pt ) +
1201  bases_val_m * ( advection_term + source_term ) ;
1202 
1203  elem_res[m] += res.val();
1204 
1205  scalar_type * const mat = elem_mat[m] ;
1206  for( unsigned n = 0; n < FunctionCount; n++) {
1207  mat[n] += res.fastAccessDx(0) * bases_vals[n] +
1208  res.fastAccessDx(1) * dpsidx[n] +
1209  res.fastAccessDx(2) * dpsidy[n] +
1210  res.fastAccessDx(3) * dpsidz[n];
1211  }
1212  }
1213  }
1214  }
1215 
1217  void operator()( const unsigned ielem ) const
1218  {
1219  double x[ FunctionCount ] ;
1220  double y[ FunctionCount ] ;
1221  double z[ FunctionCount ] ;
1222  unsigned node_index[ ElemNodeCount ];
1223 
1224  scalar_type val[ FunctionCount ] ;
1225  scalar_type elem_res[ FunctionCount ] ;
1226  scalar_type elem_mat[ FunctionCount ][ FunctionCount ] ;
1227 
1228  // Gather nodal coordinates and solution vector:
1229  this->gatherSolution( ielem, val, node_index, x, y, z, elem_res, elem_mat );
1230 
1231  // Compute nodal element residual vector and Jacobian matrix:
1232  computeElementResidualJacobian( val, x, y, z, elem_res, elem_mat );
1233 
1234  // Scatter nodal element residual and Jacobian in global vector and matrix:
1235  this->scatterResidual( ielem, node_index, elem_res, elem_mat );
1236  }
1237 }; /* ElementComputation */
1238 
1239 } /* namespace FENL */
1240 } /* namespace Example */
1241 } /* namespace Kokkos */
1242 
1243 //----------------------------------------------------------------------------
1244 
1245 #endif /* #ifndef KOKKOS_EXAMPLE_FENLFUNCTORS_HPP */
KOKKOS_INLINE_FUNCTION void operator()(const unsigned iwork) const
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 join(volatile unsigned &update, volatile const unsigned &input) const
expr val()
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
Kokkos::View< scalar_type *[FunctionCount], execution_space > elem_vectors_type
static const unsigned element_node_count
Definition: HexElement.hpp:215
CrsMatrix< ScalarType, ExecutionSpace > sparse_matrix_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
static const unsigned function_count
Definition: HexElement.hpp:217
ElementComputation< Kokkos::Example::BoxElemFixture< ExecutionSpace, Order, CoordinateMap >, CrsMatrix< ScalarType, ExecutionSpace >, Analytic > base_type
Kokkos::View< const double *[SpaceDim], Device > node_coord_type
Kokkos::View< const unsigned *[ElemNode], Device > elem_node_type
Kokkos::View< scalar_type *, Kokkos::LayoutLeft, execution_space > vector_type
KOKKOS_INLINE_FUNCTION void init(const TagFillNodeSet &, unsigned &update) const
KOKKOS_INLINE_FUNCTION void join(const TagFillNodeSet &, volatile unsigned &update, volatile const unsigned &input) const
KOKKOS_INLINE_FUNCTION void fill_graph_entries(const unsigned iset) const
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
NodeNodeGraph< elem_node_type, sparse_graph_type, ElemNodeCount >::ElemGraphType elem_graph_type
#define KOKKOS_INLINE_FUNCTION
ElemNodeIdView::execution_space execution_space
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 operator()(const TagFillNodeSet &, unsigned ielem, unsigned &count) const
ElementComputationBase(const ElementComputationBase &rhs)
KOKKOS_INLINE_FUNCTION void operator()(const unsigned irow, unsigned &update, const bool final) const
Kokkos::Example::HexElement_Data< mesh_type::ElemNode > element_data_type
CrsGraphType::row_map_type::non_const_type RowMapType
KOKKOS_INLINE_FUNCTION void sort_graph_entries(const unsigned irow) 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)
CrsMatrix(const StaticCrsGraphType &arg_graph)
Do not initialize the derivative array.
Kokkos::StaticCrsGraph< unsigned, Space, void, unsigned > StaticCrsGraphType
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
static const unsigned integration_count
Definition: HexElement.hpp:216
Kokkos::Example::BoxElemFixture< ExecutionSpace, Order, CoordinateMap > mesh_type
NodeNodeGraph(const ElemNodeIdView &arg_elem_node_id, const unsigned arg_node_count, Times &results)
KOKKOS_INLINE_FUNCTION void fill_elem_graph_map(const unsigned ielem) 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 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 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::View< scalar_type *[FunctionCount][FunctionCount], execution_space > elem_matrices_type
Kokkos::View< unsigned *[ElemNode][ElemNode], execution_space > ElemGraphType
Generate a distributed unstructured finite element mesh from a partitioned NX*NY*NZ box of elements...
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_INLINE_FUNCTION void init(unsigned &update) const
Kokkos::UnorderedMap< key_type, void, execution_space > SetType
pair< unsigned, unsigned > key_type
static const unsigned spatial_dimension
Definition: HexElement.hpp:214