DOpE
facedatacontainer_internal.h
Go to the documentation of this file.
1 
24 #ifndef FACEDATACONTAINER_INTERNAL_H_
25 #define FACEDATACONTAINER_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 fdcinternal
35  {
45  template<typename VECTOR, int dim>
47  {
48  public:
50  const std::map<std::string, const dealii::Vector<double>*> &param_values,
51  const std::map<std::string, const VECTOR*> &domain_values,
52  bool need_neighbour);
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 dealii::FEFaceValuesBase<dim>&
78  GetFEFaceValuesState() const =0;
79  virtual const dealii::FEFaceValuesBase<dim>&
80  GetFEFaceValuesControl() const = 0;
81 
82  virtual const dealii::FEFaceValuesBase<dim>&
83  GetNbrFEFaceValuesState() const = 0;
84  virtual const dealii::FEFaceValuesBase<dim>&
85  GetNbrFEFaceValuesControl() const = 0;
86 
87  /*********************************************************************/
91  const typename Triangulation<dim>::cell_iterator
92  GetElementState() const;
93 
94  /********************************************************************/
102  void
103  GetFaceValuesState(std::string name,
104  std::vector<double>& values) const;
105 
106  /*********************************************/
107  /*
108  * Same as above for the Vector valued case.
109  */
110  void
111  GetFaceValuesState(std::string name,
112  std::vector<dealii::Vector<double> >& values) const;
113 
114  /*********************************************/
115  /*
116  * Writes the values of the control variable at the quadrature points into values
117  */
118  void
119  GetFaceValuesControl(std::string name,
120  std::vector<double>& values) const;
121 
122  /*********************************************/
123  /*
124  * Same as above for the Vector valued case.
125  */
126  void
127  GetFaceValuesControl(std::string name,
128  std::vector<dealii::Vector<double> >& values) const;
129  /*********************************************/
130  /*
131  * Writes the values of the state gradient at the quadrature points into values.
132  */
133 
134  template<int targetdim>
135  void
136  GetFaceGradsState(std::string name,
137  std::vector<dealii::Tensor<1, targetdim> >& values) const;
138 
139  /*********************************************/
140  /*
141  * Same as above for the Vector valued case.
142  */
143  template<int targetdim>
144  void
145  GetFaceGradsState(std::string name,
146  std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values) const;
147 
148  /*********************************************/
149  /*
150  * Writes the values of the control gradient at the quadrature points into values.
151  */
152  template<int targetdim>
153  void
154  GetFaceGradsControl(std::string name,
155  std::vector<dealii::Tensor<1, targetdim> >& values) const;
156 
157  /*********************************************/
158  /*
159  * Same as above for the Vector valued case.
160  */
161  template<int targetdim>
162  void
163  GetFaceGradsControl(std::string name,
164  std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values) const;
165 
166  /*
167  * Writes the values of the state variable at the quadrature points into values.
168  */
169  inline void
170  GetNbrFaceValuesState(std::string name,
171  std::vector<double>& values) const;
172  /*********************************************/
173  /*
174  * Same as above for the Vector valued case.
175  */
176  inline void
177  GetNbrFaceValuesState(std::string name,
178  std::vector<Vector<double> >& values) const;
179 
180  /*********************************************/
181 
182  /*
183  * Writes the values of the control variable at the quadrature points into values
184  */
185  inline void
186  GetNbrFaceValuesControl(std::string name,
187  std::vector<double>& values) const;
188  /*********************************************/
189 
190  /*
191  * Same as above for the Vector valued case.
192  */
193  inline void
194  GetNbrFaceValuesControl(std::string name,
195  std::vector<Vector<double> >& values) const;
196  /*********************************************/
197 
198  /*
199  * Writes the values of the state gradient at the quadrature points into values.
200  */
201 
202  template<int targetdim>
203  inline void
204  GetNbrFaceGradsState(std::string name,
205  std::vector<dealii::Tensor<1, targetdim> >& values) const;
206 
207  /*********************************************/
208 
209  /*
210  * Same as avoe for the Vector valued case.
211  */
212  template<int targetdim>
213  inline void
214  GetNbrFaceGradsState(std::string name,
215  std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values) const;
216 
217  /*********************************************/
218 
219  /*
220  * Writes the values of the control gradient at the quadrature points into values.
221  */
222 
223  template<int targetdim>
224  inline void
225  GetNbrFaceGradsControl(std::string name,
226  std::vector<dealii::Tensor<1, targetdim> >& values) const;
227 
228  /*********************************************/
229  /*
230  * Same as above for the Vector valued case.
231  */
232  template<int targetdim>
233  inline void
234  GetNbrFaceGradsControl(std::string name,
235  std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values) const;
236 
237  protected:
238  void
239  SetFace(unsigned int face)
240  {
241  face_ = face;
242  }
243  unsigned int
244  GetFace() const
245  {
246  return face_;
247  }
248  void
249  SetSubFace(unsigned int subface)
250  {
251  subface_ = subface;
252  }
253  unsigned int
254  GetSubFace() const
255  {
256  return subface_;
257  }
258  bool
260  {
261  return need_neighbour_;
262  }
263 
264  private:
265  /***********************************************************/
269  void
270  GetValues(const dealii::FEFaceValuesBase<dim>& fe_values,
271  std::string name, std::vector<double>& values) const;
272  /***********************************************************/
276  void
277  GetValues(const dealii::FEFaceValuesBase<dim>& fe_values,
278  std::string name,
279  std::vector<dealii::Vector<double> >& values) const;
280  /***********************************************************/
284  template<int targetdim>
285  void
286  GetGrads(const dealii::FEFaceValuesBase<dim>& fe_values,
287  std::string name,
288  std::vector<dealii::Tensor<1, targetdim> >& values) const;
289  /***********************************************************/
293  template<int targetdim>
294  void
295  GetGrads(const dealii::FEFaceValuesBase<dim>& fe_values,
296  std::string name,
297  std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values) const;
298 
299  const std::map<std::string, const dealii::Vector<double>*> & param_values_;
300  const std::map<std::string, const VECTOR*> & domain_values_;
301 
302  unsigned int face_;
303  unsigned int subface_;
304  bool need_neighbour_;
305  };
306 
307  /**********************************************************************/
308  template<typename VECTOR, int dim>
310  const std::map<std::string, const dealii::Vector<double>*> &param_values,
311  const std::map<std::string, const VECTOR*> &domain_values,
312  bool need_neighbour)
313  : param_values_(param_values), domain_values_(domain_values), need_neighbour_(
314  need_neighbour)
315  {
316  }
317 
318  template<typename VECTOR, int dim>
319  void
321  dealii::Vector<double>& value) const
322  {
323  typename std::map<std::string, const dealii::Vector<double>*>::const_iterator it =
324  param_values_.find(name);
325  if (it == param_values_.end())
326  {
327  throw DOpEException("Did not find " + name,
328  "FaceDataContainerInternal::GetParamValues");
329  }
330  value = *(it->second);
331  }
332 
333  /*********************************************/
334  template<typename VECTOR, int dim>
335  const typename Triangulation<dim>::cell_iterator
337  {
338  return this->GetFEFaceValuesState().get_element();
339  }
340 
341  /*********************************************/
342  template<typename VECTOR, int dim>
343  void
345  std::string name, std::vector<double>& values) const
346  {
347  this->GetValues(this->GetFEFaceValuesState(), name, values);
348  }
349  /*********************************************/
350  template<typename VECTOR, int dim>
351  void
353  std::string name, std::vector<dealii::Vector<double> >& values) const
354  {
355  this->GetValues(this->GetFEFaceValuesState(), name, values);
356 
357  }
358 
359  /*********************************************/
360  template<typename VECTOR, int dim>
361  void
363  std::string name, std::vector<double>& values) const
364  {
365  this->GetValues(this->GetFEFaceValuesControl(), name, values);
366  }
367 
368  /*********************************************/
369  template<typename VECTOR, int dim>
370  void
372  std::string name, std::vector<dealii::Vector<double> >& values) const
373  {
374  this->GetValues(this->GetFEFaceValuesControl(), name, values);
375  }
376 
377  /*********************************************/
378  template<typename VECTOR, int dim>
379  template<int targetdim>
380  void
382  std::string name,
383  std::vector<dealii::Tensor<1, targetdim> >& values) const
384  {
385  this->GetGrads<targetdim>(this->GetFEFaceValuesState(), name, values);
386  }
387 
388  /*********************************************/
389  template<typename VECTOR, int dim>
390  template<int targetdim>
391  void
393  std::string name,
394  std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values) const
395  {
396  this->GetGrads<targetdim>(this->GetFEFaceValuesState(), name, values);
397  }
398 
399  /***********************************************************************/
400 
401  template<typename VECTOR, int dim>
402  template<int targetdim>
403  void
405  std::string name,
406  std::vector<dealii::Tensor<1, targetdim> >& values) const
407  {
408  this->GetGrads<targetdim>(this->GetFEFaceValuesControl(), name,
409  values);
410  }
411  /***********************************************************************/
412 
413  template<typename VECTOR, int dim>
414  template<int targetdim>
415  void
417  std::string name,
418  std::vector<std::vector<dealii::Tensor<1, targetdim> > >& values) const
419  {
420  this->GetGrads<targetdim>(this->GetFEFaceValuesControl(), name,
421  values);
422  }
423 
424  /*********************************************/
425  template<typename VECTOR, int dim>
426  void
428  std::string name, std::vector<double>& values) const
429  {
430  this->GetValues(this->GetNbrFEFaceValuesState(), name, values);
431  }
432  /*********************************************/
433  template<typename VECTOR, int dim>
434  void
436  std::string name, std::vector<Vector<double> >& values) const
437  {
438  this->GetValues(this->GetNbrFEFaceValuesState(), name, values);
439 
440  }
441 
442  /*********************************************/
443  template<typename VECTOR, int dim>
444  void
446  std::string name, std::vector<double>& values) const
447  {
448  this->GetValues(this->GetNbrFEFaceValuesControl(), name, values);
449  }
450 
451  /*********************************************/
452  template<typename VECTOR, int dim>
453  void
455  std::string name, std::vector<Vector<double> >& values) const
456  {
457  this->GetValues(this->GetNbrFEFaceValuesControl(), name, values);
458  }
459 
460  /*********************************************/
461  template<typename VECTOR, int dim>
462  template<int targetdim>
463  void
465  std::string name, std::vector<Tensor<1, targetdim> >& values) const
466  {
467  this->GetGrads<targetdim>(this->GetNbrFEFaceValuesState(), name,
468  values);
469  }
470 
471  /*********************************************/
472  template<typename VECTOR, int dim>
473  template<int targetdim>
474  void
476  std::string name,
477  std::vector<std::vector<Tensor<1, targetdim> > >& values) const
478  {
479  this->GetGrads<targetdim>(this->GetNbrFEFaceValuesState(), name,
480  values);
481  }
482 
483  /***********************************************************************/
484 
485  template<typename VECTOR, int dim>
486  template<int targetdim>
487  void
489  std::string name, std::vector<Tensor<1, targetdim> >& values) const
490  {
491  this->GetGrads<targetdim>(this->GetNbrFEFaceValuesControl(), name,
492  values);
493  }
494  /***********************************************************************/
495 
496  template<typename VECTOR, int dim>
497  template<int targetdim>
498  void
500  std::string name,
501  std::vector<std::vector<Tensor<1, targetdim> > >& values) const
502  {
503  this->GetGrads<targetdim>(this->GetNbrFEFaceValuesControl(), name,
504  values);
505  }
506 
507  /***********************************************************************/
508 
509  /***********************************************************************/
510  template<typename VECTOR, int dim>
511  void
513  const dealii::FEFaceValuesBase<dim>& fe_values, std::string name,
514  std::vector<double>& values) const
515  {
516  typename std::map<std::string, const VECTOR*>::const_iterator it =
517  this->GetDomainValues().find(name);
518  if (it == this->GetDomainValues().end())
519  {
520  throw DOpEException("Did not find " + name,
521  "ElementDataContainer::GetValues");
522  }
523  fe_values.get_function_values(*(it->second), values);
524  }
525 
526  /***********************************************************************/
527  template<typename VECTOR, int dim>
528  void
529  FaceDataContainerInternal<VECTOR, dim>::GetValues(
530  const dealii::FEFaceValuesBase<dim>& fe_values, std::string name,
531  std::vector<dealii::Vector<double> >& values) const
532  {
533  typename std::map<std::string, const VECTOR*>::const_iterator it =
534  this->GetDomainValues().find(name);
535  if (it == this->GetDomainValues().end())
536  {
537  throw DOpEException("Did not find " + name,
538  "ElementDataContainer::GetValues");
539  }
540  fe_values.get_function_values(*(it->second), values);
541  }
542 
543  /***********************************************************************/
544 
545  template<typename VECTOR, int dim>
546  template<int targetdim>
547  void
548  FaceDataContainerInternal<VECTOR, dim>::GetGrads(
549  const dealii::FEFaceValuesBase<dim>& fe_values, std::string name,
550  std::vector<dealii::Tensor<1, targetdim> >& values) const
551  {
552  typename std::map<std::string, const VECTOR*>::const_iterator it =
553  this->GetDomainValues().find(name);
554  if (it == this->GetDomainValues().end())
555  {
556  throw DOpEException("Did not find " + name,
557  "FaceDataContainerInternal::GetGrads");
558  }
559  fe_values.get_function_gradients(*(it->second), values);
560  }
561 
562  /***********************************************************************/
563 
564  template<typename VECTOR, int dim>
565  template<int targetdim>
566  void
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
570  {
571  typename std::map<std::string, const VECTOR*>::const_iterator it =
572  this->GetDomainValues().find(name);
573  if (it == this->GetDomainValues().end())
574  {
575  throw DOpEException("Did not find " + name,
576  "FaceDataContainerInternal::GetGrads");
577  }
578  fe_values.get_function_gradients(*(it->second), values);
579  }
580 
581  /***********************************************************************/
582  }
583 }
584 
585 #endif /* FACEDATACONTAINER_INTERNAL_H_ */
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 > * > &param_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