24 #ifndef MULTIMESH_FACEDATACONTAINER_H_
25 #define MULTIMESH_FACEDATACONTAINER_H_
32 #include <deal.II/dofs/dof_handler.h>
33 #include <deal.II/hp/dof_handler.h>
35 using namespace dealii;
52 template<
template<
int,
int>
class DH,
typename VECTOR,
int dim>
59 "Dummy class, this constructor should never get called.",
60 "FaceDataContainer<dealii::DoFHandler , VECTOR, dim>::FaceDataContainer"));
74 template<
typename VECTOR,
int dim>
98 template<
template<
int,
int>
class FE,
typename SPARSITYPATTERN,
int dopedim,
101 const Quadrature<dim - 1>& quad,
102 UpdateFlags update_flags,
104 VECTOR, dopedim, dealdim>& sth,
105 const typename std::vector<
typename dealii::DoFHandler<dim>::cell_iterator>& element,
106 const typename std::vector<
typename dealii::Triangulation<dim>::cell_iterator>& tria_element,
107 const std::map<std::string,
const Vector<double>*> ¶m_values,
108 const std::map<std::string, const VECTOR*> &domain_values,
bool) :
109 param_values_(param_values),
110 domain_values_(domain_values),
112 tria_element_(tria_element),
113 state_fe_values_(sth.GetFESystem(
"state"), quad,
115 control_fe_values_(sth.GetFESystem(
"control"), quad,
118 state_index_ = sth.GetStateIndex();
119 if (state_index_ == 1)
123 n_q_points_per_element_ = quad.size();
124 n_dofs_per_element_ = element[0]->get_fe().dofs_per_cell;
125 control_prolongation_ = IdentityMatrix(element_[this->GetControlIndex()]->get_fe().dofs_per_cell);
126 state_prolongation_ = IdentityMatrix(element_[this->GetStateIndex()]->get_fe().dofs_per_cell);
141 ReInit(
unsigned int coarse_index,
unsigned int fine_index,
const FullMatrix<double>& prolongation_matrix,
unsigned int face_no);
156 GetNDoFsPerElement()
const;
158 GetNbrNDoFsPerElement()
const;
162 GetMaterialId()
const;
164 GetNbrMaterialId()
const;
166 GetNbrMaterialId(
unsigned int face)
const;
168 GetIsAtBoundary()
const;
170 GetElementDiameter()
const;
172 GetBoundaryIndicator()
const;
174 GetFEFaceValuesState()
const;
176 GetFEFaceValuesControl()
const;
184 GetParamValues(std::string name, Vector<double>& value)
const;
195 GetFaceValuesState(std::string name, std::vector<double>& values)
const;
201 GetFaceValuesState(std::string name, std::vector<Vector<double> >& values)
const;
208 GetFaceValuesControl(std::string name, std::vector<double>& values)
const;
214 GetFaceValuesControl(std::string name, std::vector<Vector<double> >& values)
const;
220 template<
int targetdim>
222 GetFaceGradsState(std::string name,
223 std::vector<Tensor<1, targetdim> >& values)
const;
229 template<
int targetdim>
231 GetFaceGradsState(std::string name,
232 std::vector<std::vector<Tensor<1, targetdim> > >& values)
const;
238 template<
int targetdim>
240 GetFaceGradsControl(std::string name,
241 std::vector<Tensor<1, targetdim> >& values)
const;
247 template<
int targetdim>
249 GetFaceGradsControl(std::string name,
250 std::vector<std::vector<Tensor<1, targetdim> > >& values)
const;
257 GetStateIndex()
const;
259 GetControlIndex()
const;
261 GetCoarseIndex()
const {
return coarse_index_; }
263 GetFineIndex()
const {
return fine_index_; }
269 GetValues(
typename dealii::DoFHandler<dim>::cell_iterator element,
270 const FullMatrix<double>& prolongation,
272 std::vector<double>& values)
const;
278 GetValues(
typename dealii::DoFHandler<dim>::cell_iterator element,
279 const FullMatrix<double>& prolongation,
281 std::vector<Vector<double> >& values)
const;
286 template<
int targetdim>
288 GetGrads(
typename dealii::DoFHandler<dim>::cell_iterator element,
289 const FullMatrix<double>& prolongation,
291 std::string name, std::vector<Tensor<1, targetdim> >& values)
const;
296 template<
int targetdim>
298 GetGrads(
typename dealii::DoFHandler<dim>::cell_iterator element,
299 const FullMatrix<double>& prolongation,
301 std::string name, std::vector<std::vector<Tensor<1, targetdim> > >& values)
const;
303 inline const std::map<std::string, const VECTOR*> &
304 GetDomainValues()
const
306 return domain_values_;
310 const std::map<std::string, const Vector<double>*> ¶m_values_;
311 const std::map<std::string, const VECTOR*> &domain_values_;
312 unsigned int state_index_;
313 unsigned int control_index_;
315 const typename std::vector<typename dealii::DoFHandler<dim>::cell_iterator>& element_;
316 const typename std::vector<typename dealii::Triangulation<dim>::cell_iterator>& tria_element_;
320 FullMatrix<double> control_prolongation_;
321 FullMatrix<double> state_prolongation_;
322 unsigned int coarse_index_, fine_index_;
324 unsigned int n_q_points_per_element_;
325 unsigned int n_dofs_per_element_;
335 template<
typename VECTOR,
int dim>
338 unsigned int coarse_index,
339 unsigned int fine_index,
340 const FullMatrix<double>& prolongation_matrix,
341 unsigned int face_no)
344 coarse_index_ = coarse_index;
345 fine_index_ = fine_index;
346 assert(this->GetControlIndex() < element_.size());
348 if(coarse_index == this->GetStateIndex())
350 state_prolongation_ = prolongation_matrix;
351 control_prolongation_ = IdentityMatrix(element_[this->GetControlIndex()]->get_fe().dofs_per_cell);
355 if(coarse_index == this->GetControlIndex())
357 control_prolongation_ = prolongation_matrix;
358 state_prolongation_ = IdentityMatrix(element_[this->GetStateIndex()]->get_fe().dofs_per_cell);
362 control_prolongation_ = IdentityMatrix(element_[this->GetControlIndex()]->get_fe().dofs_per_cell);
363 state_prolongation_ = IdentityMatrix(element_[this->GetStateIndex()]->get_fe().dofs_per_cell);
368 state_fe_values_.reinit(tria_element_[GetFineIndex()], face_no);
369 control_fe_values_.reinit(tria_element_[GetFineIndex()], face_no);
372 template<
typename VECTOR,
int dim>
376 return n_dofs_per_element_;
381 template<
typename VECTOR,
int dim>
385 throw DOpEException(
"This function has not been written since we do not know what the right neigbour is!",
386 "Multimesh_FaceDataContainer::GetNbrNDoFsPerElement");
390 template<
typename VECTOR,
int dim>
394 return n_q_points_per_element_;
398 template<
typename VECTOR,
int dim>
402 return tria_element_[GetFineIndex()]->material_id();
406 template<
typename VECTOR,
int dim>
410 return this->GetNbrMaterialId(face_);
413 template<
typename VECTOR,
int dim>
418 throw DOpEException(
"This function has not been written since we do not know what the right neigbour is!",
419 "Multimesh_FaceDataContainer::GetNbrMaterialId");
423 template<
typename VECTOR,
int dim>
427 return tria_element_[GetFineIndex()]->face(face_)->at_boundary();
431 template<
typename VECTOR,
int dim>
435 throw DOpEException(
"This function has not been written since we do not know what the right Diameter!",
436 "Multimesh_FaceDataContainer::GetElementDiameter");
441 template<
typename VECTOR,
int dim>
445 return tria_element_[GetFineIndex()]->face(face_)->boundary_indicator();
449 template<
typename VECTOR,
int dim>
453 return state_fe_values_;
457 template<
typename VECTOR,
int dim>
461 return control_fe_values_;
466 template<
typename VECTOR,
int dim>
469 std::string name, Vector<double>& value)
const
471 typename std::map<std::string, const Vector<double>*>::const_iterator it =
472 param_values_.find(name);
473 if (it == param_values_.end())
476 "Multimesh_FaceDataContainer::GetParamValues");
478 value = *(it->second);
482 template<
typename VECTOR,
int dim>
485 std::string name, std::vector<double>& values)
const
487 this->GetValues(element_[this->GetStateIndex()],state_prolongation_,this->GetFEFaceValuesState(), name, values);
490 template<
typename VECTOR,
int dim>
493 std::string name, std::vector<Vector<double> >& values)
const
495 this->GetValues(element_[this->GetStateIndex()],state_prolongation_,this->GetFEFaceValuesState(), name, values);
500 template<
typename VECTOR,
int dim>
503 std::string name, std::vector<double>& values)
const
505 this->GetValues(element_[this->GetControlIndex()],control_prolongation_,this->GetFEFaceValuesControl(), name, values);
509 template<
typename VECTOR,
int dim>
512 std::string name, std::vector<Vector<double> >& values)
const
514 this->GetValues(element_[this->GetControlIndex()],control_prolongation_,this->GetFEFaceValuesControl(), name, values);
518 template<
typename VECTOR,
int dim>
519 template<
int targetdim>
522 std::string name, std::vector<Tensor<1, targetdim> >& values)
const
524 this->GetGrads<targetdim> (element_[this->GetStateIndex()],state_prolongation_,this->GetFEFaceValuesState(), name, values);
528 template<
typename VECTOR,
int dim>
529 template<
int targetdim>
532 std::string name, std::vector<std::vector<Tensor<1, targetdim> > >& values)
const
534 this->GetGrads<targetdim> (element_[this->GetStateIndex()],state_prolongation_,this->GetFEFaceValuesState(), name, values);
539 template<
typename VECTOR,
int dim>
540 template<
int targetdim>
543 std::string name, std::vector<Tensor<1, targetdim> >& values)
const
545 this->GetGrads<targetdim> (element_[this->GetControlIndex()],control_prolongation_,this->GetFEFaceValuesControl(), name, values);
549 template<
typename VECTOR,
int dim>
550 template<
int targetdim>
553 std::string name, std::vector<std::vector<Tensor<1, targetdim> > >& values)
const
555 this->GetGrads<targetdim> (element_[this->GetControlIndex()],control_prolongation_,this->GetFEFaceValuesControl(), name, values);
560 template<
typename VECTOR,
int dim>
569 template<
typename VECTOR,
int dim>
571 Multimesh_FaceDataContainer<dealii::DoFHandler, VECTOR, dim>::GetControlIndex()
const
573 return control_index_;
577 template<
typename VECTOR,
int dim>
579 Multimesh_FaceDataContainer<dealii::DoFHandler, VECTOR, dim>::GetValues(
580 typename dealii::DoFHandler<dim>::cell_iterator element,
581 const FullMatrix<double>& prolongation,
583 std::vector<double>& values)
const
585 typename std::map<std::string, const VECTOR*>::const_iterator it =
586 this->GetDomainValues().find(name);
587 if (it == this->GetDomainValues().end())
589 throw DOpEException(
"Did not find " + name,
590 "Multimesh_FaceDataContainer::GetValues");
592 unsigned int dofs_per_element = element->get_fe().dofs_per_cell;
594 dealii::Vector<double> dof_values(dofs_per_element);
595 dealii::Vector<double> dof_values_transformed(dofs_per_element);
596 element->get_dof_values (*(it->second), dof_values);
598 prolongation.vmult(dof_values_transformed,dof_values);
602 unsigned int n_quadrature_points = GetNQPoints();
603 std::fill_n (values.begin(), n_quadrature_points, 0);
604 for (
unsigned int shape_func=0; shape_func<dofs_per_element; ++shape_func)
606 const double value = dof_values_transformed(shape_func);
610 const double *shape_value_ptr = &(fe_values.shape_value(shape_func, 0));
611 for (
unsigned int point=0; point<n_quadrature_points; ++point)
612 values[point] += value * *shape_value_ptr++;
617 template<
typename VECTOR,
int dim>
619 Multimesh_FaceDataContainer<dealii::DoFHandler, VECTOR, dim>::GetValues(
620 typename dealii::DoFHandler<dim>::cell_iterator element,
621 const FullMatrix<double>& prolongation,
623 std::vector<Vector<double> >& values)
const
625 typename std::map<std::string, const VECTOR*>::const_iterator it =
626 this->GetDomainValues().find(name);
627 if (it == this->GetDomainValues().end())
629 throw DOpEException(
"Did not find " + name,
630 "Multimesh_FaceDataContainer::GetValues");
633 unsigned int dofs_per_element = element->get_fe().dofs_per_cell;
635 dealii::Vector<double> dof_values(dofs_per_element);
636 dealii::Vector<double> dof_values_transformed(dofs_per_element);
637 element->get_dof_values (*(it->second), dof_values);
639 prolongation.vmult(dof_values_transformed,dof_values);
643 const unsigned int n_components = element->get_fe().n_components();
644 unsigned int n_quadrature_points = GetNQPoints();
645 for (
unsigned i=0;i<values.size();++i)
646 std::fill_n (values[i].begin(), values[i].size(), 0);
648 for (
unsigned int shape_func=0; shape_func<dofs_per_element; ++shape_func)
650 const double value = dof_values_transformed(shape_func);
654 if (element->get_fe().is_primitive(shape_func))
656 const unsigned int comp = element->get_fe().system_to_component_index(shape_func).first;
657 for (
unsigned int point=0; point<n_quadrature_points; ++point)
658 values[point](comp) += value * fe_values.shape_value(shape_func,point);
661 for (
unsigned int c=0; c<n_components; ++c)
663 if (element->get_fe().get_nonzero_components(shape_func)[c] ==
false)
665 for (
unsigned int point=0; point<n_quadrature_points; ++point)
666 values[point](c) += value * fe_values.shape_value_component(shape_func,point,c);
673 template<
typename VECTOR,
int dim>
674 template<
int targetdim>
676 Multimesh_FaceDataContainer<dealii::DoFHandler, VECTOR, dim>::GetGrads(
677 typename dealii::DoFHandler<dim>::cell_iterator element,
678 const FullMatrix<double>& prolongation,
680 std::vector<Tensor<1, targetdim> >& values)
const
682 typename std::map<std::string, const VECTOR*>::const_iterator it =
683 this->GetDomainValues().find(name);
684 if (it == this->GetDomainValues().end())
686 throw DOpEException(
"Did not find " + name,
687 "Multimesh_FaceDataContainerBase::GetGrads");
689 unsigned int dofs_per_element = element->get_fe().dofs_per_cell;
691 dealii::Vector<double> dof_values(dofs_per_element);
692 dealii::Vector<double> dof_values_transformed(dofs_per_element);
693 element->get_dof_values (*(it->second), dof_values);
695 prolongation.vmult(dof_values_transformed,dof_values);
698 unsigned int n_quadrature_points = GetNQPoints();
699 std::fill_n (values.begin(), n_quadrature_points, Tensor<1,targetdim>());
701 for (
unsigned int shape_func=0; shape_func<dofs_per_element; ++shape_func)
703 const double value = dof_values_transformed(shape_func);
707 const Tensor<1,targetdim> *shape_gradient_ptr
708 = &(fe_values.shape_grad(shape_func,0));
709 for (
unsigned int point=0; point<n_quadrature_points; ++point)
710 values[point] += value * *shape_gradient_ptr++;
716 template<
typename VECTOR,
int dim>
717 template<
int targetdim>
719 Multimesh_FaceDataContainer<dealii::DoFHandler, VECTOR, dim>::GetGrads(
720 typename dealii::DoFHandler<dim>::cell_iterator element,
721 const FullMatrix<double>& prolongation,
723 std::vector<std::vector<Tensor<1, targetdim> > >& values)
const
725 typename std::map<std::string, const VECTOR*>::const_iterator it =
726 this->GetDomainValues().find(name);
727 if (it == this->GetDomainValues().end())
729 throw DOpEException(
"Did not find " + name,
730 "Multimesh_FaceDataContainerBase::GetGrads");
733 unsigned int dofs_per_element = element->get_fe().dofs_per_cell;
735 dealii::Vector<double> dof_values(dofs_per_element);
736 dealii::Vector<double> dof_values_transformed(dofs_per_element);
737 element->get_dof_values (*(it->second), dof_values);
739 prolongation.vmult(dof_values_transformed,dof_values);
742 const unsigned int n_components = element->get_fe().n_components();
743 unsigned int n_quadrature_points = GetNQPoints();
744 for (
unsigned i=0;i<values.size();++i)
745 std::fill_n (values[i].begin(), values[i].size(), Tensor<1,dim>());
747 for (
unsigned int shape_func=0; shape_func<dofs_per_element; ++shape_func)
749 const double value = dof_values_transformed(shape_func);
753 if (element->get_fe().is_primitive(shape_func))
755 const unsigned int comp = element->get_fe().system_to_component_index(shape_func).first;
756 for (
unsigned int point=0; point<n_quadrature_points; ++point)
757 values[point][comp] += value * fe_values.shape_grad(shape_func,point);
760 for (
unsigned int c=0; c<n_components; ++c)
762 if (element->get_fe().get_nonzero_components(shape_func)[c] ==
false)
764 for (
unsigned int point=0; point<n_quadrature_points; ++point)
765 values[point][c] += value * fe_values.shape_grad_component(shape_func,point,c);
~Multimesh_FaceDataContainer()
Definition: multimesh_facedatacontainer.h:129
Multimesh_FaceDataContainer(const Quadrature< dim-1 > &quad, UpdateFlags update_flags, SpaceTimeHandler< FE, dealii::DoFHandler, SPARSITYPATTERN, VECTOR, dopedim, dealdim > &sth, const typename std::vector< typename dealii::DoFHandler< dim >::cell_iterator > &element, const typename std::vector< typename dealii::Triangulation< dim >::cell_iterator > &tria_element, const std::map< std::string, const Vector< double > * > ¶m_values, const std::map< std::string, const VECTOR * > &domain_values, bool)
Definition: multimesh_facedatacontainer.h:100
void ReInitNbr()
Definition: multimesh_facedatacontainer.h:148
Multimesh_FaceDataContainer()
Definition: multimesh_facedatacontainer.h:56
Definition: spacetimehandler.h:71
Definition: multimesh_facedatacontainer.h:53
Definition: fevalues_wrapper.h:88
Definition: dopeexception.h:35