void IntegrationValues2<Scalar>::
  getCubatureCV(const PHX::MDField<Scalar,Cell,NODE,Dim>& in_node_coordinates)
  {
    int num_space_dim = int_rule->topology->getDimension();
    if (int_rule->isSide() && num_space_dim==1) {
       std::cout << "WARNING: 0-D quadrature rule infrastructure does not exist!!! Will not be able to do "
                 << "non-natural integration rules.";
       return; 
    }
     {
      size_type num_cells = in_node_coordinates.dimension(0);
      size_type num_nodes = in_node_coordinates.dimension(1);
      size_type num_dims = in_node_coordinates.dimension(2);

      for (size_type cell = 0; cell < num_cells;  ++cell) {
        for (size_type node = 0; node < num_nodes; ++node) {
          for (size_type dim = 0; dim < num_dims; ++dim) {
            node_coordinates(cell,node,dim) =
                in_node_coordinates(cell,node,dim);
            dyn_node_coordinates(cell,node,dim) =
                Sacado::ScalarValue<Scalar>::eval(in_node_coordinates(cell,node,dim));
          }
        }
      }
     }

    if (int_rule->cv_type == "side")
      intrepid_cubature->getCubature(dyn_phys_cub_points.get_view(),dyn_phys_cub_norms.get_view(),dyn_node_coordinates.get_view());
    else 
      intrepid_cubature->getCubature(dyn_phys_cub_points.get_view(),dyn_phys_cub_weights.get_view(),dyn_node_coordinates.get_view());

    size_type num_cells = dyn_phys_cub_points.dimension(0);
    size_type num_ip =dyn_phys_cub_points.dimension(1); 
    size_type num_dims = dyn_phys_cub_points.dimension(2);

    for (size_type cell = 0; cell < num_cells;  ++cell) {
        for (size_type ip = 0; ip < num_ip;  ++ip) {
           if (int_rule->cv_type != "side")
               weighted_measure(cell,ip) = dyn_phys_cub_weights(cell,ip);
           for (size_type dim = 0; dim < num_dims; ++dim) {
               ip_coordinates(cell,ip,dim) = dyn_phys_cub_points(cell,ip,dim);
               if (int_rule->cv_type == "side")
                  weighted_normals(cell,ip,dim) = dyn_phys_cub_norms(cell,ip,dim);
           }
        }
     }

  }
Ejemplo n.º 2
0
void IntegrationValues2<Scalar>::
evaluateValuesCV(const PHX::MDField<Scalar,Cell,NODE,Dim> & in_node_coordinates)
{

    Intrepid2::CellTools<Scalar> cell_tools;

    {
        size_type num_cells = in_node_coordinates.dimension(0);
        size_type num_nodes = in_node_coordinates.dimension(1);
        size_type num_dims = in_node_coordinates.dimension(2);

        for (size_type cell = 0; cell < num_cells;  ++cell) {
            for (size_type node = 0; node < num_nodes; ++node) {
                for (size_type dim = 0; dim < num_dims; ++dim) {
                    node_coordinates(cell,node,dim) =
                        in_node_coordinates(cell,node,dim);
                    dyn_node_coordinates(cell,node,dim) =
                        Sacado::ScalarValue<Scalar>::eval(in_node_coordinates(cell,node,dim));
                }
            }
        }
    }

    if (int_rule->cv_type == "volume")
        intrepid_cubature->getCubature(dyn_phys_cub_points,dyn_phys_cub_weights,dyn_node_coordinates);

    else if (int_rule->cv_type == "side")
        intrepid_cubature->getCubature(dyn_phys_cub_points,dyn_phys_cub_norms,dyn_node_coordinates);

    if (int_rule->cv_type == "volume")
    {
        size_type num_cells = dyn_phys_cub_points.dimension(0);
        size_type num_ip = dyn_phys_cub_points.dimension(1);
        size_type num_dims = dyn_phys_cub_points.dimension(2);

        for (size_type cell = 0; cell < num_cells;  ++cell) {
            for (size_type ip = 0; ip < num_ip;  ++ip) {
                weighted_measure(cell,ip) = dyn_phys_cub_weights(cell,ip);
                for (size_type dim = 0; dim < num_dims; ++dim)
                    ip_coordinates(cell,ip,dim) = dyn_phys_cub_points(cell,ip,dim);
            }
        }
        cell_tools.mapToReferenceFrame(ref_ip_coordinates, ip_coordinates, node_coordinates,
                                       *(int_rule->topology),-1);

        cell_tools.setJacobian(jac, ref_ip_coordinates, node_coordinates,
                               *(int_rule->topology));

    }
    else if (int_rule->cv_type == "side")
    {
        size_type num_cells = dyn_phys_cub_points.dimension(0);
        size_type num_ip = dyn_phys_cub_points.dimension(1);
        size_type num_dims = dyn_phys_cub_points.dimension(2);

        for (size_type cell = 0; cell < num_cells;  ++cell) {
            for (size_type ip = 0; ip < num_ip;  ++ip) {
                for (size_type dim = 0; dim < num_dims; ++dim) {
                    ip_coordinates(cell,ip,dim) = dyn_phys_cub_points(cell,ip,dim);
                    weighted_normals(cell,ip,dim) = dyn_phys_cub_norms(cell,ip,dim);
                }
            }
        }

        cell_tools.mapToReferenceFrame(ref_ip_coordinates, ip_coordinates, node_coordinates,
                                       *(int_rule->topology),-1);
        cell_tools.setJacobian(jac, ref_ip_coordinates, node_coordinates,
                               *(int_rule->topology));
    }

    cell_tools.setJacobianInv(jac_inv, jac);

    cell_tools.setJacobianDet(jac_det, jac);

}
Ejemplo n.º 3
0
    int FunctionSpaceTools_Test02(const bool verbose) {
      typedef ValueType value_type;

      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);
      
      *outStream \
        << "===============================================================================\n" \
        << "|                                                                             |\n" \
        << "|                      Unit Test (FunctionSpaceTools)                         |\n" \
        << "|                                                                             |\n" \
        << "|     1) basic operator transformations and integration in HCURL              |\n" \
        << "|                                                                             |\n" \
        << "|  Questions? Contact  Pavel Bochev ([email protected]) or                   |\n" \
        << "|                      Denis Ridzal ([email protected]).                     |\n" \
        << "|                                                                             |\n" \
        << "|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
        << "|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
        << "|                                                                             |\n" \
        << "===============================================================================\n";

      typedef FunctionSpaceTools<DeviceSpaceType> fst;
      typedef Kokkos::DynRankView<value_type,DeviceSpaceType> DynRankView;
#define ConstructWithLabel(obj, ...) obj(#obj, __VA_ARGS__)
      
      int errorFlag = 0;

      *outStream \
        << "\n"
        << "===============================================================================\n"\
        << "| TEST 1: correctness of math operations                                      |\n"\
        << "===============================================================================\n";

      outStream->precision(20);

      try {
        shards::CellTopology cellType = shards::getCellTopologyData< shards::Hexahedron<> >();    // cell type: hex

        /* Related to cubature. */
        DefaultCubatureFactory<double> cubFactory;                                                // create cubature factory
        int cubDegree = 20;                                                                       // cubature degree
        Teuchos::RCP<Cubature<double> > myCub = cubFactory.create(cellType, cubDegree);           // create default cubature
        int spaceDim = myCub->getDimension();                                                     // get spatial dimension 
        int numCubPoints = myCub->getNumPoints();                                                 // get number of cubature points
        /* Related to basis. */
        Basis_HCURL_HEX_I1_FEM<double, FieldContainer<double> > hexBasis;                         // create H-curl basis on a hex
        int numFields = hexBasis.getCardinality();                                                // get basis cardinality
 
        /* Cell geometries and orientations. */
        int numCells    = 4;
        int numNodes    = 8;
        int numCellData = numCells*numNodes*spaceDim;
        int numSignData = numCells*numFields;

        double 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
        };

        double edgesigns[] = {
          1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
          1, -1, 1, -1, 1, -1, 1, -1, 1, -1, 1, -1,
          -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. */
        FieldContainer<double> cub_points(numCubPoints, spaceDim);
        FieldContainer<double> cub_weights(numCubPoints);
        FieldContainer<double> cell_nodes(numCells, numNodes, spaceDim);
        FieldContainer<double>  field_signs(numCells, numFields);
        FieldContainer<double> jacobian(numCells, numCubPoints, spaceDim, spaceDim);
        FieldContainer<double> jacobian_inv(numCells, numCubPoints, spaceDim, spaceDim);
        FieldContainer<double> jacobian_det(numCells, numCubPoints);
        FieldContainer<double> weighted_measure(numCells, numCubPoints);

        FieldContainer<double> curl_of_basis_at_cub_points(numFields, numCubPoints, spaceDim);
        FieldContainer<double> transformed_curl_of_basis_at_cub_points(numCells, numFields, numCubPoints, spaceDim);
        FieldContainer<double> weighted_transformed_curl_of_basis_at_cub_points(numCells, numFields, numCubPoints, spaceDim);
        FieldContainer<double> stiffness_matrices(numCells, numFields, numFields);

        FieldContainer<double> value_of_basis_at_cub_points(numFields, numCubPoints, spaceDim);
        FieldContainer<double> transformed_value_of_basis_at_cub_points(numCells, numFields, numCubPoints, spaceDim);
        FieldContainer<double> weighted_transformed_value_of_basis_at_cub_points(numCells, numFields, numCubPoints, spaceDim);
        FieldContainer<double> mass_matrices(numCells, numFields, numFields);

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

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

        // fill cell vertex array
        cell_nodes.setValues(hexnodes, numCellData);
        // set basis function signs, for each cell
        field_signs.setValues(edgesigns, numSignData);

        // compute geometric cell information
        CellTools<double>::setJacobian(jacobian, cub_points, cell_nodes, cellType);
        CellTools<double>::setJacobianInv(jacobian_inv, jacobian);
        CellTools<double>::setJacobianDet(jacobian_det, jacobian);
        // compute weighted measure
        fst::computeCellMeasure<double>(weighted_measure, jacobian_det, cub_weights);

        // Computing stiffness matrices:
        // tabulate curls of basis functions at (reference) cubature points
        hexBasis.getValues(curl_of_basis_at_cub_points, cub_points, OPERATOR_CURL);

        // transform curls of basis functions 
        fst::HCURLtransformCURL<double>(transformed_curl_of_basis_at_cub_points,
                                        jacobian,
                                        jacobian_det,
                                        curl_of_basis_at_cub_points);
        // multiply with weighted measure
        fst::multiplyMeasure<double>(weighted_transformed_curl_of_basis_at_cub_points,
                                     weighted_measure,
                                     transformed_curl_of_basis_at_cub_points);

        // we can apply the field signs to the basis function arrays, or after the fact, see below
        fst::applyFieldSigns<double>(transformed_curl_of_basis_at_cub_points, field_signs);
        fst::applyFieldSigns<double>(weighted_transformed_curl_of_basis_at_cub_points, field_signs);

        // compute stiffness matrices
        fst::integrate<double>(stiffness_matrices,
                               transformed_curl_of_basis_at_cub_points,
                               weighted_transformed_curl_of_basis_at_cub_points,
                               COMP_CPP);
        // Computing mass matrices:
        // tabulate values of basis functions at (reference) cubature points
        hexBasis.getValues(value_of_basis_at_cub_points, cub_points, OPERATOR_VALUE);

        // transform values of basis functions 
        fst::HCURLtransformVALUE<double>(transformed_value_of_basis_at_cub_points,
                                         jacobian_inv,
                                         value_of_basis_at_cub_points);

        // multiply with weighted measure
        fst::multiplyMeasure<double>(weighted_transformed_value_of_basis_at_cub_points,
                                     weighted_measure,
                                     transformed_value_of_basis_at_cub_points);

        // compute mass matrices
        fst::integrate<double>(mass_matrices,
                               transformed_value_of_basis_at_cub_points,
                               weighted_transformed_value_of_basis_at_cub_points,
                               COMP_CPP);

        // apply field signs (after the fact, as a post-processing step)
        fst::applyLeftFieldSigns<double>(mass_matrices, field_signs);
        fst::applyRightFieldSigns<double>(mass_matrices, field_signs);
        /*******************  STOP COMPUTATION ***********************/


        /******************* START COMPARISON ***********************/
        string basedir = "./testdata";
        for (int cell_id = 0; cell_id < numCells-1; cell_id++) {

          stringstream namestream;
          string filename;
          namestream <<  basedir << "/mass_HCURL_HEX_I1_FEM" << "_" << "0" << cell_id+1 << ".dat";
          namestream >> filename;

          ifstream massfile(&filename[0]);
          if (massfile.is_open()) {
            if (compareToAnalytic<double>(&mass_matrices(cell_id, 0, 0), massfile, 1e-10, iprint) > 0)
              errorFlag++;
            massfile.close();
          }
          else {
            errorFlag = -1;
            std::cout << "End Result: TEST FAILED\n";
            return errorFlag;
          }

          namestream.clear();
          namestream << basedir << "/stiff_HCURL_HEX_I1_FEM" << "_" << "0" << cell_id+1 << ".dat";
          namestream >> filename;

          ifstream stifffile(&filename[0]);
          if (stifffile.is_open())
            {
              if (compareToAnalytic<double>(&stiffness_matrices(cell_id, 0, 0), stifffile, 1e-10, iprint) > 0)
                errorFlag++;
              stifffile.close();
            }
          else {
            errorFlag = -1;
            std::cout << "End Result: TEST FAILED\n";
            return errorFlag;
          }

        }

        for (int cell_id = 3; cell_id < numCells; cell_id++) {

          stringstream namestream;
          string filename;
          namestream <<  basedir << "/mass_fp_HCURL_HEX_I1_FEM" << "_" << "0" << cell_id+1 << ".dat";
          namestream >> filename;

          ifstream massfile(&filename[0]);
          if (massfile.is_open()) {
            if (compareToAnalytic<double>(&mass_matrices(cell_id, 0, 0), massfile, 1e-4, iprint, INTREPID2_UTILS_SCALAR) > 0)
              errorFlag++;
            massfile.close();
          }
          else {
            errorFlag = -1;
            std::cout << "End Result: TEST FAILED\n";
            return errorFlag;
          }

          namestream.clear();
          namestream << basedir << "/stiff_fp_HCURL_HEX_I1_FEM" << "_" << "0" << cell_id+1 << ".dat";
          namestream >> filename;

          ifstream stifffile(&filename[0]);
          if (stifffile.is_open())
            {
              if (compareToAnalytic<double>(&stiffness_matrices(cell_id, 0, 0), stifffile, 1e-4, iprint, INTREPID2_UTILS_SCALAR) > 0)
                errorFlag++;
              stifffile.close();
            }
          else {
            errorFlag = -1;
            std::cout << "End Result: TEST FAILED\n";
            return errorFlag;
          }

        }

        /******************* 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);
      Kokkos::finalize();
      return errorFlag;
    }
  void IntegrationValues2<Scalar>::
  evaluateValues(const PHX::MDField<Scalar,Cell,NODE,Dim>& in_node_coordinates,
                 const PHX::MDField<Scalar,Cell,IP,Dim>& other_ip_coordinates)
  {
    if (int_rule->cv_type == "none") {

       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_static_view(), ip_coordinates.get_static_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);
    }

    else {

       getCubatureCV(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_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_static_view(), ip_coordinates.get_static_view());
        Array_CellIPDim old_weighted_normals = af.template buildStaticArray<Scalar,Cell,IP,Dim>(
            "old_weighted_normals", num_cells, num_ip, num_dim);
        Array_CellIP old_weighted_measure = af.template buildStaticArray<Scalar,Cell,IP>(
            "old_weighted_measure", num_cells, num_ip);
        if (int_rule->cv_type == "side") 
           Kokkos::deep_copy(old_weighted_normals.get_static_view(), weighted_normals.get_static_view());
        else 
           Kokkos::deep_copy(old_weighted_measure.get_static_view(), weighted_measure.get_static_view());
        for (size_type cell = 0; cell < num_cells; ++cell)
        {
          for (size_type ip = 0; ip < num_ip; ++ip)
          {
            if (ip != permutation[ip]) {
              if (int_rule->cv_type == "boundary" || int_rule->cv_type == "volume") 
                  weighted_measure(cell, ip) = old_weighted_measure(cell, permutation[ip]);
              for (size_type dim = 0; dim < num_dim; ++dim)
              {
                ip_coordinates(cell, ip, dim) = old_ip_coordinates(cell, permutation[ip], dim);
                if (int_rule->cv_type == "side") 
                  weighted_normals(cell, ip, dim) = old_weighted_normals(cell, permutation[ip], dim);
                
              }
            }
          }
        }
      }

      evaluateValuesCV(in_node_coordinates);
    }
  }