/**
* Compute the vector required to
* renumber the dofs of a cell by
- * block. Furthermore,
- * compute the vector storing the
- * start indices of each local block
- * vector.
+ * block. Furthermore, compute
+ * the vector storing either the
+ * start indices or the size of
+ * each local block vector.
+ *
+ * If the @p bool parameter is
+ * true, @p block_data is filled
+ * with the start indices of each
+ * local block. If it is false,
+ * then the block sizes are
+ * returned.
*
* @todo Which way does this
* vector map the numbers?
static void compute_block_renumbering (
const FiniteElement<dim,spacedim>& fe,
std::vector<unsigned int>& renumbering,
- std::vector<unsigned int>& start_indices);
+ std::vector<unsigned int>& block_data,
+ bool return_start_indices = true);
/**
* @name Generation of local matrices
void FETools::compute_block_renumbering (
const FiniteElement<dim,spacedim>& element,
std::vector<unsigned int>& renumbering,
- std::vector<unsigned int>& start_indices)
+ std::vector<unsigned int>& block_data,
+ bool return_start_indices)
{
Assert(renumbering.size() == element.dofs_per_cell,
ExcDimensionMismatch(renumbering.size(),
element.dofs_per_cell));
- Assert(start_indices.size() == element.n_blocks(),
- ExcDimensionMismatch(start_indices.size(),
+ Assert(block_data.size() == element.n_blocks(),
+ ExcDimensionMismatch(block_data.size(),
element.n_blocks()));
unsigned int k=0;
for (unsigned int b=0;b<element.n_base_elements();++b)
for (unsigned int m=0;m<element.element_multiplicity(b);++m)
{
- start_indices[i++] = k;
+ block_data[i++] = (return_start_indices)
+ ? k
+ : (element.base_element(b).n_dofs_per_cell());
k += element.base_element(b).n_dofs_per_cell();
}
Assert (i == element.n_blocks(), ExcInternalError());
+
+ std::vector<unsigned int> start_indices(block_data.size());
+ k = 0;
+ for (unsigned int i=0;i<block_data.size();++i)
+ if (return_start_indices)
+ start_indices[i] = block_data[i];
+ else
+ {
+ start_indices[i] = k;
+ k += block_data[i];
+ }
+
//TODO:[GK] This does not work for a single RT
for (unsigned int i=0;i<element.dofs_per_cell;++i)
{
template
void FETools::compute_block_renumbering (
const FiniteElement<deal_II_dimension>& element,
- std::vector<unsigned int>&, std::vector<unsigned int>&_indices);
+ std::vector<unsigned int>&, std::vector<unsigned int>&_indices, bool);
template
void FETools::get_interpolation_matrix<deal_II_dimension>
(const FiniteElement<deal_II_dimension> &,