CubatureControlVolume<SpT,PT,WT>::
CubatureControlVolume(const shards::CellTopology cellTopology) {

    // define primary cell topology with given one
    primaryCellTopo_ = cellTopology;

    // subcell is defined either hex or quad according to dimension
    switch (primaryCellTopo_.getDimension()) {
    case 2:
        subcvCellTopo_ = shards::CellTopology(shards::getCellTopologyData<shards::Quadrilateral<4> >());
        break;
    case 3:
        subcvCellTopo_ = shards::CellTopology(shards::getCellTopologyData<shards::Hexahedron<8> >());
        break;
    }

    // computation order is always one;
    degree_ = 1;

    // create subcell cubature points and weights and cache them
    const ordinal_type subcvDegree = 2;
    auto subcvCubature = DefaultCubatureFactory::create<SpT,PT,WT>(subcvCellTopo_, subcvDegree);

    const auto numSubcvPoints = subcvCubature->getNumPoints();
    const auto subcvDim       = subcvCubature->getDimension();

    subcvCubaturePoints_  = Kokkos::DynRankView<PT,SpT>("CubatureControlVolume::subcvCubaturePoints_",
                            numSubcvPoints, subcvDim);
    subcvCubatureWeights_ = Kokkos::DynRankView<WT,SpT>("CubatureControlVolume::subcvCubatureWeights_",
                            numSubcvPoints);

    subcvCubature->getCubature(subcvCubaturePoints_,
                               subcvCubatureWeights_);
}
void IntegrationValues2<Scalar>::
evaluateValues(const PHX::MDField<Scalar,Cell,NODE,Dim>& in_node_coordinates,
               const PHX::MDField<Scalar,Cell,IP,Dim>& other_ip_coordinates)
{
    getCubature(in_node_coordinates);

    {
        // Determine the permutation.
        std::vector<size_type> permutation(other_ip_coordinates.dimension(1));
        permuteToOther(ip_coordinates, other_ip_coordinates, permutation);
        // Apply the permutation to the cubature arrays.
        MDFieldArrayFactory af(prefix, alloc_arrays);
        const size_type num_ip = dyn_cub_points.dimension(0);
        {
            const size_type num_dim = dyn_side_cub_points.dimension(1);
            DblArrayDynamic old_dyn_side_cub_points = af.template buildArray<double,IP,Dim>(
                        "old_dyn_side_cub_points", num_ip, num_dim);
            old_dyn_side_cub_points.deep_copy(dyn_side_cub_points);
            for (size_type ip = 0; ip < num_ip; ++ip)
                if (ip != permutation[ip])
                    for (size_type dim = 0; dim < num_dim; ++dim)
                        dyn_side_cub_points(ip, dim) = old_dyn_side_cub_points(permutation[ip], dim);
        }
        {
            const size_type num_dim = dyn_cub_points.dimension(1);
            DblArrayDynamic old_dyn_cub_points = af.template buildArray<double,IP,Dim>(
                    "old_dyn_cub_points", num_ip, num_dim);
            old_dyn_cub_points.deep_copy(dyn_cub_points);
            for (size_type ip = 0; ip < num_ip; ++ip)
                if (ip != permutation[ip])
                    for (size_type dim = 0; dim < num_dim; ++dim)
                        dyn_cub_points(ip, dim) = old_dyn_cub_points(permutation[ip], dim);
        }
        {
            DblArrayDynamic old_dyn_cub_weights = af.template buildArray<double,IP>(
                    "old_dyn_cub_weights", num_ip);
            old_dyn_cub_weights.deep_copy(dyn_cub_weights);
            for (size_type ip = 0; ip < dyn_cub_weights.dimension(0); ++ip)
                if (ip != permutation[ip])
                    dyn_cub_weights(ip) = old_dyn_cub_weights(permutation[ip]);
        }
        {
            const size_type num_cells = ip_coordinates.dimension(0), num_ip = ip_coordinates.dimension(1),
                            num_dim = ip_coordinates.dimension(2);
            Array_CellIPDim old_ip_coordinates = af.template buildStaticArray<Scalar,Cell,IP,Dim>(
                    "old_ip_coordinates", num_cells, num_ip, num_dim);
            Kokkos::deep_copy(old_ip_coordinates.get_kokkos_view(), ip_coordinates.get_kokkos_view());
            for (size_type cell = 0; cell < num_cells; ++cell)
                for (size_type ip = 0; ip < num_ip; ++ip)
                    if (ip != permutation[ip])
                        for (size_type dim = 0; dim < num_dim; ++dim)
                            ip_coordinates(cell, ip, dim) = old_ip_coordinates(cell, permutation[ip], dim);
        }
        // All subsequent calculations inherit the permutation.
    }

    evaluateRemainingValues(in_node_coordinates);
}
void IntegrationValues2<Scalar>::
evaluateValues(const PHX::MDField<Scalar,Cell,NODE,Dim> & in_node_coordinates)
{
    if (int_rule->cv_type != "none") {
        evaluateValuesCV(in_node_coordinates);
    }
    else {
        getCubature(in_node_coordinates);
        evaluateRemainingValues(in_node_coordinates);
    }
}
Exemple #4
0
    int FunctionSpaceTools_Test03(const bool verbose) {

      Teuchos::RCP<std::ostream> outStream;
      Teuchos::oblackholestream bhs; // outputs nothing

      if (verbose)
        outStream = Teuchos::rcp(&std::cout, false);
      else
        outStream = Teuchos::rcp(&bhs, false);

      Teuchos::oblackholestream oldFormatState;
      oldFormatState.copyfmt(std::cout);

      typedef typename
        Kokkos::Impl::is_space<DeviceSpaceType>::host_mirror_space::execution_space HostSpaceType ;

      *outStream << "DeviceSpace::  "; DeviceSpaceType::print_configuration(*outStream, false);
      *outStream << "HostSpace::    ";   HostSpaceType::print_configuration(*outStream, false);
      
      *outStream
        << "===============================================================================\n"
        << "|                                                                             |\n"
        << "|                      Unit Test (FunctionSpaceTools)                         |\n"
        << "|                                                                             |\n"
        << "|     1) Basic operator transformations and integration in HDIV:              |\n"
        << "|                                                                             |\n"
        << "|  Questions? Contact  Pavel Bochev ([email protected]) or                   |\n"
        << "|                      Denis Ridzal ([email protected]) or                   |\n"
        << "|                      Kyungjoo Kim ([email protected]).                      |\n"
        << "|                                                                             |\n"
        << "|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n"
        << "|  Trilinos website:   http://trilinos.sandia.gov                             |\n"
        << "|                                                                             |\n"
        << "===============================================================================\n";

      typedef CellTools<DeviceSpaceType> ct;      
      typedef FunctionSpaceTools<DeviceSpaceType> fst;      
      typedef Kokkos::DynRankView<ValueType,DeviceSpaceType> DynRankView;

      int errorFlag = 0;

      *outStream
        << "\n"
        << "===============================================================================\n"
        << "| TEST 1: correctness of math operations                                      |\n"
        << "===============================================================================\n";
      
      outStream->precision(20);
      
      try {
        DefaultCubatureFactory cub_factory;                            

        shards::CellTopology cell_topo = shards::getCellTopologyData< shards::Hexahedron<8> >();    

        const auto cub_degree = 20;                                                               
        auto cub = cub_factory.create<DeviceSpaceType,ValueType,ValueType>(cell_topo, cub_degree); 

        const auto space_dim = cub->getDimension();                                              
        const auto num_cub_points = cub->getNumPoints();                                          

        Basis_HDIV_HEX_I1_FEM<DeviceSpaceType> basis;                        
        const auto num_fields = basis.getCardinality();                                       
        
        /* Cell geometries and orientations. */
        const auto num_cells = 4;
        const auto num_nodes = 8;
        
        const ValueType hexnodes[] = {
          // hex 0  -- affine
          -1.0, -1.0, -1.0,
          1.0, -1.0, -1.0,
          1.0, 1.0, -1.0,
          -1.0, 1.0, -1.0,
          -1.0, -1.0, 1.0,
          1.0, -1.0, 1.0,
          1.0, 1.0, 1.0,
          -1.0, 1.0, 1.0,
          // hex 1  -- affine
          -3.0, -3.0, 1.0,
          6.0, 3.0, 1.0,
          7.0, 8.0, 0.0,
          -2.0, 2.0, 0.0,
          -3.0, -3.0, 4.0,
          6.0, 3.0, 4.0,
          7.0, 8.0, 3.0,
          -2.0, 2.0, 3.0,
          // hex 2  -- affine
          -3.0, -3.0, 0.0,
          9.0, 3.0, 0.0,
          15.0, 6.1, 0.0,
          3.0, 0.1, 0.0,
          9.0, 3.0, 0.1,
          21.0, 9.0, 0.1,
          27.0, 12.1, 0.1,
          15.0, 6.1, 0.1,
          // hex 3  -- nonaffine
          -2.0, -2.0, 0.0,
          2.0, -1.0, 0.0,
          1.0, 6.0, 0.0,
          -1.0, 1.0, 0.0,
          0.0, 0.0, 1.0,
          1.0, 0.0, 1.0,
          1.0, 1.0, 1.0,
          0.0, 1.0, 1.0
        };
        
        const ValueType facesigns[] = {
          1, 1, 1, 1, 1, 1,
          1, -1, 1, -1, 1, -1,
          -1, -1, 1, 1, -1, 1,
          -1, -1, 1, 1, -1, -1
        };
        
        /* Computational arrays. */
        DynRankView ConstructWithLabel( cub_points,  num_cub_points, space_dim);
        DynRankView ConstructWithLabel( cub_weights, num_cub_points);

        DynRankView ConstructWithLabel( cell_nodes,       num_cells, num_nodes, space_dim);
        DynRankView ConstructWithLabel( field_signs,      num_cells, num_fields);

        DynRankView ConstructWithLabel( jacobian,         num_cells, num_cub_points, space_dim, space_dim);
        DynRankView ConstructWithLabel( jacobian_det,     num_cells, num_cub_points);
        DynRankView ConstructWithLabel( weighted_measure, num_cells, num_cub_points);
        
        DynRankView ConstructWithLabel( div_of_basis_at_cub_points,                                 num_fields, num_cub_points);
        DynRankView ConstructWithLabel( transformed_div_of_basis_at_cub_points,          num_cells, num_fields, num_cub_points);
        DynRankView ConstructWithLabel( weighted_transformed_div_of_basis_at_cub_points, num_cells, num_fields, num_cub_points);
        DynRankView ConstructWithLabel( stiffness_matrices,                              num_cells, num_fields, num_fields);
        
        DynRankView ConstructWithLabel( value_of_basis_at_cub_points,                                 num_fields, num_cub_points, space_dim);
        DynRankView ConstructWithLabel( transformed_value_of_basis_at_cub_points,          num_cells, num_fields, num_cub_points, space_dim);
        DynRankView ConstructWithLabel( weighted_transformed_value_of_basis_at_cub_points, num_cells, num_fields, num_cub_points, space_dim);
        DynRankView ConstructWithLabel( mass_matrices,                                     num_cells, num_fields, num_fields);

        /******************* START COMPUTATION ***********************/

        // get cubature points and weights
        cub->getCubature(cub_points, cub_weights);

        const Kokkos::DynRankView<const ValueType,Kokkos::LayoutRight,Kokkos::HostSpace> cell_nodes_host (&hexnodes[0],  num_cells, num_nodes, space_dim); 
        const Kokkos::DynRankView<const ValueType,Kokkos::LayoutRight,Kokkos::HostSpace> field_signs_host(&facesigns[0], num_cells, num_fields); 
        
        Kokkos::deep_copy( cell_nodes,  cell_nodes_host  );
        Kokkos::deep_copy( field_signs, field_signs_host );
        
        // compute geometric cell information
        ct::setJacobian(jacobian, cub_points, cell_nodes, cell_topo);
        ct::setJacobianDet(jacobian_det, jacobian);

        // compute weighted measure
        fst::computeCellMeasure(weighted_measure, jacobian_det, cub_weights);

        // **Computing stiffness matrices:
        basis.getValues(div_of_basis_at_cub_points, cub_points, OPERATOR_DIV);
        
        // transform divergences of basis functions
        fst::HDIVtransformDIV(transformed_div_of_basis_at_cub_points,
                              jacobian_det,
                              div_of_basis_at_cub_points);
        
        // multiply with weighted measure
        fst::multiplyMeasure(weighted_transformed_div_of_basis_at_cub_points,
                             weighted_measure,
                             transformed_div_of_basis_at_cub_points);
        
        // we can apply the field signs to the basis function arrays, or after the fact, see below
        fst::applyFieldSigns(transformed_div_of_basis_at_cub_points, field_signs);
        fst::applyFieldSigns(weighted_transformed_div_of_basis_at_cub_points, field_signs);
        
        // compute stiffness matrices
        fst::integrate(stiffness_matrices,
                       transformed_div_of_basis_at_cub_points,
                       weighted_transformed_div_of_basis_at_cub_points);
        

        // **Computing mass matrices:
        basis.getValues(value_of_basis_at_cub_points, cub_points, OPERATOR_VALUE);

            
        
        // transform values of basis functions
        fst::HDIVtransformVALUE(transformed_value_of_basis_at_cub_points,
                                jacobian,
                                jacobian_det,
                                value_of_basis_at_cub_points);

        // multiply with weighted measure
        fst::multiplyMeasure(weighted_transformed_value_of_basis_at_cub_points,
                             weighted_measure,
                             transformed_value_of_basis_at_cub_points);
        
        // compute mass matrices
        fst::integrate(mass_matrices,
                       transformed_value_of_basis_at_cub_points,
                       weighted_transformed_value_of_basis_at_cub_points);
        
        // apply field signs
        fst::applyLeftFieldSigns(mass_matrices, field_signs);
        fst::applyRightFieldSigns(mass_matrices, field_signs);
        
        /*******************  STOP COMPUTATION ***********************/
        
        
        /******************* START COMPARISON ***********************/
        std::string basedir = "../testdata";
        for (auto cid=0;cid<num_cells-1;++cid) {
          std::stringstream namestream;
          std::string filename;
          namestream <<  basedir << "/mass_HDIV_HEX_I1_FEM" << "_" << "0" << cid+1 << ".dat";
          namestream >> filename;

          *outStream << "\nCell ID : " << cid << "  mass matrix comparing with " << filename << "\n\n";
          
          std::ifstream massfile(&filename[0]);
          if (massfile.is_open()) {
            const auto mass_matrix_cell = Kokkos::subdynrankview(mass_matrices, cid, Kokkos::ALL(), Kokkos::ALL());
            errorFlag += compareToAnalytic(massfile, 
                                           mass_matrix_cell, 
                                           1e-10,
                                           verbose);
            massfile.close();
          } else {
            errorFlag = -1;
            INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error,
                                          "Failed to open a file" );
          }
          
          namestream.clear();
          namestream << basedir << "/stiff_HDIV_HEX_I1_FEM" << "_" << "0" << cid+1 << ".dat";
          namestream >> filename;

          *outStream << "\nCell ID : " << cid << "  stiffness matrix comparing with " << filename << "\n\n";
          
          std::ifstream stifffile(&filename[0]);
          if (stifffile.is_open()) {
            const auto stiffness_matrix_cell = Kokkos::subdynrankview(stiffness_matrices, cid, Kokkos::ALL(), Kokkos::ALL());
            errorFlag += compareToAnalytic(stifffile,
                                           stiffness_matrix_cell,
                                           1e-10,
                                           verbose);
            stifffile.close();
          } else {
            errorFlag = -1;
            INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error,
                                          "Failed to open a file" );
          }
        }

        for (auto cid=3;cid<num_cells;++cid) {
          std::stringstream namestream;
          std::string filename;
          namestream <<  basedir << "/mass_fp_HDIV_HEX_I1_FEM" << "_" << "0" << cid+1 << ".dat";
          namestream >> filename;

          *outStream << "\nCell ID : " << cid << "  mass matrix comparing with " << filename << "\n\n";
          
          std::ifstream massfile(&filename[0]);
          if (massfile.is_open()) {
            const auto mass_matrix_cell = Kokkos::subdynrankview(mass_matrices, cid, Kokkos::ALL(), Kokkos::ALL());
            errorFlag += compareToAnalytic(massfile, 
                                           mass_matrix_cell, 
                                           1e-4, 
                                           verbose,
                                           INTREPID2_UTILS_SCALAR);
            massfile.close();
          } else {
            errorFlag = -1;
            INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error,
                                          "Failed to open a file" );
          }
          
          namestream.clear();
          namestream << basedir << "/stiff_fp_HDIV_HEX_I1_FEM" << "_" << "0" << cid+1 << ".dat";
          namestream >> filename;
          
          *outStream << "\nCell ID : " << cid << "  stiffness matrix comparing with " << filename << "\n\n";

          std::ifstream stifffile(&filename[0]);
          if (stifffile.is_open()) {
            const auto stiffness_matrix_cell = Kokkos::subdynrankview(stiffness_matrices, cid, Kokkos::ALL(), Kokkos::ALL());
            errorFlag += compareToAnalytic(stifffile,
                                           stiffness_matrix_cell,
                                           1e-4, 
                                           verbose,
                                           INTREPID2_UTILS_SCALAR);
            stifffile.close();
          } else {
            errorFlag = -1;
            INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error,
                                          "Failed to open a file" );
          }
        }
        /******************* STOP COMPARISON ***********************/

        *outStream << "\n";
      } catch (std::logic_error err) {
        *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
        *outStream << err.what() << '\n';
        *outStream << "-------------------------------------------------------------------------------" << "\n\n";
        errorFlag = -1000;
      }

      if (errorFlag != 0)
        std::cout << "End Result: TEST FAILED\n";
      else
        std::cout << "End Result: TEST PASSED\n";
      
      // reset format state of std::cout
      std::cout.copyfmt(oldFormatState);

      return errorFlag;
    }
 void IntegrationValues2<Scalar>::
 evaluateValues(const PHX::MDField<Scalar,Cell,NODE,Dim> & in_node_coordinates)
 {
   getCubature(in_node_coordinates);
   evaluateRemainingValues(in_node_coordinates);
 }
    int ComputeBasis_HGRAD_Vector(const ordinal_type nworkset,
                                  const ordinal_type C,
                                  const ordinal_type order,
                                  const bool verbose) {
      typedef Vector<VectorTagType> VectorType;

      typedef typename VectorTagType::value_type ValueType;
      constexpr int VectorLength = VectorTagType::length;

      Teuchos::RCP<std::ostream> verboseStream;
      Teuchos::oblackholestream bhs; // outputs nothing

      if (verbose) 
        verboseStream = Teuchos::rcp(&std::cout, false);
      else
        verboseStream = Teuchos::rcp(&bhs,       false);

      Teuchos::oblackholestream oldFormatState;
      oldFormatState.copyfmt(std::cout);

      typedef typename
        Kokkos::Impl::is_space<DeviceSpaceType>::host_mirror_space::execution_space HostSpaceType ;

      *verboseStream << "DeviceSpace::  "; DeviceSpaceType::print_configuration(*verboseStream, false);
      *verboseStream << "HostSpace::    ";   HostSpaceType::print_configuration(*verboseStream, false);
      *verboseStream << "VectorLength::  " << (VectorLength) << "\n";

      using BasisTypeHost = Basis_HGRAD_HEX_C1_FEM<HostSpaceType,ValueType,ValueType>;
      using ImplBasisType = Impl::Basis_HGRAD_HEX_C1_FEM;
      using range_type = Kokkos::pair<ordinal_type,ordinal_type>;

      constexpr size_t LLC_CAPACITY = 32*1024*1024;
      Intrepid2::Test::Flush<LLC_CAPACITY,DeviceSpaceType> flush;
      
      Kokkos::Impl::Timer timer;
      double t_vectorize = 0;
      int errorFlag = 0;

      BasisTypeHost hostBasis;
      const auto cellTopo = hostBasis.getBaseCellTopology();

      auto cubature = DefaultCubatureFactory::create<DeviceSpaceType,ValueType,ValueType>(cellTopo, order);

      const ordinal_type 
        numCells = C,
        numCellsAdjusted = C/VectorLength + (C%VectorLength > 0),
        numVerts = cellTopo.getVertexCount(),
        numDofs = hostBasis.getCardinality(),
        numPoints = cubature->getNumPoints(), 
        spaceDim = cubature->getDimension();

      Kokkos::DynRankView<ValueType,HostSpaceType> dofCoordsHost("dofCoordsHost", numDofs, spaceDim);
      hostBasis.getDofCoords(dofCoordsHost);
      const auto refNodesHost = Kokkos::subview(dofCoordsHost, range_type(0, numVerts), Kokkos::ALL());
      
      // pertub nodes
      Kokkos::DynRankView<VectorType,HostSpaceType> worksetCellsHost("worksetCellsHost", numCellsAdjusted, numVerts, spaceDim);
      for (ordinal_type cell=0;cell<numCells;++cell) {
        for (ordinal_type i=0;i<numVerts;++i)
          for (ordinal_type j=0;j<spaceDim;++j) {
            ValueType val = (rand()/(RAND_MAX + 1.0))*0.2 -0.1;
            worksetCellsHost(cell/VectorLength, i, j)[cell%VectorLength] = refNodesHost(i, j) + val; 
          }
      }

      auto worksetCells = Kokkos::create_mirror_view(typename DeviceSpaceType::memory_space(), worksetCellsHost);
      Kokkos::deep_copy(worksetCells, worksetCellsHost);

      Kokkos::DynRankView<ValueType,DeviceSpaceType> refPoints("refPoints", numPoints, spaceDim), refWeights("refWeights", numPoints);
      cubature->getCubature(refPoints, refWeights);

      std::cout
        << "===============================================================================\n" 
        << " Performance Test evaluating ComputeBasis \n"
        << " # of workset = " << nworkset << "\n" 
        << " Test Array Structure (C,F,P,D) = " << numCells << ", " << numDofs << ", " << numPoints << ", " << spaceDim << "\n"
        << "===============================================================================\n";

      *verboseStream
        << "\n"
        << "===============================================================================\n"
        << "TEST 1: evaluateFields vector version\n"
        << "===============================================================================\n";
      
      try {
        Kokkos::DynRankView<ValueType,DeviceSpaceType> 
          refBasisValues("refBasisValues", numDofs, numPoints),
          refBasisGrads ("refBasisGrads",  numDofs, numPoints, spaceDim);
        
        ImplBasisType::getValues<DeviceSpaceType>(refBasisValues, refPoints, OPERATOR_VALUE);
        ImplBasisType::getValues<DeviceSpaceType>(refBasisGrads,  refPoints, OPERATOR_GRAD);
        
        const ordinal_type ibegin = -3;

        // testing vertical approach
        {          
          Kokkos::DynRankView<VectorType,DeviceSpaceType> 
            weightedBasisValues("weightedBasisValues", numCellsAdjusted, numDofs, numPoints),
            weightedBasisGrads ("weightedBasisGrads",  numCellsAdjusted, numDofs, numPoints, spaceDim);

          typedef F_hgrad_eval<VectorType,ValueType,DeviceSpaceType> FunctorType;

          using range_policy_type = Kokkos::Experimental::MDRangePolicy
            < DeviceSpaceType, Kokkos::Experimental::Rank<2>, Kokkos::IndexType<ordinal_type> >;
          range_policy_type policy( {                0,         0 },
                                    { numCellsAdjusted, numPoints } );

          FunctorType functor(weightedBasisValues,
                              weightedBasisGrads,
                              refBasisGrads,
                              worksetCells,
                              refWeights,
                              refBasisValues,
                              refBasisGrads);
          
          for (ordinal_type iwork=ibegin;iwork<nworkset;++iwork) {
            flush.run();

            DeviceSpaceType::fence();
            timer.reset();
            
            Kokkos::parallel_for(policy, functor);

            DeviceSpaceType::fence();
            t_vectorize += (iwork >= 0)*timer.seconds();
          }

        }

      } catch (std::exception err) {
        *verboseStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
        *verboseStream << err.what() << '\n';
        *verboseStream << "-------------------------------------------------------------------------------" << "\n\n";
        errorFlag = -1000;
      }

      std::cout 
        << "TEST HGRAD " 
        << " t_vectorize = " << (t_vectorize/nworkset)
        << std::endl;
      
      if (errorFlag != 0)
        std::cout << "End Result: TEST FAILED\n";
      else
        std::cout << "End Result: TEST PASSED\n";
      
      // reset format state of std::cout
      std::cout.copyfmt(oldFormatState);
      
      return errorFlag;
    }
Exemple #7
0
    int FunctionSpaceTools_Test04(const bool verbose) {

      Teuchos::RCP<std::ostream> outStream;
      Teuchos::oblackholestream bhs; // outputs nothing

      if (verbose)
        outStream = Teuchos::rcp(&std::cout, false);
      else
        outStream = Teuchos::rcp(&bhs, false);

      Teuchos::oblackholestream oldFormatState;
      oldFormatState.copyfmt(std::cout);

      typedef typename
        Kokkos::Impl::is_space<DeviceSpaceType>::host_mirror_space::execution_space HostSpaceType ;

      *outStream << "DeviceSpace::  "; DeviceSpaceType::print_configuration(*outStream, false);
      *outStream << "HostSpace::    ";   HostSpaceType::print_configuration(*outStream, false);

      *outStream                                                       
        << "===============================================================================\n"
        << "|                                                                             |\n"
        << "|                      Unit Test (FunctionSpaceTools)                         |\n"
        << "|                                                                             |\n"
        << "|     1) volume integration on tetrahedra, testing dataIntegral               |\n"
        << "|                                                                             |\n"
        << "|  Questions? Contact  Pavel Bochev ([email protected]) or                   |\n"
        << "|                      Denis Ridzal ([email protected]) or                   |\n"
        << "|                      Kyungjoo Kim ([email protected]).                      |\n"
        << "|                                                                             |\n"
        << "|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n"
        << "|  Trilinos website:   http://trilinos.sandia.gov                             |\n"
        << "|                                                                             |\n"
        << "===============================================================================\n";

      typedef CellTools<DeviceSpaceType> ct;
      typedef FunctionSpaceTools<DeviceSpaceType> fst;
      typedef Kokkos::DynRankView<ValueType,DeviceSpaceType> DynRankView;

      const auto tol = tolerence<ValueType>();

      int errorFlag = 0;

      *outStream
        << "\n"
        << "===============================================================================\n"
        << "| TEST 1: correctness of cell volumes                                         |\n"
        << "===============================================================================\n";

      outStream->precision(20);

      try {
        DefaultCubatureFactory cub_factory;

        shards::CellTopology cell_topo = shards::getCellTopologyData< shards::Tetrahedron<4> >();

        const auto cub_degree = 0;
        auto cub = cub_factory.create<DeviceSpaceType,ValueType,ValueType>(cell_topo, cub_degree);

        const auto space_dim = cub->getDimension();
        const auto num_cub_points = cub->getNumPoints();

        /* Cell geometries. */
        const auto num_cells = 4;
        const auto num_nodes = 4;

        const ValueType tetnodes[] = {
          // tet 0
          0.0, 0.0, 0.0,
          1.0, 0.0, 0.0,
          0.0, 1.0, 0.0,
          0.0, 0.0, 1.0,
          // tet 1
          4.0, 5.0, 1.0,
          -6.0, 2.0, 0.0,
          4.0, -3.0, -1.0,
          0.0, 2.0, 5.0,
          // tet 2
          -6.0, -3.0, 1.0,
          9.0, 2.0, 1.0,
          8.9, 2.1, 0.9,
          8.9, 2.1, 1.1,
          // tet 3
          -6.0, -3.0, 1.0,
          12.0, 3.0, 1.0,
          2.9, 0.1, 0.9,
          2.9, 0.1, 1.1
        };
        
        /* Analytic volumes. */
        const ValueType tetvols[] = {1.0/6.0, 194.0/3.0, 1.0/15.0, 2.0/25.0};

      /* Computational arrays. */
        DynRankView ConstructWithLabel( cub_points,  num_cub_points, space_dim);
        DynRankView ConstructWithLabel( cub_weights, num_cub_points);

        DynRankView ConstructWithLabel( cell_nodes,       num_cells, num_nodes, space_dim);

        DynRankView ConstructWithLabel( jacobian,         num_cells, num_cub_points, space_dim, space_dim);
        DynRankView ConstructWithLabel( jacobian_det,     num_cells, num_cub_points);
        DynRankView ConstructWithLabel( weighted_measure, num_cells, num_cub_points);

        DynRankView ConstructWithLabel( data_one, num_cells, num_cub_points);
        DynRankView ConstructWithLabel( volumes, num_cells);

        /******************* START COMPUTATION ***********************/

        // get cubature points and weights
        cub->getCubature(cub_points, cub_weights);

        const Kokkos::DynRankView<const ValueType,Kokkos::LayoutRight,HostSpaceType> cell_nodes_host (&tetnodes[0],  num_cells, num_nodes, space_dim);
        Kokkos::deep_copy( cell_nodes,  cell_nodes_host  );
        
        // compute geometric cell information
        ct::setJacobian(jacobian, cub_points, cell_nodes, cell_topo);
        ct::setJacobianDet(jacobian_det, jacobian);

        // compute weighted measure
        fst::computeCellMeasure(weighted_measure, jacobian_det, cub_weights);

        Kokkos::deep_copy(data_one, 1.0);

        // compute volumes
        fst::integrate(volumes, data_one, weighted_measure);
        
        /******************* STOP COMPUTATION ***********************/
        
        // memcpy do implicit sync
        const auto volumes_host = Kokkos::create_mirror_view(typename HostSpaceType::memory_space(), volumes);
        Kokkos::deep_copy(volumes_host, volumes);

        /******************* START COMPARISON ***********************/
        for (auto cid=0;cid<num_cells-1;++cid) {
          *outStream << "Volume of cell " << cid 
                     << " = " << std::setw(24) << volumes_host(cid) 
                     << "    vs.    Analytic value =  " << std::setw(24) << tetvols[cid] << "\n";
          errorFlag += (std::fabs(volumes_host(cid)-tetvols[cid]) > tol);
          
        }
        /******************* STOP COMPARISON ***********************/
        
        *outStream << "\n";
      } catch (std::logic_error err) {
        *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
        *outStream << err.what() << '\n';
        *outStream << "-------------------------------------------------------------------------------" << "\n\n";
        errorFlag = -1000;
      }

      if (errorFlag != 0)
        std::cout << "End Result: TEST FAILED\n";
      else
        std::cout << "End Result: TEST PASSED\n";

      // reset format state of std::cout
      std::cout.copyfmt(oldFormatState);

      return errorFlag;
    }