Flow123d  release_2.2.0-41-g0958a8d
assembly_lmh.hh
Go to the documentation of this file.
1 /*
2  * assembly_lmh.hh
3  *
4  * Created on: Apr 30, 2016
5  * Author: jb
6  */
7 
8 #ifndef SRC_FLOW_ASSEMBLY_LMH_HH_
9 #define SRC_FLOW_ASSEMBLY_LMH_HH_
10 
12 #include "soil_models.hh"
13 
14 
15 
16 
17 /**
18  * Prove of concept for general assembly scheme.
19  * Ideas:
20  * - Local assembly class for given DIM provides assembly methods for:
21  * - local assembly on elements of dimension DIM
22  * - local assembly on edges/faces of neighbouring elemnts of dimension DIM
23  * - local assembly of interactiong elements DIM and OTHER_DIM with OTHER_DIM < DIM
24  * - Such class template would be necessary to create particular instance of Assembler class, that
25  * takes job of pass through the mesh, precomputing necessary fields, etc.
26  * This allows to manage order of assembled entities in arbitrary way in order to use caches and vectorisation efficiently.
27  * - Separate local assembly template would be necessary for every pass through the mesh, but it may delegate its actions
28  * to a single general local assembly class.
29  * - The local assembly class gets an Accessor object for particular domain of the assembly, i.e. ElementAccessor for element assembly,
30  * FaceAccessor for interface intergrals assembly, and InteractionAccessor for assembly of interaction of two elements, e.g. for XFEM
31  * there may be an interaction between every 3D element within the enrichment with appropriate 1D elemnt on the well.
32  * These accessors provides methods to access fields values as well as the DOF indices on the element(s).
33  * - Results of the assembly are passed to the global linear algebra objects collected in the RichardSystem class, global indices
34  * (which are still local indicies of the single process) provided by accessors are used.
35  *
36  *
37  * TODO:
38  * - finish functional Richards model
39  * - move whole internals of the assembly_mh_matrix into local assembly classes
40  * - mean while: imporve accessors, imporve mesh and DOF handler interface; possibly create new mesh implementation used in Darcy first and then
41  * apply it to other equations
42  */
43 
44 
45 template<int dim>
46 class AssemblyLMH : public AssemblyMH<dim> {
47 public:
48 
49  typedef std::shared_ptr<RichardsLMH::EqData> AssemblyDataPtr;
50 
51  AssemblyLMH(AssemblyDataPtr data)
52  : AssemblyMH<dim>(data),
53  ad_(data),
54  system_(data->system_),
55  genuchten_on(false),
56  cross_section(1.0),
57  soil_model(data->soil_model_)
58  {}
59 
61  genuchten_on = (this->ad_->genuchten_p_head_scale.field_result({ele.region()}) != result_zeros);
62  if (genuchten_on) {
63  SoilData soil_data;
64  soil_data.n = this->ad_->genuchten_n_exponent.value(ele.centre(), ele.element_accessor());
65  soil_data.alpha = this->ad_->genuchten_p_head_scale.value(ele.centre(), ele.element_accessor());
66  soil_data.Qr = this->ad_->water_content_residual.value(ele.centre(), ele.element_accessor());
67  soil_data.Qs = this->ad_->water_content_saturated.value(ele.centre(), ele.element_accessor());
68  soil_data.Ks = this->ad_->conductivity.value(ele.centre(), ele.element_accessor());
69  //soil_data.cut_fraction = 0.99; // set by model
70 
71  soil_model->reset(soil_data);
72  }
73 
74  }
75 
77  reset_soil_model(ele);
78  double storativity = this->ad_->storativity.value(ele.centre(), ele.element_accessor());
79  FOR_ELEMENT_SIDES(ele.full_iter(), i) {
80  double capacity = 0;
81  double water_content = 0;
82  double phead = ad_->phead_edge_[ele.edge_local_idx(i)];
83  if (genuchten_on) {
84 
85  fadbad::B<double> x_phead(phead);
86  fadbad::B<double> evaluated( soil_model->water_content_diff(x_phead) );
87  evaluated.diff(0,1);
88  water_content = evaluated.val();
89  capacity = x_phead.d(0);
90  }
91  ad_->capacity[ele.side_local_idx(i)] = capacity + storativity;
92  ad_->water_content_previous_it[ele.side_local_idx(i)] = water_content + storativity * phead;
93  }
94  }
95 
97  {
98  reset_soil_model(ele);
99  cross_section = this->ad_->cross_section.value(ele.centre(), ele.element_accessor());
100 
101 
102  double conductivity, head;
103  if (genuchten_on) {
104  conductivity=0;
105  head=0;
106  FOR_ELEMENT_SIDES(ele.full_iter(), i)
107  {
108  uint local_edge = ele.edge_local_idx(i);
109  double phead = ad_->phead_edge_[local_edge];
110  conductivity += soil_model->conductivity(phead);
111  head += ad_->phead_edge_[local_edge];
112  }
113  conductivity /= ele.n_sides();
114  head /= ele.n_sides();
115  } else {
116  conductivity = this->ad_->conductivity.value(ele.centre(), ele.element_accessor());
117  }
118 
119  double scale = 1 / cross_section / conductivity;
121 
123  }
124 
125  /***
126  * Called from assembly_local_matrix, assumes precomputed:
127  * cross_section, genuchten_on, soil_model
128  */
130 
131  // set lumped source
132  double diagonal_coef = ele.measure() * cross_section / ele.n_sides();
133 
134 
135  double source_diagonal = diagonal_coef * this->ad_->water_source_density.value(ele.centre(), ele.element_accessor());
136 
138  FOR_ELEMENT_SIDES(ele.full_iter(), i)
139  {
140 
141  uint local_edge = ele.edge_local_idx(i);
142  uint local_side = ele.side_local_idx(i);
143  uint edge_row = ele.edge_row(i);
144  if (ad_->system_.dirichlet_edge[i] == 0) {
145 
146  double capacity = this->ad_->capacity[local_side];
147  double water_content_diff = -ad_->water_content_previous_it[local_side] + ad_->water_content_previous_time[local_side];
148  double mass_diagonal = diagonal_coef * capacity;
149 
150  /*
151  DebugOut().fmt("w diff: {:g} mass: {:g} w prev: {:g} w time: {:g} c: {:g} p: {:g} z: {:g}",
152  water_content_diff,
153  mass_diagonal * ad_->phead_edge_[local_edge],
154  -ad_->water_content_previous_it[local_side],
155  ad_->water_content_previous_time[local_side],
156  capacity,
157  ad_->phead_edge_[local_edge],
158  ele.centre()[0] );
159  */
160 
161 
162  double mass_rhs = mass_diagonal * ad_->phead_edge_[local_edge] / this->ad_->time_step_
163  + diagonal_coef * water_content_diff / this->ad_->time_step_;
164 
165 
166  system_.lin_sys->mat_set_value(edge_row, edge_row, -mass_diagonal/this->ad_->time_step_ );
167  system_.lin_sys->rhs_set_value(edge_row, -source_diagonal - mass_rhs);
168  }
169 
170  if (system_.balance != nullptr) {
171  system_.balance->add_mass_vec_value(ad_->water_balance_idx_, ele.region().bulk_idx(),
172  diagonal_coef*ad_->water_content_previous_it[local_side]);
173  system_.balance->add_source_vec_values(ad_->water_balance_idx_, ele.region().bulk_idx(), {(int)edge_row}, {source_diagonal});
174  }
175  }
176 
177  }
178 
179  AssemblyDataPtr ad_;
181 
184  std::shared_ptr<SoilModelBase> soil_model;
185 };
186 
187 
188 
189 #endif /* SRC_FLOW_ASSEMBLY_LMH_HH_ */
const arma::vec3 centre() const
AssemblyLMH(AssemblyDataPtr data)
Definition: assembly_lmh.hh:51
unsigned int uint
double cross_section
ElementAccessor< 3 > element_accessor()
#define FOR_ELEMENT_SIDES(i, j)
Definition: elements.h:188
RichardsSystem system_
void reset_soil_model(LocalElementAccessorBase< 3 > ele)
Definition: assembly_lmh.hh:60
LinSys * lin_sys
double alpha
Definition: soil_models.hh:38
std::shared_ptr< SoilModelBase > soil_model
AssemblyDataPtr ad_
void mat_set_value(int row, int col, double val)
Definition: linsys.hh:325
std::shared_ptr< arma::mat > local_matrix
double Qr
Definition: soil_models.hh:39
double Qs
Definition: soil_models.hh:40
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
std::shared_ptr< RichardsLMH::EqData > AssemblyDataPtr
Definition: assembly_lmh.hh:49
double Ks
Definition: soil_models.hh:42
std::shared_ptr< Balance > balance
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definition: region.hh:90
void assembly_source_term(LocalElementAccessorBase< 3 > ele)
void assembly_local_matrix(LocalElementAccessorBase< 3 > ele) override
Definition: assembly_lmh.hh:96
void update_water_content(LocalElementAccessorBase< 3 > ele) override
Definition: assembly_lmh.hh:76
arma::mat::fixed< dim+1, dim+1 > assembly_local_geometry_matrix(ElementFullIter ele)
double n
Definition: soil_models.hh:37
ElementFullIter full_iter()
void rhs_set_value(int row, double val)
Definition: linsys.hh:337