Flow123d  release_2.2.0-37-g336ee74
first_order_reaction_base.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file first_order_reaction_base.cc
15  * @brief
16  */
17 
20 
24 
25 
26 #include "system/global_defs.h"
27 #include "system/sys_profiler.hh"
28 
29 #include "mesh/mesh.h"
30 #include "la/distribution.hh"
31 #include "input/accessors.hh"
32 
33 // FLOW123D_FORCE_LINK_IN_PARENT(padeApproximant);
34 FLOW123D_FORCE_LINK_IN_PARENT(linearODEAnalytic);
35 
36 
37 using namespace Input::Type;
38 
39 
41  : ReactionTerm(init_mesh, in_rec)
42 {
44  linear_ode_solver_ = std::make_shared<PadeApproximant>(*num_it);
45 }
46 
48 {
49 }
50 
52 {
53  OLD_ASSERT(distribution_ != nullptr, "Distribution has not been set yet.\n");
54  OLD_ASSERT(time_ != nullptr, "Time governor has not been set yet.\n");
56 
59 
60  // allocation
61  prev_conc_.resize(n_substances_);
62  reaction_matrix_.resize(n_substances_, n_substances_);
63  molar_matrix_.resize(n_substances_, n_substances_);
64  molar_mat_inverse_.resize(n_substances_, n_substances_);
65 
66  // initialize diagonal matrices with molar masses
67  molar_matrix_.zeros();
68  molar_mat_inverse_.zeros();
69  for (unsigned int i=0; i<n_substances_; ++i)
70  {
71  molar_matrix_(i,i) = substances_[i].molar_mass();
72  molar_mat_inverse_(i,i) = 1./substances_[i].molar_mass();
73  }
74 }
75 
76 
78 {
79  OLD_ASSERT(distribution_ != nullptr, "Distribution has not been set yet.\n");
80  OLD_ASSERT(time_ != nullptr, "Time governor has not been set yet.\n");
82 
84  // make scaling that takes into account different molar masses of substances
86 
87  linear_ode_solver_->set_system_matrix(reaction_matrix_);
88 }
89 
90 
91 double **FirstOrderReactionBase::compute_reaction(double **concentrations, int loc_el) //multiplication of concentrations array by reaction matrix
92 {
93  unsigned int rows; // row in the concentration matrix, regards the substance index
94  arma::vec new_conc;
95 
96  // save previous concentrations to column vector
97  for(rows = 0; rows < n_substances_; rows++)
98  prev_conc_(rows) = concentrations[rows][loc_el];
99 
100  // compute new concetrations R*c
101  linear_ode_solver_->update_solution(prev_conc_, new_conc);
102 
103  // save new concentrations to the concentration matrix
104  for(rows = 0; rows < n_substances_; rows++)
105  concentrations[rows][loc_el] = new_conc(rows);
106 
107  return concentrations;
108 }
109 
111 {
112  //DebugOut() << "FirstOrderReactionBases - update solution\n";
113  if(time_->is_changed_dt())
114  {
115  linear_ode_solver_->set_step(time_->dt());
116  }
117 
118  START_TIMER("linear reaction step");
119 
120  for (unsigned int loc_el = 0; loc_el < distribution_->lsize(); loc_el++)
122 
123  END_TIMER("linear reaction step");
124 }
125 
126 
127 unsigned int FirstOrderReactionBase::find_subst_name(const string &name)
128 {
129  unsigned int k=0;
130  for(; k < n_substances_; k++)
131  if (name == substances_[k].name()) return k;
132 
133  return k;
134 }
135 
136 
138 {
139  if (!linear_ode_solver_->evaluate_time_constraint(time_constraint)) return false;
140 
141  DebugOut().fmt("CFL constraint(first order reaction): {}.\n", time_constraint);
142 
143  return true;
144 }
virtual void initialize_from_input()=0
Initializes private members of sorption from the input record.
arma::mat molar_matrix_
Diagonal matrix with molar masses of substances.
double ** concentration_matrix_
void initialize() override
Prepares the object to usage.
void zero_time_step() override
Moves the model to zero time.
SubstanceList substances_
Names belonging to substances.
Definition: mesh.h:97
Iterator< Ret > find(const string &key) const
unsigned int find_subst_name(const std::string &name)
virtual void assemble_ode_matrix(void)=0
Assembles the matrix of the ODEs.
unsigned int size() const
Definition: substance.hh:87
Class ReactionTerm is an abstract class representing reaction term in transport.
#define OLD_ASSERT(...)
Definition: global_defs.h:131
void update_solution(void) override
Updates the solution.
Global macros to enhance readability and debugging, general constants.
arma::mat reaction_matrix_
Reaction matrix.
Distribution * distribution_
Pointer to reference to distribution of elements between processors.
FirstOrderReactionBase(Mesh &init_mesh, Input::Record in_rec)
Constructor.
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
bool is_changed_dt() const
#define START_TIMER(tag)
Starts a timer with specified tag.
#define OLD_ASSERT_LESS(a, b)
Definition: global_defs.h:134
Support classes for parallel programing.
arma::vec prev_conc_
Column vector storing previous concetrations on an element.
Input::Record input_record_
Definition: equation.hh:225
double dt() const
arma::mat molar_mat_inverse_
Inverse of molar_matrix_.
FLOW123D_FORCE_LINK_IN_PARENT(linearODEAnalytic)
#define END_TIMER(tag)
Ends a timer with specified tag.
unsigned int n_substances_
Number of all transported substances. It is the dimension of the reaction matrix. ...
virtual double ** compute_reaction(double **concentrations, int loc_el) override
Computes the reaction on a specified element.
#define DebugOut()
Macro defining &#39;debug&#39; record of log.
Definition: logger.hh:242
~FirstOrderReactionBase(void)
Destructor.
std::shared_ptr< PadeApproximant > linear_ode_solver_
bool evaluate_time_constraint(double &time_constraint) override
Computes a constraint for time step.
TimeGovernor * time_
Definition: equation.hh:224
unsigned int lsize(int proc) const
get local size