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) 2009-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 //
16 // ----------------------------------------------------------------------------
17 
18 #ifndef _FUEL_CELL_APPLICATION_CORE_MESH_LOOP_INFO_OBJECTS_H_
19 #define _FUEL_CELL_APPLICATION_CORE_MESH_LOOP_INFO_OBJECTS_H_
20 
21 //-- deal.II
22 #include <deal.II/base/geometry_info.h>
23 #include <deal.II/base/quadrature_lib.h>
24 #include <deal.II/base/vector_slice.h>
25 #include <deal.II/lac/block_indices.h>
26 #include <deal.II/fe/fe.h>
27 #include <deal.II/fe/fe_tools.h>
28 #include <deal.II/fe/fe_values.h>
29 #include <deal.II/dofs/dof_tools.h>
30 
31 //-- OpenFCST
33 
34 //-- boost libraries:
35 #include <boost/shared_ptr.hpp>
36 
37 //----
38 namespace dealii
39 {
40  template<int,int> class DoFHandler;
41 }
42 
43 using namespace dealii;
44 
63 namespace FuelCell
64 {
65 namespace ApplicationCore
66 {
67 
91 // +++ IMPLEMENTATION +++++++++++++++++++++++++++++++++++++++++++++++++++++++++
92 
93 template<int dim, typename TYPE>
94 inline
95 void
96 fill_data(const FEValuesBase<dim>& fe_values,
97  const FEVector& fe_vector,
98  const std::vector<unsigned int>& local_dof_indices,
99  unsigned int first_index,
100  unsigned int n_indices,
101  std::vector< std::vector<TYPE> >& result)
102 {
103  Assert(false, ExcNotImplemented());
104 }
105 
106 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
107 
108 template<int dim>
109 inline
110 void
111 fill_data(const FEValuesBase<dim>& fe_values,
112  const FEVector& fe_vector,
113  const std::vector<unsigned int>& local_dof_indices,
114  unsigned int first_index,
115  unsigned int n_indices,
116  std::vector< std::vector<double> >& result)
117 {
118 
119  VectorSlice<std::vector< std::vector<double> > > aux(result);
120  fe_values.get_function_values(fe_vector,
121  make_slice(local_dof_indices, first_index, n_indices),
122  aux,
123  true);
124 }
125 
126 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
127 
128 template<int dim>
129 inline
130 void
131 fill_data(const FEValuesBase<dim>& fe_values,
132  const FEVector& fe_vector,
133  const std::vector<unsigned int>& local_dof_indices,
134  unsigned int first_index,
135  unsigned int n_indices,
136  std::vector< std::vector< Tensor<1,dim> > >& result)
137 {
138  VectorSlice< std::vector< std::vector< Tensor<1,dim> > > > aux(result);
139  fe_values.get_function_gradients(fe_vector,
140  make_slice(local_dof_indices, first_index, n_indices),
141  aux,
142  true);
143 }
144 
145 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
146 
147 template<int dim>
148 inline
149 void
150 fill_data(const FEValuesBase<dim>& fe_values,
151  const FEVector& fe_vector,
152  const std::vector<unsigned int>& local_dof_indices,
153  unsigned int first_index,
154  unsigned int n_indices,
155  std::vector< std::vector< Tensor<2,dim> > >& result)
156 {
157  fe_values.get_function_hessians(fe_vector,
158  make_slice(local_dof_indices, first_index, n_indices),
159  result,
160  true);
161 }
162 
163 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
164 
173 struct BlockInfo : public Subscriptor
174 {
178  BlockIndices global;
179 
183  BlockIndices local;
184 
188  std::vector<BlockIndices> levels;
189 
195  std::vector<unsigned int> base_element;
196 
206  std::vector<unsigned int> local_renumbering;
207 
211  template<int dim, int spacedim>
212  void initialize(const DoFHandler<dim, spacedim>& dof_handler)
213  {
214  const FiniteElement<dim, spacedim>& fe = dof_handler.get_fe();
215  std::vector<unsigned int> sizes(fe.n_blocks());
216  DoFTools::count_dofs_per_block(dof_handler, sizes);
217  global.reinit(sizes);
218  }
219 
223  template<int dim, int spacedim>
225  {
226  const FiniteElement<dim, spacedim>& fe = dof_handler.get_fe();
227  std::vector<unsigned int> sizes(fe.n_blocks());
228 
229  base_element.resize(fe.n_blocks());
230 
231  for(unsigned int i = 0; i < base_element.size(); ++i)
232  base_element[i] = fe.block_to_base_index(i).first;
233 
234  local_renumbering.resize(fe.n_dofs_per_cell());
235  FETools::compute_block_renumbering(fe,
236  local_renumbering,
237  sizes,
238  false);
239  local.reinit(sizes);
240  }
241 };
242 
243 
267 template<int dim, int spacedim = dim>
268 class DoFInfo
269 {
270 public:
271 
273 
274 
278  DoFInfo(const BlockInfo& block_info);
279 
284  DoFInfo(const FEVectors&,
285  const BlockInfo& block_info);
286 
288 
290 
291 
296  template<typename DHCellIterator>
297  void reinit(const DHCellIterator& c);
298 
303  template<typename DHCellIterator, typename DHFaceIterator>
304  void reinit(const DHCellIterator& c,
305  const DHFaceIterator& f,
306  const unsigned int fn);
307 
312  template<typename DHCellIterator, typename DHFaceIterator>
313  void reinit(const DHCellIterator& c,
314  const DHFaceIterator& f,
315  const unsigned int fn,
316  const unsigned int sn);
317 
319 
321 
322 
327 
332 
337 
342 
347 
357  unsigned int face_number;
358 
369  unsigned int sub_number;
370 
372 
374 
375 
380  std::vector<unsigned int> indices;
381 
385  SmartPointer<const BlockInfo> block_info;
386 
388 
389 private:
390 
394  void get_indices(const typename DoFHandler<dim, spacedim>::cell_iterator c);
395 
399  std::vector<unsigned int> indices_org;
400 };
401 
402 // +++ IMPLEMENTATION +++++++++++++++++++++++++++++++++++++++++++++++++++++++++
403 
404 template<int dim, int spacedim>
405 inline
407 :
408 block_info(&block_info,
409  typeid(*this).name())
410 { }
411 
412 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
413 
414 template<int dim, int spacedim>
415 inline
417  const BlockInfo& block_info)
418 :
419 block_info(&block_info,
420  typeid(*this).name())
421 { }
422 
423 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
424 
425 template<int dim, int spacedim>
426 template<typename DHCellIterator>
427 inline
428 void
429 DoFInfo<dim, spacedim>::reinit(const DHCellIterator& c)
430 {
431  get_indices(c);
432 
433  dof_cell = c;
434  dof_active_cell = c;
435  cell = static_cast<typename Triangulation<dim>::cell_iterator> (c);
436 
437  face_number = static_cast<unsigned int>(-1);
438  sub_number = static_cast<unsigned int>(-1);
439 }
440 
441 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
442 
443 template<int dim, int spacedim>
444 template<typename DHCellIterator, typename DHFaceIterator>
445 inline
446 void
447 DoFInfo<dim, spacedim>::reinit(const DHCellIterator& c,
448  const DHFaceIterator& f,
449  const unsigned int fn)
450 {
451  if( cell.state() != IteratorState::valid || cell != static_cast<typename Triangulation<dim>::cell_iterator> (c) )
452  {
453  get_indices(c);
454  }
455 
456  dof_cell = c;
457  dof_active_cell = c;
458  cell = static_cast<typename Triangulation<dim>::cell_iterator> (c);
459 
460  dof_face = f;
461  face = static_cast<typename Triangulation<dim>::face_iterator> (f);
462 
463  face_number = fn;
464  sub_number = static_cast<unsigned int>(-1);
465 }
466 
467 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
468 
469 template<int dim, int spacedim>
470 template<typename DHCellIterator, typename DHFaceIterator>
471 inline
472 void
473 DoFInfo<dim, spacedim>::reinit(const DHCellIterator& c,
474  const DHFaceIterator& f,
475  const unsigned int fn,
476  const unsigned int sn)
477 {
478  if( cell.state() != IteratorState::valid || cell != static_cast<typename Triangulation<dim>::cell_iterator> (c) )
479  {
480  get_indices(c);
481  }
482 
483  dof_cell = c;
484  dof_active_cell = c;
485  cell = static_cast<typename Triangulation<dim>::cell_iterator> (c);
486 
487  dof_face = f;
488  face = static_cast<typename Triangulation<dim>::face_iterator> (f);
489 
490  face_number = fn;
491  sub_number = sn;
492 }
493 
494 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
495 
496 template<int dim, int spacedim>
497 inline
498 void
500 {
501  indices.resize(c->get_fe().dofs_per_cell);
502 
503  if(block_info->local_renumbering.size() == 0)
504  {
505  c->get_dof_indices(indices);
506  }
507  else
508  {
509  indices_org.resize(c->get_fe().dofs_per_cell);
510  c->get_dof_indices(indices_org);
511 
512  for(unsigned int i = 0; i < indices.size(); ++i)
513  indices[this->block_info->local_renumbering[i]] = indices_org[i];
514  }
515 }
516 
517 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
518 
543 template<int dim, typename FEVALUESBASE>
544 class IntegrationInfo : public DoFInfo<dim>
545 {
546 public:
547 
549 
550 
555 
559  IntegrationInfo(const FEVectors& data,
560  const BlockInfo& block_info);
561 
586  template<typename FEVALUES>
587  void initialize(const FEVALUES* FE_VALUES,
588  const FiniteElement<dim>& fe,
589  const Mapping<dim>& mapping,
590  const Quadrature<FEVALUES::integral_dimension>& quadrature,
591  const UpdateFlags flags);
592 
596  void initialize_data(const FEVectors& data);
597 
599 
601 
602 
606  template<typename DHCellIterator>
607  void reinit(const DHCellIterator& c);
608 
612  template<typename DHCellIterator, typename DHFaceIterator>
613  void reinit(const DHCellIterator& c,
614  const DHFaceIterator& f,
615  const unsigned int fn);
616 
620  template<typename DHCellIterator, typename DHFaceIterator>
621  void reinit(const DHCellIterator& c,
622  const DHFaceIterator& f,
623  const unsigned int fn,
624  const unsigned int sn);
625 
632  template<typename TYPE>
633  void fill_local_data(std::vector< std::vector< std::vector<TYPE> > >& data,
634  bool split_fevalues) const;
635 
637 
639 
640 
644  const FEVALUESBASE& fe() const;
645 
649  const FEVALUESBASE& fe(unsigned int i) const;
650 
654  const FEVALUESBASE& get_fe_val_unsplit() const;
655 
659  void clear();
660 
662 
664 
665 
669  bool multigrid;
670 
674  SmartPointer<const FEVectors> global_data;
675 
682  std::vector< std::vector< std::vector<double> > > values;
683 
690  std::vector< std::vector< std::vector< Tensor<1,dim> > > > gradients;
691 
695  std::vector< std::vector< std::vector< Tensor<1,dim> > > >& derivatives;
696 
703  std::vector< std::vector< std::vector< Tensor<2,dim> > > > hessians;
704 
706 
707 private:
708 
710 
711 
715  std::vector< boost::shared_ptr<FEVALUESBASE> > fevalv;
716 
723  boost::shared_ptr<FEVALUESBASE> fe_val_unsplit;
724 
726 };
727 
728 // +++ IMPLEMENTATION +++++++++++++++++++++++++++++++++++++++++++++++++++++++++
729 
730 template<int dim, typename FEVALUESBASE>
731 inline
733 :
734 DoFInfo<dim>(block_info),
735 multigrid(false),
736 global_data(0, typeid(*this).name()),
737 derivatives(gradients),
738 fevalv(0)
739 { }
740 
741 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
742 
743 template<int dim, typename FEVALUESBASE>
744 inline
746  const BlockInfo& block_info)
747 :
748 DoFInfo<dim>(block_info),
749 multigrid(false),
750 global_data(&data),
751 derivatives(gradients),
752 fevalv(0)
753 { }
754 
755 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
756 
757 template<int dim, typename FEVALUESBASE>
758 template<typename FEVALUES>
759 inline
760 void
762  const FiniteElement<dim>& fe,
763  const Mapping<dim>& mapping,
764  const Quadrature<FEVALUES::integral_dimension>& quadrature,
765  const UpdateFlags flags)
766 {
767  if( this->block_info->local_renumbering.size() != 0 ) // MODE1
768  {
769  fevalv.resize(fe.n_base_elements());
770  for(unsigned int i = 0; i < fevalv.size(); ++i)
771  {
772  fevalv[i] = boost::shared_ptr<FEVALUESBASE>( new FEVALUES(mapping,
773  fe.base_element(i),
774  quadrature,
775  flags)
776  );
777  }
778 
779  // --- fe_val_unsplit initialization ---
780 
781  fe_val_unsplit = boost::shared_ptr<FEVALUESBASE>( new FEVALUES(mapping,
782  fe,
783  quadrature,
784  flags)
785  );
786  }
787  else // MODE2
788  {
789  fevalv.resize(1);
790  fevalv[0] = boost::shared_ptr<FEVALUESBASE>( new FEVALUES(mapping,
791  fe,
792  quadrature,
793  flags)
794  );
795  }
796 
797  values.resize(global_data->n_vectors(),
798  std::vector< std::vector<double> >(fe.n_components(),
799  std::vector<double>(quadrature.size())
800  )
801  );
802  gradients.resize(global_data->n_vectors(),
803  std::vector< std::vector< Tensor<1,dim> > >(fe.n_components(),
804  std::vector< Tensor<1,dim> >(quadrature.size())
805  )
806  );
807  hessians.resize(global_data->n_vectors(),
808  std::vector< std::vector< Tensor<2,dim> > >(fe.n_components(),
809  std::vector< Tensor<2,dim> >(quadrature.size())
810  )
811  );
812 }
813 
814 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
815 
816 template<int dim, typename FEVALUESBASE>
817 inline
818 void
820 {
821  global_data = &data;
822 }
823 
824 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
825 
826 template<int dim, typename FEVALUESBASE>
827 template<typename DHCellIterator>
828 inline
829 void
831 {
833 
834  for(unsigned int i = 0; i < fevalv.size(); ++i)
835  {
836  FEVALUESBASE& fe_values_base = *fevalv[i];
837  FEValues<dim>& fe_values = dynamic_cast<FEValues<dim>&>(fe_values_base);
838  fe_values.reinit(this->cell);
839  }
840 
841  FEVALUESBASE& fe_values_base_unsplit = *fe_val_unsplit;
842  FEValues<dim>& fe_values_unsplit = dynamic_cast<FEValues<dim>&>(fe_values_base_unsplit);
843  fe_values_unsplit.reinit(this->dof_active_cell);
844 
845  const bool split_fevalues = this->block_info->local_renumbering.size() != 0;
846  fill_local_data(values,
847  split_fevalues);
848  fill_local_data(gradients,
849  split_fevalues);
850  //fill_local_data(hessians,
851  // split_fevalues);
852 }
853 
854 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
855 
856 template<int dim, typename FEVALUESBASE>
857 template<typename DHCellIterator, typename DHFaceIterator>
858 inline
859 void
861  const DHFaceIterator& f,
862  const unsigned int fn)
863 {
865  f,
866  fn);
867 
868  for(unsigned int i = 0; i < fevalv.size(); ++i)
869  {
870  FEVALUESBASE& fe_values_base = *fevalv[i];
871  FEFaceValues<dim>& fe_face_values = dynamic_cast<FEFaceValues<dim>&>(fe_values_base);
872  fe_face_values.reinit(this->cell,
873  fn);
874  }
875 
876  FEVALUESBASE& fe_values_base_unsplit = *fe_val_unsplit;
877  FEFaceValues<dim>& fe_face_values_unsplit = dynamic_cast<FEFaceValues<dim>&>(fe_values_base_unsplit);
878  fe_face_values_unsplit.reinit(this->dof_active_cell,
879  fn);
880 
881  const bool split_fevalues = this->block_info->local_renumbering.size() != 0;
882  fill_local_data(values,
883  split_fevalues);
884  fill_local_data(gradients,
885  split_fevalues);
886  //fill_local_data(hessians,
887  // split_fevalues);
888 }
889 
890 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
891 
892 template<int dim, typename FEVALUESBASE>
893 template<typename DHCellIterator, typename DHFaceIterator>
894 inline
895 void
897  const DHFaceIterator& f,
898  const unsigned int fn,
899  const unsigned int sn)
900 {
901  DoFInfo<dim>::reinit(c, f, fn, sn);
902 
903  for(unsigned int i = 0; i < fevalv.size(); ++i)
904  {
905  FEVALUESBASE& fe_values_base = *fevalv[i];
906  FESubfaceValues<dim>& fe_subface_values = dynamic_cast<FESubfaceValues<dim>&>(fe_values_base);
907  fe_subface_values.reinit(this->cell,
908  fn,
909  sn);
910  }
911 
912  FEVALUESBASE& fe_values_base_unsplit = *fe_val_unsplit;
913  FESubfaceValues<dim>& fe_subface_values_unsplit = dynamic_cast<FESubfaceValues<dim>&>(fe_values_base_unsplit);
914  fe_subface_values_unsplit.reinit(this->dof_active_cell,
915  fn,
916  sn);
917 
918  const bool split_fevalues = this->block_info->local_renumbering.size() != 0;
919  fill_local_data(values,
920  split_fevalues);
921  fill_local_data(gradients,
922  split_fevalues);
923  //fill_local_data(hessians,
924  // split_fevalues);
925 }
926 
927 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
928 
929 template<int dim, typename FEVALUESBASE>
930 template<typename TYPE>
931 inline
932 void
933 IntegrationInfo<dim, FEVALUESBASE>::fill_local_data(std::vector< std::vector< std::vector<TYPE> > >& data,
934  bool split_fevalues) const
935 {
936  if(split_fevalues)
937  {
938  std::vector< std::vector<TYPE> > local_data;
939 
940  unsigned int comp = 0;
941  for(unsigned int b = 0; b < this->block_info->base_element.size(); ++b)
942  {
943  const unsigned int no_fe = this->block_info->base_element[b];
944  const FEValuesBase<dim>& fe_values_base = this->fe(no_fe);
945 
946  const unsigned int n_comp = fe_values_base.get_fe().n_components();
947  local_data.resize(n_comp,
948  std::vector<TYPE>(fe_values_base.n_quadrature_points));
949 
950  for(unsigned int i = 0; i < this->global_data->n_vectors(); ++i)
951  {
952  const FEVector& src = this->global_data->vector(i);
953  fill_data(fe_values_base,
954  src,
955  this->indices,
956  this->block_info->local.block_start(b),
957  fe_values_base.dofs_per_cell,
958  local_data);
959 
960  for(unsigned int c = 0; c < local_data.size(); ++c)
961  for(unsigned int k = 0; k < local_data[c].size(); ++k)
962  data[i][comp+c][k] = local_data[c][k];
963  }
964 
965  comp += n_comp;
966  }
967  }
968  else
969  {
970  for(unsigned int i = 0; i < this->global_data->n_vectors(); ++i)
971  {
972  const FEVector& src = this->global_data->vector(i);
973  const FEValuesBase<dim>& fe_values_base = this->fe();
974 
975  fill_data(fe_values_base,
976  src,
977  this->indices,
978  0,
979  fe_values_base.get_fe().dofs_per_cell,
980  data[i]);
981  }
982  }
983 }
984 
985 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
986 
987 template<int dim, typename FEVALUESBASE>
988 inline
989 const FEVALUESBASE&
991 {
992  AssertDimension(fevalv.size(), 1);
993  return *fevalv[0];
994 }
995 
996 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
997 
998 template<int dim, typename FEVALUESBASE>
999 inline
1000 const FEVALUESBASE&
1002 {
1003  Assert( i < fevalv.size(), ExcIndexRange(i, 0, fevalv.size()) );
1004  return *fevalv[i];
1005 }
1006 
1007 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1008 
1009 template<int dim, typename FEVALUESBASE>
1010 inline
1011 const FEVALUESBASE&
1013 {
1014  return *fe_val_unsplit;
1015 }
1016 
1017 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1018 
1019 template<int dim, typename FEVALUESBASE>
1020 inline
1021 void
1023 {
1024  fevalv.clear();
1025  fevalv.resize(0);
1026 }
1027 
1028 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1029 
1030 } // ApplicationCore
1031 
1032 } // FuelCell
1033 
1034 #endif
Definition: dof_application.h:67
void reinit(const DHCellIterator &c)
Set the current cell and fill indices.
Definition: mesh_loop_info_objects.h:429
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:761
const unsigned int dim
Definition: fcst_constants.h:23
unsigned int sub_number
The number of the current subface on the current face of the current cell.
Definition: mesh_loop_info_objects.h:369
Definition: dof_application.h:68
A small structure collecting the different BlockIndices of FEVector vectors (for instance, solution) involved in the computations.
Definition: mesh_loop_info_objects.h:173
BlockIndices global
The block structure of the global FEVector solution.
Definition: mesh_loop_info_objects.h:178
SmartPointer< const FEVectors > global_data
The smart pointer to the FEVectors object called global_data.
Definition: mesh_loop_info_objects.h:674
void initialize_data(const FEVectors &data)
Initialize global_data.
Definition: mesh_loop_info_objects.h:819
void clear()
Resize fevalv to 0.
Definition: mesh_loop_info_objects.h:1022
void reinit(const DHCellIterator &c)
Reinitialize internal data for use on a cell.
Definition: mesh_loop_info_objects.h:830
unsigned int face_number
The number of the current face on the current cell.
Definition: mesh_loop_info_objects.h:357
Definition: dof_application.h:69
bool multigrid
true if we assemble for multigrid.
Definition: mesh_loop_info_objects.h:669
std::vector< std::vector< std::vector< Tensor< 1, dim > > > > & derivatives
Definition: mesh_loop_info_objects.h:695
Very basic info class only containing information on geometry and degrees of freedom on a mesh entity...
Definition: mesh_loop_info_objects.h:268
SmartPointer< const BlockInfo > block_info
The block structure of the problem at hand.
Definition: mesh_loop_info_objects.h:385
IntegrationInfo(const BlockInfo &block_info)
Constructor.
Definition: mesh_loop_info_objects.h:732
This class is created for the objects handed to the mesh loops.
Definition: mesh_loop_info_objects.h:544
DoFInfo(const BlockInfo &block_info)
Constructor.
Definition: mesh_loop_info_objects.h:406
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:723
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:690
std::vector< BlockIndices > levels
The multilevel block structure of the global FEVector solution.
Definition: mesh_loop_info_objects.h:188
void initialize(const DoFHandler< dim, spacedim > &dof_handler)
Fill the structure globally with values describing the DoFHandler.
Definition: mesh_loop_info_objects.h:212
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:96
void initialize_local(const DoFHandler< dim, spacedim > &dof_handler)
Fill the structure locally with values describing the DoFHandler.
Definition: mesh_loop_info_objects.h:224
Triangulation< dim >::face_iterator face
The current Triangulation&lt;dim&gt; face.
Definition: mesh_loop_info_objects.h:346
std::vector< unsigned int > indices
The local dof indices associated with the current cell.
Definition: mesh_loop_info_objects.h:380
BlockIndices local
The block structure of a local FEVector solution.
Definition: mesh_loop_info_objects.h:183
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:703
void get_indices(const typename DoFHandler< dim, spacedim >::cell_iterator c)
Fill indices.
Definition: mesh_loop_info_objects.h:499
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:682
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:206
Definition: dof_application.h:70
const FEVALUESBASE & get_fe_val_unsplit() const
Access to fe_val_unsplit.
Definition: mesh_loop_info_objects.h:1012
DoFHandler< dim >::active_cell_iterator dof_active_cell
The current DoFHandler&lt;dim&gt; active cell.
Definition: mesh_loop_info_objects.h:331
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:933
BlockVector< double > FEVector
The vector class used by applications.
Definition: application_data.h:46
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:715
DoFHandler< dim >::face_iterator dof_face
The current DoFHandler&lt;dim&gt; face.
Definition: mesh_loop_info_objects.h:336
DoFHandler< dim >::cell_iterator dof_cell
The current DoFHandler&lt;dim&gt; cell.
Definition: mesh_loop_info_objects.h:326
std::vector< unsigned int > base_element
A vector of base elements.
Definition: mesh_loop_info_objects.h:195
Triangulation< dim >::cell_iterator cell
The current Triangulation&lt;dim&gt; cell.
Definition: mesh_loop_info_objects.h:341
std::vector< unsigned int > indices_org
Auxiliary vector.
Definition: mesh_loop_info_objects.h:399
const FEVALUESBASE & fe() const
Access to a single actual FEVALUES object.
Definition: mesh_loop_info_objects.h:990