OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Private Attributes | List of all members
FuelCell::Application::AppLaplace< dim > Class Template Reference

This class is used to develop and test new applications. More...

#include <app_laplace.h>

Inheritance diagram for FuelCell::Application::AppLaplace< dim >:
Inheritance graph
[legend]
Collaboration diagram for FuelCell::Application::AppLaplace< dim >:
Collaboration graph
[legend]

Public Member Functions

 AppLaplace (boost::shared_ptr< AppFrame::ApplicationData > data=boost::shared_ptr< AppFrame::ApplicationData >())
 Constructor.
 
 ~AppLaplace ()
 Destructor.
 
virtual void declare_parameters (ParameterHandler &param)
 Declare all parameters that are needed for:
 
virtual void set_parameters (const std::vector< std::string > &name_dvar, const std::vector< double > &value_dvar, ParameterHandler &param)
 Function called by optimization loop in order to set the values in the ParameterHandler to the new design parameters.
 
void _initialize (ParameterHandler &param)
 Set up how many equations are needed and read in parameters for the parameter handler in order to initialize data.
 
virtual void initialize (ParameterHandler &param)
 Call the other initialize routines from the inherited classes.
 
void init_solution (AppFrame::FEVector &src)
 Initialize nonlinear solution Note that this routine is extremely important since it is used to resize the solution vector.
 
virtual void cell_matrix (MatrixVector &cell_matrices, const typename DoFApplication< dim >::CellInfo &cell)
 Integration of local bilinear form.
 
virtual void cell_residual (AppFrame::FEVector &cell_vector, const typename DoFApplication< dim >::CellInfo &cell)
 Integration of the rhs of the equations.
 
virtual void dirichlet_bc (std::map< unsigned int, double > &boundary_values) const
 Integration of local bilinear form.
 
virtual void solve (AppFrame::FEVector &start, const AppFrame::FEVectors &rhs)
 Solve the linear system of equations.
 
virtual double estimate (const AppFrame::FEVectors &sol)
 Estimate error per cell.
 
virtual void cell_responses (std::vector< double > &resp, const typename DoFApplication< dim >::CellInfo &info, const AppFrame::FEVector &sol)
 This class is called by responses to make sure that all responses requested are implemented in either cell_responses, global_responses or face_responses.
 
virtual void global_responses (std::vector< double > &resp, const AppFrame::FEVector &sol)
 This class is used to evaluate all responses that do not require looping over cells.
 
virtual double evaluate (const AppFrame::FEVectors &src)
 This class is used to evaluate the derivative of all the functionals that require looping over cells with respect to the design variables.
 
- Public Member Functions inherited from AppFrame::OptimizationBlockMatrixApplication< dim >
 OptimizationBlockMatrixApplication (boost::shared_ptr< AppFrame::ApplicationData > data=boost::shared_ptr< AppFrame::ApplicationData >())
 Constructor for an object, owning its own mesh and dof handler.
 
 OptimizationBlockMatrixApplication (AppFrame::DoFApplication< dim > &, bool triangulation_only)
 Constructor for an object, borrowing mesh and dof handler from another object.
 
 ~OptimizationBlockMatrixApplication ()
 Destructor.
 
virtual void responses (std::vector< double > &f, const AppFrame::FEVectors &vectors)
 Post-processing.
 
virtual void check_responses ()
 This class is called by responses to make sure that all responses requested are implemented in either cell_responses, global_responses or face_responses.
 
virtual void print_responses (std::vector< double > &resp)
 This function is used to print the responses to screen (deallog)
 
virtual void dresponses_dl (std::vector< std::vector< double > > &df_dl, const AppFrame::FEVectors &src)
 Post-processing.
 
virtual void cell_dresponses_dl (std::vector< std::vector< double > > &cell_df_dl, const typename DoFApplication< dim >::CellInfo &info, const AppFrame::FEVector &src)
 This class is used to evaluate the derivative of all the functionals that require looping over cells with respect to the design variables.
 
virtual void global_dresponses_dl (std::vector< std::vector< double > > &df_dl, const AppFrame::FEVector &src)
 This class is used to evaluate the sensitivities of all responses that do not require looping over cells with respect of the design variables.
 
virtual void dresponses_du (std::vector< AppFrame::FEVector > &dst, const AppFrame::FEVectors &src)
 Loop over all cells and assemble the vector

\[ \frac{\partial f_m}{\partial u_i} \]

that is used to solve the sensitivity equations by using the local matrix integration functions cell_dfunctional_du(), bdry_dfunctional_dlu() and edge_dfunctinal_du() provided by the derived class.

 
virtual void cell_dresponses_du (std::vector< AppFrame::FEVector > &df_du, const typename DoFApplication< dim >::CellInfo &info, std::vector< std::vector< double > > &src)
 Integration of local bilinear form.
 
virtual void global_dresponses_du (std::vector< AppFrame::FEVector > &df_du, const AppFrame::FEVector &src)
 This class is used to evaluate the sensitivities of all responses that do not require looping over cells with respect of the design variables.
 
virtual void dresidual_dlambda (std::vector< AppFrame::FEVector > &dst, const AppFrame::FEVectors &src)
 Loop over all cells and assemble the vector

\[ \frac{\partial R_n}{\partial \lambda_k} \]

that is used to solve the sensitivity equations by using the local matrix integration functions cell_dresidual_dlambda(), bdry_dresidual_dlambda() and edge_dresidual_dlambda() provided by the derived class.

 
virtual void cell_dresidual_dlambda (std::vector< AppFrame::FEVector > &cell_vector, const typename DoFApplication< dim >::CellInfo &cell, std::vector< std::vector< double > > &src)
 Integration of local bilinear form.
 
void solve_direct (std::vector< std::vector< double > > &df_dl, const AppFrame::FEVectors &sol)
 Solver in order to obtain the analytical sensitivities for the given objective and constraints using the direct method.
 
void solve_adjoint (std::vector< std::vector< double > > &df_dl, const AppFrame::FEVector &sol)
 Solver in order to obtain the analytical sensitivities for the given objective and constraints using the adjoint method.
 
unsigned int get_n_resp () const
 Member function that returns the number of responses.
 
unsigned int get_n_dvar () const
 Member function that returns the number of design variables.
 
std::vector< std::string > get_name_dvar () const
 Member function that returns the name of the design variables.
 
std::vector< std::string > get_name_responses () const
 Member function that returns the name of the responses.
 
void print_default_parameter_file ()
 Print the default parameter handler file.
 
bool get_bool_transfer_solution ()
 Function to see if we are transferring a solution on a refined grid to the initial coarse grid.
 
virtual void set_read_initial_solution (bool &read)
 
void set_optimization_parameters (unsigned int &n_dvar, unsigned int &n_resp, std::vector< std::string > &name_design_var, std::vector< std::string > &name_responses)
 This routine assigns the value of n_dvar, n_resp, name_design_var and name_responses to the internal variables n_dvar, n_resp, name_design_var and name_response.
 
void set_output_variables (std::vector< std::string > &dakota_name_responses)
 
- Public Member Functions inherited from AppFrame::BlockMatrixApplication< dim >
void do_assemble (const FEVectors &, typename DoFHandler< dim >::active_cell_iterator begin, typename DoFHandler< dim >::active_cell_iterator end, Threads::Mutex &mutex)
 Multithreaded assemble function called by assemble().
 
virtual void post_cell_assemble ()
 A call back function for assemble(), called after all cell matrices have been entered, but before the face matrices are computed.
 
 BlockMatrixApplication (boost::shared_ptr< ApplicationData > data=boost::shared_ptr< ApplicationData >())
 Constructor for an object, owning its own mesh and dof handler.
 
 BlockMatrixApplication (DoFApplication< dim > &, bool triangulation_only)
 Constructor for an object, borrowing mesh and dof handler from another object.
 
void _initialize (ParameterHandler &param)
 Initialize data of this class.
 
void add_vector_for_cell_matrix (std::string name, unsigned int block, unsigned int components, bool values, bool derivatives)
 Add a named vector from data to be handed to the cell_matrix() function via CellInfo::values and CellInfo::derivatives.
 
void add_vector_for_flux_matrix (std::string name, unsigned int block, unsigned int components, bool values, bool derivatives)
 Add a named vector from data to be handed to the bdry_matrix() and face_matrix() functions via FaceInfo::values and FaceInfo::derivatives.
 
void remesh_matrices ()
 Initialize sparsity patterns and matrices for the new mesh.
 
virtual void remesh ()
 Refine grid accordingly to the Refinement options under "Grid Generation" in the parameter file.
 
void assemble (const FEVectors &)
 Loop over all cells and assemble the system #matrices by using the local matrix integration functions cell_matrix(), bdry_matrix() and face_matrix() provided by the derived class.
 
void assemble_numerically (const FEVectors &src, const double delta=1e-6)
 Compute the Jacobian of the system of equations,

\[ J(i,j) = \frac{\partial R_i}{\partial u_j} \]

by using forward differences.

 
void residual_constraints (FEVector &dst) const
 Redefinition of residual_constraints() in DoFHandler.
 
virtual void bdry_matrix (MatrixVector &face_matrices, const typename DoFApplication< dim >::FaceInfo &face)
 Integration of local bilinear form.
 
virtual void face_matrix (MatrixVector &matrices11, MatrixVector &matrices12, MatrixVector &matrices21, MatrixVector &matrices22, const typename DoFApplication< dim >::FaceInfo &face1, const typename DoFApplication< dim >::FaceInfo &face2)
 Integration of local bilinear form.
 
- Public Member Functions inherited from AppFrame::DoFApplication< dim >
 DoFApplication (boost::shared_ptr< ApplicationData > data=boost::shared_ptr< ApplicationData >(), bool multigrid=false)
 Constructor for an object, owning its own mesh and dof handler.
 
 DoFApplication (bool multigrid)
 Constructor for an object owning its own mesh and dof handler and creating new ApplicationData.
 
 DoFApplication (DoFApplication< dim > &, bool triangulation_only, bool multigrid=false)
 Constructor for an object, borrowing mesh and dof handler from another object.
 
 ~DoFApplication ()
 Destructor which deletes owned objects.
 
void _initialize (ParameterHandler &param)
 Initialize from parameter values.
 
void add_vector_for_transfer (FEVector *)
 Add the vector to be transfered from one mesh to the next.
 
void delete_vector_for_transfer ()
 Delete the vector to be transfered from one mesh to the next.
 
virtual void remesh_dofs ()
 Initialize dof handler, count the dofs in each block and renumber the dofs.
 
virtual void init_vector (FEVector &dst) const
 Reinitialize the BlockVector dst such that it contains block_size.size() blocks.
 
virtual double residual (FEVector &dst, const FEVectors &src, bool apply_boundaries=true)
 Loop over all cells to compute the residual.
 
unsigned int memory_consumption () const
 Compute the amount of memory needed by this object.
 
virtual double cell_estimate (const CellInfo &)
 Integration of estimate inside a cell.
 
virtual double bdry_estimate (const FaceInfo &)
 Integration of estimate on a boundary face.
 
virtual double face_estimate (const FaceInfo &, const FaceInfo &)
 Integration of estimate on an interior face.
 
virtual void cell_residual (FEVector &cell_vector, const CellInfo &cell)
 Integration of local linear form.
 
virtual void bdry_residual (FEVector &face_vector, const FaceInfo &face)
 Integration of local linear form.
 
virtual void face_residual (FEVector &face_vector1, FEVector &face_vector2, const FaceInfo &face1, const FaceInfo &face2)
 Integration of local linear form.
 
template<class BOX , class WORKER >
void integrate_1form (BOX &box, WORKER &local_forms, FEVector &dst, const FEVectors &src)
 The function integrating 1-forms as for instance used by residual().
 
void integrate_functional (LocalEstimate< dim > &local_forms, FEVectors &dst, const FEVectors &src)
 The function integrating functionals using local integrators.
 
void store_triangulation (Triangulation< dim > &new_tr)
 Function to copy a triangulation object for use after refinement.
 
void transfer_solution_to_coarse_mesh (Triangulation< dim > &tr_coarse, AppFrame::FEVector &coarse_solution, AppFrame::FEVector &refined_solution)
 Function to perform the transfer of a solution on a refined grid to the initial coarse grid.
 
void read_init_solution (AppFrame::FEVector &dst, bool &good_solution) const
 Function to transfer a solution on a refined grid to the initial coarse grid.
 
template<class INFO , typename TYPE >
void fill_local_data (const INFO &info, const std::vector< VectorSelector > &data_vector, std::vector< std::vector< std::vector< TYPE > > > &data) const
 Fill local data vector with values of the stored finite element function vectors.
 
virtual void grid_out (const std::string &basename)
 Output the grid used to solve the problem.
 
virtual void data_out (const std::string &basename, const FEVectors &src)
 
virtual void data_out (const std::string &basename, const FEVectors &src, const std::vector< std::string > &solution_names)
 
virtual void data_out (const std::string &basename, const FEVector &solution, const std::vector< std::string > &solution_names, const FEVector &postprocessing=FEVector(), const std::vector< std::string > &postprocessing_names=std::vector< std::string >())
 This function outputs the results of a computation.
 
- Public Member Functions inherited from AppFrame::ApplicationBase
 ApplicationBase (boost::shared_ptr< ApplicationData > ex_data=boost::shared_ptr< ApplicationData >())
 Constructor for an application.
 
 ApplicationBase (const ApplicationBase &other)
 Copy constructor.
 
virtual ~ApplicationBase ()
 Virtual destructor.
 
void print_parameters_to_file (ParameterHandler &param, const std::string &file_name, const ParameterHandler::OutputStyle &style)
 Print default parameters for the application to a file.
 
virtual void start_vector (FEVector &dst, std::string caller) const
 Initialize vector to problem size.
 
virtual void Tsolve (FEVector &start, const FEVectors &rhs)
 Solve the dual system assembled with right hand side rhs and return the result in start.
 
virtual void grid_out (const std::string &filename) const
 Write the mesh in the format specified by the ParameterHandler.
 
virtual void data_out (const std::string &filename, const FEVectors &src, const std::vector< std::string >)
 
boost::shared_ptr
< ApplicationData
get_data ()
 Get access to the protected variable data.
 
const boost::shared_ptr
< ApplicationData
get_data () const
 Get read-only access to the protected variable data.
 
virtual std::string id () const
 Return a unique identification string for this application.
 
virtual void notify (const Event &reason)
 Add a reason for assembling.
 
virtual void clear ()
 Reset the application class such that a call to initialize() is possible again and will produce the same result as if the object was fresh.
 
void clear_events ()
 Clear all notifications.
 

Private Attributes

FuelCell::SystemManagement system_management
 Reimplementation of the routine in the base class BaseApplication in namespace AppFrame so that the right labels are outputed and so that I can compute and output the source term.
 
boost::shared_ptr
< FuelCellShop::Geometry::GridBase
< dim > > 
grid
 Create an object to generate the mesh.
 
FuelCell::OperatingConditions OC
 Operating conditions class object.
 
boost::shared_ptr
< FuelCellShop::Layer::GasDiffusionLayer
< dim > > 
CGDL
 Create an object to generate the mesh.
 
FuelCellShop::Material::Oxygen oxygen
 Create an oxygen gas object.
 
FuelCellShop::Material::Nitrogen nitrogen
 Create an oxygen gas object.
 
FuelCellShop::Equation::FicksTransportEquation
< dim
FicksLaw
 Create equations.
 

Additional Inherited Members

- Public Types inherited from AppFrame::DoFApplication< dim >
typedef
MeshWorker::InfoObjects::IntegrationInfo
< dim, FEValuesBase< dim > > 
CellInfo
 This is a redefinition of CellInfo in order to call MeshWorker objects instead of the previous implementation in AppFrame.
 
typedef
MeshWorker::InfoObjects::IntegrationInfo
< dim, FEFaceValuesBase< dim > > 
FaceInfo
 This is a redefinition of CellInfo in order to call MeshWorker objects instead of the previous implementation in AppFrame.
 
- Public Attributes inherited from AppFrame::DoFApplication< dim >
unsigned int verbosity
 Controls verbosity of certain functions.
 
std::vector
< DataComponentInterpretation::DataComponentInterpretation > 
solution_interpretations
 solution_interpretations identifies whether some solution_names are scalars or parts of a vector.
 
std::vector
< DataComponentInterpretation::DataComponentInterpretation > 
postprocessing_interpretations
 postprocessing_interpretations identifies whether some postprocessing_names are scalars or parts of a vector.
 
bool output_materials_and_levels
 output_materials_and_levels if true then visualized, otherwise suppressed.
 
Vector< double > output_materials
 Vector that will be used by data_out to store the material ids.
 
Vector< double > output_levels
 Vector that will be used by data_out to store the refinement levels at each cell.
 
- Static Public Attributes inherited from AppFrame::OptimizationBlockMatrixApplication< dim >
static const AppFrame::Event sensitivity_analysis
 The Event set by OptimizationMatrixApplication if if we want to compute the sensitivities and a new matrix should be assembled.
 
- Protected Types inherited from AppFrame::OptimizationBlockMatrixApplication< dim >
typedef void(* CELL_Dvalues )(std::vector< AppFrame::FEVector > &, const typename DoFApplication< dim >::CellInfo &, std::vector< std::vector< double > > &)
 Definition of a pointer to a function that is called in a loop over cells to compute the derivatives of either residual or responses with respect to design variables.
 
- Protected Member Functions inherited from AppFrame::OptimizationBlockMatrixApplication< dim >
void print_dresponses_dl (std::vector< std::vector< double > > pdf_pdl)
 Auxiliary routine to print the values df_dl This routine should be called once df_df is assembled.
 
void print_dresponses_du (std::vector< AppFrame::FEVector > df_du)
 Auxiliary routine to print the values of df_du This routine should be called once df_du is assembled.
 
void dfunction (std::vector< AppFrame::FEVector > &dst, const AppFrame::FEVectors &src, bool dfunctional_du, bool dresidual_dlambda)
 This is an auxiliary function called by dresidual_dlambda and dfunctional_du.
 
- Protected Attributes inherited from AppFrame::OptimizationBlockMatrixApplication< dim >
unsigned int n_dvar
 Number of design variables.
 
unsigned int n_obj
 Number of objective functions.
 
unsigned int n_resp
 Number of responses = n_obj + n_con.
 
std::vector< std::string > name_design_var
 Member that stores the name of the design variables.
 
std::vector< std::string > name_responses
 Member that stores the name of the responses, i.e.
 
std::vector< std::string > name_output_var
 Member that stores the name of the output variables, These names will be written to name_responses if optimization is not being used and ignored otherwise.
 
bool output_coarse_solution
 Decision variable for whether the solution is to be output for transfer.
 
bool read_in_initial_solution
 Decision variable for whether the initial solution should be read from file.
 
bool optimization
 Decision variable for whether the application is being used for optimization.
 

Detailed Description

template<int dim>
class FuelCell::Application::AppLaplace< dim >

This class is used to develop and test new applications.

All functions are blank so all will have to be implemented.

Constructor & Destructor Documentation

template<int dim>
FuelCell::Application::AppLaplace< dim >::AppLaplace ( boost::shared_ptr< AppFrame::ApplicationData data = boost::shared_ptr< AppFrame::ApplicationData >())

Constructor.

Note
the pointer data is initialized to boost::shared_ptr<> (), this means that the pointer is empty and when we do data.get() it will return 0. This is good because at ApplicationBase constructor an ApplicationData will be constructed.
template<int dim>
FuelCell::Application::AppLaplace< dim >::~AppLaplace ( )

Destructor.

Member Function Documentation

template<int dim>
void FuelCell::Application::AppLaplace< dim >::_initialize ( ParameterHandler &  param)

Set up how many equations are needed and read in parameters for the parameter handler in order to initialize data.

template<int dim>
virtual void FuelCell::Application::AppLaplace< dim >::cell_matrix ( MatrixVector cell_matrices,
const typename DoFApplication< dim >::CellInfo cell 
)
virtual

Integration of local bilinear form.

Here we loop over the quadrature points and over degrees of freedom in order to compute the matrix for the cell This routine depends on the problem at hand and is called by assemble() in DoF_Handler class The matrix to be assembled is:

\[ \begin{array}{l} M(i,j).block(0) = \int_{\Omega} a \nabla \phi_i \nabla \phi_j d\Omega + \int_{\Omega} \phi_i \frac{\partial f}{\partial u_0}\Big|_n \phi_j d\Omega \\ \end{array} \]

Reimplemented from AppFrame::BlockMatrixApplication< dim >.

template<int dim>
virtual void FuelCell::Application::AppLaplace< dim >::cell_residual ( AppFrame::FEVector cell_vector,
const typename DoFApplication< dim >::CellInfo cell 
)
virtual

Integration of the rhs of the equations.

Here we loop over the quadrature points and over degrees of freedom in order to compute the right hand side for each cell This routine depends on the problem at hand and is called by residual() in DoF_Handler class

Note
This function is called residual because in the case of nonlinear systems the rhs is equivalent to the residual
template<int dim>
virtual void FuelCell::Application::AppLaplace< dim >::cell_responses ( std::vector< double > &  resp,
const typename DoFApplication< dim >::CellInfo info,
const AppFrame::FEVector sol 
)
inlinevirtual

This class is called by responses to make sure that all responses requested are implemented in either cell_responses, global_responses or face_responses.

Note
Every time we add a new response that we can calculate we need to update this file. Compute the value of all objective function and constraints

Reimplemented from AppFrame::OptimizationBlockMatrixApplication< dim >.

template<int dim>
virtual void FuelCell::Application::AppLaplace< dim >::declare_parameters ( ParameterHandler &  param)
virtual

Declare all parameters that are needed for:

  • the computation of the equation coefficients
  • the control of the linear system solution
  • ...

Reimplemented from AppFrame::OptimizationBlockMatrixApplication< dim >.

template<int dim>
virtual void FuelCell::Application::AppLaplace< dim >::dirichlet_bc ( std::map< unsigned int, double > &  boundary_values) const
virtual

Integration of local bilinear form.

Member function used to set dirichlet boundary conditions. This function is application specific and it only computes the boundary_value values that are used to constraint the linear system of equations that is being solved

Reimplemented from AppFrame::BlockMatrixApplication< dim >.

template<int dim>
virtual double FuelCell::Application::AppLaplace< dim >::estimate ( const AppFrame::FEVectors sol)
virtual

Estimate error per cell.

Reimplemented from AppFrame::DoFApplication< dim >.

template<int dim>
virtual double FuelCell::Application::AppLaplace< dim >::evaluate ( const AppFrame::FEVectors src)
virtual

This class is used to evaluate the derivative of all the functionals that require looping over cells with respect to the design variables.

This class is called by responses to evaluate the response at each cell. This class is used to evaluate the sensitivities of all responses that do not require looping over cells with respect of the design variables. An example of one of this types of constraints is the solid volume fraction. This class is used to evaluate the derivative of all the functionals that require looping over cells with respect of the unknowns of the system of governing equations. This class is called by responses to evaluate the response at each cell. This class is used to evaluate the sensitivities of all responses that do not require looping over cells with respecto of the unknowns of the system of governing equations. An example of one of this types of constraints is the solid volume fraction. Post-processing. Evaluate a functional such as the objective function of an optimization problem

Reimplemented from AppFrame::DoFApplication< dim >.

template<int dim>
virtual void FuelCell::Application::AppLaplace< dim >::global_responses ( std::vector< double > &  resp,
const AppFrame::FEVector sol 
)
inlinevirtual

This class is used to evaluate all responses that do not require looping over cells.

An example of one of this types of constraints is the solid volume fraction.

Reimplemented from AppFrame::OptimizationBlockMatrixApplication< dim >.

template<int dim>
void FuelCell::Application::AppLaplace< dim >::init_solution ( AppFrame::FEVector src)
virtual

Initialize nonlinear solution Note that this routine is extremely important since it is used to resize the solution vector.

Therefore, it always has to exist. Otherwise you will receive an error such as: An error occurred in line <306> of file </home/secanell/Programs/fcst/contrib/deal.II/deal.II/include/deal.II/lac/block_indices.h> in function std::pair<unsigned int, unsigned int> dealii::BlockIndices::global_to_local(unsigned int) const The violated condition was: i<total_size() The name and call sequence of the exception was: ExcIndexRange(i, 0, total_size()) Additional Information: Index 0 is not in [0,0[ To initialize the solution vector at least you need the following: src.reinit(this->block_info.global);

Implements AppFrame::OptimizationBlockMatrixApplication< dim >.

template<int dim>
virtual void FuelCell::Application::AppLaplace< dim >::initialize ( ParameterHandler &  param)
virtual

Call the other initialize routines from the inherited classes.

Reimplemented from AppFrame::OptimizationBlockMatrixApplication< dim >.

template<int dim>
virtual void FuelCell::Application::AppLaplace< dim >::set_parameters ( const std::vector< std::string > &  name_dvar,
const std::vector< double > &  value_dvar,
ParameterHandler &  param 
)
virtual

Function called by optimization loop in order to set the values in the ParameterHandler to the new design parameters.

Since ParameterHandler depends on the problem we are solving, set_parameters() is set at the most inner loop of the application.

Reimplemented from AppFrame::OptimizationBlockMatrixApplication< dim >.

template<int dim>
virtual void FuelCell::Application::AppLaplace< dim >::solve ( AppFrame::FEVector start,
const AppFrame::FEVectors rhs 
)
virtual

Solve the linear system of equations.

Reimplemented from AppFrame::ApplicationBase.

Member Data Documentation

template<int dim>
boost::shared_ptr<FuelCellShop::Layer::GasDiffusionLayer<dim> > FuelCell::Application::AppLaplace< dim >::CGDL
private

Create an object to generate the mesh.

Create equations.

template<int dim>
boost::shared_ptr< FuelCellShop::Geometry::GridBase<dim> > FuelCell::Application::AppLaplace< dim >::grid
private

Create an object to generate the mesh.

Create an oxygen gas object.

Operating conditions class object.

Create an oxygen gas object.

template<int dim>
FuelCell::SystemManagement FuelCell::Application::AppLaplace< dim >::system_management
private

Reimplementation of the routine in the base class BaseApplication in namespace AppFrame so that the right labels are outputed and so that I can compute and output the source term.

The object that knows everything about FCST equations, variables, couplings, etc.


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