* The function above adjusted with the variable cell_normal_vectors for the
* case of codimension 1
*/
- virtual void
+ virtual
+ CellSimilarity::Similarity
fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const Quadrature<dim> &quadrature,
InternalDataBase &internal,
std::vector<DerivativeForm<2,dim,spacedim> > &jacobian_grads,
std::vector<DerivativeForm<1,spacedim,dim> > &inverse_jacobians,
std::vector<Point<spacedim> > &cell_normal_vectors,
- CellSimilarity::Similarity &cell_similarity
+ const CellSimilarity::Similarity cell_similarity
) const=0;
get_subface_data (const UpdateFlags flags,
const Quadrature<dim-1>& quadrature) const;
- virtual void
+ virtual
+ CellSimilarity::Similarity
fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const Quadrature<dim> &quadrature,
typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
std::vector<DerivativeForm<2,dim,spacedim> > &jacobian_grads,
std::vector<DerivativeForm<1,spacedim,dim> > &inverse_jacobians,
std::vector<Point<spacedim> > &,
- CellSimilarity::Similarity &cell_similarity) const ;
+ const CellSimilarity::Similarity cell_similarity) const ;
virtual void
/**
* Implementation of the interface in Mapping.
*/
- virtual void
+ virtual
+ CellSimilarity::Similarity
fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const Quadrature<dim> &quadrature,
typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
std::vector<DerivativeForm<2,dim,spacedim> > &jacobian_grads,
std::vector<DerivativeForm<1,spacedim,dim> > &inverse_jacobians,
std::vector<Point<spacedim> > &cell_normal_vectors,
- CellSimilarity::Similarity &cell_similarity) const ;
+ const CellSimilarity::Similarity cell_similarity) const ;
/**
* Implementation of the interface in Mapping.
/**
* Implementation of the interface in Mapping.
*/
- virtual void
+ virtual
+ CellSimilarity::Similarity
fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const Quadrature<dim> &quadrature,
typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
std::vector<DerivativeForm<2,dim,spacedim> > &jacobian_grads,
std::vector<DerivativeForm<1,spacedim,dim> > &inverse_jacobians,
std::vector<Point<spacedim> > &cell_normal_vectors,
- CellSimilarity::Similarity &cell_similarity) const ;
+ const CellSimilarity::Similarity cell_similarity) const ;
/**
* Implementation of the interface in Mapping.
/**
* Implementation of the interface in Mapping.
*/
- virtual void
+ virtual
+ CellSimilarity::Similarity
fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const Quadrature<dim> &quadrature,
typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
std::vector<DerivativeForm<2,dim,spacedim> > &jacobian_grads,
std::vector<DerivativeForm<1,spacedim,dim> > &inverse_jacobians,
std::vector<Point<spacedim> > &cell_normal_vectors,
- CellSimilarity::Similarity &cell_similarity) const;
+ const CellSimilarity::Similarity cell_similarity) const;
/**
* Implementation of the interface in Mapping.
* Implementation of the interface in MappingQ1. Overrides the function in
* the base class, since we cannot use any cell similarity for this class.
*/
- virtual void
+ virtual
+ CellSimilarity::Similarity
fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const Quadrature<dim> &quadrature,
typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
std::vector<DerivativeForm<2,dim,spacedim> > &jacobian_grads,
std::vector<DerivativeForm<1,spacedim,dim> > &inverse_jacobians,
std::vector<Point<spacedim> > &cell_normal_vectors,
- CellSimilarity::Similarity &cell_similarity) const;
+ const CellSimilarity::Similarity cell_similarity) const;
/**
* Reference to the vector of shifts.
* Implementation of the interface in MappingQ. Overrides the function in
* the base class, since we cannot use any cell similarity for this class.
*/
- virtual void
+ virtual
+ CellSimilarity::Similarity
fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const Quadrature<dim> &quadrature,
typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
std::vector<DerivativeForm<2,dim,spacedim> > &jacobian_grads,
std::vector<DerivativeForm<1,spacedim,dim> > &inverse_jacobians,
std::vector<Point<spacedim> > &cell_normal_vectors,
- CellSimilarity::Similarity &cell_similarity) const;
+ const CellSimilarity::Similarity cell_similarity) const;
/**
* Reference to the vector of shifts.
template <int dim, int spacedim>
void FEValues<dim,spacedim>::do_reinit ()
{
- this->get_mapping().fill_fe_values(*this->present_cell,
- quadrature,
- *this->mapping_data,
- this->quadrature_points,
- this->JxW_values,
- this->jacobians,
- this->jacobian_grads,
- this->inverse_jacobians,
- this->normal_vectors,
- this->cell_similarity);
-
+ // first call the mapping and let it generate the data
+ // specific to the mapping. also let it inspect the
+ // cell similarity flag and, if necessary, update
+ // it
+ this->cell_similarity
+ = this->get_mapping().fill_fe_values(*this->present_cell,
+ quadrature,
+ *this->mapping_data,
+ this->quadrature_points,
+ this->JxW_values,
+ this->jacobians,
+ this->jacobian_grads,
+ this->inverse_jacobians,
+ this->normal_vectors,
+ this->cell_similarity);
+
+ // then call the finite element and, with the data
+ // already filled by the mapping, let it compute the
+ // data for the mapped shape function values, gradients,
+ // etc.
this->get_fe().fill_fe_values(this->get_mapping(),
*this->present_cell,
quadrature,
template<int dim, int spacedim>
-void
+CellSimilarity::Similarity
MappingCartesian<dim, spacedim>::
fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const Quadrature<dim> &q,
std::vector<DerivativeForm<2,dim,spacedim> > &jacobian_grads,
std::vector<DerivativeForm<1,spacedim,dim> > &inverse_jacobians,
std::vector<Point<spacedim> > &,
- CellSimilarity::Similarity &cell_similarity) const
+ const CellSimilarity::Similarity cell_similarity) const
{
// convert data object to internal
// data for this class. fails with
for (unsigned int j=0; j<dim; ++j)
inverse_jacobians[j][j]=1./data.length[j];
}
+
+ return cell_similarity;
}
// Note that the CellSimilarity flag is modifyable, since MappingFEField can need to
// recalculate data even when cells are similar.
template<int dim, int spacedim, class VECTOR, class DH>
-void
+CellSimilarity::Similarity
MappingFEField<dim,spacedim,VECTOR,DH>::fill_fe_values (
const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const Quadrature<dim> &q,
std::vector<DerivativeForm<2,dim,spacedim> > &jacobian_grads,
std::vector<DerivativeForm<1,spacedim,dim> > &inverse_jacobians,
std::vector<Point<spacedim> > &normal_vectors,
- CellSimilarity::Similarity &cell_similarity) const
+ const CellSimilarity::Similarity cell_similarity) const
{
AssertDimension(fe->dofs_per_cell, local_dof_values.get().size());
Assert(local_dof_values.get().size()>0, ExcMessage("Cannot do anything with zero degrees of freedom"));
Assert (dynamic_cast<InternalData *> (&mapping_data) != 0, ExcInternalError());
InternalData &data = static_cast<InternalData &> (mapping_data);
- if (get_degree() > 1)
- cell_similarity = CellSimilarity::invalid_next_cell;
-
const unsigned int n_q_points=q.size();
-
- compute_fill (cell, n_q_points, QProjector<dim>::DataSetDescriptor::cell (), cell_similarity,
+ const CellSimilarity::Similarity updated_cell_similarity
+ = (get_degree() == 1
+ ?
+ cell_similarity
+ :
+ CellSimilarity::invalid_next_cell);
+
+ compute_fill (cell, n_q_points, QProjector<dim>::DataSetDescriptor::cell (),
+ updated_cell_similarity,
data, quadrature_points);
const UpdateFlags update_flags(data.current_update_flags());
inverse_jacobians[point] = data.covariant[point].transpose();
}
+ return updated_cell_similarity;
}
// Note that the CellSimilarity flag is modifyable, since MappingQ can need to
// recalculate data even when cells are similar.
template<int dim, int spacedim>
-void
+CellSimilarity::Similarity
MappingQ<dim,spacedim>::fill_fe_values (
const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const Quadrature<dim> &q,
std::vector<DerivativeForm<2,dim,spacedim> > &jacobian_grads,
std::vector<DerivativeForm<1,spacedim,dim> > &inverse_jacobians,
std::vector<Point<spacedim> > &normal_vectors,
- CellSimilarity::Similarity &cell_similarity) const
+ const CellSimilarity::Similarity cell_similarity) const
{
// convert data object to internal data for this class. fails with an
// exception if that is not possible
|| cell->has_boundary_lines());
// depending on this result, use this or the other data object for the
- // mapping. furthermore, we need to ensure that the flag indicating whether
+ // mapping.
+ typename MappingQ1<dim,spacedim>::InternalData *
+ p_data = (data.use_mapping_q1_on_current_cell
+ ?
+ &data.mapping_q1_data
+ :
+ &data);
+
+ // call the base class. we need to ensure that the flag indicating whether
// we can use some similarity has to be modified - for a general MappingQ,
// the data needs to be recomputed anyway since then the mapping changes the
// data. this needs to be known also for later operations, so modify the
// variable here. this also affects the calculation of the next cell -- if
// we use Q1 data on the next cell, the data will still be invalid.
- typename MappingQ1<dim,spacedim>::InternalData *p_data=0;
- if (data.use_mapping_q1_on_current_cell)
- p_data=&data.mapping_q1_data;
- else
- {
- p_data=&data;
- if (get_degree() > 1)
- cell_similarity = CellSimilarity::invalid_next_cell;
- }
-
+ const CellSimilarity::Similarity updated_cell_similarity
+ = ((data.use_mapping_q1_on_current_cell == false)
+ &&
+ (get_degree() > 1)
+ ?
+ CellSimilarity::invalid_next_cell
+ :
+ cell_similarity);
MappingQ1<dim,spacedim>::fill_fe_values(cell, q, *p_data,
quadrature_points, JxW_values,
jacobians, jacobian_grads, inverse_jacobians,
- normal_vectors, cell_similarity);
+ normal_vectors,
+ updated_cell_similarity);
+ return updated_cell_similarity;
}
template<int dim, int spacedim>
-void
+CellSimilarity::Similarity
MappingQ1<dim,spacedim>::fill_fe_values (
const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const Quadrature<dim> &q,
std::vector<DerivativeForm<2,dim,spacedim> > &jacobian_grads,
std::vector<DerivativeForm<1,spacedim,dim> > &inverse_jacobians,
std::vector<Point<spacedim> > &normal_vectors,
- CellSimilarity::Similarity &cell_similarity) const
+ const CellSimilarity::Similarity cell_similarity) const
{
// ensure that the following static_cast is really correct:
Assert (dynamic_cast<InternalData *>(&mapping_data) != 0,
inverse_jacobians[point] = data.covariant[point].transpose();
}
+ return cell_similarity;
}
template<int dim, class EulerVectorType, int spacedim>
-void
+CellSimilarity::Similarity
MappingQ1Eulerian<dim,EulerVectorType,spacedim>::fill_fe_values (
const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const Quadrature<dim> &q,
std::vector<DerivativeForm<2,dim,spacedim> > &jacobian_grads,
std::vector<DerivativeForm<1,spacedim,dim> > &inverse_jacobians,
std::vector<Point<spacedim> > &normal_vectors,
- CellSimilarity::Similarity &cell_similarity) const
+ const CellSimilarity::Similarity) const
{
- // disable any previously detected
- // similarity and then enter the
- // respective function of the base class.
- cell_similarity = CellSimilarity::invalid_next_cell;
+ // call the function of the base class, but ignoring
+ // any potentially detected cell similarity between
+ // the current and the previous cell
MappingQ1<dim,spacedim>::fill_fe_values (cell, q, mapping_data,
quadrature_points, JxW_values, jacobians,
jacobian_grads, inverse_jacobians,
- normal_vectors, cell_similarity);
+ normal_vectors,
+ CellSimilarity::invalid_next_cell);
+ // also return the updated flag since any detected
+ // similarity wasn't based on the mapped field, but
+ // the original vertices which are meaningless
+ return CellSimilarity::invalid_next_cell;
}
template<int dim, class EulerVectorType, int spacedim>
-void
+CellSimilarity::Similarity
MappingQEulerian<dim,EulerVectorType,spacedim>::fill_fe_values (
const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const Quadrature<dim> &q,
std::vector<DerivativeForm<2,dim,spacedim> > &jacobian_grads,
std::vector<DerivativeForm<1,spacedim,dim> > &inverse_jacobians,
std::vector<Point<spacedim> > &normal_vectors,
- CellSimilarity::Similarity &cell_similarity) const
+ const CellSimilarity::Similarity ) const
{
- // disable any previously detected similarity and hand on to the respective
- // function of the base class.
- cell_similarity = CellSimilarity::invalid_next_cell;
+ // call the function of the base class, but ignoring
+ // any potentially detected cell similarity between
+ // the current and the previous cell
MappingQ<dim,spacedim>::fill_fe_values (cell, q, mapping_data,
quadrature_points, JxW_values, jacobians,
jacobian_grads, inverse_jacobians,
- normal_vectors, cell_similarity);
+ normal_vectors,
+ CellSimilarity::invalid_next_cell);
+ // also return the updated flag since any detected
+ // similarity wasn't based on the mapped field, but
+ // the original vertices which are meaningless
+ return CellSimilarity::invalid_next_cell;
}