49 #ifndef __INTREPID2_FUNCTIONSPACETOOLS_DEF_HPP__ 
   50 #define __INTREPID2_FUNCTIONSPACETOOLS_DEF_HPP__ 
   55   template<
typename SpT>
 
   56   template<
typename outputValueType, 
class ...outputProperties,
 
   57            typename inputValueType,  
class ...inputProperties>
 
   61                        const Kokkos::DynRankView<inputValueType, inputProperties...>  input ) {
 
   62     if(output.rank() == input.rank()) {
 
   63 #ifdef HAVE_INTREPID2_DEBUG 
   65       for (size_type i=0;i< input.rank();++i) {
 
   66         INTREPID2_TEST_FOR_EXCEPTION( (input.extent(i) != output.extent(i)), std::invalid_argument,
 
   67                                         ">>> ERROR (FunctionSpaceTools::HGRADtransformVALUE): Dimensions of input and output fields containers do not match.");
 
   79   namespace FunctorFunctionSpaceTools {
 
   83     template <
typename OutputViewType,
 
   84               typename jacInverseViewType,
 
   85               typename inputViewType,
 
   86               ordinal_type spaceDim>
 
   88             OutputViewType     _output;
 
   89       const jacInverseViewType  _jacInverse;
 
   90       const inputViewType _input;
 
   93       KOKKOS_INLINE_FUNCTION
 
   95                            jacInverseViewType  jacInverse_,
 
   98           _jacInverse(jacInverse_), 
 
  101       KOKKOS_INLINE_FUNCTION
 
  102       void operator()(
const ordinal_type cl,
 
  103                       const ordinal_type bf,
 
  104                       const ordinal_type pt)
 const {
 
  105          auto y = Kokkos::subview(_output,     cl, bf, pt, Kokkos::ALL());
 
  106         const auto A = Kokkos::subview(_jacInverse, cl,     pt, Kokkos::ALL(), Kokkos::ALL());
 
  107         const auto x = Kokkos::subview(_input,      bf,     pt, Kokkos::ALL());
 
  110           Kernels::Serial::matvec_trans_product_d2( y, A, x );
 
  112           Kernels::Serial::matvec_trans_product_d3( y, A, x );
 
  118   template<
typename SpT>
 
  119   template<
typename outputValValueType,       
class ...outputValProperties,
 
  120            typename jacobianInverseValueType, 
class ...jacobianInverseProperties,
 
  121            typename inputValValueType,        
class ...inputValProperties>
 
  125                       const Kokkos::DynRankView<jacobianInverseValueType,jacobianInverseProperties...> jacobianInverse,
 
  126                       const Kokkos::DynRankView<inputValValueType,       inputValProperties...>        inputVals ) {
 
  127     return HCURLtransformVALUE(outputVals, jacobianInverse, inputVals);
 
  168   template<
typename SpT>
 
  169   template<
typename outputValValueType,       
class ...outputValProperties,
 
  170            typename jacobianInverseValueType, 
class ...jacobianInverseProperties,
 
  171            typename inputValValueType,        
class ...inputValProperties>
 
  175                        const Kokkos::DynRankView<jacobianInverseValueType,jacobianInverseProperties...> jacobianInverse,
 
  176                        const Kokkos::DynRankView<inputValValueType,       inputValProperties...>        inputVals ) {
 
  182   template<
typename SpT>
 
  183   template<
typename outputValValueType,   
class ...outputValProperties,
 
  184            typename jacobianValueType,    
class ...jacobianProperties,
 
  185            typename jacobianDetValueType, 
class ...jacobianDetProperties,
 
  186            typename inputValValueType,    
class ...inputValProperties>
 
  190                       const Kokkos::DynRankView<jacobianValueType,   jacobianProperties...>    jacobian,
 
  191                       const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
 
  192                       const Kokkos::DynRankView<inputValValueType,   inputValProperties...>    inputVals ) {
 
  193     if(jacobian.data()==NULL || jacobian.extent(2)==2) 
 
  194       return HVOLtransformVALUE(outputVals, jacobianDet, inputVals);
 
  196       return HDIVtransformVALUE(outputVals, jacobian, jacobianDet, inputVals);
 
  201   template<
typename SpT>
 
  202   template<
typename outputValValueType,   
class ...outputValProperties,
 
  203            typename jacobianDetValueType, 
class ...jacobianDetProperties,
 
  204            typename inputValValueType,    
class ...inputValProperties>
 
  208                       const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
 
  209                       const Kokkos::DynRankView<inputValValueType,   inputValProperties...>    inputVals ) {
 
  210 #ifdef HAVE_INTREPID2_DEBUG 
  212       INTREPID2_TEST_FOR_EXCEPTION( outputVals.rank() == 4, std::invalid_argument,
 
  213                                     ">>> ERROR (FunctionSpaceTools::HCURLtransformCURL): Output rank must have rank 3.\n If these are 3D fields, then use the appropriate overload of this function.");
 
  216     return HVOLtransformVALUE(outputVals, jacobianDet, inputVals);
 
  221   template<
typename SpT>
 
  222   template<
typename outputValValueType,   
class ...outputValProperties,
 
  223            typename jacobianValueType,    
class ...jacobianProperties,
 
  224            typename jacobianDetValueType, 
class ...jacobianDetProperties,
 
  225            typename inputValValueType,    
class ...inputValProperties>
 
  229                       const Kokkos::DynRankView<jacobianValueType,   jacobianProperties...>    jacobian,
 
  230                       const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
 
  231                       const Kokkos::DynRankView<inputValValueType,   inputValProperties...>    inputVals ) {
 
  232 #ifdef HAVE_INTREPID2_DEBUG 
  234       INTREPID2_TEST_FOR_EXCEPTION( outputVals.extent(3)!=2, std::invalid_argument,
 
  235                                     ">>> ERROR (FunctionSpaceTools::HGRADtransformCURL):\n output field is 3D by the function is meant for 2D fields");
 
  238     return HDIVtransformVALUE(outputVals, jacobian, jacobianDet, inputVals);
 
  243   template<
typename SpT>
 
  244   template<
typename outputValValueType,   
class ...outputValProperties,
 
  245            typename jacobianValueType,    
class ...jacobianProperties,
 
  246            typename jacobianDetValueType, 
class ...jacobianDetProperties,
 
  247            typename inputValValueType,    
class ...inputValProperties>
 
  251                       const Kokkos::DynRankView<jacobianValueType,   jacobianProperties...>    jacobian,
 
  252                       const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
 
  253                       const Kokkos::DynRankView<inputValValueType,   inputValProperties...>    inputVals ) {
 
  260   template<
typename SpT>
 
  261   template<
typename outputValValueType,   
class ...outputValProperties,
 
  262            typename jacobianDetValueType, 
class ...jacobianDetProperties,
 
  263            typename inputValValueType,    
class ...inputValProperties>
 
  266   HDIVtransformDIV(       Kokkos::DynRankView<outputValValueType,  outputValProperties...>   outputVals,
 
  267                     const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
 
  268                     const Kokkos::DynRankView<inputValValueType,   inputValProperties...>    inputVals ) {
 
  269     return HVOLtransformVALUE(outputVals, jacobianDet, inputVals);
 
  274   template<
typename SpT>
 
  275   template<
typename outputValValueType,   
class ...outputValProperties,
 
  276            typename jacobianDetValueType, 
class ...jacobianDetProperties,
 
  277            typename inputValValueType,    
class ...inputValProperties>
 
  281                       const Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
 
  282                       const Kokkos::DynRankView<inputValValueType,   inputValProperties...>    inputVals ) {
 
  288   template<
typename SpT>  
 
  289   template<
typename outputValueValueType, 
class ...outputValueProperties,
 
  290            typename leftValueValueType,   
class ...leftValueProperties,
 
  291            typename rightValueValueType,  
class ...rightValueProperties>
 
  294   integrate(       Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
 
  295              const Kokkos::DynRankView<leftValueValueType,  leftValueProperties...>   leftValues,
 
  296              const Kokkos::DynRankView<rightValueValueType, rightValueProperties...>  rightValues,
 
  297              const bool sumInto ) {
 
  299 #ifdef HAVE_INTREPID2_DEBUG 
  301       INTREPID2_TEST_FOR_EXCEPTION( leftValues.rank() < 2 ||
 
  302                                     leftValues.rank() > 4, std::invalid_argument,
 
  303                                     ">>> ERROR (FunctionSpaceTools::integrate): Left data must have rank 2, 3 or 4.");
 
  304       INTREPID2_TEST_FOR_EXCEPTION( outputValues.rank() < 1 ||
 
  305                                     outputValues.rank() > 3, std::invalid_argument,
 
  306                                     ">>> ERROR (FunctionSpaceTools::integrate): Output values must have rank 1, 2 or 3.");
 
  310     const ordinal_type outRank  = outputValues.rank();
 
  311     const ordinal_type leftRank = leftValues.rank();
 
  312     const ordinal_type mode = outRank*10 + leftRank;
 
  375       INTREPID2_TEST_FOR_EXCEPTION( outRank < 1 || outRank > 3, std::runtime_error,
 
  376                                     ">>> ERROR (FunctionSpaceTools::integrate): outRank must be 1,2, or 3.");
 
  377       INTREPID2_TEST_FOR_EXCEPTION( leftRank < 2 || leftRank > 4, std::runtime_error,
 
  378                                     ">>> ERROR (FunctionSpaceTools::integrate): leftRank must be 1,2, 3 or 4.");
 
  384   namespace FunctorFunctionSpaceTools {
 
  388     template<
typename outputValViewType,
 
  389              typename inputDetViewType,
 
  390              typename inputWeightViewType>
 
  392             outputValViewType   _outputVals;
 
  393       const inputDetViewType    _inputDet;
 
  394       const inputWeightViewType    _inputWeight;
 
  396       KOKKOS_INLINE_FUNCTION
 
  398                             inputDetViewType inputDet_,
 
  399                             inputWeightViewType inputWeight_)
 
  400         : _outputVals(outputVals_), 
 
  401           _inputDet(inputDet_),
 
  402           _inputWeight(inputWeight_) {}
 
  404       typedef ordinal_type value_type;
 
  417       KOKKOS_INLINE_FUNCTION
 
  418       void operator()(
const size_type cl, value_type &dst)
 const {
 
  420         const bool hasNegativeDet = (_inputDet(cl, 0) < 0.0);
 
  421         dst |= hasNegativeDet;
 
  424         const auto sign = (hasNegativeDet ? -1.0 : 1.0);
 
  425         const ordinal_type pt_end = _outputVals.extent(1);
 
  426         for (ordinal_type pt=0;pt<pt_end;++pt) 
 
  427           _outputVals(cl, pt) = sign*_inputDet(cl, pt)*_inputWeight(pt);
 
  432   template<
typename SpT>  
 
  433   template<
typename outputValValueType,   
class ...outputValProperties,
 
  434            typename inputDetValueType,    
class ...inputDetProperties,
 
  435            typename inputWeightValueType, 
class ...inputWeightProperties>
 
  439                       const Kokkos::DynRankView<inputDetValueType,   inputDetProperties...>    inputDet,
 
  440                       const Kokkos::DynRankView<inputWeightValueType,inputWeightProperties...> inputWeights ) {
 
  441 #ifdef HAVE_INTREPID2_DEBUG 
  443       INTREPID2_TEST_FOR_EXCEPTION( inputDet.rank()     != 2 || 
 
  444                                     inputWeights.rank() != 1 || 
 
  445                                     outputVals.rank()   != 2, std::invalid_argument,
 
  446                                     ">>> ERROR (FunctionSpaceTools::computeCellMeasure): Ranks are not compatible.");
 
  447       INTREPID2_TEST_FOR_EXCEPTION( outputVals.extent(0) != inputDet.extent(0), std::invalid_argument,
 
  448                                     ">>> ERROR (FunctionSpaceTools::computeCellMeasure): Cell dimension does not match.");
 
  449       INTREPID2_TEST_FOR_EXCEPTION( inputDet.extent(1)      !=  outputVals.extent(1) ||
 
  450                                     inputWeights.extent(0) !=  outputVals.extent(1), std::invalid_argument,
 
  451                                     ">>> ERROR (FunctionSpaceTools::computeCellMeasure): Point dimension does not match.");
 
  454     typedef          Kokkos::DynRankView<outputValValueType,  outputValProperties...>         outputValViewType;
 
  455     typedef          Kokkos::DynRankView<inputDetValueType,   inputDetProperties...>          inputDetViewType;
 
  456     typedef          Kokkos::DynRankView<inputWeightValueType,inputWeightProperties...>       inputWeightViewType;
 
  458                      <outputValViewType,inputDetViewType,inputWeightViewType> FunctorType;
 
  460     const ordinal_type C = inputDet.extent(0);
 
  461     Kokkos::RangePolicy<SpT,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
 
  463     typename FunctorType::value_type hasNegativeDet = 
false;
 
  464     Kokkos::parallel_reduce( policy, FunctorType(outputVals, inputDet, inputWeights), hasNegativeDet );
 
  466     return hasNegativeDet;
 
  471   template<
typename SpT>
 
  472   template<
typename outputValValueType,   
class ...outputValProperties,
 
  473            typename inputJacValueType,    
class ...inputJacProperties,
 
  474            typename inputWeightValueType, 
class ...inputWeightProperties,
 
  475            typename scratchValueType,     
class ...scratchProperties>
 
  479                       const Kokkos::DynRankView<inputJacValueType,  inputJacProperties...>     inputJac,
 
  480                       const Kokkos::DynRankView<inputWeightValueType,inputWeightProperties...> inputWeights,
 
  481                       const ordinal_type                   whichFace,
 
  482                       const shards::CellTopology  parentCell,
 
  483                       const Kokkos::DynRankView<scratchValueType,    scratchProperties...>     scratch ) {
 
  484 #ifdef HAVE_INTREPID2_DEBUG 
  485     INTREPID2_TEST_FOR_EXCEPTION( inputJac.rank() != 4, std::invalid_argument,
 
  486                                   ">>> ERROR (FunctionSpaceTools::computeFaceMeasure): Input Jacobian container must have rank 4.");
 
  487     INTREPID2_TEST_FOR_EXCEPTION( scratch.rank() != 1, std::invalid_argument,
 
  488                                   ">>> ERROR (FunctionSpaceTools::computeFaceMeasure): Scratch view imust have rank 1.");
 
  489     INTREPID2_TEST_FOR_EXCEPTION( scratch.span() < inputJac.span(), std::invalid_argument,
 
  490                                   ">>> ERROR (FunctionSpaceTools::computeFaceMeasure): Scratch storage must be greater than or equal to inputJac's one.");
 
  498     auto vcprop = Kokkos::common_view_alloc_prop(scratch);
 
  500     typedef Kokkos::DynRankView<scratchValueType, SpT> viewType;
 
  501     viewType faceNormals(Kokkos::view_wrap(scratch.data(), vcprop),
 
  518   template<
typename SpT>
 
  519   template<
typename outputValValueType,   
class ...outputValProperties,
 
  520            typename inputJacValueType,    
class ...inputJacProperties,
 
  521            typename inputWeightValueType, 
class ...inputWeightProperties,
 
  522            typename scratchValueType,     
class ...scratchProperties>
 
  526                       const Kokkos::DynRankView<inputJacValueType,  inputJacProperties...>   inputJac,
 
  527                       const Kokkos::DynRankView<inputWeightValueType,inputWeightProperties...> inputWeights,
 
  528                       const ordinal_type                   whichEdge,
 
  529                       const shards::CellTopology  parentCell,
 
  530                       const Kokkos::DynRankView<scratchValueType,    scratchProperties...>    scratch ) {
 
  531 #ifdef HAVE_INTREPID2_DEBUG 
  532     INTREPID2_TEST_FOR_EXCEPTION( (inputJac.rank() != 4), std::invalid_argument,
 
  533                                   ">>> ERROR (FunctionSpaceTools::computeEdgeMeasure): Input Jacobian container must have rank 4.");
 
  534     INTREPID2_TEST_FOR_EXCEPTION( scratch.rank() != 1, std::invalid_argument,
 
  535                                   ">>> ERROR (FunctionSpaceTools::computeEdgeMeasure): Scratch view must have a rank 1.");
 
  536     INTREPID2_TEST_FOR_EXCEPTION( scratch.span() < inputJac.span(), std::invalid_argument,
 
  537                                   ">>> ERROR (FunctionSpaceTools::computeEdgeMeasure): Scratch storage must be greater than or equal to inputJac'one.");
 
  545     auto vcprop = Kokkos::common_view_alloc_prop(scratch);
 
  547     typedef Kokkos::DynRankView<scratchValueType, SpT> viewType;
 
  548     viewType edgeTangents(Kokkos::view_wrap(scratch.data(), vcprop),
 
  565   template<
typename SpT>
 
  566   template<
typename outputValValueType,    
class ...outputValProperties,
 
  567            typename inputMeasureValueType, 
class ...inputMeasureProperties,
 
  568            typename inputValValueType,     
class ...inputValProperteis>
 
  571   multiplyMeasure(       Kokkos::DynRankView<outputValValueType,   outputValProperties...>    outputVals,
 
  572                    const Kokkos::DynRankView<inputMeasureValueType,inputMeasureProperties...> inputMeasure,
 
  573                    const Kokkos::DynRankView<inputValValueType,    inputValProperteis...>     inputVals ) {
 
  574     scalarMultiplyDataField( outputVals, 
 
  581   template<
typename SpT>
 
  582   template<
typename outputFieldValueType, 
class ...outputFieldProperties,
 
  583            typename inputDataValueType,   
class ...inputDataProperties,
 
  584            typename inputFieldValueType,  
class ...inputFieldProperties>
 
  588                            const Kokkos::DynRankView<inputDataValueType,  inputDataProperties...>   inputData,
 
  589                            const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...>  inputFields,
 
  590                            const bool reciprocal ) {
 
  599   template<
typename SpT>
 
  600   template<
typename outputDataValuetype,     
class ...outputDataProperties,
 
  601            typename inputDataLeftValueType,  
class ...inputDataLeftProperties,
 
  602            typename inputDataRightValueType, 
class ...inputDataRightProperties>
 
  606                           const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...>  inputDataLeft,
 
  607                           const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight,
 
  608                           const bool reciprocal ) {
 
  617   template<
typename SpT>
 
  618   template<
typename outputFieldValueType, 
class ...outputFieldProperties,
 
  619            typename inputDataValueType,   
class ...inputDataProperties,
 
  620            typename inputFieldValueType,  
class ...inputFieldProperties>
 
  624                         const Kokkos::DynRankView<inputDataValueType,  inputDataProperties...>   inputData,
 
  625                         const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...>  inputFields ) {
 
  633   template<
typename SpT>
 
  634   template<
typename outputDataValueType,     
class ...outputDataProperties,
 
  635            typename inputDataLeftValueType,  
class ...inputDataLeftProperties,
 
  636            typename inputDataRightValueType, 
class ...inputDataRightProperties>
 
  640                        const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...>  inputDataLeft,
 
  641                        const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight ) {
 
  649   template<
typename SpT>
 
  650   template<
typename outputFieldValueType, 
class ...outputFieldProperties,
 
  651            typename inputDataValueType,   
class ...inputDataProperties,
 
  652            typename inputFieldValueType,  
class ...inputFieldProperties>
 
  656                            const Kokkos::DynRankView<inputDataValueType,  inputDataProperties...>   inputData,
 
  657                            const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...>  inputFields ) {
 
  658     const auto outRank = outputFields.rank();
 
  672       INTREPID2_TEST_FOR_EXCEPTION( outRank < 3 && outRank > 5, std::runtime_error,
 
  673                                     ">>> ERROR (FunctionSpaceTools::vectorMultiplyDataField): Output container must have rank 3, 4 or 5.");
 
  680   template<
typename SpT>
 
  681   template<
typename outputDataValueType,     
class ...outputDataProperties,
 
  682            typename inputDataLeftValueType,  
class ...inputDataLeftProperties,
 
  683            typename inputDataRightValueType, 
class ...inputDataRightProperties>
 
  687                           const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...>  inputDataLeft,
 
  688                           const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight ) {
 
  689     const auto outRank = outputData.rank();
 
  703       INTREPID2_TEST_FOR_EXCEPTION( outRank < 2 && outRank > 4, std::runtime_error,
 
  704                                     ">>> ERROR (FunctionSpaceTools::vectorMultiplyDataData): Output container must have rank 2, 3 or 4.");
 
  711   template<
typename SpT>
 
  712   template<
typename outputFieldValueType, 
class ...outputFieldProperties,
 
  713            typename inputDataValueType,   
class ...inputDataProperties,
 
  714            typename inputFieldValueType,  
class ...inputFieldProperties>
 
  718                            const Kokkos::DynRankView<inputDataValueType,  inputDataProperties...>   inputData,
 
  719                            const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...>  inputFields,
 
  720                            const char transpose ) {
 
  722     const auto outRank = outputFields.rank();
 
  737       INTREPID2_TEST_FOR_EXCEPTION( outRank < 4 && outRank > 5, std::runtime_error,
 
  738                                     ">>> ERROR (FunctionSpaceTools::tensorMultiplyDataField): Output container must have rank 4 or 5.");
 
  745   template<
typename SpT>
 
  746   template<
typename outputDataValueType,     
class ...outputDataProperties,
 
  747            typename inputDataLeftValueType,  
class ...inputDataLeftProperties,
 
  748            typename inputDataRightValueType, 
class ...inputDataRightProperties>
 
  752                           const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...>  inputDataLeft,
 
  753                           const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight,
 
  754                           const char transpose ) {
 
  755     const auto outRank = outputData.rank();
 
  770       INTREPID2_TEST_FOR_EXCEPTION( outRank < 4 && outRank > 5, std::runtime_error,
 
  771                                     ">>> ERROR (FunctionSpaceTools::tensorMultiplyDataField): Output container must have rank 4 or 5.");
 
  778   namespace FunctorFunctionSpaceTools {
 
  783     template<
typename inoutOperatorViewType,
 
  784              typename fieldSignViewType>
 
  786             inoutOperatorViewType _inoutOperator;
 
  787       const fieldSignViewType     _fieldSigns;
 
  789       KOKKOS_INLINE_FUNCTION
 
  791                              fieldSignViewType     fieldSigns_ )
 
  792         : _inoutOperator(inoutOperator_), _fieldSigns(fieldSigns_) {}
 
  794       KOKKOS_INLINE_FUNCTION
 
  795       void operator()(
const ordinal_type cl)
 const {
 
  796         const ordinal_type nlbf = _inoutOperator.extent(1);
 
  797         const ordinal_type nrbf = _inoutOperator.extent(2);
 
  799         for (ordinal_type lbf=0;lbf<nlbf;++lbf)
 
  800           for (ordinal_type rbf=0;rbf<nrbf;++rbf)
 
  801             _inoutOperator(cl, lbf, rbf) *= _fieldSigns(cl, lbf);
 
  806   template<
typename SpT>  
 
  807   template<
typename inoutOperatorValueType, 
class ...inoutOperatorProperties,
 
  808            typename fieldSignValueType,     
class ...fieldSignProperties>
 
  811   applyLeftFieldSigns(       Kokkos::DynRankView<inoutOperatorValueType,inoutOperatorProperties...> inoutOperator,
 
  812                        const Kokkos::DynRankView<fieldSignValueType,    fieldSignProperties...>     fieldSigns ) {
 
  814 #ifdef HAVE_INTREPID2_DEBUG 
  815     INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.rank() != 3, std::invalid_argument,
 
  816                                   ">>> ERROR (FunctionSpaceTools::applyLeftFieldSigns): Input operator container must have rank 3.");
 
  817     INTREPID2_TEST_FOR_EXCEPTION( fieldSigns.rank() != 2, std::invalid_argument,
 
  818                                   ">>> ERROR (FunctionSpaceTools::applyLeftFieldSigns): Input field signs container must have rank 2.");
 
  819     INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.extent(0) != fieldSigns.extent(0), std::invalid_argument,
 
  820                                   ">>> ERROR (FunctionSpaceTools::applyLeftFieldSigns): Zeroth dimensions (number of cells) of the operator and field signs containers must agree!");
 
  821     INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.extent(1) != fieldSigns.extent(1), std::invalid_argument,
 
  822                                   ">>> ERROR (FunctionSpaceTools::applyLeftFieldSigns): First dimensions (number of left fields) of the operator and field signs containers must agree!");
 
  825     typedef          Kokkos::DynRankView<inoutOperatorValueType,inoutOperatorProperties...>        inoutOperatorViewType;
 
  826     typedef          Kokkos::DynRankView<fieldSignValueType,    fieldSignProperties...>            fieldSignViewType;
 
  828                  <inoutOperatorViewType,fieldSignViewType>                                     FunctorType;
 
  829     typedef typename ExecSpace<typename inoutOperatorViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
 
  831     const ordinal_type C = inoutOperator.extent(0);
 
  832     Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
 
  833     Kokkos::parallel_for( policy, FunctorType(inoutOperator, fieldSigns) );
 
  838   namespace FunctorFunctionSpaceTools {
 
  842     template<
typename inoutOperatorViewType,
 
  843              typename fieldSignViewType>
 
  845         inoutOperatorViewType _inoutOperator;
 
  846       const fieldSignViewType     _fieldSigns;
 
  848       KOKKOS_INLINE_FUNCTION
 
  850                               fieldSignViewType     fieldSigns_ )
 
  851         : _inoutOperator(inoutOperator_), _fieldSigns(fieldSigns_) {}
 
  853       KOKKOS_INLINE_FUNCTION
 
  854       void operator()(
const ordinal_type cl)
 const {
 
  855         const ordinal_type nlbf = _inoutOperator.extent(1);
 
  856         const ordinal_type nrbf = _inoutOperator.extent(2);
 
  858         for (ordinal_type lbf=0;lbf<nlbf;++lbf)
 
  859           for (ordinal_type rbf=0;rbf<nrbf;++rbf)
 
  860             _inoutOperator(cl, lbf, rbf) *= _fieldSigns(cl, rbf);
 
  865   template<
typename SpT>  
 
  866   template<
typename inoutOperatorValueType, 
class ...inoutOperatorProperties,
 
  867            typename fieldSignValueType,     
class ...fieldSignProperties>
 
  871                         const Kokkos::DynRankView<fieldSignValueType,    fieldSignProperties...>     fieldSigns ) {
 
  873 #ifdef HAVE_INTREPID2_DEBUG 
  874     INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.rank() != 3, std::invalid_argument,
 
  875                                   ">>> ERROR (FunctionSpaceTools::applyRightFieldSigns): Input operator container must have rank 3.");
 
  876     INTREPID2_TEST_FOR_EXCEPTION( fieldSigns.rank() != 2, std::invalid_argument,
 
  877                                   ">>> ERROR (FunctionSpaceTools::applyRightFieldSigns): Input field signs container must have rank 2.");
 
  878     INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.extent(0) != fieldSigns.extent(0), std::invalid_argument,
 
  879                                   ">>> ERROR (FunctionSpaceTools::applyRightFieldSigns): Zeroth dimensions (number of cells) of the operator and field signs containers must agree!");
 
  880     INTREPID2_TEST_FOR_EXCEPTION( inoutOperator.extent(2) != fieldSigns.extent(1), std::invalid_argument,
 
  881                                   ">>> ERROR (FunctionSpaceTools::applyRightFieldSigns): Second dimension of the operator container and first dimension of the field signs container (number of right fields) must agree!");
 
  884     typedef          Kokkos::DynRankView<inoutOperatorValueType,inoutOperatorProperties...>        inoutOperatorViewType;
 
  885     typedef          Kokkos::DynRankView<fieldSignValueType,    fieldSignProperties...>            fieldSignViewType;
 
  887                      <inoutOperatorViewType,fieldSignViewType>                                     FunctorType;
 
  888     typedef typename ExecSpace<typename inoutOperatorViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
 
  890     const ordinal_type C = inoutOperator.extent(0);
 
  891     Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
 
  892     Kokkos::parallel_for( policy, FunctorType(inoutOperator, fieldSigns) );
 
  897   namespace FunctorFunctionSpaceTools {
 
  901     template<
typename inoutFunctionViewType,
 
  902              typename fieldSignViewType>
 
  904             inoutFunctionViewType _inoutFunction;
 
  905       const fieldSignViewType     _fieldSigns;
 
  907       KOKKOS_INLINE_FUNCTION
 
  909                          fieldSignViewType     fieldSigns_)
 
  910         : _inoutFunction(inoutFunction_), _fieldSigns(fieldSigns_) {}
 
  912       KOKKOS_INLINE_FUNCTION
 
  913       void operator()(
const ordinal_type cl)
 const {
 
  914         const ordinal_type nbfs = _inoutFunction.extent(1);
 
  915         const ordinal_type npts = _inoutFunction.extent(2);
 
  916         const ordinal_type iend = _inoutFunction.extent(3);
 
  917         const ordinal_type jend = _inoutFunction.extent(4);
 
  919         for (ordinal_type bf=0;bf<nbfs;++bf) 
 
  920           for (ordinal_type pt=0;pt<npts;++pt)
 
  921             for (ordinal_type i=0;i<iend;++i) 
 
  922               for (ordinal_type j=0;j<jend;++j) 
 
  923                 _inoutFunction(cl, bf, pt, i, j) *= _fieldSigns(cl, bf);
 
  928   template<
typename SpT>  
 
  929   template<
typename inoutFunctionValueType, 
class ...inoutFunctionProperties,
 
  930            typename fieldSignValueType,     
class ...fieldSignProperties>
 
  933   applyFieldSigns(       Kokkos::DynRankView<inoutFunctionValueType,inoutFunctionProperties...> inoutFunction,
 
  934                    const Kokkos::DynRankView<fieldSignValueType,    fieldSignProperties...>     fieldSigns ) {
 
  936 #ifdef HAVE_INTREPID2_DEBUG 
  937     INTREPID2_TEST_FOR_EXCEPTION( inoutFunction.rank() < 2 || inoutFunction.rank() > 5, std::invalid_argument,
 
  938                                   ">>> ERROR (FunctionSpaceTools::applyFieldSigns): Input function container must have rank 2, 3, 4, or 5.");
 
  939     INTREPID2_TEST_FOR_EXCEPTION( fieldSigns.rank() != 2, std::invalid_argument,
 
  940                                   ">>> ERROR (FunctionSpaceTools::applyFieldSigns): Input field signs container must have rank 2.");
 
  941     INTREPID2_TEST_FOR_EXCEPTION( inoutFunction.extent(0) != fieldSigns.extent(0), std::invalid_argument,
 
  942                                   ">>> ERROR (FunctionSpaceTools::applyFieldSigns): Zeroth dimensions (number of integration domains) of the function and field signs containers must agree!");
 
  943     INTREPID2_TEST_FOR_EXCEPTION( inoutFunction.extent(1) != fieldSigns.extent(1), std::invalid_argument,
 
  944                                   ">>> ERROR (FunctionSpaceTools::applyFieldSigns): First dimensions (number of fields) of the function and field signs containers must agree!");
 
  948     typedef          Kokkos::DynRankView<inoutFunctionValueType,inoutFunctionProperties...>        inoutFunctionViewType;
 
  949     typedef          Kokkos::DynRankView<fieldSignValueType,    fieldSignProperties...>            fieldSignViewType;
 
  951                      <inoutFunctionViewType,fieldSignViewType>                                     FunctorType;
 
  952     typedef typename ExecSpace<typename inoutFunctionViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
 
  954     const ordinal_type C = inoutFunction.extent(0);
 
  955     Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
 
  956     Kokkos::parallel_for( policy, FunctorType(inoutFunction, fieldSigns) );
 
  961   namespace FunctorFunctionSpaceTools {
 
  966     template<
typename outputPointViewType,
 
  967              typename inputCoeffViewType,
 
  968              typename inputFieldViewType>
 
  970             outputPointViewType _outputPointVals;
 
  971       const inputCoeffViewType  _inputCoeffs;
 
  972       const inputFieldViewType  _inputFields;
 
  974       KOKKOS_INLINE_FUNCTION
 
  975       F_evaluate( outputPointViewType outputPointVals_,
 
  976                   inputCoeffViewType  inputCoeffs_,
 
  977                   inputFieldViewType  inputFields_ )
 
  978         : _outputPointVals(outputPointVals_), _inputCoeffs(inputCoeffs_), _inputFields(inputFields_) {}
 
  980       KOKKOS_INLINE_FUNCTION
 
  981       void operator()(
const ordinal_type cl)
 const {
 
  982         const ordinal_type nbfs = _inputFields.extent(1);
 
  983         const ordinal_type npts = _inputFields.extent(2);
 
  985         const ordinal_type iend = _inputFields.extent(3);
 
  986         const ordinal_type jend = _inputFields.extent(4);
 
  988         for (ordinal_type bf=0;bf<nbfs;++bf) 
 
  989           for (ordinal_type pt=0;pt<npts;++pt)
 
  990             for (ordinal_type i=0;i<iend;++i) 
 
  991               for (ordinal_type j=0;j<jend;++j) 
 
  992                 _outputPointVals(cl, pt, i, j) += _inputCoeffs(cl, bf) * _inputFields(cl, bf, pt, i, j);
 
  997   template<
typename SpT>    
 
  998   template<
typename outputPointValueType, 
class ...outputPointProperties,
 
  999            typename inputCoeffValueType,  
class ...inputCoeffProperties,
 
 1000            typename inputFieldValueType,  
class ...inputFieldProperties>
 
 1003   evaluate(       Kokkos::DynRankView<outputPointValueType,outputPointProperties...> outputPointVals,
 
 1004             const Kokkos::DynRankView<inputCoeffValueType, inputCoeffProperties...>  inputCoeffs,
 
 1005             const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...>  inputFields ) {
 
 1007 #ifdef HAVE_INTREPID2_DEBUG 
 1008     INTREPID2_TEST_FOR_EXCEPTION( inputFields.rank() < 3 || inputFields.rank() > 5, std::invalid_argument,
 
 1009                                   ">>> ERROR (FunctionSpaceTools::evaluate): Input fields container must have rank 3, 4, or 5.");
 
 1010     INTREPID2_TEST_FOR_EXCEPTION( inputCoeffs.rank() != 2, std::invalid_argument,
 
 1011                                   ">>> ERROR (FunctionSpaceTools::evaluate): Input coefficient container must have rank 2.");
 
 1012     INTREPID2_TEST_FOR_EXCEPTION( outputPointVals.rank() != (inputFields.rank()-1), std::invalid_argument,
 
 1013                                   ">>> ERROR (FunctionSpaceTools::evaluate): Output values container must have rank one less than the rank of the input fields container.");
 
 1014     INTREPID2_TEST_FOR_EXCEPTION( inputCoeffs.extent(0) != inputFields.extent(0), std::invalid_argument,
 
 1015                                   ">>> ERROR (FunctionSpaceTools::evaluate): Zeroth dimensions (number of cells) of the coefficient and fields input containers must agree!");
 
 1016     INTREPID2_TEST_FOR_EXCEPTION( inputCoeffs.extent(1) != inputFields.extent(1), std::invalid_argument,
 
 1017                                   ">>> ERROR (FunctionSpaceTools::evaluate): First dimensions (number of fields) of the coefficient and fields input containers must agree!");
 
 1018     INTREPID2_TEST_FOR_EXCEPTION( outputPointVals.extent(0) != inputFields.extent(0), std::invalid_argument,
 
 1019                                   ">>> ERROR (FunctionSpaceTools::evaluate): Zeroth dimensions (number of cells) of the input fields container and the output values container must agree!");
 
 1020     for (size_type i=1;i<outputPointVals.rank();++i)
 
 1021       INTREPID2_TEST_FOR_EXCEPTION( outputPointVals.extent(i) != inputFields.extent(i+1), std::invalid_argument, 
 
 1022                                     ">>> ERROR (FunctionSpaceTools::evaluate): outputPointVals dimension(i) does not match to inputFields dimension(i+1).");
 
 1025     typedef          Kokkos::DynRankView<outputPointValueType,outputPointProperties...>         outputPointValViewType;
 
 1026     typedef          Kokkos::DynRankView<inputCoeffValueType, inputCoeffProperties...>          inputCoeffViewType;
 
 1027     typedef          Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...>          inputFieldViewType; 
 
 1029                      <outputPointValViewType,inputCoeffViewType,inputFieldViewType>             FunctorType;
 
 1030     typedef typename ExecSpace<typename inputCoeffViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
 
 1032     const ordinal_type C = inputFields.extent(0);
 
 1033     Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, C);
 
 1034     Kokkos::parallel_for( policy, FunctorType(outputPointVals, inputCoeffs, inputFields) );