OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
mesh_loop_info_objects.h
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 //
3 // FCST: Fuel Cell Simulation Toolbox
4 //
5 // Copyright (C) 2006-2009 by Guido Kanschat
6 // Copyright (C) 2006-2014 by Energy Systems Design Laboratory, University of Alberta
7 //
8 // This software is distributed under the MIT License
9 // For more information, see the README file in /doc/LICENSE
10 //
11 // - Class: mesh_loop_info_objects.h
12 // - Description: Objects used for looping over mesh
13 // - Developers: Guido Kanschat, Texas A&M University
14 // Valentin N. Zingan, University of Alberta
15 // - Id: $Id: mesh_loop_info_objects.h 2605 2014-08-15 03:36:44Z secanell $
16 //
17 // ----------------------------------------------------------------------------
18 
19 #ifndef _FUEL_CELL_APPLICATION_CORE_MESH_LOOP_INFO_OBJECTS_H_
20 #define _FUEL_CELL_APPLICATION_CORE_MESH_LOOP_INFO_OBJECTS_H_
21 
22 #include <base/geometry_info.h>
23 #include <base/quadrature_lib.h>
24 #include <base/mg_level_object.h>
25 
26 #include <lac/block_indices.h>
27 
28 #include <fe/fe.h>
29 #include <fe/fe_tools.h>
30 #include <fe/fe_values.h>
31 
32 #include <multigrid/mg_tools.h>
33 
34 #include "fe_vectors.h"
35 
36 #include <boost/shared_ptr.hpp>
37 
38 namespace dealii
39 {
40  template<int,int> class DoFHandler;
41  template<int,int> class MGDoFHandler;
42 }
43 
44 using namespace dealii;
45 
64 namespace FuelCell
65 {
66 namespace ApplicationCore
67 {
68 
102 // +++ IMPLEMENTATION +++++++++++++++++++++++++++++++++++++++++++++++++++++++++
103 
104 template<int dim, typename TYPE>
105 inline
106 void
107 fill_data(const FEValuesBase<dim>& fe_values,
108  const FEVector& fe_vector,
109  const std::vector<unsigned int>& local_dof_indices,
110  unsigned int first_index,
111  unsigned int n_indices,
112  std::vector< std::vector<TYPE> >& result)
113 {
114  Assert(false, ExcNotImplemented());
115 }
116 
117 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
118 
119 template<int dim>
120 inline
121 void
122 fill_data(const FEValuesBase<dim>& fe_values,
123  const FEVector& fe_vector,
124  const std::vector<unsigned int>& local_dof_indices,
125  unsigned int first_index,
126  unsigned int n_indices,
127  std::vector< std::vector<double> >& result)
128 {
129  fe_values.get_function_values(fe_vector,
130  make_slice(local_dof_indices, first_index, n_indices),
131  result,
132  true);
133 }
134 
135 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
136 
137 template<int dim>
138 inline
139 void
140 fill_data(const FEValuesBase<dim>& fe_values,
141  const FEVector& fe_vector,
142  const std::vector<unsigned int>& local_dof_indices,
143  unsigned int first_index,
144  unsigned int n_indices,
145  std::vector< std::vector< Tensor<1,dim> > >& result)
146 {
147  fe_values.get_function_grads(fe_vector,
148  make_slice(local_dof_indices, first_index, n_indices),
149  result,
150  true);
151 }
152 
153 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
154 
155 template<int dim>
156 inline
157 void
158 fill_data(const FEValuesBase<dim>& fe_values,
159  const FEVector& fe_vector,
160  const std::vector<unsigned int>& local_dof_indices,
161  unsigned int first_index,
162  unsigned int n_indices,
163  std::vector< std::vector< Tensor<2,dim> > >& result)
164 {
165  fe_values.get_function_hessians(fe_vector,
166  make_slice(local_dof_indices, first_index, n_indices),
167  result,
168  true);
169 }
170 
171 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
172 
182 struct BlockInfo : public Subscriptor
183 {
187  BlockIndices global;
188 
192  BlockIndices local;
193 
197  std::vector<BlockIndices> levels;
198 
204  std::vector<unsigned int> base_element;
205 
217  std::vector<unsigned int> local_renumbering;
218 
222  template<int dim, int spacedim>
223  void initialize(const DoFHandler<dim, spacedim>& dof_handler);
224 
228  template<int dim, int spacedim>
229  void initialize_local(const DoFHandler<dim, spacedim>& dof_handler);
230 
234  template<int dim, int spacedim>
235  void initialize(const MGDoFHandler<dim, spacedim>& mg_dof_handler);
236 };
237 
238 // +++ IMPLEMENTATION +++++++++++++++++++++++++++++++++++++++++++++++++++++++++
239 
240 template<int dim, int spacedim>
241 inline
242 void
243 BlockInfo::initialize(const DoFHandler<dim, spacedim>& dof_handler)
244 {
245  const FiniteElement<dim, spacedim>& fe = dof_handler.get_fe();
246  std::vector<unsigned int> sizes(fe.n_blocks());
247  DoFTools::count_dofs_per_block(dof_handler, sizes);
248  global.reinit(sizes);
249 }
250 
251 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
252 
253 template<int dim, int spacedim>
254 inline
255 void
256 BlockInfo::initialize_local(const DoFHandler<dim, spacedim>& dof_handler)
257 {
258  const FiniteElement<dim, spacedim>& fe = dof_handler.get_fe();
259  std::vector<unsigned int> sizes(fe.n_blocks());
260 
261  base_element.resize(fe.n_blocks());
262 
263  for(unsigned int i = 0; i < base_element.size(); ++i)
264  base_element[i] = fe.block_to_base_index(i).first;
265 
266  local_renumbering.resize(fe.n_dofs_per_cell());
267  FETools::compute_block_renumbering(fe,
268  local_renumbering,
269  sizes,
270  false);
271  local.reinit(sizes);
272 }
273 
274 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
275 
276 template<int dim, int spacedim>
277 inline
278 void
279 BlockInfo::initialize(const MGDoFHandler<dim, spacedim>& mg_dof_handler)
280 {
281  initialize((DoFHandler<dim, spacedim>&) mg_dof_handler);
282 
283  std::vector< std::vector<unsigned int> > sizes(mg_dof_handler.get_tria().n_levels(),
284  std::vector<unsigned int>(mg_dof_handler.get_fe().n_blocks())
285  );
286 
287  MGTools::count_dofs_per_block(mg_dof_handler, sizes);
288 
289  levels.resize(sizes.size());
290  for(unsigned int i = 0; i < sizes.size(); ++i)
291  levels[i].reinit(sizes[i]);
292 }
293 
294 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
295 
320 template<int dim, int spacedim = dim>
321 class DoFInfo
322 {
323 public:
324 
326 
327 
331  DoFInfo(const BlockInfo& block_info);
332 
337  DoFInfo(const FEVectors&,
338  const BlockInfo& block_info);
339 
341 
343 
344 
349  template<typename DHCellIterator>
350  void reinit(const DHCellIterator& c);
351 
356  template<typename DHCellIterator, typename DHFaceIterator>
357  void reinit(const DHCellIterator& c,
358  const DHFaceIterator& f,
359  const unsigned int fn);
360 
365  template<typename DHCellIterator, typename DHFaceIterator>
366  void reinit(const DHCellIterator& c,
367  const DHFaceIterator& f,
368  const unsigned int fn,
369  const unsigned int sn);
370 
372 
374 
375 
380 
385 
390 
395 
400 
410  unsigned int face_number;
411 
422  unsigned int sub_number;
423 
425 
427 
428 
433  std::vector<unsigned int> indices;
434 
438  SmartPointer<const BlockInfo> block_info;
439 
441 
442 private:
443 
447  void get_indices(const typename DoFHandler<dim, spacedim>::cell_iterator c);
448 
452  void get_indices(const typename MGDoFHandler<dim, spacedim>::cell_iterator c);
453 
457  std::vector<unsigned int> indices_org;
458 };
459 
460 // +++ IMPLEMENTATION +++++++++++++++++++++++++++++++++++++++++++++++++++++++++
461 
462 template<int dim, int spacedim>
463 inline
465 :
466 block_info(&block_info,
467  typeid(*this).name())
468 { }
469 
470 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
471 
472 template<int dim, int spacedim>
473 inline
475  const BlockInfo& block_info)
476 :
477 block_info(&block_info,
478  typeid(*this).name())
479 { }
480 
481 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
482 
483 template<int dim, int spacedim>
484 template<typename DHCellIterator>
485 inline
486 void
487 DoFInfo<dim, spacedim>::reinit(const DHCellIterator& c)
488 {
489  get_indices(c);
490 
491  dof_cell = c;
492  dof_active_cell = c;
493  cell = static_cast<typename Triangulation<dim>::cell_iterator> (c);
494 
495  face_number = deal_II_numbers::invalid_unsigned_int;
496  sub_number = deal_II_numbers::invalid_unsigned_int;
497 }
498 
499 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
500 
501 template<int dim, int spacedim>
502 template<typename DHCellIterator, typename DHFaceIterator>
503 inline
504 void
505 DoFInfo<dim, spacedim>::reinit(const DHCellIterator& c,
506  const DHFaceIterator& f,
507  const unsigned int fn)
508 {
509  if( cell.state() != IteratorState::valid || cell != static_cast<typename Triangulation<dim>::cell_iterator> (c) )
510  {
511  get_indices(c);
512  }
513 
514  dof_cell = c;
515  dof_active_cell = c;
516  cell = static_cast<typename Triangulation<dim>::cell_iterator> (c);
517 
518  dof_face = f;
519  face = static_cast<typename Triangulation<dim>::face_iterator> (f);
520 
521  face_number = fn;
522  sub_number = deal_II_numbers::invalid_unsigned_int;
523 }
524 
525 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
526 
527 template<int dim, int spacedim>
528 template<typename DHCellIterator, typename DHFaceIterator>
529 inline
530 void
531 DoFInfo<dim, spacedim>::reinit(const DHCellIterator& c,
532  const DHFaceIterator& f,
533  const unsigned int fn,
534  const unsigned int sn)
535 {
536  if( cell.state() != IteratorState::valid || cell != static_cast<typename Triangulation<dim>::cell_iterator> (c) )
537  {
538  get_indices(c);
539  }
540 
541  dof_cell = c;
542  dof_active_cell = c;
543  cell = static_cast<typename Triangulation<dim>::cell_iterator> (c);
544 
545  dof_face = f;
546  face = static_cast<typename Triangulation<dim>::face_iterator> (f);
547 
548  face_number = fn;
549  sub_number = sn;
550 }
551 
552 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
553 
554 template<int dim, int spacedim>
555 inline
556 void
558 {
559  indices.resize(c->get_fe().dofs_per_cell);
560 
561  if(block_info->local_renumbering.size() == 0)
562  {
563  c->get_dof_indices(indices);
564  }
565  else
566  {
567  indices_org.resize(c->get_fe().dofs_per_cell);
568  c->get_dof_indices(indices_org);
569 
570  for(unsigned int i = 0; i < indices.size(); ++i)
571  indices[this->block_info->local_renumbering[i]] = indices_org[i];
572  }
573 }
574 
575 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
576 
577 template<int dim, int spacedim>
578 inline
579 void
581 {
582  indices.resize(c->get_fe().dofs_per_cell);
583 
584  if(block_info->local_renumbering.size() == 0)
585  {
586  c->get_mg_dof_indices(indices);
587  }
588  else
589  {
590  indices_org.resize(c->get_fe().dofs_per_cell);
591  c->get_mg_dof_indices(indices_org);
592 
593  for(unsigned int i = 0; i < indices.size(); ++i)
594  indices[this->block_info->local_renumbering[i]] = indices_org[i];
595  }
596 }
597 
598 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
599 
624 template<int dim, typename FEVALUESBASE>
625 class IntegrationInfo : public DoFInfo<dim>
626 {
627 public:
628 
630 
631 
636 
640  IntegrationInfo(const FEVectors& data,
641  const BlockInfo& block_info);
642 
685  template<typename FEVALUES>
686  void initialize(const FEVALUES* FE_VALUES,
687  const FiniteElement<dim>& fe,
688  const Mapping<dim>& mapping,
689  const Quadrature<FEVALUES::integral_dimension>& quadrature,
690  const UpdateFlags flags);
691 
695  void initialize_data(const FEVectors& data);
696 
698 
700 
701 
706  template<typename DHCellIterator>
707  void reinit(const DHCellIterator& c);
708 
713  template<typename DHCellIterator, typename DHFaceIterator>
714  void reinit(const DHCellIterator& c,
715  const DHFaceIterator& f,
716  const unsigned int fn);
717 
722  template<typename DHCellIterator, typename DHFaceIterator>
723  void reinit(const DHCellIterator& c,
724  const DHFaceIterator& f,
725  const unsigned int fn,
726  const unsigned int sn);
727 
736  template<typename TYPE>
737  void fill_local_data(std::vector< std::vector< std::vector<TYPE> > >& data,
738  bool split_fevalues) const;
739 
741 
743 
744 
748  const FEVALUESBASE& fe() const;
749 
753  const FEVALUESBASE& fe(unsigned int i) const;
754 
758  const FEVALUESBASE& get_fe_val_unsplit() const;
759 
763  void clear();
764 
766 
768 
769 
773  bool multigrid;
774 
780  SmartPointer<const FEVectors> global_data;
781 
790  std::vector< std::vector< std::vector<double> > > values;
791 
800  std::vector< std::vector< std::vector< Tensor<1,dim> > > > gradients;
801 
805  std::vector< std::vector< std::vector< Tensor<1,dim> > > >& derivatives;
806 
815  std::vector< std::vector< std::vector< Tensor<2,dim> > > > hessians;
816 
818 
819 private:
820 
822 
823 
827  std::vector< boost::shared_ptr<FEVALUESBASE> > fevalv;
828 
835  boost::shared_ptr<FEVALUESBASE> fe_val_unsplit;
836 
838 };
839 
840 // +++ IMPLEMENTATION +++++++++++++++++++++++++++++++++++++++++++++++++++++++++
841 
842 template<int dim, typename FEVALUESBASE>
843 inline
845 :
846 DoFInfo<dim>(block_info),
847 multigrid(false),
848 global_data(0, typeid(*this).name()),
849 derivatives(gradients),
850 fevalv(0)
851 { }
852 
853 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
854 
855 template<int dim, typename FEVALUESBASE>
856 inline
858  const BlockInfo& block_info)
859 :
860 DoFInfo<dim>(block_info),
861 multigrid(false),
862 global_data(&data),
863 derivatives(gradients),
864 fevalv(0)
865 { }
866 
867 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
868 
869 template<int dim, typename FEVALUESBASE>
870 template<typename FEVALUES>
871 inline
872 void
874  const FiniteElement<dim>& fe,
875  const Mapping<dim>& mapping,
876  const Quadrature<FEVALUES::integral_dimension>& quadrature,
877  const UpdateFlags flags)
878 {
879  if( this->block_info->local_renumbering.size() != 0 ) // MODE1
880  {
881  fevalv.resize(fe.n_base_elements());
882  for(unsigned int i = 0; i < fevalv.size(); ++i)
883  {
884  fevalv[i] = boost::shared_ptr<FEVALUESBASE>( new FEVALUES(mapping,
885  fe.base_element(i),
886  quadrature,
887  flags)
888  );
889  }
890 
891  // --- fe_val_unsplit initialization ---
892 
893  fe_val_unsplit = boost::shared_ptr<FEVALUESBASE>( new FEVALUES(mapping,
894  fe,
895  quadrature,
896  flags)
897  );
898  }
899  else // MODE2
900  {
901  fevalv.resize(1);
902  fevalv[0] = boost::shared_ptr<FEVALUESBASE>( new FEVALUES(mapping,
903  fe,
904  quadrature,
905  flags)
906  );
907  }
908 
909  values.resize(global_data->n_vectors(),
910  std::vector< std::vector<double> >(fe.n_components(),
911  std::vector<double>(quadrature.size())
912  )
913  );
914  gradients.resize(global_data->n_vectors(),
915  std::vector< std::vector< Tensor<1,dim> > >(fe.n_components(),
916  std::vector< Tensor<1,dim> >(quadrature.size())
917  )
918  );
919  hessians.resize(global_data->n_vectors(),
920  std::vector< std::vector< Tensor<2,dim> > >(fe.n_components(),
921  std::vector< Tensor<2,dim> >(quadrature.size())
922  )
923  );
924 }
925 
926 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
927 
928 template<int dim, typename FEVALUESBASE>
929 inline
930 void
932 {
933  global_data = &data;
934 }
935 
936 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
937 
938 template<int dim, typename FEVALUESBASE>
939 template<typename DHCellIterator>
940 inline
941 void
943 {
945 
946  for(unsigned int i = 0; i < fevalv.size(); ++i)
947  {
948  FEVALUESBASE& fe_values_base = *fevalv[i];
949  FEValues<dim>& fe_values = dynamic_cast<FEValues<dim>&>(fe_values_base);
950  fe_values.reinit(this->cell);
951  }
952 
953  FEVALUESBASE& fe_values_base_unsplit = *fe_val_unsplit;
954  FEValues<dim>& fe_values_unsplit = dynamic_cast<FEValues<dim>&>(fe_values_base_unsplit);
955  fe_values_unsplit.reinit(this->dof_active_cell);
956 
957  const bool split_fevalues = this->block_info->local_renumbering.size() != 0;
958  fill_local_data(values,
959  split_fevalues);
960  fill_local_data(gradients,
961  split_fevalues);
962  //fill_local_data(hessians,
963  // split_fevalues);
964 }
965 
966 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
967 
968 template<int dim, typename FEVALUESBASE>
969 template<typename DHCellIterator, typename DHFaceIterator>
970 inline
971 void
973  const DHFaceIterator& f,
974  const unsigned int fn)
975 {
977  f,
978  fn);
979 
980  for(unsigned int i = 0; i < fevalv.size(); ++i)
981  {
982  FEVALUESBASE& fe_values_base = *fevalv[i];
983  FEFaceValues<dim>& fe_face_values = dynamic_cast<FEFaceValues<dim>&>(fe_values_base);
984  fe_face_values.reinit(this->cell,
985  fn);
986  }
987 
988  FEVALUESBASE& fe_values_base_unsplit = *fe_val_unsplit;
989  FEFaceValues<dim>& fe_face_values_unsplit = dynamic_cast<FEFaceValues<dim>&>(fe_values_base_unsplit);
990  fe_face_values_unsplit.reinit(this->dof_active_cell,
991  fn);
992 
993  const bool split_fevalues = this->block_info->local_renumbering.size() != 0;
994  fill_local_data(values,
995  split_fevalues);
996  fill_local_data(gradients,
997  split_fevalues);
998  //fill_local_data(hessians,
999  // split_fevalues);
1000 }
1001 
1002 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1003 
1004 template<int dim, typename FEVALUESBASE>
1005 template<typename DHCellIterator, typename DHFaceIterator>
1006 inline
1007 void
1009  const DHFaceIterator& f,
1010  const unsigned int fn,
1011  const unsigned int sn)
1012 {
1014  f,
1015  fn,
1016  sn);
1017 
1018  for(unsigned int i = 0; i < fevalv.size(); ++i)
1019  {
1020  FEVALUESBASE& fe_values_base = *fevalv[i];
1021  FESubfaceValues<dim>& fe_subface_values = dynamic_cast<FESubfaceValues<dim>&>(fe_values_base);
1022  fe_subface_values.reinit(this->cell,
1023  fn,
1024  sn);
1025  }
1026 
1027  FEVALUESBASE& fe_values_base_unsplit = *fe_val_unsplit;
1028  FESubfaceValues<dim>& fe_subface_values_unsplit = dynamic_cast<FESubfaceValues<dim>&>(fe_values_base_unsplit);
1029  fe_subface_values_unsplit.reinit(this->dof_active_cell,
1030  fn,
1031  sn);
1032 
1033  const bool split_fevalues = this->block_info->local_renumbering.size() != 0;
1034  fill_local_data(values,
1035  split_fevalues);
1036  fill_local_data(gradients,
1037  split_fevalues);
1038  //fill_local_data(hessians,
1039  // split_fevalues);
1040 }
1041 
1042 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1043 
1044 template<int dim, typename FEVALUESBASE>
1045 template<typename TYPE>
1046 inline
1047 void
1048 IntegrationInfo<dim, FEVALUESBASE>::fill_local_data(std::vector< std::vector< std::vector<TYPE> > >& data,
1049  bool split_fevalues) const
1050 {
1051  if(split_fevalues)
1052  {
1053  std::vector< std::vector<TYPE> > local_data;
1054 
1055  unsigned int comp = 0;
1056  for(unsigned int b = 0; b < this->block_info->base_element.size(); ++b)
1057  {
1058  const unsigned int no_fe = this->block_info->base_element[b];
1059  const FEValuesBase<dim>& fe_values_base = this->fe(no_fe);
1060 
1061  const unsigned int n_comp = fe_values_base.get_fe().n_components();
1062  local_data.resize(n_comp,
1063  std::vector<TYPE>(fe_values_base.n_quadrature_points));
1064 
1065  for(unsigned int i = 0; i < this->global_data->n_vectors(); ++i)
1066  {
1067  const FEVector& src = this->global_data->vector(i);
1068  fill_data(fe_values_base,
1069  src,
1070  this->indices,
1071  this->block_info->local.block_start(b),
1072  fe_values_base.dofs_per_cell,
1073  local_data);
1074 
1075  for(unsigned int c = 0; c < local_data.size(); ++c)
1076  for(unsigned int k = 0; k < local_data[c].size(); ++k)
1077  data[i][comp+c][k] = local_data[c][k];
1078  }
1079 
1080  comp += n_comp;
1081  }
1082  }
1083  else
1084  {
1085  for(unsigned int i = 0; i < this->global_data->n_vectors(); ++i)
1086  {
1087  const FEVector& src = this->global_data->vector(i);
1088  const FEValuesBase<dim>& fe_values_base = this->fe();
1089 
1090  fill_data(fe_values_base,
1091  src,
1092  this->indices,
1093  0,
1094  fe_values_base.get_fe().dofs_per_cell,
1095  data[i]);
1096  }
1097  }
1098 }
1099 
1100 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1101 
1102 template<int dim, typename FEVALUESBASE>
1103 inline
1104 const FEVALUESBASE&
1106 {
1107  AssertDimension(fevalv.size(), 1);
1108  return *fevalv[0];
1109 }
1110 
1111 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1112 
1113 template<int dim, typename FEVALUESBASE>
1114 inline
1115 const FEVALUESBASE&
1117 {
1118  Assert( i < fevalv.size(), ExcIndexRange(i, 0, fevalv.size()) );
1119  return *fevalv[i];
1120 }
1121 
1122 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1123 
1124 template<int dim, typename FEVALUESBASE>
1125 inline
1126 const FEVALUESBASE&
1128 {
1129  return *fe_val_unsplit;
1130 }
1131 
1132 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1133 
1134 template<int dim, typename FEVALUESBASE>
1135 inline
1136 void
1138 {
1139  fevalv.clear();
1140  fevalv.resize(0);
1141 }
1142 
1143 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1144 
1145 } // ApplicationCore
1146 
1147 } // FuelCell
1148 
1149 #endif
Definition: dof_application.h:82
void reinit(const DHCellIterator &c)
Set the current cell and fill indices.
Definition: mesh_loop_info_objects.h:487
void initialize(const FEVALUES *FE_VALUES, const FiniteElement< dim > &fe, const Mapping< dim > &mapping, const Quadrature< FEVALUES::integral_dimension > &quadrature, const UpdateFlags flags)
Build internal structures fevalv and allocate memory for values, gradients, hessians.
Definition: mesh_loop_info_objects.h:873
const unsigned int dim
Definition: fcst_constants.h:24
unsigned int sub_number
The number of the current subface on the current face of the current cell.
Definition: mesh_loop_info_objects.h:422
Definition: dof_application.h:83
A small structure collecting the different BlockIndices of FEVector vectors (for instance, solution) involved in the computations.
Definition: mesh_loop_info_objects.h:182
Definition: mesh_loop_info_objects.h:41
BlockIndices global
The block structure of the global FEVector solution.
Definition: mesh_loop_info_objects.h:187
SmartPointer< const FEVectors > global_data
The smart pointer to the FEVectors object called global_data.
Definition: mesh_loop_info_objects.h:780
void initialize_data(const FEVectors &data)
Initialize global_data.
Definition: mesh_loop_info_objects.h:931
void clear()
Resize fevalv to 0.
Definition: mesh_loop_info_objects.h:1137
void reinit(const DHCellIterator &c)
Reinitialize internal data for use on a cell.
Definition: mesh_loop_info_objects.h:942
unsigned int face_number
The number of the current face on the current cell.
Definition: mesh_loop_info_objects.h:410
Definition: dof_application.h:84
bool multigrid
true if we assemble for multigrid.
Definition: mesh_loop_info_objects.h:773
std::vector< std::vector< std::vector< Tensor< 1, dim > > > > & derivatives
Definition: mesh_loop_info_objects.h:805
Very basic info class only containing information on geometry and degrees of freedom on a mesh entity...
Definition: mesh_loop_info_objects.h:321
SmartPointer< const BlockInfo > block_info
The block structure of the problem at hand.
Definition: mesh_loop_info_objects.h:438
IntegrationInfo(const BlockInfo &block_info)
Constructor.
Definition: mesh_loop_info_objects.h:844
This class is created for the objects handed to the mesh loops.
Definition: mesh_loop_info_objects.h:625
DoFInfo(const BlockInfo &block_info)
Constructor.
Definition: mesh_loop_info_objects.h:464
boost::shared_ptr< FEVALUESBASE > fe_val_unsplit
In the MODE1 (or in the main mode), the class splits an FEVALUESBASE object into several sub-objects ...
Definition: mesh_loop_info_objects.h:835
std::vector< std::vector< std::vector< Tensor< 1, dim > > > > gradients
The vector containing the gradients of finite element functions in the quadrature points...
Definition: mesh_loop_info_objects.h:800
std::vector< BlockIndices > levels
The multilevel block structure of the global FEVector solution.
Definition: mesh_loop_info_objects.h:197
void fill_data(const FEValuesBase< dim > &fe_values, const FEVector &fe_vector, const std::vector< unsigned int > &local_dof_indices, unsigned int first_index, unsigned int n_indices, std::vector< std::vector< TYPE > > &result)
Helper functions computing the desired data in each quadrature point of a mesh entity by calling FEVa...
Definition: mesh_loop_info_objects.h:107
Triangulation< dim >::face_iterator face
The current Triangulation&lt;dim&gt; face.
Definition: mesh_loop_info_objects.h:399
std::vector< unsigned int > indices
The local dof indices associated with the current cell.
Definition: mesh_loop_info_objects.h:433
BlockIndices local
The block structure of a local FEVector solution.
Definition: mesh_loop_info_objects.h:192
std::vector< std::vector< std::vector< Tensor< 2, dim > > > > hessians
The vector containing the hessians of finite element functions in the quadrature points.
Definition: mesh_loop_info_objects.h:815
void get_indices(const typename DoFHandler< dim, spacedim >::cell_iterator c)
Fill indices.
Definition: mesh_loop_info_objects.h:557
std::vector< std::vector< std::vector< double > > > values
The vector containing the values of finite element functions in the quadrature points.
Definition: mesh_loop_info_objects.h:790
std::vector< unsigned int > local_renumbering
A vector containing the internal renumbering of degrees of freedom from the standard ordering on a ce...
Definition: mesh_loop_info_objects.h:217
Definition: dof_application.h:85
const FEVALUESBASE & get_fe_val_unsplit() const
Access to fe_val_unsplit.
Definition: mesh_loop_info_objects.h:1127
DoFHandler< dim >::active_cell_iterator dof_active_cell
The current DoFHandler&lt;dim&gt; active cell.
Definition: mesh_loop_info_objects.h:384
void fill_local_data(std::vector< std::vector< std::vector< TYPE > > > &data, bool split_fevalues) const
Use the finite element functions in global_data and fill the vectors values, gradients, hessians.
Definition: mesh_loop_info_objects.h:1048
BlockVector< double > FEVector
The vector class used by applications.
Definition: application_data.h:39
The data type used in function calls of Application.
Definition: fe_vectors.h:59
std::vector< boost::shared_ptr< FEVALUESBASE > > fevalv
The vector of smart pointers to FEVALUESBASE objects.
Definition: mesh_loop_info_objects.h:827
DoFHandler< dim >::face_iterator dof_face
The current DoFHandler&lt;dim&gt; face.
Definition: mesh_loop_info_objects.h:389
DoFHandler< dim >::cell_iterator dof_cell
The current DoFHandler&lt;dim&gt; cell.
Definition: mesh_loop_info_objects.h:379
std::vector< unsigned int > base_element
A vector of base elements.
Definition: mesh_loop_info_objects.h:204
Triangulation< dim >::cell_iterator cell
The current Triangulation&lt;dim&gt; cell.
Definition: mesh_loop_info_objects.h:394
std::vector< unsigned int > indices_org
Auxiliary vector.
Definition: mesh_loop_info_objects.h:457
const FEVALUESBASE & fe() const
Access to a single actual FEVALUES object.
Definition: mesh_loop_info_objects.h:1105