DOpE
elementdatacontainer_internal.h
Go to the documentation of this file.
1 
24 #ifndef ELEMENTDATACONTAINER_INTERNAL_H_
25 #define ELEMENTDATACONTAINER_INTERNAL_H_
26 
27 #include <deal.II/lac/vector.h>
28 
29 #include "fevalues_wrapper.h"
30 #include "dopeexception.h"
31 
32 namespace DOpE
33 {
34  namespace edcinternal
35  {
45  template<typename VECTOR, int dim>
47  {
48  public:
50  const std::map<std::string, const dealii::Vector<double>*> &param_values
51  ,
52  const std::map<std::string, const VECTOR*> &domain_values);
53 
54  virtual
56  {
57  }
58  ;
59 
64  void
65  GetParamValues(std::string name, dealii::Vector<double>& value) const;
66 
70  const std::map<std::string, const VECTOR*> &
72  {
73  return domain_values_;
74  }
75  ;
76 
77  virtual const DOpEWrapper::FEValues<dim>&
78  GetFEValuesState() const = 0;
79 
80  virtual const DOpEWrapper::FEValues<dim>&
81  GetFEValuesControl() const = 0;
82 
83  /*********************************************************************/
87  const typename Triangulation<dim>::cell_iterator
88  GetElementState() const;
89 
90 
91  /********************************************************************/
99  void
100  GetValuesState(std::string name, std::vector<double>& values) const;
101 
102  /*********************************************/
103  /*
104  * Same as above for the Vector valued case.
105  */
106  void
107  GetValuesState(std::string name,
108  std::vector<dealii::Vector<double> >& values) const;
109 
110  /*********************************************/
111  /*
112  * Writes the values of the control variable at the quadrature points into values
113  */
114  void
115  GetValuesControl(std::string name, std::vector<double>& values) const;
116 
117  /*********************************************/
118  /*
119  * Same as above for the Vector valued case.
120  */
121  void
122  GetValuesControl(std::string name,
123  std::vector<dealii::Vector<double> >& values) const;
124  /*********************************************/
125  /*
126  * Writes the values of the state gradient at the quadrature points into values.
127  */
128 
129  template<int targetdim>
130  void
131  GetGradsState(std::string name,
132  std::vector<dealii::Tensor<1, targetdim> >& values) const;
133 
134  /*********************************************/
135  /*
136  * Same as above for the Vector valued case.
137  */
138  template<int targetdim>
139  void
141  std::string name,
142  std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values) const;
143 
144  /*********************************************/
145  /*
146  * Writes the values of the control gradient at the quadrature points into values.
147  */
148  template<int targetdim>
149  void
150  GetGradsControl(std::string name,
151  std::vector<dealii::Tensor<1, targetdim> >& values) const;
152 
153  /*********************************************/
154  /*
155  * Same as above for the Vector valued case.
156  */
157  template<int targetdim>
158  void
160  std::string name,
161  std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values) const;
162  /*********************************************/
163  /*
164  * Writes the values of the state hessian at the quadrature points into values.
165  */
166  template<int targetdim>
167  void
168  GetHessiansState(std::string name,
169  std::vector<dealii::Tensor<2, targetdim> >& values) const;
170 
171  /*********************************************/
172  /*
173  * Same as above for the Vector valued case.
174  */
175  template<int targetdim>
176  void
178  std::string name,
179  std::vector<std::vector<dealii::Tensor<2, targetdim> > >& values) const;
180 
181  /*********************************************/
182  /*
183  * Writes the values of the control hessian at the quadrature points into values.
184  */
185  template<int targetdim>
186  void
187  GetHessiansControl(std::string name,
188  std::vector<dealii::Tensor<2, targetdim> >& values) const;
189 
190  /*********************************************/
191  /*
192  * Same as above for the Vector valued case.
193  */
194  template<int targetdim>
195  void
197  std::string name,
198  std::vector<std::vector<dealii::Tensor<2, targetdim> > >& values) const;
199 
200  /*********************************************/
201  /*
202  * Writes the values of the state laplacian
203  * at the quadrature points into values.
204  */
205 
206  void
207  GetLaplaciansState(std::string name,
208  std::vector<double> & values) const;
209 
210  /*********************************************/
211  /*
212  * Same as above for the Vector valued case.
213  */
214 
215  void
216  GetLaplaciansState(std::string name,
217  std::vector<dealii::Vector<double> >& values) const;
218 
219  /*********************************************/
220  /*
221  * Writes the values of the control laplacian
222  * at the quadrature points into values.
223  */
224 
225  void
226  GetLaplaciansControl(std::string name,
227  std::vector<double> & values) const;
228 
229  /*********************************************/
230  /*
231  * Same as above for the Vector valued case.
232  */
233 
234  void
235  GetLaplaciansControl(std::string name,
236  std::vector<dealii::Vector<double> >& values) const;
237 
238  private:
239  /***********************************************************/
243  void
244  GetValues(const DOpEWrapper::FEValues<dim>& fe_values,
245  std::string name, std::vector<double>& values) const;
246  /***********************************************************/
250  void
251  GetValues(const DOpEWrapper::FEValues<dim>& fe_values,
252  std::string name,
253  std::vector<dealii::Vector<double> >& values) const;
254  /***********************************************************/
258  template<int targetdim>
259  void
260  GetGrads(const DOpEWrapper::FEValues<dim>& fe_values,
261  std::string name,
262  std::vector<dealii::Tensor<1, targetdim> >& values) const;
263  /***********************************************************/
267  template<int targetdim>
268  void
269  GetGrads(
270  const DOpEWrapper::FEValues<dim>& fe_values,
271  std::string name,
272  std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values) const;
273  /***********************************************************/
277  void
278  GetLaplacians(const DOpEWrapper::FEValues<dim>& fe_values,
279  std::string name, std::vector<double>& values) const;
280 
281  /***********************************************************/
285  void
286  GetLaplacians(const DOpEWrapper::FEValues<dim>& fe_values,
287  std::string name,
288  std::vector<dealii::Vector<double> >& values) const;
289 
290  /***********************************************************/
294  template<int targetdim>
295  void
296  GetHessians(const DOpEWrapper::FEValues<dim>& fe_values,
297  std::string name,
298  std::vector<dealii::Tensor<2, targetdim> > & values) const;
299 
300  /***********************************************************/
304  template<int targetdim>
305  void
306  GetHessians(
307  const DOpEWrapper::FEValues<dim>& fe_values,
308  std::string name,
309  std::vector<std::vector<dealii::Tensor<2, targetdim> > >& values) const;
310 
311  const std::map<std::string, const dealii::Vector<double>*> &param_values_;
312  const std::map<std::string, const VECTOR*> &domain_values_;
313  };
314 
315  /**********************************************************************/
316  template<typename VECTOR, int dim>
318  const std::map<std::string, const dealii::Vector<double>*> &param_values,
319  const std::map<std::string, const VECTOR*> &domain_values)
320  : param_values_(param_values), domain_values_(domain_values)
321  {
322  }
323 
324  template<typename VECTOR, int dim>
325  void
327  dealii::Vector<double>& value) const
328  {
329  typename std::map<std::string, const dealii::Vector<double>*>::const_iterator it =
330  param_values_.find(name);
331  if (it == param_values_.end())
332  {
333  throw DOpEException("Did not find " + name,
334  "ElementDataContainerInternal::GetParamValues");
335  }
336  value = *(it->second);
337  }
338 
339  /*********************************************/
340  template<typename VECTOR, int dim>
341  void
343  std::vector<double>& values) const
344  {
345  this->GetValues(this->GetFEValuesState(), name, values);
346  }
347  /*********************************************/
348  template<typename VECTOR, int dim>
349  void
351  std::vector<dealii::Vector<double> >& values) const
352  {
353  this->GetValues(this->GetFEValuesState(), name, values);
354 
355  }
356 
357 
358  /*********************************************/
359  template<typename VECTOR, int dim>
360  const typename Triangulation<dim>::cell_iterator
362  {
363  return this->GetFEValuesState().get_element();
364  }
365 
366  /*********************************************/
367  template<typename VECTOR, int dim>
368  void
370  std::vector<double>& values) const
371  {
372  this->GetValues(this->GetFEValuesControl(), name, values);
373  }
374 
375  /*********************************************/
376  template<typename VECTOR, int dim>
377  void
379  std::vector<dealii::Vector<double> >& values) const
380  {
381  this->GetValues(this->GetFEValuesControl(), name, values);
382  }
383 
384  /*********************************************/
385  template<typename VECTOR, int dim>
386  template<int targetdim>
387  void
389  std::vector<dealii::Tensor<1, targetdim> >& values) const
390  {
391  this->GetGrads<targetdim>(this->GetFEValuesState(), name, values);
392  }
393 
394  /*********************************************/
395  template<typename VECTOR, int dim>
396  template<int targetdim>
397  void
399  std::string name,
400  std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values) const
401  {
402  this->GetGrads<targetdim>(this->GetFEValuesState(), name, values);
403  }
404 
405  /***********************************************************************/
406 
407  template<typename VECTOR, int dim>
408  template<int targetdim>
409  void
411  std::string name,
412  std::vector<dealii::Tensor<1, targetdim> >& values) const
413  {
414  this->GetGrads<targetdim>(this->GetFEValuesControl(), name, values);
415  }
416 
417  /***********************************************************************/
418 
419  template<typename VECTOR, int dim>
420  template<int targetdim>
421  void
423  std::string name,
424  std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values) const
425  {
426  this->GetGrads<targetdim>(this->GetFEValuesControl(), name, values);
427  }
428 
429  /***********************************************************************/
430 
431  template<typename VECTOR, int dim>
432  template<int targetdim>
433  void
435  std::string name,
436  std::vector<std::vector<dealii::Tensor<2, targetdim> > >& values) const
437  {
438  this->GetHessians<targetdim>(this->GetFEValuesState(), name, values);
439  }
440 
441  /***********************************************************************/
442 
443  template<typename VECTOR, int dim>
444  template<int targetdim>
445  void
447  std::string name,
448  std::vector<dealii::Tensor<2, targetdim> >& values) const
449  {
450  this->GetHessians<targetdim>(this->GetFEValuesState(), name, values);
451  }
452 
453  /***********************************************************************/
454 
455  template<typename VECTOR, int dim>
456  template<int targetdim>
457  void
459  std::string name,
460  std::vector<std::vector<dealii::Tensor<2, targetdim> > >& values) const
461  {
462  this->GetHessians<targetdim>(this->GetFEValuesControl(), name, values);
463  }
464 
465  /***********************************************************************/
466 
467  template<typename VECTOR, int dim>
468  template<int targetdim>
469  void
471  std::string name,
472  std::vector<dealii::Tensor<2, targetdim> >& values) const
473  {
474  this->GetHessians<targetdim>(this->GetFEValuesControl(), name,
475  values);
476  }
477 
478  /***********************************************************************/
479  template<typename VECTOR, int dim>
480  void
482  std::string name, std::vector<double> & values) const
483  {
484  this->GetLaplacians(this->GetFEValuesState(), name, values);
485  }
486 
487  /***********************************************************************/
488  template<typename VECTOR, int dim>
489  void
491  std::string name, std::vector<dealii::Vector<double> >& values) const
492  {
493  this->GetLaplacians(this->GetFEValuesState(), name, values);
494  }
495  /***********************************************************************/
496  template<typename VECTOR, int dim>
497  void
499  std::string name, std::vector<double> & values) const
500  {
501  this->GetLaplacians(this->GetFEValuesControl(), name, values);
502  }
503 
504  /***********************************************************************/
505  template<typename VECTOR, int dim>
506  void
508  std::string name, std::vector<dealii::Vector<double> >& values) const
509  {
510  this->GetLaplacians(this->GetFEValuesControl(), name, values);
511  }
512  /***********************************************************************/
513  template<typename VECTOR, int dim>
514  void
516  const DOpEWrapper::FEValues<dim>& fe_values, std::string name,
517  std::vector<double>& values) const
518  {
519  typename std::map<std::string, const VECTOR*>::const_iterator it =
520  this->GetDomainValues().find(name);
521  if (it == this->GetDomainValues().end())
522  {
523  throw DOpEException("Did not find " + name,
524  "ElementDataContainer::GetValues");
525  }
526  fe_values.get_function_values(*(it->second), values);
527  }
528 
529  /***********************************************************************/
530  template<typename VECTOR, int dim>
531  void
532  ElementDataContainerInternal<VECTOR, dim>::GetValues(
533  const DOpEWrapper::FEValues<dim>& fe_values, std::string name,
534  std::vector<dealii::Vector<double> >& values) const
535  {
536  typename std::map<std::string, const VECTOR*>::const_iterator it =
537  this->GetDomainValues().find(name);
538  if (it == this->GetDomainValues().end())
539  {
540  throw DOpEException("Did not find " + name,
541  "ElementDataContainer::GetValues");
542  }
543  fe_values.get_function_values(*(it->second), values);
544  }
545 
546  /***********************************************************************/
547 
548  template<typename VECTOR, int dim>
549  template<int targetdim>
550  void
551  ElementDataContainerInternal<VECTOR, dim>::GetGrads(
552  const DOpEWrapper::FEValues<dim>& fe_values, std::string name,
553  std::vector<dealii::Tensor<1, targetdim> >& values) const
554  {
555  typename std::map<std::string, const VECTOR*>::const_iterator it =
556  this->GetDomainValues().find(name);
557  if (it == this->GetDomainValues().end())
558  {
559  throw DOpEException("Did not find " + name,
560  "ElementDataContainerInternal::GetGrads");
561  }
562  fe_values.get_function_gradients(*(it->second), values);
563  }
564 
565  /***********************************************************************/
566 
567  template<typename VECTOR, int dim>
568  template<int targetdim>
569  void
570  ElementDataContainerInternal<VECTOR, dim>::GetGrads(
571  const DOpEWrapper::FEValues<dim>& fe_values,
572  std::string name,
573  std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values) const
574  {
575  typename std::map<std::string, const VECTOR*>::const_iterator it =
576  this->GetDomainValues().find(name);
577  if (it == this->GetDomainValues().end())
578  {
579  throw DOpEException("Did not find " + name,
580  "ElementDataContainerInternal::GetGrads");
581  }
582  fe_values.get_function_gradients(*(it->second), values);
583  }
584 
585  /***********************************************************************/
586 
587  template<typename VECTOR, int dim>
588  void
589  ElementDataContainerInternal<VECTOR, dim>::GetLaplacians(
590  const DOpEWrapper::FEValues<dim>& fe_values, std::string name,
591  std::vector<double> & values) const
592  {
593  typename std::map<std::string, const VECTOR*>::const_iterator it =
594  this->GetDomainValues().find(name);
595  if (it == this->GetDomainValues().end())
596  {
597  throw DOpEException("Did not find " + name,
598  "ElementDataContainerInternal::GetLaplacians");
599  }
600  fe_values.get_function_laplacians(*(it->second), values);
601  }
602 
603  /***********************************************************************/
604 
605  template<typename VECTOR, int dim>
606  void
607  ElementDataContainerInternal<VECTOR, dim>::GetLaplacians(
608  const DOpEWrapper::FEValues<dim>& fe_values, std::string name,
609  std::vector<dealii::Vector<double> >& values) const
610  {
611  typename std::map<std::string, const VECTOR*>::const_iterator it =
612  this->GetDomainValues().find(name);
613  if (it == this->GetDomainValues().end())
614  {
615  throw DOpEException("Did not find " + name,
616  "ElementDataContainerInternal::GetLaplacians");
617  }
618  fe_values.get_function_laplacians(*(it->second), values);
619  }
620 
621  /***********************************************************************/
622 
623  template<typename VECTOR, int dim>
624  template<int targetdim>
625  void
626  ElementDataContainerInternal<VECTOR, dim>::GetHessians(
627  const DOpEWrapper::FEValues<dim>& fe_values,
628  std::string name,
629  std::vector<std::vector<dealii::Tensor<2, targetdim> > >& values) const
630  {
631  typename std::map<std::string, const VECTOR*>::const_iterator it =
632  this->GetDomainValues().find(name);
633  if (it == this->GetDomainValues().end())
634  {
635  throw DOpEException("Did not find " + name,
636  "ElementDataContainerInternal::GetGrads");
637  }
638  fe_values.get_function_hessians(*(it->second), values);
639  }
640 
641  /***********************************************************************/
642 
643  template<typename VECTOR, int dim>
644  template<int targetdim>
645  void
646  ElementDataContainerInternal<VECTOR, dim>::GetHessians(
647  const DOpEWrapper::FEValues<dim>& fe_values, std::string name,
648  std::vector<dealii::Tensor<2, targetdim> >& values) const
649  {
650  typename std::map<std::string, const VECTOR*>::const_iterator it =
651  this->GetDomainValues().find(name);
652  if (it == this->GetDomainValues().end())
653  {
654  throw DOpEException("Did not find " + name,
655  "ElementDataContainerInternal::GetGrads");
656  }
657  fe_values.get_function_hessians(*(it->second), values);
658  }
659 
660  } //end of namespace edcinternal
661 }
662 
663 #endif /* ElementDataContainer_INTERNAL_H_ */
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 > * > &param_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