24 #ifndef FACEDATACONTAINER_INTERNAL_H_
25 #define FACEDATACONTAINER_INTERNAL_H_
27 #include <deal.II/lac/vector.h>
45 template<
typename VECTOR,
int dim>
50 const std::map<std::string,
const dealii::Vector<double>*> ¶m_values,
51 const std::map<std::string, const VECTOR*> &domain_values,
65 GetParamValues(std::string name, dealii::Vector<double>& value)
const;
70 const std::map<std::string, const VECTOR*> &
73 return domain_values_;
77 virtual const dealii::FEFaceValuesBase<dim>&
79 virtual const dealii::FEFaceValuesBase<dim>&
82 virtual const dealii::FEFaceValuesBase<dim>&
84 virtual const dealii::FEFaceValuesBase<dim>&
91 const typename Triangulation<dim>::cell_iterator
104 std::vector<double>& values)
const;
112 std::vector<dealii::Vector<double> >& values)
const;
120 std::vector<double>& values)
const;
128 std::vector<dealii::Vector<double> >& values)
const;
134 template<
int targetdim>
137 std::vector<dealii::Tensor<1, targetdim> >& values)
const;
143 template<
int targetdim>
146 std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values)
const;
152 template<
int targetdim>
155 std::vector<dealii::Tensor<1, targetdim> >& values)
const;
161 template<
int targetdim>
164 std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values)
const;
171 std::vector<double>& values)
const;
178 std::vector<Vector<double> >& values)
const;
187 std::vector<double>& values)
const;
195 std::vector<Vector<double> >& values)
const;
202 template<
int targetdim>
205 std::vector<dealii::Tensor<1, targetdim> >& values)
const;
212 template<
int targetdim>
215 std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values)
const;
223 template<
int targetdim>
226 std::vector<dealii::Tensor<1, targetdim> >& values)
const;
232 template<
int targetdim>
235 std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values)
const;
261 return need_neighbour_;
270 GetValues(
const dealii::FEFaceValuesBase<dim>& fe_values,
271 std::string name, std::vector<double>& values)
const;
277 GetValues(
const dealii::FEFaceValuesBase<dim>& fe_values,
279 std::vector<dealii::Vector<double> >& values)
const;
284 template<
int targetdim>
286 GetGrads(
const dealii::FEFaceValuesBase<dim>& fe_values,
288 std::vector<dealii::Tensor<1, targetdim> >& values)
const;
293 template<
int targetdim>
295 GetGrads(
const dealii::FEFaceValuesBase<dim>& fe_values,
297 std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values)
const;
299 const std::map<std::string, const dealii::Vector<double>*> & param_values_;
300 const std::map<std::string, const VECTOR*> & domain_values_;
303 unsigned int subface_;
304 bool need_neighbour_;
308 template<
typename VECTOR,
int dim>
310 const std::map<std::string,
const dealii::Vector<double>*> ¶m_values,
311 const std::map<std::string, const VECTOR*> &domain_values,
313 : param_values_(param_values), domain_values_(domain_values), need_neighbour_(
318 template<
typename VECTOR,
int dim>
321 dealii::Vector<double>& value)
const
323 typename std::map<std::string, const dealii::Vector<double>*>::const_iterator it =
324 param_values_.find(name);
325 if (it == param_values_.end())
328 "FaceDataContainerInternal::GetParamValues");
330 value = *(it->second);
334 template<
typename VECTOR,
int dim>
335 const typename Triangulation<dim>::cell_iterator
338 return this->GetFEFaceValuesState().get_element();
342 template<
typename VECTOR,
int dim>
345 std::string name, std::vector<double>& values)
const
347 this->GetValues(this->GetFEFaceValuesState(), name, values);
350 template<
typename VECTOR,
int dim>
353 std::string name, std::vector<dealii::Vector<double> >& values)
const
355 this->GetValues(this->GetFEFaceValuesState(), name, values);
360 template<
typename VECTOR,
int dim>
363 std::string name, std::vector<double>& values)
const
365 this->GetValues(this->GetFEFaceValuesControl(), name, values);
369 template<
typename VECTOR,
int dim>
372 std::string name, std::vector<dealii::Vector<double> >& values)
const
374 this->GetValues(this->GetFEFaceValuesControl(), name, values);
378 template<
typename VECTOR,
int dim>
379 template<
int targetdim>
383 std::vector<dealii::Tensor<1, targetdim> >& values)
const
385 this->GetGrads<targetdim>(this->GetFEFaceValuesState(), name, values);
389 template<
typename VECTOR,
int dim>
390 template<
int targetdim>
394 std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values)
const
396 this->GetGrads<targetdim>(this->GetFEFaceValuesState(), name, values);
401 template<
typename VECTOR,
int dim>
402 template<
int targetdim>
406 std::vector<dealii::Tensor<1, targetdim> >& values)
const
408 this->GetGrads<targetdim>(this->GetFEFaceValuesControl(), name,
413 template<
typename VECTOR,
int dim>
414 template<
int targetdim>
418 std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values)
const
420 this->GetGrads<targetdim>(this->GetFEFaceValuesControl(), name,
425 template<
typename VECTOR,
int dim>
428 std::string name, std::vector<double>& values)
const
430 this->GetValues(this->GetNbrFEFaceValuesState(), name, values);
433 template<
typename VECTOR,
int dim>
436 std::string name, std::vector<Vector<double> >& values)
const
438 this->GetValues(this->GetNbrFEFaceValuesState(), name, values);
443 template<
typename VECTOR,
int dim>
446 std::string name, std::vector<double>& values)
const
448 this->GetValues(this->GetNbrFEFaceValuesControl(), name, values);
452 template<
typename VECTOR,
int dim>
455 std::string name, std::vector<Vector<double> >& values)
const
457 this->GetValues(this->GetNbrFEFaceValuesControl(), name, values);
461 template<
typename VECTOR,
int dim>
462 template<
int targetdim>
465 std::string name, std::vector<Tensor<1, targetdim> >& values)
const
467 this->GetGrads<targetdim>(this->GetNbrFEFaceValuesState(), name,
472 template<
typename VECTOR,
int dim>
473 template<
int targetdim>
477 std::vector<std::vector<Tensor<1, targetdim> > >& values)
const
479 this->GetGrads<targetdim>(this->GetNbrFEFaceValuesState(), name,
485 template<
typename VECTOR,
int dim>
486 template<
int targetdim>
489 std::string name, std::vector<Tensor<1, targetdim> >& values)
const
491 this->GetGrads<targetdim>(this->GetNbrFEFaceValuesControl(), name,
496 template<
typename VECTOR,
int dim>
497 template<
int targetdim>
501 std::vector<std::vector<Tensor<1, targetdim> > >& values)
const
503 this->GetGrads<targetdim>(this->GetNbrFEFaceValuesControl(), name,
510 template<
typename VECTOR,
int dim>
513 const dealii::FEFaceValuesBase<dim>& fe_values, std::string name,
514 std::vector<double>& values)
const
516 typename std::map<std::string, const VECTOR*>::const_iterator it =
517 this->GetDomainValues().find(name);
518 if (it == this->GetDomainValues().end())
521 "ElementDataContainer::GetValues");
523 fe_values.get_function_values(*(it->second), values);
527 template<
typename VECTOR,
int dim>
529 FaceDataContainerInternal<VECTOR, dim>::GetValues(
530 const dealii::FEFaceValuesBase<dim>& fe_values, std::string name,
531 std::vector<dealii::Vector<double> >& values)
const
533 typename std::map<std::string, const VECTOR*>::const_iterator it =
534 this->GetDomainValues().find(name);
535 if (it == this->GetDomainValues().end())
537 throw DOpEException(
"Did not find " + name,
538 "ElementDataContainer::GetValues");
540 fe_values.get_function_values(*(it->second), values);
545 template<
typename VECTOR,
int dim>
546 template<
int targetdim>
548 FaceDataContainerInternal<VECTOR, dim>::GetGrads(
549 const dealii::FEFaceValuesBase<dim>& fe_values, std::string name,
550 std::vector<dealii::Tensor<1, targetdim> >& values)
const
552 typename std::map<std::string, const VECTOR*>::const_iterator it =
553 this->GetDomainValues().find(name);
554 if (it == this->GetDomainValues().end())
556 throw DOpEException(
"Did not find " + name,
557 "FaceDataContainerInternal::GetGrads");
559 fe_values.get_function_gradients(*(it->second), values);
564 template<
typename VECTOR,
int dim>
565 template<
int targetdim>
567 FaceDataContainerInternal<VECTOR, dim>::GetGrads(
568 const dealii::FEFaceValuesBase<dim>& fe_values, std::string name,
569 std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values)
const
571 typename std::map<std::string, const VECTOR*>::const_iterator it =
572 this->GetDomainValues().find(name);
573 if (it == this->GetDomainValues().end())
575 throw DOpEException(
"Did not find " + name,
576 "FaceDataContainerInternal::GetGrads");
578 fe_values.get_function_gradients(*(it->second), values);
void GetFaceValuesControl(std::string name, std::vector< double > &values) const
Definition: facedatacontainer_internal.h:362
void SetSubFace(unsigned int subface)
Definition: facedatacontainer_internal.h:249
Definition: facedatacontainer_internal.h:46
FaceDataContainerInternal(const std::map< std::string, const dealii::Vector< double > * > ¶m_values, const std::map< std::string, const VECTOR * > &domain_values, bool need_neighbour)
Definition: facedatacontainer_internal.h:309
bool NeedNeighbour() const
Definition: facedatacontainer_internal.h:259
unsigned int GetFace() const
Definition: facedatacontainer_internal.h:244
void GetNbrFaceValuesState(std::string name, std::vector< double > &values) const
Definition: facedatacontainer_internal.h:427
void GetNbrFaceGradsControl(std::string name, std::vector< dealii::Tensor< 1, targetdim > > &values) const
const Triangulation< dim >::cell_iterator GetElementState() const
Definition: facedatacontainer_internal.h:336
void GetNbrFaceGradsState(std::string name, std::vector< dealii::Tensor< 1, targetdim > > &values) const
void GetParamValues(std::string name, dealii::Vector< double > &value) const
Definition: facedatacontainer_internal.h:320
virtual const dealii::FEFaceValuesBase< dim > & GetNbrFEFaceValuesState() const =0
virtual const dealii::FEFaceValuesBase< dim > & GetNbrFEFaceValuesControl() const =0
void GetFaceGradsState(std::string name, std::vector< dealii::Tensor< 1, targetdim > > &values) const
Definition: facedatacontainer_internal.h:381
unsigned int GetSubFace() const
Definition: facedatacontainer_internal.h:254
const std::map< std::string, const VECTOR * > & GetDomainValues() const
Definition: facedatacontainer_internal.h:71
void GetFaceGradsControl(std::string name, std::vector< dealii::Tensor< 1, targetdim > > &values) const
Definition: facedatacontainer_internal.h:404
virtual const dealii::FEFaceValuesBase< dim > & GetFEFaceValuesControl() const =0
virtual const dealii::FEFaceValuesBase< dim > & GetFEFaceValuesState() const =0
virtual ~FaceDataContainerInternal()
Definition: facedatacontainer_internal.h:55
void SetFace(unsigned int face)
Definition: facedatacontainer_internal.h:239
void GetFaceValuesState(std::string name, std::vector< double > &values) const
Definition: facedatacontainer_internal.h:344
void GetNbrFaceValuesControl(std::string name, std::vector< double > &values) const
Definition: facedatacontainer_internal.h:445
Definition: dopeexception.h:35