24 #ifndef ELEMENTDATACONTAINER_INTERNAL_H_
25 #define ELEMENTDATACONTAINER_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
52 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_;
87 const typename Triangulation<dim>::cell_iterator
100 GetValuesState(std::string name, std::vector<double>& values)
const;
108 std::vector<dealii::Vector<double> >& values)
const;
123 std::vector<dealii::Vector<double> >& values)
const;
129 template<
int targetdim>
132 std::vector<dealii::Tensor<1, targetdim> >& values)
const;
138 template<
int targetdim>
142 std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values)
const;
148 template<
int targetdim>
151 std::vector<dealii::Tensor<1, targetdim> >& values)
const;
157 template<
int targetdim>
161 std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values)
const;
166 template<
int targetdim>
169 std::vector<dealii::Tensor<2, targetdim> >& values)
const;
175 template<
int targetdim>
179 std::vector<std::vector<dealii::Tensor<2, targetdim> > >& values)
const;
185 template<
int targetdim>
188 std::vector<dealii::Tensor<2, targetdim> >& values)
const;
194 template<
int targetdim>
198 std::vector<std::vector<dealii::Tensor<2, targetdim> > >& values)
const;
208 std::vector<double> & values)
const;
217 std::vector<dealii::Vector<double> >& values)
const;
227 std::vector<double> & values)
const;
236 std::vector<dealii::Vector<double> >& values)
const;
245 std::string name, std::vector<double>& values)
const;
253 std::vector<dealii::Vector<double> >& values)
const;
258 template<
int targetdim>
262 std::vector<dealii::Tensor<1, targetdim> >& values)
const;
267 template<
int targetdim>
272 std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values)
const;
279 std::string name, std::vector<double>& values)
const;
288 std::vector<dealii::Vector<double> >& values)
const;
294 template<
int targetdim>
298 std::vector<dealii::Tensor<2, targetdim> > & values)
const;
304 template<
int targetdim>
309 std::vector<std::vector<dealii::Tensor<2, targetdim> > >& values)
const;
311 const std::map<std::string, const dealii::Vector<double>*> ¶m_values_;
312 const std::map<std::string, const VECTOR*> &domain_values_;
316 template<
typename VECTOR,
int dim>
318 const std::map<std::string,
const dealii::Vector<double>*> ¶m_values,
319 const std::map<std::string, const VECTOR*> &domain_values)
320 : param_values_(param_values), domain_values_(domain_values)
324 template<
typename VECTOR,
int dim>
327 dealii::Vector<double>& value)
const
329 typename std::map<std::string, const dealii::Vector<double>*>::const_iterator it =
330 param_values_.find(name);
331 if (it == param_values_.end())
334 "ElementDataContainerInternal::GetParamValues");
336 value = *(it->second);
340 template<
typename VECTOR,
int dim>
343 std::vector<double>& values)
const
345 this->GetValues(this->GetFEValuesState(), name, values);
348 template<
typename VECTOR,
int dim>
351 std::vector<dealii::Vector<double> >& values)
const
353 this->GetValues(this->GetFEValuesState(), name, values);
359 template<
typename VECTOR,
int dim>
360 const typename Triangulation<dim>::cell_iterator
363 return this->GetFEValuesState().get_element();
367 template<
typename VECTOR,
int dim>
370 std::vector<double>& values)
const
372 this->GetValues(this->GetFEValuesControl(), name, values);
376 template<
typename VECTOR,
int dim>
379 std::vector<dealii::Vector<double> >& values)
const
381 this->GetValues(this->GetFEValuesControl(), name, values);
385 template<
typename VECTOR,
int dim>
386 template<
int targetdim>
389 std::vector<dealii::Tensor<1, targetdim> >& values)
const
391 this->GetGrads<targetdim>(this->GetFEValuesState(), name, values);
395 template<
typename VECTOR,
int dim>
396 template<
int targetdim>
400 std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values)
const
402 this->GetGrads<targetdim>(this->GetFEValuesState(), name, values);
407 template<
typename VECTOR,
int dim>
408 template<
int targetdim>
412 std::vector<dealii::Tensor<1, targetdim> >& values)
const
414 this->GetGrads<targetdim>(this->GetFEValuesControl(), name, values);
419 template<
typename VECTOR,
int dim>
420 template<
int targetdim>
424 std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values)
const
426 this->GetGrads<targetdim>(this->GetFEValuesControl(), name, values);
431 template<
typename VECTOR,
int dim>
432 template<
int targetdim>
436 std::vector<std::vector<dealii::Tensor<2, targetdim> > >& values)
const
438 this->GetHessians<targetdim>(this->GetFEValuesState(), name, values);
443 template<
typename VECTOR,
int dim>
444 template<
int targetdim>
448 std::vector<dealii::Tensor<2, targetdim> >& values)
const
450 this->GetHessians<targetdim>(this->GetFEValuesState(), name, values);
455 template<
typename VECTOR,
int dim>
456 template<
int targetdim>
460 std::vector<std::vector<dealii::Tensor<2, targetdim> > >& values)
const
462 this->GetHessians<targetdim>(this->GetFEValuesControl(), name, values);
467 template<
typename VECTOR,
int dim>
468 template<
int targetdim>
472 std::vector<dealii::Tensor<2, targetdim> >& values)
const
474 this->GetHessians<targetdim>(this->GetFEValuesControl(), name,
479 template<
typename VECTOR,
int dim>
482 std::string name, std::vector<double> & values)
const
484 this->GetLaplacians(this->GetFEValuesState(), name, values);
488 template<
typename VECTOR,
int dim>
491 std::string name, std::vector<dealii::Vector<double> >& values)
const
493 this->GetLaplacians(this->GetFEValuesState(), name, values);
496 template<
typename VECTOR,
int dim>
499 std::string name, std::vector<double> & values)
const
501 this->GetLaplacians(this->GetFEValuesControl(), name, values);
505 template<
typename VECTOR,
int dim>
508 std::string name, std::vector<dealii::Vector<double> >& values)
const
510 this->GetLaplacians(this->GetFEValuesControl(), name, values);
513 template<
typename VECTOR,
int dim>
517 std::vector<double>& values)
const
519 typename std::map<std::string, const VECTOR*>::const_iterator it =
520 this->GetDomainValues().find(name);
521 if (it == this->GetDomainValues().end())
524 "ElementDataContainer::GetValues");
526 fe_values.get_function_values(*(it->second), values);
530 template<
typename VECTOR,
int dim>
532 ElementDataContainerInternal<VECTOR, dim>::GetValues(
534 std::vector<dealii::Vector<double> >& values)
const
536 typename std::map<std::string, const VECTOR*>::const_iterator it =
537 this->GetDomainValues().find(name);
538 if (it == this->GetDomainValues().end())
540 throw DOpEException(
"Did not find " + name,
541 "ElementDataContainer::GetValues");
543 fe_values.get_function_values(*(it->second), values);
548 template<
typename VECTOR,
int dim>
549 template<
int targetdim>
551 ElementDataContainerInternal<VECTOR, dim>::GetGrads(
553 std::vector<dealii::Tensor<1, targetdim> >& values)
const
555 typename std::map<std::string, const VECTOR*>::const_iterator it =
556 this->GetDomainValues().find(name);
557 if (it == this->GetDomainValues().end())
559 throw DOpEException(
"Did not find " + name,
560 "ElementDataContainerInternal::GetGrads");
562 fe_values.get_function_gradients(*(it->second), values);
567 template<
typename VECTOR,
int dim>
568 template<
int targetdim>
570 ElementDataContainerInternal<VECTOR, dim>::GetGrads(
573 std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values)
const
575 typename std::map<std::string, const VECTOR*>::const_iterator it =
576 this->GetDomainValues().find(name);
577 if (it == this->GetDomainValues().end())
579 throw DOpEException(
"Did not find " + name,
580 "ElementDataContainerInternal::GetGrads");
582 fe_values.get_function_gradients(*(it->second), values);
587 template<
typename VECTOR,
int dim>
589 ElementDataContainerInternal<VECTOR, dim>::GetLaplacians(
591 std::vector<double> & values)
const
593 typename std::map<std::string, const VECTOR*>::const_iterator it =
594 this->GetDomainValues().find(name);
595 if (it == this->GetDomainValues().end())
597 throw DOpEException(
"Did not find " + name,
598 "ElementDataContainerInternal::GetLaplacians");
600 fe_values.get_function_laplacians(*(it->second), values);
605 template<
typename VECTOR,
int dim>
607 ElementDataContainerInternal<VECTOR, dim>::GetLaplacians(
609 std::vector<dealii::Vector<double> >& values)
const
611 typename std::map<std::string, const VECTOR*>::const_iterator it =
612 this->GetDomainValues().find(name);
613 if (it == this->GetDomainValues().end())
615 throw DOpEException(
"Did not find " + name,
616 "ElementDataContainerInternal::GetLaplacians");
618 fe_values.get_function_laplacians(*(it->second), values);
623 template<
typename VECTOR,
int dim>
624 template<
int targetdim>
626 ElementDataContainerInternal<VECTOR, dim>::GetHessians(
629 std::vector<std::vector<dealii::Tensor<2, targetdim> > >& values)
const
631 typename std::map<std::string, const VECTOR*>::const_iterator it =
632 this->GetDomainValues().find(name);
633 if (it == this->GetDomainValues().end())
635 throw DOpEException(
"Did not find " + name,
636 "ElementDataContainerInternal::GetGrads");
638 fe_values.get_function_hessians(*(it->second), values);
643 template<
typename VECTOR,
int dim>
644 template<
int targetdim>
646 ElementDataContainerInternal<VECTOR, dim>::GetHessians(
648 std::vector<dealii::Tensor<2, targetdim> >& values)
const
650 typename std::map<std::string, const VECTOR*>::const_iterator it =
651 this->GetDomainValues().find(name);
652 if (it == this->GetDomainValues().end())
654 throw DOpEException(
"Did not find " + name,
655 "ElementDataContainerInternal::GetGrads");
657 fe_values.get_function_hessians(*(it->second), values);
virtual const DOpEWrapper::FEValues< dim > & GetFEValuesControl() const =0
void GetLaplaciansState(std::string name, std::vector< double > &values) const
Definition: elementdatacontainer_internal.h:481
ElementDataContainerInternal(const std::map< std::string, const dealii::Vector< double > * > ¶m_values, const std::map< std::string, const VECTOR * > &domain_values)
Definition: elementdatacontainer_internal.h:317
void GetHessiansControl(std::string name, std::vector< dealii::Tensor< 2, targetdim > > &values) const
Definition: elementdatacontainer_internal.h:470
virtual const DOpEWrapper::FEValues< dim > & GetFEValuesState() const =0
virtual ~ElementDataContainerInternal()
Definition: elementdatacontainer_internal.h:55
const Triangulation< dim >::cell_iterator GetElementState() const
Definition: elementdatacontainer_internal.h:361
void GetLaplaciansControl(std::string name, std::vector< double > &values) const
Definition: elementdatacontainer_internal.h:498
void GetGradsControl(std::string name, std::vector< dealii::Tensor< 1, targetdim > > &values) const
Definition: elementdatacontainer_internal.h:410
Definition: elementdatacontainer_internal.h:46
void GetValuesState(std::string name, std::vector< double > &values) const
Definition: elementdatacontainer_internal.h:342
void GetGradsState(std::string name, std::vector< dealii::Tensor< 1, targetdim > > &values) const
Definition: elementdatacontainer_internal.h:388
void GetParamValues(std::string name, dealii::Vector< double > &value) const
Definition: elementdatacontainer_internal.h:326
Definition: fevalues_wrapper.h:43
const std::map< std::string, const VECTOR * > & GetDomainValues() const
Definition: elementdatacontainer_internal.h:71
Definition: dopeexception.h:35
void GetValuesControl(std::string name, std::vector< double > &values) const
Definition: elementdatacontainer_internal.h:369
void GetHessiansState(std::string name, std::vector< dealii::Tensor< 2, targetdim > > &values) const
Definition: elementdatacontainer_internal.h:446