OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes | List of all members
FuelCell::ApplicationCore::ThetaScheme Class Reference

#include <theta_scheme.h>

Inheritance diagram for FuelCell::ApplicationCore::ThetaScheme:
Inheritance graph
[legend]
Collaboration diagram for FuelCell::ApplicationCore::ThetaScheme:
Collaboration graph
[legend]

Public Member Functions

 ThetaScheme (ApplicationBase &app)
 Constructor, receiving an application with transient model and a data object. More...
 
virtual void declare_parameters (ParameterHandler &param)
 Declare input parameters. More...
 
virtual void initialize (ParameterHandler &param)
 Read parameters. More...
 
virtual void solve (FuelCell::ApplicationCore::FEVector &u, const FuelCell::ApplicationCore::FEVectors &in_vectors)
 Actual transient solver. More...
 
- Public Member Functions inherited from FuelCell::ApplicationCore::TransientBase
 TransientBase (ApplicationBase &app)
 Constructor, receiving an application with transient model and a data object. More...
 
- Public Member Functions inherited from FuelCell::ApplicationCore::ApplicationWrapper
 ApplicationWrapper (ApplicationBase &app)
 Constructor for a derived application. More...
 
 ~ApplicationWrapper ()
 Destructor. More...
 
virtual void remesh ()
 Generate the next mesh depending on the mesh generation parameters. More...
 
virtual void init_vector (FEVector &dst) const
 Initialize vector to problem size. More...
 
virtual double residual (FEVector &dst, const FEVectors &src, bool apply_boundaries=true)
 Compute residual of src and store it into dst. More...
 
virtual void Tsolve (FEVector &dst, const FEVectors &src)
 Solve the dual system assembled with right hand side rhs and return the result in start. More...
 
virtual double estimate (const FEVectors &src)
 Estimate cell-wise errors. More...
 
virtual double evaluate (const FEVectors &src)
 Evaluate a functional. More...
 
virtual void grid_out (const std::string &filename) const
 
virtual void data_out (const std::string &filename, const FEVectors &src)
 Write data in the format specified by the ParameterHandler. More...
 
virtual std::string id () const
 Return a unique identification string for this application. More...
 
virtual void notify (const Event &reason)
 Add a reason for assembling. More...
 
- Public Member Functions inherited from FuelCell::ApplicationCore::ApplicationBase
 ApplicationBase (boost::shared_ptr< ApplicationData > data=boost::shared_ptr< ApplicationData >())
 Constructor for an application. More...
 
 ApplicationBase (const ApplicationBase &other)
 Copy constructor. More...
 
virtual ~ApplicationBase ()
 Virtual destructor. More...
 
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. More...
 
virtual void start_vector (FEVector &dst, std::string) const
 Initialize vector to problem size. More...
 
virtual void grid_out (const std::string &)
 Write the mesh in the format specified by the ParameterHandler. More...
 
boost::shared_ptr
< ApplicationData
get_data ()
 Get access to the protected variable data. More...
 
const boost::shared_ptr
< ApplicationData
get_data () const
 Get read-only access to the protected variable data. More...
 
virtual void clear ()
 All true in notifications. More...
 
virtual void clear_events ()
 All false in notifications. More...
 
virtual unsigned int get_solution_index ()
 Returns solution index. More...
 

Static Public Attributes

static const
FuelCell::ApplicationCore::Event 
new_assembly
 Event used to make application reassemble its system matrix and rhs. More...
 

Private Member Functions

void solve_operations (FuelCell::ApplicationCore::FEVector &u, const FuelCell::ApplicationCore::FEVectors &in_vectors, bool print_solution)
 Performs one transient iteration with the given time step tau. More...
 

Private Attributes

double theta
 Parameter $ \theta \in [0,\,1] $. More...
 
double current_time
 
float current_layer
 Variable that contains current time layer. More...
 
unsigned short int scheme_order
 

Additional Inherited Members

- Protected Member Functions inherited from FuelCell::ApplicationCore::TransientBase
void debug_output (const FEVector &sol)
 Function used to output any necessary information at each time layer for debugging purposes. More...
 
void responses_output ()
 Function used to output boundary responses computed by application. More...
 
void error_output (const double &err)
 Function used to output the solution error on each layer. More...
 
void tau_output ()
 Function used to output time step on each time layer. More...
 
void initialize_u_old (const FEVector u, const double n_ref, const FEVector old_solution, const bool force=true)
 If we are inside the adaptive refinement loop, then we receive the old solution from class since it was transferred to a finer mesh and has now a different size. More...
 
void ats_Richardson (FEVector &u, const FEVector u1, const FEVector u2, int scheme_order, bool &solution_accepted)
 This routine estimates solution error and changes the time step using Richardson extrapolation. More...
 
bool file_exists (const std::string &file)
 
- Protected Member Functions inherited from FuelCell::ApplicationCore::ApplicationWrapper
SmartPointer< ApplicationBaseget_wrapped_application ()
 Gain access to the inner application. More...
 
- Protected Member Functions inherited from FuelCell::ApplicationCore::ApplicationBase
void print_caller_name (const std::string &caller_name) const
 Print caller name. More...
 
- Protected Attributes inherited from FuelCell::ApplicationCore::TransientBase
bool print_last_layer_only
 If this variable is set to true, then only the last time layer data will be printed. More...
 
bool debug_solution
 Print updated solution after each step into file solution_NNNNN.dat. More...
 
unsigned int debug
 Write debug output to FcstUtilities::log; the higher the number, the more output. More...
 
unsigned int step
 Number of the current time layer. More...
 
std::string ats_method
 Decision varible for selecting the adaptive time stepping method. More...
 
double tau
 Initial time step size. More...
 
double max_tau
 Maximal time step size. More...
 
double min_tau
 Minimal time step size. More...
 
int max_ats_it
 Maximal number of adaptive time step iterations used to resolve each time layer. More...
 
FEVector u_old
 A vector to store the solution from the previous time layer in. More...
 
int ar_freq
 If adaptive refinement is selected, then it will be performed at each ar_freq-th layer. More...
 
double total_time
 Total time of simulation in seconds. More...
 
double time
 Current time of simulation. More...
 
double abs_err_estimate
 
double rel_err_estimate
 
double abs_error_tolerance
 
double rel_error_tolerance
 
bool output_error
 
bool output_tau
 
std::ofstream responses_outfile
 Variable used in printing responses data into a file. More...
 
std::ofstream error_outfile
 Variable used in printing solution error data into a file. More...
 
std::ofstream tau_outfile
 Variable used in printing time step into a file. More...
 
const std::string fname_err ="Solution_error.dat"
 File name for solution error output. More...
 
const std::string fname_tau ="Time_step.dat"
 File name for solution error output. More...
 
std::vector< double > all_time_values
 A vector of all time values. More...
 
bool everything_is_fine =true
 
unsigned int n_ref
 
- Protected Attributes inherited from FuelCell::ApplicationCore::ApplicationWrapper
SmartPointer< ApplicationBaseapp
 Pointer to the application this one depends upon. More...
 
- Protected Attributes inherited from FuelCell::ApplicationCore::ApplicationBase
boost::shared_ptr
< ApplicationData
data
 Object for auxiliary data. More...
 
Event notifications
 Accumulate reasons for assembling here. More...
 

Detailed Description

Theta scheme

Theta scheme is a time discretization method for PDEs arising from finite difference approximation of temporal derivative and implicit-explicit treatment of spatial operators and source terms.

In general, a theta scheme for a PDE

\[ u_{,t}+L(u)=f(u,\,t), \]

where $ L(u) $ is a spatial differential operator and $ f(u,\,t) $ is a source term, will look like

\[ \frac{u^{n+1}-u^n}{\tau}+\theta L(u^{n+1}) + (1-\theta) L(u^n) = \theta f(u^{n+1},\,t) + (1-\theta) f(u^n,\,t). \]

Here we used first order approximation of $ u_{,t} $. After finite element method is applied, resulting equation is of the form

\[ [ M + \tau \theta A ] u^{n+1} = [ M + \tau (1-\theta) A ] u^n + \tau \theta f(u^{n+1},\,t) + \tau (1-\theta) f(u^n,\,t), \]

where $ M $ is a mass matrix and $ A $ is a matrix corresponding to operator $ L(u) $. Depending on parameter $ \theta \in [0,\,1]$, one can obtain explicit Euler scheme ( $ \theta = 0 $), Crank-Nicolson scheme ( $ \theta = 0.5 $), implicit Euler scheme ( $ \theta = 1 $), or anything inbetween.

Parameters

Parameters that control theta scheme solver are defined in Theta_Scheme subsection in the data input file. Default parameters are as follows:

subsection Transient
set Total time of simulation, [s] = 10
subsection Theta Scheme
set Theta = 0.5
end
end

Usage

In order to use the class, first create an object passing in a linear application and data application

boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData > data = FuelCell::ApplicationCore::ApplicationData ();
FuelCell::ApplicationCore::ThetaScheme> transient_solver (app_lin,data);

Then, read data and initialize code using a ParameterHandler object:

transient_solver.declare_parameters(param);
// Use ParameterHandler param to read from file ...
transient_solver.initialize(param);

Finally, execute:

transient_solver.solve(u_out, in_vectors)
Author
Aslan Kosakian and Marc Secanell, 2015

Constructor & Destructor Documentation

FuelCell::ApplicationCore::ThetaScheme::ThetaScheme ( ApplicationBase app)

Constructor, receiving an application with transient model and a data object.

Member Function Documentation

virtual void FuelCell::ApplicationCore::ThetaScheme::declare_parameters ( ParameterHandler &  param)
virtual

Declare input parameters.

The following are the parameters used and their default values in section "Transient".

  • Total time of simulation, [s] = 10;
  • Theta = 0.5.

These values are set in parameter file.

Reimplemented from FuelCell::ApplicationCore::TransientBase.

virtual void FuelCell::ApplicationCore::ThetaScheme::initialize ( ParameterHandler &  param)
virtual

Read parameters.

Reimplemented from FuelCell::ApplicationCore::TransientBase.

virtual void FuelCell::ApplicationCore::ThetaScheme::solve ( FuelCell::ApplicationCore::FEVector u,
const FuelCell::ApplicationCore::FEVectors in_vectors 
)
virtual

Actual transient solver.

Performs loop over time layers.

Parameters:

  • u : Output solution at the next time step
  • in_vectors : collection of vectors with information at each DOF.

Implements FuelCell::ApplicationCore::TransientBase.

void FuelCell::ApplicationCore::ThetaScheme::solve_operations ( FuelCell::ApplicationCore::FEVector u,
const FuelCell::ApplicationCore::FEVectors in_vectors,
bool  print_solution 
)
private

Performs one transient iteration with the given time step tau.

Member Data Documentation

float FuelCell::ApplicationCore::ThetaScheme::current_layer
private

Variable that contains current time layer.

Note
The type is float so that partial steps can be shown, e.g. step $ n+1/2 $.
double FuelCell::ApplicationCore::ThetaScheme::current_time
private
const FuelCell::ApplicationCore::Event FuelCell::ApplicationCore::ThetaScheme::new_assembly
static

Event used to make application reassemble its system matrix and rhs.

unsigned short int FuelCell::ApplicationCore::ThetaScheme::scheme_order
private
double FuelCell::ApplicationCore::ThetaScheme::theta
private

Parameter $ \theta \in [0,\,1] $.

Depending on parameter $ \theta \in [0,\,1]$, one can obtain explicit Euler scheme ( $ \theta = 0 $), Crank-Nicolson scheme ( $ \theta = 0.5 $), implicit Euler scheme ( $ \theta = 1 $), or anything inbetween.


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