OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
response_water_sorption.h
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 //
3 // FCST: Fuel Cell Simulation Toolbox
4 //
5 // Copyright (C) 2014 by Energy Systems Design Laboratory, University of Alberta
6 //
7 // This software is distributed under the MIT License
8 // For more information, see the README file in /doc/LICENSE
9 //
10 // - Class: response_water_sorption.h
11 // - Description: This file defines response class computing amount of water sorbed in the catalyst layers.
12 // - Developers: Marc Secanell Gallart and Madhur Bhaiya
13 // - $Id:
14 //
15 // ----------------------------------------------------------------------------
16 
17 #ifndef _FUELCELLSHOP__RESPONSE_WATER_SORPTION_H
18 #define _FUELCELLSHOP__RESPONSE_WATER_SORPTION_H
19 
20 // Include deal.II classes
21 #include <base/parameter_handler.h>
22 #include <base/point.h>
23 #include <base/function.h>
24 #include <lac/vector.h>
25 #include <fe/fe_values.h>
26 
27 //Include STL
28 #include <cmath>
29 #include <iostream>
30 #include <exception> // std::exception
31 
32 // Include OpenFCST routines:
35 
38 
39 #include "layers/base_layer.h"
40 #include "layers/catalyst_layer.h"
41 
43 
44 using namespace dealii;
45 
46 namespace FuelCellShop
47 {
48 
49  namespace PostProcessing
50  {
51 
119  template <int dim>
121  {
122  public:
124 
127  :
128  BaseResponse<dim>(sm),
129  sorption_source(sst)
130  {}
131 
133 
137  void initialize(ParameterHandler& param)
138  {
139  if ( this->system_management->solution_in_userlist("membrane_water_content") )
140  {
141  lambda.solution_index = this->system_management->solution_name_to_index("membrane_water_content");
142  lambda.fetype_index = this->system_management->block_info->base_element[lambda.solution_index];
143  lambda.indices_exist = true;
144  }
145  else
146  throw std::runtime_error("membrane_water_content variable required for WaterSorptionResponse");
147 
148  if ( this->system_management->solution_in_userlist("water_molar_fraction") )
149  {
150  x_w.solution_index = this->system_management->solution_name_to_index("water_molar_fraction");
151  x_w.fetype_index = this->system_management->block_info->base_element[x_w.solution_index];
152  x_w.indices_exist = true;
153  }
154  else
155  throw std::runtime_error("water_molar_fraction variable required for WaterSorptionResponse");
156 
157  if ( this->system_management->solution_in_userlist("temperature_of_REV") )
158  {
159  tRev.solution_index = this->system_management->solution_name_to_index("temperature_of_REV");
160  tRev.fetype_index = this->system_management->block_info->base_element[tRev.solution_index];
161  tRev.indices_exist = true;
162  }
163 
164  time_k = sorption_source->get_time_constant();
165  }
167 
169 
170 
186  std::map<FuelCellShop::PostProcessing::ResponsesNames, double>& respMap) const
187  {
188  respMap.clear();
189 
190  // Make sure you are in a CL
191  const std::type_info& base_layer = layer->get_base_type();
192  const std::type_info& CatalystLayer = typeid(FuelCellShop::Layer::CatalystLayer<dim>);
193  AssertThrow(base_layer == CatalystLayer,
194  ExcMessage("WaterSorptionResponse can only be used with a CatalystLayer object"));
195 
196  // Get number of quadrature points in the cell:
197  unsigned int n_q_points_cell = (info.fe(lambda.fetype_index)).n_quadrature_points;
198 
199  // Find where the solution is stored in info:
200  unsigned int solIndex = info.global_data->find_vector("Solution");
201 
202  // Create CL:
204 
205  double rhoDry = catalyst_layer->get_electrolyte()->get_density();
206  double EW = catalyst_layer->get_electrolyte()->get_EW();
207  catalyst_layer->get_electrolyte()->set_water_molar_fraction( FuelCellShop::SolutionVariable(&info.values[solIndex][x_w.solution_index], water_molar_fraction) );
208  if (tRev.indices_exist)
209  catalyst_layer->get_electrolyte()->set_temperature( FuelCellShop::SolutionVariable(&info.values[solIndex][tRev.solution_index], temperature_of_REV) );
210 
211  std::vector<double> lambdaValue( info.values[solIndex][lambda.solution_index] );
212 
213  std::vector<double> lambdaEq(n_q_points_cell, 0.0);
214  catalyst_layer->get_electrolyte()->sorption_isotherm(lambdaEq);
215 
216  for (unsigned int q = 0; q < n_q_points_cell; ++q){
217  double JxW = info.fe(lambda.fetype_index).JxW(q);
218  respMap[FuelCellShop::PostProcessing::ResponsesNames::sorbed_water] += ( ( ((time_k*rhoDry)/EW)*(lambdaEq[q] - lambdaValue[q]) ) * JxW );
219  }
220  }
226  void compute_responses(std::vector< FuelCellShop::SolutionVariable > solution_variables,
227  const typename DoFApplication<dim>::CellInfo& info,
229  std::map<FuelCellShop::PostProcessing::ResponsesNames, double>& resp) const
230  {
231  throw std::runtime_error("WaterSorptionResponse::compute_responses(solution_variables, info, layer, resp) not implemented");
232  }
233  //
234  private:
239 
244  double time_k;
245 
251 
257 
263  };
264 
265  }
266 }
267 
268 #endif
WaterSorptionResponse(const FuelCell::SystemManagement &sm, const FuelCellShop::Equation::SorptionSourceTerms< dim > *sst)
Definition: response_water_sorption.h:125
const unsigned int dim
Definition: fcst_constants.h:24
This class assembles source terms corresponding to sorption/desorption of water inside the catalyst l...
Definition: sorption_source_terms.h:103
double get_EW() const
Get Equivalent Weight (grams of dry polymer electrolyte per moles of ) of the polymer electrolyte mat...
Definition: polymer_electrolyte_material_base.h:448
FuelCellShop::Equation::VariableInfo x_w
VariableInfo structure corresponding to the &quot;water_molar_fraction&quot;.
Definition: response_water_sorption.h:256
SmartPointer< const FEVectors > global_data
The smart pointer to the FEVectors object called global_data.
Definition: mesh_loop_info_objects.h:780
This structure is used to encapsulate data from constant values and variable solution data that is us...
Definition: fcst_variables.h:86
void set_water_molar_fraction(const FuelCellShop::SolutionVariable &x_in)
Set the solution variable, water vapor molar fraction .
Definition: polymer_electrolyte_material_base.h:526
Virtual class used to develop a common interface to a set of functions used to evaluate functionals t...
Definition: base_response.h:129
virtual void sorption_isotherm(std::vector< double > &) const
Compute the equilibrium water content, , inside the polymer electrolyte for vapor-equilibriated case...
Definition: polymer_electrolyte_material_base.h:142
Definition: system_management.h:75
double get_density() const
Get the density [gm/cm^3] of the dry polymer electrolyte material.
Definition: polymer_electrolyte_material_base.h:430
Class used to calculate the amount of water sorbed inside the catalyst layer.
Definition: response_water_sorption.h:120
This class is created for the objects handed to the mesh loops.
Definition: mesh_loop_info_objects.h:625
void initialize(ParameterHandler &param)
Initialize class parameters.
Definition: response_water_sorption.h:137
virtual const std::type_info & get_base_type() const
This member function return the name of the type of layer, i.e.
Definition: base_layer.h:177
void compute_responses(const typename DoFApplication< dim >::CellInfo &info, FuelCellShop::Layer::BaseLayer< dim > *const layer, std::map< FuelCellShop::PostProcessing::ResponsesNames, double > &respMap) const
This member function computes the water adsorbed/desorbed from the electrolyte in the catalyst layer...
Definition: response_water_sorption.h:184
virtual FuelCellShop::Material::PolymerElectrolyteBase * get_electrolyte() const
Method to provide access to pointer of the electrolyte object of the catalyst layer.
Definition: catalyst_layer.h:809
Definition: system_management.h:67
void compute_responses(std::vector< FuelCellShop::SolutionVariable > solution_variables, const typename DoFApplication< dim >::CellInfo &info, FuelCellShop::Layer::BaseLayer< dim > *const layer, std::map< FuelCellShop::PostProcessing::ResponsesNames, double > &resp) const
Routine used in order to compute the response with a modified solution (not the one stored in info) ...
Definition: response_water_sorption.h:226
double time_k
Rate of sorption/desorption.
Definition: response_water_sorption.h:244
const FuelCellShop::Equation::SorptionSourceTerms< dim > * sorption_source
Pointer to SorptionSourceTerms object.
Definition: response_water_sorption.h:238
Definition: base_response.h:60
~WaterSorptionResponse()
Definition: response_water_sorption.h:132
IMPORTANT: Add all new solution variables and equations here !
Definition: system_management.h:271
This simple structure stores certain information regarding a particular variable for the equation (al...
Definition: equation_base.h:121
FuelCellShop::Equation::VariableInfo lambda
VariableInfo structure corresponding to the &quot;membrane_water_content&quot;.
Definition: response_water_sorption.h:250
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
void set_temperature(const FuelCellShop::SolutionVariable &T_in)
Set the solution variable, temperature [Kelvin].
Definition: polymer_electrolyte_material_base.h:509
Virtual class used to characterize a generic layer interface.
Definition: base_layer.h:58
FuelCellShop::Equation::VariableInfo tRev
VariableInfo structure corresponding to the &quot;temperature_of_REV&quot;.
Definition: response_water_sorption.h:262
const FEVALUESBASE & fe() const
Access to a single actual FEVALUES object.
Definition: mesh_loop_info_objects.h:1105
Virtual class used to provide the interface for all CatalystLayer children.
Definition: catalyst_layer.h:130