19 #ifndef _FUEL_CELL_APPLICATION_CORE_MESH_LOOP_INFO_OBJECTS_H_
20 #define _FUEL_CELL_APPLICATION_CORE_MESH_LOOP_INFO_OBJECTS_H_
22 #include <base/geometry_info.h>
23 #include <base/quadrature_lib.h>
24 #include <base/mg_level_object.h>
26 #include <lac/block_indices.h>
29 #include <fe/fe_tools.h>
30 #include <fe/fe_values.h>
32 #include <multigrid/mg_tools.h>
36 #include <boost/shared_ptr.hpp>
44 using namespace dealii;
66 namespace ApplicationCore
104 template<
int dim,
typename TYPE>
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)
114 Assert(
false, ExcNotImplemented());
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)
129 fe_values.get_function_values(fe_vector,
130 make_slice(local_dof_indices, first_index, n_indices),
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)
147 fe_values.get_function_grads(fe_vector,
148 make_slice(local_dof_indices, first_index, n_indices),
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)
165 fe_values.get_function_hessians(fe_vector,
166 make_slice(local_dof_indices, first_index, n_indices),
222 template<
int dim,
int spacedim>
228 template<
int dim,
int spacedim>
234 template<
int dim,
int spacedim>
240 template<
int dim,
int spacedim>
246 std::vector<unsigned int> sizes(fe.n_blocks());
247 DoFTools::count_dofs_per_block(dof_handler, sizes);
248 global.reinit(sizes);
253 template<
int dim,
int spacedim>
259 std::vector<unsigned int> sizes(fe.n_blocks());
261 base_element.resize(fe.n_blocks());
263 for(
unsigned int i = 0; i < base_element.size(); ++i)
264 base_element[i] = fe.block_to_base_index(i).first;
266 local_renumbering.resize(fe.n_dofs_per_cell());
267 FETools::compute_block_renumbering(fe,
276 template<
int dim,
int spacedim>
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())
287 MGTools::count_dofs_per_block(mg_dof_handler, sizes);
289 levels.resize(sizes.size());
290 for(
unsigned int i = 0; i < sizes.size(); ++i)
291 levels[i].reinit(sizes[i]);
320 template<
int dim,
int spacedim = dim>
349 template<
typename DHCellIterator>
350 void reinit(
const DHCellIterator& c);
356 template<
typename DHCellIterator,
typename DHFaceIterator>
357 void reinit(
const DHCellIterator& c,
358 const DHFaceIterator& f,
359 const unsigned int fn);
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);
462 template<
int dim,
int spacedim>
466 block_info(&block_info,
467 typeid(*this).name())
472 template<
int dim,
int spacedim>
477 block_info(&block_info,
478 typeid(*this).name())
483 template<
int dim,
int spacedim>
484 template<
typename DHCellIterator>
495 face_number = deal_II_numbers::invalid_unsigned_int;
496 sub_number = deal_II_numbers::invalid_unsigned_int;
501 template<
int dim,
int spacedim>
502 template<
typename DHCellIterator,
typename DHFaceIterator>
506 const DHFaceIterator& f,
507 const unsigned int fn)
522 sub_number = deal_II_numbers::invalid_unsigned_int;
527 template<
int dim,
int spacedim>
528 template<
typename DHCellIterator,
typename DHFaceIterator>
532 const DHFaceIterator& f,
533 const unsigned int fn,
534 const unsigned int sn)
554 template<
int dim,
int spacedim>
559 indices.resize(c->get_fe().dofs_per_cell);
561 if(block_info->local_renumbering.size() == 0)
563 c->get_dof_indices(indices);
567 indices_org.resize(c->get_fe().dofs_per_cell);
568 c->get_dof_indices(indices_org);
570 for(
unsigned int i = 0; i < indices.size(); ++i)
571 indices[this->block_info->local_renumbering[i]] = indices_org[i];
577 template<
int dim,
int spacedim>
582 indices.resize(c->get_fe().dofs_per_cell);
584 if(block_info->local_renumbering.size() == 0)
586 c->get_mg_dof_indices(indices);
590 indices_org.resize(c->get_fe().dofs_per_cell);
591 c->get_mg_dof_indices(indices_org);
593 for(
unsigned int i = 0; i < indices.size(); ++i)
594 indices[this->block_info->local_renumbering[i]] = indices_org[i];
624 template<
int dim,
typename FEVALUESBASE>
685 template<
typename FEVALUES>
689 const Quadrature<FEVALUES::integral_dimension>& quadrature,
690 const UpdateFlags flags);
706 template<
typename DHCellIterator>
707 void reinit(
const DHCellIterator& c);
713 template<
typename DHCellIterator,
typename DHFaceIterator>
714 void reinit(
const DHCellIterator& c,
715 const DHFaceIterator& f,
716 const unsigned int fn);
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);
736 template<
typename TYPE>
737 void fill_local_data(std::vector< std::vector< std::vector<TYPE> > >& data,
738 bool split_fevalues)
const;
748 const FEVALUESBASE&
fe()
const;
753 const FEVALUESBASE&
fe(
unsigned int i)
const;
790 std::vector< std::vector< std::vector<double> > >
values;
800 std::vector< std::vector< std::vector< Tensor<1,dim> > > >
gradients;
805 std::vector< std::vector< std::vector< Tensor<1,dim> > > >&
derivatives;
815 std::vector< std::vector< std::vector< Tensor<2,dim> > > >
hessians;
827 std::vector< boost::shared_ptr<FEVALUESBASE> >
fevalv;
842 template<
int dim,
typename FEVALUESBASE>
848 global_data(0, typeid(*this).name()),
849 derivatives(gradients),
855 template<
int dim,
typename FEVALUESBASE>
863 derivatives(gradients),
869 template<
int dim,
typename FEVALUESBASE>
870 template<
typename FEVALUES>
876 const Quadrature<FEVALUES::integral_dimension>& quadrature,
877 const UpdateFlags flags)
879 if( this->block_info->local_renumbering.size() != 0 )
881 fevalv.resize(fe.n_base_elements());
882 for(
unsigned int i = 0; i < fevalv.size(); ++i)
884 fevalv[i] = boost::shared_ptr<FEVALUESBASE>(
new FEVALUES(mapping,
893 fe_val_unsplit = boost::shared_ptr<FEVALUESBASE>(
new FEVALUES(mapping,
902 fevalv[0] = boost::shared_ptr<FEVALUESBASE>(
new FEVALUES(mapping,
909 values.resize(global_data->n_vectors(),
910 std::vector< std::vector<double> >(fe.n_components(),
911 std::vector<double>(quadrature.size())
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())
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())
928 template<
int dim,
typename FEVALUESBASE>
938 template<
int dim,
typename FEVALUESBASE>
939 template<
typename DHCellIterator>
946 for(
unsigned int i = 0; i < fevalv.size(); ++i)
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);
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);
957 const bool split_fevalues = this->block_info->local_renumbering.size() != 0;
958 fill_local_data(values,
960 fill_local_data(gradients,
968 template<
int dim,
typename FEVALUESBASE>
969 template<
typename DHCellIterator,
typename DHFaceIterator>
973 const DHFaceIterator& f,
974 const unsigned int fn)
980 for(
unsigned int i = 0; i < fevalv.size(); ++i)
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,
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,
993 const bool split_fevalues = this->block_info->local_renumbering.size() != 0;
994 fill_local_data(values,
996 fill_local_data(gradients,
1004 template<
int dim,
typename FEVALUESBASE>
1005 template<
typename DHCellIterator,
typename DHFaceIterator>
1009 const DHFaceIterator& f,
1010 const unsigned int fn,
1011 const unsigned int sn)
1018 for(
unsigned int i = 0; i < fevalv.size(); ++i)
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,
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,
1033 const bool split_fevalues = this->block_info->local_renumbering.size() != 0;
1034 fill_local_data(values,
1036 fill_local_data(gradients,
1044 template<
int dim,
typename FEVALUESBASE>
1045 template<
typename TYPE>
1049 bool split_fevalues)
const
1053 std::vector< std::vector<TYPE> > local_data;
1055 unsigned int comp = 0;
1056 for(
unsigned int b = 0; b < this->block_info->base_element.size(); ++b)
1058 const unsigned int no_fe = this->block_info->base_element[b];
1059 const FEValuesBase<dim>& fe_values_base = this->fe(no_fe);
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));
1065 for(
unsigned int i = 0; i < this->global_data->n_vectors(); ++i)
1067 const FEVector& src = this->global_data->vector(i);
1071 this->block_info->local.block_start(b),
1072 fe_values_base.dofs_per_cell,
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];
1085 for(
unsigned int i = 0; i < this->global_data->n_vectors(); ++i)
1087 const FEVector& src = this->global_data->vector(i);
1088 const FEValuesBase<dim>& fe_values_base = this->fe();
1094 fe_values_base.get_fe().dofs_per_cell,
1102 template<
int dim,
typename FEVALUESBASE>
1107 AssertDimension(fevalv.size(), 1);
1113 template<
int dim,
typename FEVALUESBASE>
1118 Assert( i < fevalv.size(), ExcIndexRange(i, 0, fevalv.size()) );
1124 template<
int dim,
typename FEVALUESBASE>
1129 return *fe_val_unsplit;
1134 template<
int dim,
typename FEVALUESBASE>
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<dim> 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<dim> 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<dim> face.
Definition: mesh_loop_info_objects.h:389
DoFHandler< dim >::cell_iterator dof_cell
The current DoFHandler<dim> 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<dim> 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