KOKKOS_INLINE_FUNCTION
    void
    Basis_HGRAD_LINE_Cn_FEM::Serial<opType>::
    getValues(       outputViewType output,
               const inputViewType  input,
                     workViewType   work,
               const vinvViewType   vinv,
               const ordinal_type   operatorDn ) {    
      ordinal_type opDn = operatorDn;

      const ordinal_type card = vinv.extent(0);
      const ordinal_type npts = input.extent(0);

      const ordinal_type order = card - 1;
      const double alpha = 0.0, beta = 0.0;

      typedef typename Kokkos::DynRankView<typename workViewType::value_type, typename workViewType::memory_space> viewType;
      auto vcprop = Kokkos::common_view_alloc_prop(work);
      
      switch (opType) {
      case OPERATOR_VALUE: {
        viewType phis(Kokkos::view_wrap(work.data(), vcprop), card, npts);     

        Impl::Basis_HGRAD_LINE_Cn_FEM_JACOBI::
          Serial<opType>::getValues(phis, input, order, alpha, beta);

        for (ordinal_type i=0;i<card;++i) 
          for (ordinal_type j=0;j<npts;++j) {
            output.access(i,j) = 0.0;
            for (ordinal_type k=0;k<card;++k)
              output.access(i,j) += vinv(k,i)*phis.access(k,j);
          }
        break;
      }
      case OPERATOR_GRAD:
      case OPERATOR_D1:
      case OPERATOR_D2:
      case OPERATOR_D3:
      case OPERATOR_D4:
      case OPERATOR_D5:
      case OPERATOR_D6:
      case OPERATOR_D7:
      case OPERATOR_D8:
      case OPERATOR_D9:
      case OPERATOR_D10: 
        opDn = getOperatorOrder(opType);
      case OPERATOR_Dn: {
        // dkcard is always 1 for 1D element
        const ordinal_type dkcard = 1;
        viewType phis(Kokkos::view_wrap(work.data(), vcprop), card, npts, dkcard);     
        Impl::Basis_HGRAD_LINE_Cn_FEM_JACOBI::
          Serial<opType>::getValues(phis, input, order, alpha, beta, opDn);

        for (ordinal_type i=0;i<card;++i) 
          for (ordinal_type j=0;j<npts;++j) 
            for (ordinal_type k=0;k<dkcard;++k) {
              output.access(i,j,k) = 0.0;
              for (ordinal_type l=0;l<card;++l)
                output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
            }
        break;
      }
      default: {
        INTREPID2_TEST_FOR_ABORT( true,
                                  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::Serial::getValues) operator is not supported." );
      }
      }
    }
    KOKKOS_INLINE_FUNCTION                                                                                                                                               
    void                                                                                                                                                                 
    Basis_HGRAD_TET_C2_FEM::Serial<opType>::                                                                                                                             
    getValues(       outputViewType output,                                                                                                                              
               const inputViewType input ) {   
      switch (opType) {
      case OPERATOR_VALUE: {
        const auto x = input(0);
        const auto y = input(1);
        const auto z = input(2);
        
        // output is a rank-2 array with dimensions (basisCardinality_, dim0)
        output.access(0) = (-1. + x + y + z)*(-1. + 2.*x + 2.*y + 2.*z);
        output.access(1) = x*(-1. + 2.*x);
        output.access(2) = y*(-1. + 2.*y);
        output.access(3) = z*(-1. + 2.*z);

        output.access(4) = -4.*x*(-1. + x + y + z);
        output.access(5) =  4.*x*y;
        output.access(6) = -4.*y*(-1. + x + y + z);
        output.access(7) = -4.*z*(-1. + x + y + z);
        output.access(8) =  4.*x*z;
        output.access(9) =  4.*y*z;
        break;
      }
      case OPERATOR_D1:
      case OPERATOR_GRAD: {
        const auto x = input(0);
        const auto y = input(1);
        const auto z = input(2);
        
        // output.access is a rank-3 array with dimensions (basisCardinality_, dim0, spaceDim)
        output.access(0, 0) = -3.+ 4.*x + 4.*y + 4.*z;
        output.access(0, 1) = -3.+ 4.*x + 4.*y + 4.*z;
        output.access(0, 2) = -3.+ 4.*x + 4.*y + 4.*z; 
      
        output.access(1, 0) = -1.+ 4.*x; 
        output.access(1, 1) =  0.;
        output.access(1, 2) =  0.;
      
        output.access(2, 0) =  0.;        
        output.access(2, 1) = -1.+ 4.*y;
        output.access(2, 2) =  0.;

        output.access(3, 0) =  0.;         
        output.access(3, 1) =  0.;
        output.access(3, 2) = -1.+ 4.*z;
 
      
        output.access(4, 0) = -4.*(-1.+ 2*x + y + z);         
        output.access(4, 1) = -4.*x;
        output.access(4, 2) = -4.*x;
      
        output.access(5, 0) =  4.*y;       
        output.access(5, 1) =  4.*x; 
        output.access(5, 2) =  0.;

        output.access(6, 0) = -4.*y;          
        output.access(6, 1) = -4.*(-1.+ x + 2*y + z);
        output.access(6, 2) = -4.*y; 

        output.access(7, 0) = -4.*z;          
        output.access(7, 1) = -4.*z;
        output.access(7, 2) = -4.*(-1.+ x + y + 2*z);

        output.access(8, 0) =  4.*z;     
        output.access(8, 1) =  0.;
        output.access(8, 2) =  4.*x;

        output.access(9, 0) =  0.;         
        output.access(9, 1) =  4.*z;
        output.access(9, 2) =  4.*y;
        break;
      }
      case OPERATOR_D2: {
        output.access(0, 0) =  4.;
        output.access(0, 1) =  4.;
        output.access(0, 2) =  4.;
        output.access(0, 3) =  4.;
        output.access(0, 4) =  4.;
        output.access(0, 5) =  4.;
      
        output.access(1, 0) =  4.;
        output.access(1, 1) =  0.;
        output.access(1, 2) =  0.;
        output.access(1, 3) =  0.;
        output.access(1, 4) =  0.;
        output.access(1, 5) =  0.;

        output.access(2, 0) =  0.;
        output.access(2, 1) =  0.;
        output.access(2, 2) =  0.;
        output.access(2, 3) =  4.;
        output.access(2, 4) =  0.;
        output.access(2, 5) =  0.;

        output.access(3, 0) =  0.;
        output.access(3, 1) =  0.;
        output.access(3, 2) =  0.;
        output.access(3, 3) =  0.;
        output.access(3, 4) =  0.;
        output.access(3, 5) =  4.;

        output.access(4, 0) = -8.;
        output.access(4, 1) = -4.;
        output.access(4, 2) = -4.;
        output.access(4, 3) =  0.;
        output.access(4, 4) =  0.;
        output.access(4, 5) =  0.;

        output.access(5, 0) =  0.;
        output.access(5, 1) =  4.;
        output.access(5, 2) =  0.;
        output.access(5, 3) =  0.;
        output.access(5, 4) =  0.;
        output.access(5, 5) =  0.;

        output.access(6, 0) =  0.;
        output.access(6, 1) = -4.;
        output.access(6, 2) =  0.;
        output.access(6, 3) = -8.;
        output.access(6, 4) = -4.;
        output.access(6, 5) =  0;

        output.access(7, 0) =  0.;
        output.access(7, 1) =  0.;
        output.access(7, 2) = -4.;
        output.access(7, 3) =  0.;
        output.access(7, 4) = -4.;
        output.access(7, 5) = -8.;

        output.access(8, 0) =  0.;
        output.access(8, 1) =  0.;
        output.access(8, 2) =  4.;
        output.access(8, 3) =  0.;
        output.access(8, 4) =  0.;
        output.access(8, 5) =  0.;

        output.access(9, 0) =  0.;
        output.access(9, 1) =  0.;
        output.access(9, 2) =  0.;
        output.access(9, 3) =  0.;
        output.access(9, 4) =  4.;
        output.access(9, 5) =  0.;
        break;
      }
      case OPERATOR_MAX: {
        const ordinal_type jend = output.extent(1);
        const ordinal_type iend = output.extent(0);

        for (ordinal_type j=0;j<jend;++j)
          for (ordinal_type i=0;i<iend;++i)
            output.access(i, j) = 0.0;
        break;
      }
      default: {
        INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
                                  opType != OPERATOR_GRAD &&
                                  opType != OPERATOR_D1 &&
                                  opType != OPERATOR_D2 &&
                                  opType != OPERATOR_MAX,
                                  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::Serial::getValues) operator is not supported");
      }
      }
    }