pi-DoMUS: Parallel Deal.II MUltiphysics Solver
Interface< dim, spacedim, n_components, LAC > Class Template Reference

Interface. More...

#include <old_interface.h>

Inheritance diagram for Interface< dim, spacedim, n_components, LAC >:
ConservativeInterface< dim, spacedim, n_components, Implementation, LAC > ConservativeInterface< dim, dim, 1, HeatEquation< dim, LAC >, LAC > ConservativeInterface< dim, spacedim, dim+2, FreeSwellingThreeFields< dim, spacedim >, LAC > NonConservativeInterface< dim, spacedim, n_components, Implementation, LAC > HeatEquation< dim, LAC > FreeSwellingThreeFields< dim, spacedim >

Public Member Functions

virtual ~Interface ()
 
 Interface (const std::string &name="", const std::string &default_fe="FE_Q(1)", const std::string &default_component_names="u", const std::string &default_coupling="", const std::string &default_preconditioner_coupling="", const std::string &default_differential_components="")
 
virtual void declare_parameters (ParameterHandler &prm)
 
virtual void parse_parameters_call_back ()
 
const std::vector< unsigned int > get_differential_blocks () const
 
virtual void set_time (const double &t) const
 update time of all parsed mapped functions More...
 
virtual void postprocess_newly_created_triangulation (Triangulation< dim, spacedim > &tria) const
 This function is used to modify triangulation using boundary_id or manifold_id. More...
 
virtual void apply_dirichlet_bcs (const DoFHandler< dim, spacedim > &dof_handler, ConstraintMatrix &constraints) const
 Applies Dirichlet boundary conditions. More...
 
template<typename Number >
void apply_neumann_bcs (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Scratch &scratch, CopySystem &data, std::vector< Number > &local_residual) const
 Applies Neumann boundary conditions. More...
 
template<typename Number >
void apply_forcing_terms (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Scratch &scratch, CopySystem &data, std::vector< Number > &local_residual) const
 Applies CONSERVATIVE forcing terms. More...
 
virtual void initialize_data (const typename LAC::VectorType &solution, const typename LAC::VectorType &solution_dot, const double t, const double alpha) const
 Initialize all data required for the system. More...
 
virtual void get_preconditioner_energy (const typename DoFHandler< dim, spacedim >::active_cell_iterator &, Scratch &, CopySystem &, Sdouble &) const
 Build the energy needed to get the preconditioner in the case it is required just one derivative. More...
 
virtual void get_preconditioner_energy (const typename DoFHandler< dim, spacedim >::active_cell_iterator &, Scratch &, CopyPreconditioner &, SSdouble &) const
 Build the energy needed to get the preconditioner in the case it is required two derivatives. More...
 
virtual void get_system_energy (const typename DoFHandler< dim, spacedim >::active_cell_iterator &, Scratch &, CopySystem &, Sdouble &) const
 Build the energy needed to get the system matrix in the case it is required two derivatives. More...
 
virtual void get_system_energy (const typename DoFHandler< dim, spacedim >::active_cell_iterator &, Scratch &, CopySystem &, SSdouble &) const
 Build the energy needed to get the system matrix in the case it is required two derivatives. More...
 
virtual void get_system_residual (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Scratch &scratch, CopySystem &data, std::vector< Sdouble > &local_residual) const
 Build the residual needed to get the system matrix in the case it is required two derivatives. More...
 
virtual void get_system_residual (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Scratch &scratch, CopySystem &data, std::vector< double > &local_residual) const
 Build the residual needed to get the system matrix in the case it is required just one derivative. More...
 
virtual void get_preconditioner_residual (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Scratch &scratch, CopyPreconditioner &data, std::vector< Sdouble > &local_residual) const
 Build the residual needed to get the preconditioner matrix in the case two derivatives are required. More...
 
virtual void compute_system_operators (const DoFHandler< dim, spacedim > &, const typename LAC::BlockMatrix &, const typename LAC::BlockMatrix &, const std::vector< shared_ptr< typename LAC::BlockMatrix > >, LinearOperator< typename LAC::VectorType > &, LinearOperator< typename LAC::VectorType > &) const
 Compute linear operators needed by the problem. More...
 
virtual void assemble_local_system (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Scratch &scratch, CopySystem &data) const
 
virtual void assemble_local_preconditioner (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Scratch &scratch, CopyPreconditioner &data) const
 
virtual shared_ptr< Mapping< dim, spacedim > > get_mapping (const DoFHandler< dim, spacedim > &, const typename LAC::VectorType &) const
 
virtual UpdateFlags get_jacobian_flags () const
 
virtual UpdateFlags get_residual_flags () const
 
virtual UpdateFlags get_jacobian_preconditioner_flags () const
 
virtual UpdateFlags get_face_flags () const
 
void fix_solution_dot_derivative (FEValuesCache< dim, spacedim > &, double) const
 
void fix_solution_dot_derivative (FEValuesCache< dim, spacedim > &fe_cache, Sdouble alpha) const
 
void fix_solution_dot_derivative (FEValuesCache< dim, spacedim > &fe_cache, SSdouble alpha) const
 
template<typename Number >
void reinit (const Number &alpha, const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, FEValuesCache< dim, spacedim > &fe_cache) const
 
template<typename Number >
void reinit (const Number &alpha, const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, const unsigned int face_no, FEValuesCache< dim, spacedim > &fe_cache) const
 
virtual void get_aux_matrix_residuals (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Scratch &scratch, CopyPreconditioner &data, std::vector< std::vector< double > > &local_residuals) const
 
virtual void get_aux_matrix_residuals (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Scratch &scratch, CopyPreconditioner &data, std::vector< std::vector< Sdouble > > &local_residuals) const
 
virtual void assemble_local_aux_matrices (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Scratch &scratch, CopyPreconditioner &data) const
 
UpdateFlags get_aux_matrix_flags (const unsigned int &i) const
 
const Table< 2, DoFTools::Coupling > & get_aux_matrix_coupling (const unsigned int &i) const
 
virtual unsigned int get_number_of_aux_matrices () const
 

Public Attributes

std::vector< UpdateFlagsaux_matrix_update_flags
 
std::vector< Table< 2, DoFTools::Coupling > > aux_matrix_coupling
 

Protected Attributes

ParsedMappedFunctions< spacedim, n_components > forcing_terms
 
ParsedMappedFunctions< spacedim, n_components > neumann_bcs
 
ParsedDirichletBCs< dim, spacedim, n_components > dirichlet_bcs
 
std::string str_diff_comp
 
std::vector< unsigned int > _diff_comp
 
const LAC::VectorTypesolution
 Solution vector evaluated at time t. More...
 
LAC::VectorType old_solution
 Solution vector evaluated at time t-dt. More...
 
const LAC::VectorTypesolution_dot
 Time derivative solution vector evaluated at time t. More...
 
double t
 Current time step. More...
 
double old_t
 Previous time step. More...
 
double alpha
 
unsigned int dofs_per_cell
 
unsigned int n_q_points
 
unsigned int n_face_q_points
 

Detailed Description

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
class Interface< dim, spacedim, n_components, LAC >

Interface.

This class has two child: conservative_interface.h and non_conservative_interface.h Users should not derive directly from this class, but from its specialization classes. For istance, a stokes problem should consist in a interface class derived from conservative_interface.h. (see include/interfaces/ for some examples)

Goal: provide a derivable interface to solve a particular PDEs problem (time depending, first-order, non linear).

Usage: This class requires some arguments related to the setting of the problem: finite elements, boundary conditions, and initial conditions. Moreover, it helps to write the system matrix and the preconditioner matrix. (see conservative_interface.h and non_conservative_interface.h)

Varibles:

  • General:
    • Finite Elements
    • Boundary conditions ( Dirichlet, Neumann, and Robin )
    • Initial conditions ( y(0) and d/dt y (0) )
  • System Matrix:
    • coupling (coupling variable is a matrix of zeroes and ones: if the row and the coloumn are indipendent there will be 0, otherwise 1)
  • Preconditioner:
    • preconditioner coupling (same as above)
  • default_differential_components : this variable is a list of ones and zeroes. 1 in the case the corresponding variable should be differentiable and 0 otherwise. TODO: add flags

Constructor & Destructor Documentation

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual Interface< dim, spacedim, n_components, LAC >::~Interface ( )
inlinevirtual
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
Interface< dim, spacedim, n_components, LAC >::Interface ( const std::string &  name = "",
const std::string &  default_fe = "FE_Q(1)",
const std::string &  default_component_names = "u",
const std::string &  default_coupling = "",
const std::string &  default_preconditioner_coupling = "",
const std::string &  default_differential_components = "" 
)
inline

Member Function Documentation

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::apply_dirichlet_bcs ( const DoFHandler< dim, spacedim > &  dof_handler,
ConstraintMatrix constraints 
) const
virtual

Applies Dirichlet boundary conditions.

This function is used to applies Dirichlet boundary conditions. It takes as argument a DoF handler dof_handler and a constraint matrix constraints.

template<int dim, int spacedim, int n_components, typename LAC >
template<typename Number >
void Interface< dim, spacedim, n_components, LAC >::apply_forcing_terms ( const typename DoFHandler< dim, spacedim >::active_cell_iterator cell,
Scratch &  scratch,
CopySystem &  data,
std::vector< Number > &  local_residual 
) const

Applies CONSERVATIVE forcing terms.

This function applies the conservative forcing terms, which can be defined by expressions in the parameter file.

If the problem involves NON-conservative loads, they must be included in the residual formulation.

template<int dim, int spacedim, int n_components, typename LAC >
template<typename Number >
void Interface< dim, spacedim, n_components, LAC >::apply_neumann_bcs ( const typename DoFHandler< dim, spacedim >::active_cell_iterator cell,
Scratch &  scratch,
CopySystem &  data,
std::vector< Number > &  local_residual 
) const

Applies Neumann boundary conditions.

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::assemble_local_aux_matrices ( const typename DoFHandler< dim, spacedim >::active_cell_iterator cell,
Scratch &  scratch,
CopyPreconditioner &  data 
) const
virtual
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::assemble_local_preconditioner ( const typename DoFHandler< dim, spacedim >::active_cell_iterator cell,
Scratch &  scratch,
CopyPreconditioner &  data 
) const
virtual
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::assemble_local_system ( const typename DoFHandler< dim, spacedim >::active_cell_iterator cell,
Scratch &  scratch,
CopySystem &  data 
) const
virtual
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::compute_system_operators ( const DoFHandler< dim, spacedim > &  ,
const typename LAC::BlockMatrix ,
const typename LAC::BlockMatrix ,
const std::vector< shared_ptr< typename LAC::BlockMatrix > >  ,
LinearOperator< typename LAC::VectorType > &  ,
LinearOperator< typename LAC::VectorType > &   
) const
virtual

Compute linear operators needed by the problem.

This function is used to assemble linear operators related to the problem. It takes a reference to DoF Handler, two references to block sparse matrices representing the system matrix and the preconditioner, and two references to LinearOperator.

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
void Interface< dim, spacedim, n_components, LAC >::fix_solution_dot_derivative ( FEValuesCache< dim, spacedim > &  ,
double   
) const
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
void Interface< dim, spacedim, n_components, LAC >::fix_solution_dot_derivative ( FEValuesCache< dim, spacedim > &  fe_cache,
Sdouble  alpha 
) const
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
void Interface< dim, spacedim, n_components, LAC >::fix_solution_dot_derivative ( FEValuesCache< dim, spacedim > &  fe_cache,
SSdouble  alpha 
) const
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
const Table<2,DoFTools::Coupling>& Interface< dim, spacedim, n_components, LAC >::get_aux_matrix_coupling ( const unsigned int &  i) const
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
UpdateFlags Interface< dim, spacedim, n_components, LAC >::get_aux_matrix_flags ( const unsigned int &  i) const
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::get_aux_matrix_residuals ( const typename DoFHandler< dim, spacedim >::active_cell_iterator cell,
Scratch &  scratch,
CopyPreconditioner &  data,
std::vector< std::vector< double > > &  local_residuals 
) const
virtual
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::get_aux_matrix_residuals ( const typename DoFHandler< dim, spacedim >::active_cell_iterator cell,
Scratch &  scratch,
CopyPreconditioner &  data,
std::vector< std::vector< Sdouble > > &  local_residuals 
) const
virtual
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
const std::vector<unsigned int> Interface< dim, spacedim, n_components, LAC >::get_differential_blocks ( ) const
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual UpdateFlags Interface< dim, spacedim, n_components, LAC >::get_face_flags ( ) const
virtual
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual UpdateFlags Interface< dim, spacedim, n_components, LAC >::get_jacobian_flags ( ) const
virtual
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual UpdateFlags Interface< dim, spacedim, n_components, LAC >::get_jacobian_preconditioner_flags ( ) const
virtual
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual shared_ptr<Mapping<dim,spacedim> > Interface< dim, spacedim, n_components, LAC >::get_mapping ( const DoFHandler< dim, spacedim > &  ,
const typename LAC::VectorType  
) const
virtual
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual unsigned int Interface< dim, spacedim, n_components, LAC >::get_number_of_aux_matrices ( ) const
virtual
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::get_preconditioner_energy ( const typename DoFHandler< dim, spacedim >::active_cell_iterator ,
Scratch &  ,
CopySystem &  ,
Sdouble &   
) const
virtual

Build the energy needed to get the preconditioner in the case it is required just one derivative.

This function is used to build the energy associated to the preconditioner in the case it is required just one derivative. It takes as argument a reference to the active cell (DoFHandler<dim,spacedim>::active_cell_iterator), all the informations of the system such as fe values, quadrature points, AnyData (Scratch), all the informations related to the PDE (CopySystem) and the energy (Sdouble)

Reimplemented in ConservativeInterface< dim, dim, 1, HeatEquation< dim, LAC >, LAC >, and ConservativeInterface< dim, spacedim, dim+2, FreeSwellingThreeFields< dim, spacedim >, LAC >.

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::get_preconditioner_energy ( const typename DoFHandler< dim, spacedim >::active_cell_iterator ,
Scratch &  ,
CopyPreconditioner &  ,
SSdouble &   
) const
virtual

Build the energy needed to get the preconditioner in the case it is required two derivatives.

This function is used to build the energy associated to the preconditioner in the case the Jacobian is automatically constructed using the derivative of the residual. It takes as argument a reference to the active cell (DoFHandler<dim,spacedim>::active_cell_iterator), all the informations of the system such as fe values, quadrature points, AnyData (Scratch), all the informations related to the PDE (CopySystem) and the energy (SSdouble)

Reimplemented in ConservativeInterface< dim, dim, 1, HeatEquation< dim, LAC >, LAC >, and ConservativeInterface< dim, spacedim, dim+2, FreeSwellingThreeFields< dim, spacedim >, LAC >.

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::get_preconditioner_residual ( const typename DoFHandler< dim, spacedim >::active_cell_iterator cell,
Scratch &  scratch,
CopyPreconditioner &  data,
std::vector< Sdouble > &  local_residual 
) const
virtual

Build the residual needed to get the preconditioner matrix in the case two derivatives are required.

This function is used to build the residual associated to the preconditioner in the case it is required just one derivatice. It takes as argument a reference to the active cell cell, all the informations of the system scratch ( fe values, quadrature points, AnyData ), all the informations related to the PDE data and a reference to the local residual local_residual.

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual UpdateFlags Interface< dim, spacedim, n_components, LAC >::get_residual_flags ( ) const
virtual
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::get_system_energy ( const typename DoFHandler< dim, spacedim >::active_cell_iterator ,
Scratch &  ,
CopySystem &  ,
Sdouble &   
) const
virtual

Build the energy needed to get the system matrix in the case it is required two derivatives.

This function is used to build the energy associated to the system matrix in the case two derivatives are required. It takes as argument a reference to the active cell (DoFHandler<dim,spacedim>::active_cell_iterator), all the informations of the system such as fe values, quadrature points, AnyData (Scratch), all the informations related to the PDE (CopySystem) and the energy (SSdouble)

Reimplemented in ConservativeInterface< dim, dim, 1, HeatEquation< dim, LAC >, LAC >, and ConservativeInterface< dim, spacedim, dim+2, FreeSwellingThreeFields< dim, spacedim >, LAC >.

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::get_system_energy ( const typename DoFHandler< dim, spacedim >::active_cell_iterator ,
Scratch &  ,
CopySystem &  ,
SSdouble &   
) const
virtual

Build the energy needed to get the system matrix in the case it is required two derivatives.

This function is used to build the energy associated to the system matrix in the case two derivatives are required. It takes as argument a reference to the active cell (DFHandler<dim,spacedim>::active_cell_iterator), all the informations of the system such as fe values, quadrature points, AnyData (Scratch), all the informations related to the PDE (CopySystem) and the energy (SSdouble)

Reimplemented in ConservativeInterface< dim, dim, 1, HeatEquation< dim, LAC >, LAC >, and ConservativeInterface< dim, spacedim, dim+2, FreeSwellingThreeFields< dim, spacedim >, LAC >.

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::get_system_residual ( const typename DoFHandler< dim, spacedim >::active_cell_iterator cell,
Scratch &  scratch,
CopySystem &  data,
std::vector< Sdouble > &  local_residual 
) const
virtual

Build the residual needed to get the system matrix in the case it is required two derivatives.

This function is used to build the residual associated to the system in the case two derivatives are required. It takes as argument a reference to the active cell cell, all the informations of the system scratch ( fe values, quadrature points, AnyData ), all the informations related to the PDE data and a reference to the local residual local_residual.

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::get_system_residual ( const typename DoFHandler< dim, spacedim >::active_cell_iterator cell,
Scratch &  scratch,
CopySystem &  data,
std::vector< double > &  local_residual 
) const
virtual

Build the residual needed to get the system matrix in the case it is required just one derivative.

This function is used to build the residual associated to the system in the case it is required just one derivatice. It takes as argument a reference to the active cell cell, all the informations of the system scratch ( fe values, quadrature points, AnyData ), all the informations related to the PDE data and a reference to the local residual local_residual.

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::initialize_data ( const typename LAC::VectorType solution,
const typename LAC::VectorType solution_dot,
const double  t,
const double  alpha 
) const
virtual

Initialize all data required for the system.

This function is used to initialize the varibale AnyData d that contains all data of the problem (solutions, DoF, quadrature points, solutions vector, etc ). It takes as argument the number of DoF per cell dofs_per_cell, the number of quadrature points n_q_points, the number of quadrature points per face n_face_q_points, the reference to solutions vectors sol and the reference to the AnyData d.

TODO: add current_time and current_alpha

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::postprocess_newly_created_triangulation ( Triangulation< dim, spacedim > &  tria) const
virtual

This function is used to modify triangulation using boundary_id or manifold_id.

In the case a signal is required, this is the function to modify.

template<int dim, int spacedim, int n_components, typename LAC >
template<typename Number >
void Interface< dim, spacedim, n_components, LAC >::reinit ( const Number &  alpha,
const typename DoFHandler< dim, spacedim >::active_cell_iterator cell,
FEValuesCache< dim, spacedim > &  fe_cache 
) const
template<int dim, int spacedim, int n_components, typename LAC >
template<typename Number >
void Interface< dim, spacedim, n_components, LAC >::reinit ( const Number &  alpha,
const typename DoFHandler< dim, spacedim >::active_cell_iterator cell,
const unsigned int  face_no,
FEValuesCache< dim, spacedim > &  fe_cache 
) const
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
virtual void Interface< dim, spacedim, n_components, LAC >::set_time ( const double &  t) const
virtual

update time of all parsed mapped functions

Member Data Documentation

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
std::vector<unsigned int> Interface< dim, spacedim, n_components, LAC >::_diff_comp
protected
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
double Interface< dim, spacedim, n_components, LAC >::alpha
mutableprotected
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
std::vector<Table<2,DoFTools::Coupling> > Interface< dim, spacedim, n_components, LAC >::aux_matrix_coupling
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
std::vector<UpdateFlags> Interface< dim, spacedim, n_components, LAC >::aux_matrix_update_flags
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
ParsedDirichletBCs<dim,spacedim,n_components> Interface< dim, spacedim, n_components, LAC >::dirichlet_bcs
mutableprotected
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
unsigned int Interface< dim, spacedim, n_components, LAC >::dofs_per_cell
mutableprotected
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
ParsedMappedFunctions<spacedim,n_components> Interface< dim, spacedim, n_components, LAC >::forcing_terms
mutableprotected
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
unsigned int Interface< dim, spacedim, n_components, LAC >::n_face_q_points
mutableprotected
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
unsigned int Interface< dim, spacedim, n_components, LAC >::n_q_points
mutableprotected
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
ParsedMappedFunctions<spacedim,n_components> Interface< dim, spacedim, n_components, LAC >::neumann_bcs
mutableprotected
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
LAC::VectorType Interface< dim, spacedim, n_components, LAC >::old_solution
mutableprotected

Solution vector evaluated at time t-dt.

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
double Interface< dim, spacedim, n_components, LAC >::old_t
mutableprotected

Previous time step.

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
const LAC::VectorType* Interface< dim, spacedim, n_components, LAC >::solution
mutableprotected

Solution vector evaluated at time t.

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
const LAC::VectorType* Interface< dim, spacedim, n_components, LAC >::solution_dot
mutableprotected

Time derivative solution vector evaluated at time t.

template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
std::string Interface< dim, spacedim, n_components, LAC >::str_diff_comp
protected
template<int dim, int spacedim = dim, int n_components = 1, typename LAC = LATrilinos>
double Interface< dim, spacedim, n_components, LAC >::t
mutableprotected

Current time step.


The documentation for this class was generated from the following file: