Namespace containing the basic application framework used to loop over cells and create objects passed to FuelCellShop::Equation objects.
More...
|
| typedef BlockVector< double > | FEVector |
| | The vector class used by applications. More...
|
| |
typedef std::vector
< MatrixBlock< FullMatrix
< double > > > | MatrixVector |
| | The matrix vector used in the mesh loops. More...
|
| |
typedef std::vector< Vector
< double > > | VectorVector |
| | The std::vector of dealii::Vectors used in the mesh loops. More...
|
| |
| typedef void(* | fsub_ptr )(double &, double[], double[], double[]) |
| |
| typedef void(* | dfsub_ptr )(double &, double[], double[], double[]) |
| |
| typedef void(* | gsub_ptr )(int &, double[], double &) |
| |
| typedef void(* | dgsub_ptr )(int &, double[], double[]) |
| |
| typedef void(* | guess_ptr )(double &, double[], double[], double[]) |
| |
|
| template<int dim, typename TYPE > |
| 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 FEValuesBase::get_function_values(), FEValuesBase::get_function_grads(), and FEValuesBase::get_function_hessians(). More...
|
| |
| template<int dim> |
| 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< double > > &result) |
| |
| template<int dim> |
| 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< Tensor< 1, dim > > > &result) |
| |
| template<int dim> |
| 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< Tensor< 2, dim > > > &result) |
| |
| void | DAE_dummy_guess (double &x, double z[], double y[], double df[]) |
| | A dummy guess function to be provided to COLDAE when a user wishes to provide none. More...
|
| |
| void | for_to_c_matrix (int rows, int cols, double *fmat, double **cmat) |
| | Converts a FORTRAN 2D array to a C/C++ 2D array. More...
|
| |
| void | c_to_for_matrix (int rows, int cols, double **cmat, double *fmat) |
| | Converts a C/C++ 2D array to a Fortran 2D array. More...
|
| |
Namespace containing the basic application framework used to loop over cells and create objects passed to FuelCellShop::Equation objects.
This namespace encapsulates all the routines that are used to generate an application. ApplicationBase and ApplicationWrapper are declared here. All applications inherit from the objects developed here.
- Note
- This namespace should only be modified by advanced users as all applications depend on these routines
template<int dim, typename TYPE >
| void FuelCell::ApplicationCore::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 |
|
) |
| |
|
inline |
Helper functions computing the desired data in each quadrature point of a mesh entity by calling FEValuesBase::get_function_values(), FEValuesBase::get_function_grads(), and FEValuesBase::get_function_hessians().
- Parameters
-
| fe_values,: | The FEValues object. |
| fe_vector,: | The global finite element function in the form of a global nodal FEVector. |
| local_dof_indices,: | The local DoF indices associated with the current mesh entity. |
| first_index,: | The first index in local_dof_indices to be used. |
| n_indices,: | The number of indices in local_dof_indices to be used. |
| result,: | The result. The access to ith component in qth quadrature point is result[i][q]. |
- Author
- Guido Kanschat
Referenced by FuelCell::ApplicationCore::IntegrationInfo< dim, FEVALUESBASE >::fill_local_data().