Esempio n. 1
0
void Probes::add_positions(const Array<double>& x, const FunctionSpace& V)
{
  const std::size_t gdim = V.mesh()->geometry().dim();
  const std::size_t N = x.size() / gdim;
  Array<double> _x(gdim);
  const std::size_t old_N = total_number_probes;
  const std::size_t old_local_size = local_size();  
  total_number_probes += N;

  for (std::size_t i=0; i<N; i++)
  {
    for (std::size_t j=0; j<gdim; j++)
      _x[j] = x[i*gdim + j];
    try
    {
      Probe* probe = new Probe(_x, V);
      std::pair<std::size_t, Probe*> newprobe = std::make_pair(old_N+i, &(*probe));
      _allprobes.push_back(newprobe);
    } 
    catch (std::exception &e)
    { // do-nothing
    }
  }
  //cout << local_size() - old_local_size << " of " << N  << " probes found on processor " << MPI::process_number() << endl;
}
Esempio n. 2
0
//-----------------------------------------------------------------------------
DirichletBC::LocalData::LocalData(const FunctionSpace& V)
  : w(V.dofmap()->max_element_dofs(), 0.0),
    facet_dofs(V.dofmap()->num_facet_dofs(), 0),
    coordinates(boost::extents[V.dofmap()->max_element_dofs()][V.mesh()->geometry().dim()])
{
  // Do nothing
}
Esempio n. 3
0
    void Prolongation(const FunctionSpace& V, const Array<double>& Gx, const Array<double>& Gy, const Array<double>& Gz, double * dataX, double * dataY, double * dataZ) {
            const Mesh mesh = *V.mesh();
            const GenericDofMap& dofmap_u = *V.dofmap();
            const int n = mesh.num_cells();
            const int m = mesh.num_edges();

            std::vector<int> edge2dof;
            edge2dof.reserve(m);

            Cell* dolfin_cell;

            for (int i = 0; i < n; ++i) {
                dolfin_cell = new Cell(mesh, i);
                std::vector<int> cellDOF = dofmap_u.cell_dofs(i);
                const uint* edgeVALUES = dolfin_cell->entities(1);
                for (int j = 0; j < cellDOF.size(); ++j) {
                    edge2dof[edgeVALUES[j]]=cellDOF[j];
                }
            }

            Edge* dolfin_edge;
            int k = -1;
            for (int i = 0; i < m; ++i) {
                k = k + 1;
                dolfin_edge = new Edge(mesh, i);
                dataX[k] =Gx[dolfin_edge->index()]/2;
                dataY[k] =Gy[dolfin_edge->index()]/2;
                dataZ[k] =Gz[dolfin_edge->index()]/2;
                k = k + 1;
                dataX[k] =Gx[dolfin_edge->index()]/2;
                dataY[k] =Gy[dolfin_edge->index()]/2;
                dataZ[k] =Gz[dolfin_edge->index()]/2;
            }
        }
Esempio n. 4
0
Probes::Probes(const Array<double>& x, const FunctionSpace& V)
{
  const std::size_t Nd = V.mesh()->geometry().dim();
  const std::size_t N = x.size() / Nd;
  Array<double> _x(Nd);
  total_number_probes = N;
  _value_size = 1;
  _num_evals = 0;
  for (std::size_t i = 0; i < V.element()->value_rank(); i++)
    _value_size *= V.element()->value_dimension(i);

  for (std::size_t i=0; i<N; i++)
  {
    for (std::size_t j=0; j<Nd; j++)
      _x[j] = x[i*Nd + j];
    try
    {
      Probe* probe = new Probe(_x, V);
      std::pair<std::size_t, Probe*> newprobe = std::make_pair(i, &(*probe));
      _allprobes.push_back(newprobe);
    } 
    catch (std::exception &e)
    { // do-nothing
    }
  }
  //cout << _allprobes.size() << " of " << N  << " probes found on processor " << MPI::process_number() << endl;
}
Esempio n. 5
0
//-----------------------------------------------------------------------------
void Extrapolation::build_unique_dofs(
  std::set<std::size_t>& unique_dofs,
  std::map<std::size_t, std::map<std::size_t, std::size_t>>& cell2dof2row,
  const Cell& cell0,
  const FunctionSpace& V)
{
  // Counter for matrix row index
  std::size_t row = 0;
  dolfin_assert(V.mesh());

  // Get unique set of surrounding cells (including cell0)
  std::set<std::size_t> cell_set;
  for (VertexIterator vtx(cell0); !vtx.end(); ++vtx)
  {
    for (CellIterator cell1(*vtx); !cell1.end(); ++cell1)
      cell_set.insert(cell1->index());
  }

  // Compute unique dofs on patch
  for (auto cell_it : cell_set)
  {
    Cell cell1(cell0.mesh(), cell_it);
    cell2dof2row[cell_it] = compute_unique_dofs(cell1, V, row,
                                                unique_dofs);
  }

}
Esempio n. 6
0
//-----------------------------------------------------------------------------
SubSpace::SubSpace(const FunctionSpace& V,
                   const std::vector<std::size_t>& component)
    : FunctionSpace(V.mesh(), V.element(), V.dofmap())
{
    // Extract subspace and assign
    std::shared_ptr<FunctionSpace> _function_space(V.extract_sub_space(component));
    *static_cast<FunctionSpace*>(this) = *_function_space;
}
Esempio n. 7
0
Probe::Probe(const Array<double>& x, const FunctionSpace& V) :
  _element(V.element()), _num_evals(0)
{

  const Mesh& mesh = *V.mesh();
  std::size_t gdim = mesh.geometry().dim();
    
  // Find the cell that contains probe
  const Point point(gdim, x.data());
  unsigned int cell_id = mesh.bounding_box_tree()->compute_first_entity_collision(point);

  // If the cell is on this process, then create an instance 
  // of the Probe class. Otherwise raise a dolfin_error.
  if (cell_id != std::numeric_limits<unsigned int>::max())
  {
    // Store position of probe
    for (std::size_t i = 0; i < 3; i++) 
      _x[i] = (i < gdim ? x[i] : 0.0);
    
    // Compute in tensor (one for scalar function, . . .)
    _value_size_loc = 1;
    for (std::size_t i = 0; i < _element->value_rank(); i++)
       _value_size_loc *= _element->value_dimension(i);

    _probes.resize(_value_size_loc);

    // Create cell that contains point
    dolfin_cell.reset(new Cell(mesh, cell_id));
    dolfin_cell->get_cell_data(ufc_cell);
    
    coefficients.resize(_element->space_dimension());
    
    // Cell vertices
    dolfin_cell->get_vertex_coordinates(vertex_coordinates);
        
    // Create work vector for basis
    basis_matrix.resize(_value_size_loc);
    for (std::size_t i = 0; i < _value_size_loc; ++i)
      basis_matrix[i].resize(_element->space_dimension());
        
    std::vector<double> basis(_value_size_loc);
    const int cell_orientation = 0;
    for (std::size_t i = 0; i < _element->space_dimension(); ++i)
    {
      _element->evaluate_basis(i, &basis[0], &x[0], 
                               vertex_coordinates.data(), 
                               cell_orientation);
      for (std::size_t j = 0; j < _value_size_loc; ++j)
        basis_matrix[j][i] = basis[j];
    }
  }  
  else
  {  
    dolfin_error("Probe.cpp","set probe","Probe is not found on processor");
  }
}
Esempio n. 8
0
//-----------------------------------------------------------------------------
SubSpace::SubSpace(const FunctionSpace& V,
                   std::size_t component,
                   std::size_t sub_component)
    : FunctionSpace(V.mesh(), V.element(), V.dofmap())
{
    // Create array
    std::vector<std::size_t> c = {{component, sub_component}};

    // Extract subspace and assign
    std::shared_ptr<FunctionSpace> _function_space(V.extract_sub_space(c));
    *static_cast<FunctionSpace*>(this) = *_function_space;
}
Esempio n. 9
0
//-----------------------------------------------------------------------------
void Extrapolation::build_unique_dofs(std::set<std::size_t>& unique_dofs,
                                      std::map<std::size_t, std::map<std::size_t, std::size_t> >& cell2dof2row,
                                      const Cell& cell0,
                                      const ufc::cell& c0,
                                      const FunctionSpace& V)
{
  // Counter for matrix row index
  std::size_t row = 0;
  dolfin_assert(V.mesh());
  UFCCell c1(*V.mesh());

  // Compute unique dofs on center cell
  cell2dof2row[cell0.index()] = compute_unique_dofs(cell0, c0, V, row, unique_dofs);

  // Compute unique dofs on neighbouring cells
  for (CellIterator cell1(cell0); !cell1.end(); ++cell1)
  {
    c1.update(*cell1);
    cell2dof2row[cell1->index()] = compute_unique_dofs(*cell1, c1, V, row, unique_dofs);
  }
}
Esempio n. 10
0
    void Gradient(const FunctionSpace& V, const Array<int>& Mapping, double * column, double * row, double * data) {

        const Mesh mesh = *V.mesh();
        const GenericDofMap& dofmap_u = *V.dofmap();
        const int n = mesh.num_cells();
        const int m = mesh.num_edges();

        std::vector<int> edge2dof;
        edge2dof.reserve(m);

        Cell* dolfin_cell;

        for (int i = 0; i < n; ++i) {
            dolfin_cell = new Cell(mesh, i);
            std::vector<int> cellDOF = dofmap_u.cell_dofs(i);
            const uint* edgeVALUES = dolfin_cell->entities(1);
            for (int j = 0; j < cellDOF.size(); ++j) {
                edge2dof[edgeVALUES[j]]=cellDOF[j];
            }

        }


        Edge* dolfin_edge;
        int k = -1;
        for (int i = 0; i < m; ++i) {
            dolfin_edge = new Edge(mesh, i);
            const uint* edgeVERTICES = dolfin_edge->entities(0);

            k = k+1;
            row[k]=edge2dof[dolfin_edge->index()];
            column[k]=Mapping[edgeVERTICES[0]];
            data[k]=-1;
            k = k+1;
            row[k]=edge2dof[dolfin_edge->index()];
            column[k]=Mapping[edgeVERTICES[1]];
            data[k]=1;


        }

    }
Esempio n. 11
0
StatisticsProbes::StatisticsProbes(const Array<double>& x, const FunctionSpace& V, bool segregated)
{
  const std::size_t Nd = V.mesh()->geometry().dim();
  const std::size_t N = x.size() / Nd;
  Array<double> _x(Nd);
  total_number_probes = N;
  _num_evals = 0;
  _value_size = 1;
  for (std::size_t i = 0; i < V.element()->value_rank(); i++)
    _value_size *= V.element()->value_dimension(i);
  
  if (segregated)
  {
    assert(V.element()->value_rank() == 0);
    _value_size *= V.element()->geometric_dimension();
  }
    
  // Symmetric statistics. Velocity: u, v, w, uu, vv, ww, uv, uw, vw
  _value_size = _value_size*(_value_size+3)/2.;

  for (std::size_t i=0; i<N; i++)
  {
    for (std::size_t j=0; j<Nd; j++)
      _x[j] = x[i*Nd + j];
    try
    {
      StatisticsProbe* probe = new StatisticsProbe(_x, V, segregated);
      std::pair<std::size_t, StatisticsProbe*> newprobe = std::make_pair(i, &(*probe));
      _allprobes.push_back(newprobe);
    } 
    catch (std::exception &e)
    { // do-nothing
    }
  }
  //cout << local_size() << " of " << N  << " probes found on processor " << MPI::process_number() << endl;
}
Esempio n. 12
0
 void extract_dof_component_map(std::map<std::size_t, std::size_t>& dof_component_map, 
                                const FunctionSpace& VV, int* component)
 { // Extract sub dofmaps recursively and store dof to component map
   boost::unordered_map<std::size_t, std::size_t> collapsed_map;
   boost::unordered_map<std::size_t, std::size_t>::iterator map_it;
   std::vector<std::size_t> comp(1);
   
   if (VV.element()->num_sub_elements() == 0)
   {
     boost::shared_ptr<GenericDofMap> dummy = VV.dofmap()->collapse(collapsed_map, *VV.mesh());
     (*component)++;
     for (map_it =collapsed_map.begin(); map_it!=collapsed_map.end(); ++map_it)
       dof_component_map[map_it->second] = (*component);  
   }
   else
   {
     for (std::size_t i=0; i<VV.element()->num_sub_elements(); i++)
     {
       comp[0] = i;
       boost::shared_ptr<FunctionSpace> Vs = VV.extract_sub_space(comp);
       extract_dof_component_map(dof_component_map, *Vs, component);
     }
   }
 }
Esempio n. 13
0
    void ProlongationGrad(const FunctionSpace& V, double * dataX, double * dataY, double * dataZ, double * data, double * row, double * column) {
        const Mesh mesh = *V.mesh();
        const GenericDofMap& dofmap_u = *V.dofmap();
        const int n = mesh.num_cells();
        const int m = mesh.num_edges();
        const int dim = mesh.geometry().dim();
        std::vector<double> coord = mesh.coordinates();
        /*
        even are the x coords...
        odd are the y coords...
        */
        std::vector<double> X;
        std::vector<double> Y;
        std::vector<double> Z;

        for (int i = 0; i < coord.size()/dim; ++i)
        {
            if (dim == 2) {
                // std::cout << 2*i << std::endl;
                X.push_back(coord[2*i]);
                Y.push_back(coord[2*i+1]);
            }
            else {
                // std::cout << 3*i << std::endl;
                X.push_back(coord[3*i]);
                Y.push_back(coord[3*i+1]);
                Z.push_back(coord[3*i+2]);
            }
        }


        std::vector<double> tangent;
        std::vector<double> Ntangent;
        tangent.reserve(dim);
        Ntangent.reserve(dim);
        Edge* dolfin_edge;
        int k = -1;
        for (int i = 0; i < m; ++i) {
            k = k + 1;
            dolfin_edge = new Edge(mesh, i);
            const uint* edgeVERTICES = dolfin_edge->entities(0);

            tangent[0] = (X[edgeVERTICES[1]]-X[edgeVERTICES[0]]);
            tangent[1] = (Y[edgeVERTICES[1]]-Y[edgeVERTICES[0]]);
            if (dim == 3) {
                tangent[2] = Z[edgeVERTICES[1]]-Z[edgeVERTICES[0]];
            }

            for (int j = 0; j < dim; ++j) {
                if (dim == 3) {
                    if (tangent[0] == 0 && tangent[1] == 0 && tangent[2] == 0){
                        Ntangent[j] =0;
                    }
                    else {
                        Ntangent[j] = tangent[j]/(sqrt(pow(tangent[0],2.0)+pow(tangent[1],2.0)+pow(tangent[2],2.0)));
                    }
                }
                else {

                    if (tangent[0] == 0 && tangent[1] == 0){
                        Ntangent[j] =0;
                    }
                    else {
                        Ntangent[j] = tangent[j]/(sqrt(pow(tangent[0],2.0)+pow(tangent[1],2.0)));
                    }
                }
            }

            // double len = sqrt(pow(tangent[0],2.0)+pow(tangent[1],2.0));


            // std::cout << len << " , " << dolfin_edge->length() << std::endl;
            dataX[k]=0.5*dolfin_edge->length()*Ntangent[0];
            dataY[k]=0.5*dolfin_edge->length()*Ntangent[1];

            if (dim == 3) {
                dataZ[k]=0.5*dolfin_edge->length()*Ntangent[2];
            }
            data[k] = -1;
            row[k] = i;
            column[k] = edgeVERTICES[0];
            // std::cout << dataX[k] << std::endl;
            k = k + 1;
            dataX[k]=0.5*dolfin_edge->length()*Ntangent[0];
            dataY[k]=0.5*dolfin_edge->length()*Ntangent[1];

            if (dim == 3) {
                dataZ[k]=0.5*dolfin_edge->length()*Ntangent[2];
            }
            data[k] = 1;
            row[k] = i;
            column[k] = edgeVERTICES[1];
            // std::cout << dataX[k] << std::endl;

        }


    }
Esempio n. 14
0
    void ProlongationP(const FunctionSpace& V, const Array<int>& Mapping, const Array<double>& X, const Array<double>& Y, const Array<double>& Z, double * dataX, double * dataY, double * dataZ) {
        const Mesh mesh = *V.mesh();
        const GenericDofMap& dofmap_u = *V.dofmap();
        const int n = mesh.num_cells();
        const int m = mesh.num_edges();
        const int dim = mesh.geometry().dim();
        // std::vector<double> coord = mesh.coordinates();
        // /*
        // even are the x coords...
        // odd are the y coords...
        // */
        // std::vector<double> X;
        // std::vector<double> Y;
        // std::vector<double> Z;

        // for (int i = 0; i < coord.size()/dim; ++i)
        // {
        //     if (dim == 2) {
        //         std::cout << 2*i << std::endl;
        //         X.push_back(coord[2*i]);
        //         Y.push_back(coord[2*i+1]);
        //     }
        //     else {
        //         std::cout << 3*i << std::endl;
        //         X.push_back(coord[3*i]);
        //         Y.push_back(coord[3*i+1]);
        //         Z.push_back(coord[3*i+2]);
        //     }
        // }


        std::vector<double> tangent;
        std::vector<double> Ntangent;
        tangent.reserve(dim);
        Ntangent.reserve(dim);
        Edge* dolfin_edge;
        int k = -1;
        for (int i = 0; i < m; ++i) {
            k = k + 1;
            dolfin_edge = new Edge(mesh, i);
            const uint* edgeVERTICES = dolfin_edge->entities(0);

            tangent[0] = (X[edgeVERTICES[1]]-X[edgeVERTICES[0]]);
            tangent[1] = (Y[edgeVERTICES[1]]-Y[edgeVERTICES[0]]);
            if (dim == 3) {
                tangent[2] = Z[edgeVERTICES[1]]-Z[edgeVERTICES[0]];
            }

            for (int j = 0; j < dim; ++j) {
                if (dim == 3) {
                    if (tangent[0] == 0 && tangent[1] == 0 && tangent[2] == 0){
                        Ntangent[j] =0;
                    }
                    else {
                        Ntangent[j] = tangent[j]/(sqrt(pow(tangent[0],2.0)+pow(tangent[1],2.0)+pow(tangent[2],2.0)));
                    }
                }
                else {

                    if (tangent[0] == 0 && tangent[1] == 0){
                        Ntangent[j] =0;
                    }
                    else {
                        Ntangent[j] = tangent[j]/(sqrt(pow(tangent[0],2.0)+pow(tangent[1],2.0)));
                    }
                }
            }

            double len = sqrt(pow(tangent[0],2.0)+pow(tangent[1],2.0));


            // std::cout << len << " , " << dolfin_edge->length() << std::endl;
            dataX[k]=0.5*dolfin_edge->length()*Ntangent[0];
            dataY[k]=0.5*dolfin_edge->length()*Ntangent[1];

            if (dim == 3) {
                dataZ[k]=0.5*dolfin_edge->length()*Ntangent[2];
            }

            k = k + 1;
            dataX[k]=0.5*dolfin_edge->length()*Ntangent[0];
            dataY[k]=0.5*dolfin_edge->length()*Ntangent[1];

            if (dim == 3) {
                dataZ[k]=0.5*dolfin_edge->length()*Ntangent[2];
            }

        }


    }
Esempio n. 15
0
// Base case for all divergence computations. 
// Compute divergence of vector field u.
void cr_divergence_matrix(GenericMatrix& M, GenericMatrix& A,
                       const FunctionSpace& DGscalar, 
                       const FunctionSpace& CRvector)
{
  std::shared_ptr<const GenericDofMap>
    CR1_dofmap = CRvector.dofmap(),
    DG0_dofmap = DGscalar.dofmap();

  // Figure out about the local dofs of DG0 
  std::pair<std::size_t, std::size_t>
  first_last_dof = DG0_dofmap->ownership_range();
  std::size_t first_dof = first_last_dof.first;
  std::size_t last_dof = first_last_dof.second;
  std::size_t n_local_dofs = last_dof - first_dof;

  // Get topological dimension so that we know what Facet is
  const Mesh mesh = *DGscalar.mesh();
  std::size_t tdim = mesh.topology().dim(); 
  std::size_t gdim = mesh.geometry().dim();
  
  std::vector<std::size_t> columns;
  std::vector<double> values;
  
  // Fill the values
  for(CellIterator cell(mesh); !cell.end(); ++cell)
  {
    auto dg_dofs = DG0_dofmap->cell_dofs(cell->index());
    // There is only one DG0 dof per cell
    dolfin::la_index cell_dof = dg_dofs[0];
    
    Point cell_mp = cell->midpoint();
    double cell_volume = cell->volume();
    std::size_t local_facet_index = 0;      
    
    auto cr_dofs = CR1_dofmap->cell_dofs(cell->index());
    
    for(FacetIterator facet(*cell); !facet.end(); ++facet)
    {
      double facet_measure=0;
      if(tdim == 2)
        facet_measure = Edge(mesh, facet->index()).length();
      else if(tdim == 3)
        facet_measure = Face(mesh, facet->index()).area();
      // Tdim 1 will not happen because CR is not defined there
      
      Point facet_normal = facet->normal();

      // Flip the normal if it is not outer already
      Point facet_mp = facet->midpoint();
      double sign = (facet_normal.dot(facet_mp - cell_mp) > 0.0) ? 1.0 : -1.0;
      facet_normal *= (sign*facet_measure/cell_volume);
      
      // Dofs of CR on the facet, local order
      std::vector<std::size_t> facet_dofs;
      CR1_dofmap->tabulate_facet_dofs(facet_dofs, local_facet_index);
      
      for (std::size_t j = 0 ; j < facet_dofs.size(); j++)
      {   
        columns.push_back(cr_dofs[facet_dofs[j]]);
        values.push_back(facet_normal[j]);
      }        
      local_facet_index += 1;
    }
    M.setrow(cell_dof, columns, values);
    columns.clear();
    values.clear();
  }
  M.apply("insert");
  //std::shared_ptr<GenericMatrix> Cp = MatMatMult(M, A);
  //return Cp;
}
Esempio n. 16
0
//-----------------------------------------------------------------------------
void Extrapolation::compute_coefficients(std::vector<std::vector<double> >& coefficients,
                                         const Function& v,
                                         const FunctionSpace& V,
                                         const FunctionSpace& W,
                                         const Cell& cell0,
                                         const ufc::cell& c0,
                                         const std::vector<dolfin::la_index>& dofs,
                                         std::size_t& offset)
{
  // Call recursively for mixed elements
  dolfin_assert(V.element());
  const std::size_t num_sub_spaces = V.element()->num_sub_elements();
  if (num_sub_spaces > 0)
  {
    for (std::size_t k = 0; k < num_sub_spaces; k++)
    {
      compute_coefficients(coefficients, v[k], *V[k], *W[k], cell0, c0,
                           dofs, offset);
    }
    return;
  }

  // Build data structures for keeping track of unique dofs
  std::map<std::size_t, std::map<std::size_t, std::size_t> > cell2dof2row;
  std::set<std::size_t> unique_dofs;
  build_unique_dofs(unique_dofs, cell2dof2row, cell0, c0, V);

  // Compute size of linear system
  dolfin_assert(W.element());
  const std::size_t N = W.element()->space_dimension();
  const std::size_t M = unique_dofs.size();

  // Check size of system
  if (M < N)
  {
    dolfin_error("Extrapolation.cpp",
                 "compute extrapolation",
                 "Not enough degrees of freedom on local patch to build extrapolation");
  }

  // Create matrix and vector for linear system
  arma::mat A(M, N);
  arma::vec b(M);

  // Add equations on cell and neighboring cells
  add_cell_equations(A, b, cell0, cell0, c0, c0, V, W, v, cell2dof2row[cell0.index()]);
  dolfin_assert(V.mesh());
  UFCCell c1(*V.mesh());
  for (CellIterator cell1(cell0); !cell1.end(); ++cell1)
  {
    if (cell2dof2row[cell1->index()].empty())
      continue;

    c1.update(*cell1);
    add_cell_equations(A, b, cell0, *cell1, c0, c1, V, W, v, cell2dof2row[cell1->index()]);
  }

  // Solve least squares system
  arma::Col<double> x = arma::solve(A, b);

  // Insert resulting coefficients into global coefficient vector
  dolfin_assert(W.dofmap());
  for (std::size_t i = 0; i < W.dofmap()->cell_dimension(cell0.index()); ++i)
    coefficients[dofs[i + offset]].push_back(x[i]);

  // Increase offset
  offset += W.dofmap()->cell_dimension(cell0.index());
}
Esempio n. 17
0
//-----------------------------------------------------------------------------
void Extrapolation::compute_coefficients(
  std::vector<std::vector<double>>& coefficients,
  const Function& v,
  const FunctionSpace& V,
  const FunctionSpace& W,
  const Cell& cell0,
  const std::vector<double>& coordinate_dofs0,
  const ufc::cell& c0,
  const ArrayView<const dolfin::la_index>& dofs,
  std::size_t& offset)
{
  // Call recursively for mixed elements
  dolfin_assert(V.element());
  const std::size_t num_sub_spaces = V.element()->num_sub_elements();
  if (num_sub_spaces > 0)
  {
    for (std::size_t k = 0; k < num_sub_spaces; k++)
    {
      compute_coefficients(coefficients, v[k], *V[k], *W[k], cell0,
                           coordinate_dofs0, c0, dofs, offset);
    }
    return;
  }

  // Build data structures for keeping track of unique dofs
  std::map<std::size_t, std::map<std::size_t, std::size_t>> cell2dof2row;
  std::set<std::size_t> unique_dofs;
  build_unique_dofs(unique_dofs, cell2dof2row, cell0, V);

  // Compute size of linear system
  dolfin_assert(W.element());
  const std::size_t N = W.element()->space_dimension();
  const std::size_t M = unique_dofs.size();

  // Check size of system
  if (M < N)
  {
    dolfin_error("Extrapolation.cpp",
                 "compute extrapolation",
                 "Not enough degrees of freedom on local patch to build extrapolation");
  }

  // Create matrix and vector for linear system
  Eigen::MatrixXd A(M, N);
  Eigen::VectorXd b(M);

  // Add equations on cell and neighboring cells
  dolfin_assert(V.mesh());
  ufc::cell c1;
  std::vector<double> coordinate_dofs1;

  // Get unique set of surrounding cells (including cell0)
  std::set<std::size_t> cell_set;
  for (VertexIterator vtx(cell0); !vtx.end(); ++vtx)
  {
    for (CellIterator cell1(*vtx); !cell1.end(); ++cell1)
      cell_set.insert(cell1->index());
  }

  for (auto cell_it : cell_set)
  {
    if (cell2dof2row[cell_it].empty())
      continue;

    Cell cell1(cell0.mesh(), cell_it);

    cell1.get_coordinate_dofs(coordinate_dofs1);
    cell1.get_cell_data(c1);
    add_cell_equations(A, b, cell0, cell1,
                       coordinate_dofs0, coordinate_dofs1,
                       c0, c1, V, W, v,
                       cell2dof2row[cell_it]);
  }

  // Solve least squares system
  const Eigen::VectorXd x
    = A.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);

  // Insert resulting coefficients into global coefficient vector
  dolfin_assert(W.dofmap());
  for (std::size_t i = 0; i < W.dofmap()->num_element_dofs(cell0.index()); ++i)
    coefficients[dofs[i + offset]].push_back(x[i]);

  // Increase offset
  offset += W.dofmap()->num_element_dofs(cell0.index());
}
Esempio n. 18
0
//-----------------------------------------------------------------------------
void XMLFunctionData::build_dof_map(std::vector<std::vector<uint> >& dof_map,
                                    const FunctionSpace& V)
{
  // Get mesh and dofmap
  dolfin_assert(V.mesh());
  dolfin_assert(V.dofmap());
  const Mesh& mesh = *V.mesh();
  const GenericDofMap& dofmap = *V.dofmap();

  const uint num_cells = MPI::sum(mesh.num_cells());

  std::vector<std::vector<std::vector<uint > > > gathered_dofmap;
  std::vector<std::vector<uint > > local_dofmap(mesh.num_cells());

  if (MPI::num_processes() > 1)
  {
    // Get local-to-global cell numbering
    dolfin_assert(mesh.parallel_data().have_global_entity_indices(mesh.topology().dim()));
    const MeshFunction<uint>& global_cell_indices
      = mesh.parallel_data().global_entity_indices(mesh.topology().dim());

    // Build dof map data with global cell indices
    for (CellIterator cell(mesh); !cell.end(); ++cell)
    {
      const uint local_cell_index = cell->index();
      const uint global_cell_index = global_cell_indices[*cell];
      local_dofmap[local_cell_index] = dofmap.cell_dofs(local_cell_index);
      local_dofmap[local_cell_index].push_back(global_cell_index);
    }
  }
  else
  {
    // Build dof map data
    for (CellIterator cell(mesh); !cell.end(); ++cell)
    {
      const uint local_cell_index = cell->index();
      local_dofmap[local_cell_index] = dofmap.cell_dofs(local_cell_index);
      local_dofmap[local_cell_index].push_back(local_cell_index);
    }
  }

  // Gather dof map data on root process
  MPI::gather(local_dofmap, gathered_dofmap);


  // Build global dofmap on root process
  if (MPI::process_number() == 0)
  {
    dof_map.resize(num_cells);

    // Loop of dof map from each process
    std::vector<std::vector<std::vector<uint > > > ::const_iterator proc_dofmap;
    for (proc_dofmap = gathered_dofmap.begin(); proc_dofmap != gathered_dofmap.end(); ++proc_dofmap)
    {
      std::vector<std::vector<uint> >::const_iterator cell_dofmap;
      for (cell_dofmap = proc_dofmap->begin(); cell_dofmap != proc_dofmap->end(); ++cell_dofmap)
      {
        const std::vector<uint>& cell_dofs = *cell_dofmap;
        const uint global_cell_index = cell_dofs.back();
        dolfin_assert(global_cell_index < dof_map.size());
        dof_map[global_cell_index] = *cell_dofmap;
        dof_map[global_cell_index].pop_back();
      }
    }
  }
}
Esempio n. 19
0
//-----------------------------------------------------------------------------
void XMLFunctionData::build_dof_map(std::vector<std::vector<la_index>>& dof_map,
                                    const FunctionSpace& V)
{
  // Get mesh and dofmap
  dolfin_assert(V.mesh());
  dolfin_assert(V.dofmap());
  const Mesh& mesh = *V.mesh();
  const GenericDofMap& dofmap = *V.dofmap();

  // Get local-to-global map
  std::vector<std::size_t> local_to_global_dof;
  dofmap.tabulate_local_to_global_dofs(local_to_global_dof);

  // Get global number of cells
  const std::size_t num_cells = MPI::sum(mesh.mpi_comm(), mesh.num_cells());

  std::vector<dolfin::la_index> local_dofmap;
  if (MPI::size(mesh.mpi_comm()) > 1)
  {
    // Check that local-to-global cell numbering is available
    const std::size_t D = mesh.topology().dim();
    dolfin_assert(mesh.topology().have_global_indices(D));

    // Build dof map data with global cell indices
    std::vector<la_index> cell_dofs_global;
    for (CellIterator cell(mesh); !cell.end(); ++cell)
    {
      const std::size_t local_cell_index = cell->index();
      const std::size_t global_cell_index = cell->global_index();
      auto cell_dofs = dofmap.cell_dofs(local_cell_index);
      local_dofmap.push_back(global_cell_index);
      local_dofmap.push_back(cell_dofs.size());

      cell_dofs_global.resize(cell_dofs.size());
      for(Eigen::Index i = 0; i < cell_dofs.size(); ++i)
        cell_dofs_global[i] = local_to_global_dof[cell_dofs[i]];

      // Insert global dof indices
      local_dofmap.insert(local_dofmap.end(), cell_dofs_global.data(),
                          cell_dofs_global.data() + cell_dofs_global.size());
    }
  }
  else
  {
    // Build dof map data
    for (CellIterator cell(mesh); !cell.end(); ++cell)
    {
      const std::size_t local_cell_index = cell->index();
      local_dofmap.push_back(local_cell_index);
      local_dofmap.push_back(dofmap.cell_dofs(local_cell_index).size());

      auto dmap = dofmap.cell_dofs(local_cell_index);
      local_dofmap.insert(local_dofmap.end(),
                          dmap.data(), dmap.data() + dmap.size());
    }
  }

  // Gather dof map data on root process
  std::vector<dolfin::la_index> gathered_dofmap;
  MPI::gather(mesh.mpi_comm(), local_dofmap, gathered_dofmap);

  // Build global dofmap on root process
  if (MPI::rank(mesh.mpi_comm()) == 0)
  {
    dof_map.resize(num_cells);
    for (std::size_t i = 0; i < gathered_dofmap.size(); )
    {
      const std::size_t global_cell_index = gathered_dofmap[i++];
      const std::size_t num_dofs     = gathered_dofmap[i++];
      for (std::size_t j = 0; j < num_dofs; ++j)
        dof_map[global_cell_index].push_back(gathered_dofmap[i++]);
    }
  }
}
Esempio n. 20
0
//-----------------------------------------------------------------------------
std::vector<dolfin::la_index>
dolfin::vertex_to_dof_map(const FunctionSpace& space)
{
  // Get the mesh
  dolfin_assert(space.mesh());
  dolfin_assert(space.dofmap());
  const Mesh& mesh = *space.mesh();
  const GenericDofMap& dofmap = *space.dofmap();

  if (dofmap.is_view())
  {
    dolfin_error("fem_utils.cpp",
                 "tabulate vertex to dof map",
                 "Cannot tabulate vertex_to_dof_map for a subspace");
  }

  // Initialize vertex to cell connections
  const std::size_t top_dim = mesh.topology().dim();
  mesh.init(0, top_dim);

  // Num dofs per vertex
  const std::size_t dofs_per_vertex = dofmap.num_entity_dofs(0);
  const std::size_t vert_per_cell = mesh.topology()(top_dim, 0).size(0);
  if (vert_per_cell*dofs_per_vertex != dofmap.max_element_dofs())
  {
    dolfin_error("DofMap.cpp",
                 "tabulate dof to vertex map",
                 "Can only tabulate dofs on vertices");
  }

  // Allocate data for tabulating local to local map
  std::vector<std::size_t> local_to_local_map(dofs_per_vertex);

  // Create return data structure
  std::vector<dolfin::la_index>
    return_map(dofs_per_vertex*mesh.num_entities(0));

  // Iterate over vertices
  std::size_t local_vertex_ind = 0;
  for (VertexIterator vertex(mesh, "all"); !vertex.end(); ++vertex)
  {
    // Get the first cell connected to the vertex
    const Cell cell(mesh, vertex->entities(top_dim)[0]);

    // Find local vertex number
#ifdef DEBUG
    bool vertex_found = false;
#endif
    for (std::size_t i = 0; i < cell.num_entities(0); i++)
    {
      if (cell.entities(0)[i] == vertex->index())
      {
        local_vertex_ind = i;
#ifdef DEBUG
        vertex_found = true;
#endif
        break;
      }
    }
    dolfin_assert(vertex_found);

    // Get all cell dofs
    const ArrayView<const dolfin::la_index> cell_dofs
      = dofmap.cell_dofs(cell.index());

    // Tabulate local to local map of dofs on local vertex
    dofmap.tabulate_entity_dofs(local_to_local_map, 0,
				local_vertex_ind);

    // Fill local dofs for the vertex
    for (std::size_t local_dof = 0; local_dof < dofs_per_vertex; local_dof++)
    {
      const dolfin::la_index global_dof
        = cell_dofs[local_to_local_map[local_dof]];
      return_map[dofs_per_vertex*vertex->index() + local_dof] = global_dof;
    }
  }

  // Return the map
  return return_map;
}
Esempio n. 21
0
//-----------------------------------------------------------------------------
void XMLFunctionData::build_global_to_cell_dof(
  std::vector<std::vector<std::pair<uint, uint> > >& global_dof_to_cell_dof,
  const FunctionSpace& V)
{
  // Get mesh and dofmap
  dolfin_assert(V.mesh());
  dolfin_assert(V.dofmap());
  const Mesh& mesh = *V.mesh();
  const GenericDofMap& dofmap = *V.dofmap();

  std::vector<std::vector<std::vector<uint > > > gathered_dofmap;
  std::vector<std::vector<uint > > local_dofmap(mesh.num_cells());

  if (MPI::num_processes() > 1)
  {
    // Get local-to-global cell numbering
    dolfin_assert(mesh.parallel_data().have_global_entity_indices(mesh.topology().dim()));
    const MeshFunction<uint>& global_cell_indices
      = mesh.parallel_data().global_entity_indices(mesh.topology().dim());

    // Build dof map data with global cell indices
    for (CellIterator cell(mesh); !cell.end(); ++cell)
    {
      const uint local_cell_index = cell->index();
      const uint global_cell_index = global_cell_indices[*cell];
      local_dofmap[local_cell_index] = dofmap.cell_dofs(local_cell_index);
      local_dofmap[local_cell_index].push_back(global_cell_index);
    }
  }
  else
  {
    // Build dof map data
    for (CellIterator cell(mesh); !cell.end(); ++cell)
    {
      const uint local_cell_index = cell->index();
      local_dofmap[local_cell_index] = dofmap.cell_dofs(local_cell_index);
      local_dofmap[local_cell_index].push_back(local_cell_index);
    }
  }

  // Gather dof map data on root process
  MPI::gather(local_dofmap, gathered_dofmap);

  // Build global dof - (global cell, local dof) map on root process
  if (MPI::process_number() == 0)
  {
    global_dof_to_cell_dof.resize(dofmap.global_dimension());

    std::vector<std::vector<std::vector<uint > > > ::const_iterator proc_dofmap;
    for (proc_dofmap = gathered_dofmap.begin(); proc_dofmap != gathered_dofmap.end(); ++proc_dofmap)
    {
      std::vector<std::vector<uint> >::const_iterator cell_dofmap;
      for (cell_dofmap = proc_dofmap->begin(); cell_dofmap != proc_dofmap->end(); ++cell_dofmap)
      {
        const std::vector<uint>& cell_dofs = *cell_dofmap;
        const uint global_cell_index = cell_dofs.back();
        for (uint i = 0; i < cell_dofs.size() - 1; ++i)
          global_dof_to_cell_dof[cell_dofs[i]].push_back(std::make_pair(global_cell_index, i));
      }
    }
  }
}