Flow123d  release_3.0.0-968-gc359144
mortar_assembly.cc
Go to the documentation of this file.
1 /*
2  * mortar_assembly.cc
3  *
4  * Created on: Feb 22, 2017
5  * Author: jb
6  */
7 
10 #include "la/linsys.hh"
11 #include "mesh/accessors.hh"
13 #include <armadillo>
14 
16 : MortarAssemblyBase(data),
17  tensor_average_(16),
18  col_average_(4),
19  quadrature_(*(data->mesh)),
20  slave_ac_(data->mh_dh)
21 {
22  isec_data_list.reserve(30);
23 
24 
25  for(uint row_dim=0; row_dim<4; row_dim ++)
26  for(uint col_dim=0; col_dim<4; col_dim ++) {
27  arma::vec row_avg = arma::vec(row_dim+1);
28  row_avg.fill(1.0 / (row_dim+1));
29  arma::vec col_avg = arma::vec(col_dim+1);
30  col_avg.fill(1.0 / (col_dim+1));
31  arma::mat avg = row_avg * col_avg.t(); // tensor product
32  tensor_average(row_dim, col_dim) = avg;
33  col_average_[row_dim] = row_avg;
34  }
35 }
36 
37 
39  isec_data_list.push_back(IsecData());
40  IsecData &i_data = isec_data_list.back();
41 
42  i_data.dim= ele_ac.dim();
43  i_data.delta = delta;
44  i_data.dofs.resize(ele_ac.n_sides());
45  i_data.vel_dofs.resize(ele_ac.n_sides());
46  i_data.dirichlet_dofs.resize(ele_ac.n_sides());
47  i_data.dirichlet_sol.resize(ele_ac.n_sides());
48  i_data.n_dirichlet=0;
49  i_data.ele_z_coord_=ele_ac.centre()[2];
50 
51  for(unsigned int i_side=0; i_side < ele_ac.n_sides(); i_side++ ) {
52  i_data.dofs[i_side]=ele_ac.edge_row(i_side);
53  i_data.vel_dofs[i_side] = ele_ac.side_row(i_side);
54  //i_data.z_sides[i_side]=ele_ac.side(i_side)->centre()[2];
55  //DebugOut().fmt("edge: {} {}", i_side, ele_ac.edge_row(i_side));
56  Boundary * bcd = ele_ac.element_accessor().side(i_side)->cond();
57  if (bcd) {
58  ElementAccessor<3> b_ele = bcd->element_accessor();
59  auto type = (DarcyMH::EqData::BC_Type)data_->bc_type.value(b_ele.centre(), b_ele);
60  //DebugOut().fmt("bcd id: {} sidx: {} type: {}\n", ele->id(), i_side, type);
61  if (type == DarcyMH::EqData::dirichlet) {
62  //DebugOut().fmt("Dirichlet: {}\n", ele->index());
63  double bc_pressure = data_->bc_pressure.value(b_ele.centre(), b_ele);
64  i_data.dirichlet_dofs[i_data.n_dirichlet] = i_side;
65  i_data.dirichlet_sol[i_data.n_dirichlet] = bc_pressure;
66  i_data.n_dirichlet++;
67  }
68  }
69  }
70 
71 }
72 
73 
74 
75 
76 
77 /**
78  * Works well but there is large error next to the boundary.
79  */
81 {
82 
83 
84  // on the intersection element we consider
85  // * intersection dofs for master and slave
86  // those are dofs of the space into which we interpolate
87  // base functions from individual master and slave elements
88  // For the master dofs both are usualy eqivalent.
89  // * original dofs - for master same as intersection dofs, for slave
90  // all dofs of slave elements
91 
92  // form list of intersection dofs, in this case pressures in barycenters
93  // but we do not use those form MH system in order to allow second schur somplement. We rather map these
94  // dofs to pressure traces, i.e. we use averages of traces as barycentric values.
95 
96  /**
97  * TODO:
98  * - pass through the master and all slaves and collect global dofs , bcd, solution.
99  * I.e. call Nx pressure_diff not NxNx.
100  *
101  * - Is it safe to have duplicate rows in local_system?
102  * - Is it better to have more smaller local system then single big one?
103  * - use one big or more smaller local systems to set.
104  */
105 
106 
107  if (master_ac.dim() > 2) return; // supported only for 1D and 2D master elements
108  auto &isec_list = mixed_mesh_.element_intersections_[master_ac.ele_global_idx()];
109  if (isec_list.size() == 0) return; // skip empty masters
110 
111  //slave_ac_.setup(master_ac);
112 
113  ElementAccessor<3> ele = master_ac.element_accessor();
114  arma::vec3 ele_centre = ele.centre();
115  double m_sigma = data_->sigma.value( ele_centre, ele);
116  double m_conductivity = data_->conductivity.value( ele_centre, ele);
117  double m_crossection = data_->cross_section.value( ele_centre, ele );
118 
119  double master_sigma = 2*m_sigma*m_conductivity *
120  2/ m_crossection;
121  /**
122  * ?? How to deal with anisotropy ??
123  * 3d-2d : compute nv of 2d triangle
124  * 2d-2d : interpret as 2d-1d-2d, should be symmetric master-slave
125  * 2d-1d : nv is tangent to 2d and normal to 1d
126  arma::dot(data_->anisotropy.value( ele_centre, ele->element_accessor())*nv, nv)
127  */
128 
129 
130 
131 
132 
133 
134 
135  isec_data_list.clear();
136  double cs_sqr_avg = 0.0;
137  double isec_sum = 0.0;
138  uint i = 0;
139  for(; i < isec_list.size(); ++i) {
140  bool non_zero = quadrature_.reinit(isec_list[i].second);
142  if (slave_ac_.dim() == master_ac.dim()) break;
143  if (! non_zero) continue; // skip quadratures close to zero
144 
145  double cs = data_->cross_section.value(slave_ac_.element_accessor().centre(), slave_ac_.element_accessor());
146  double isec_measure = quadrature_.measure();
147  //DebugOut() << print_var(cs) << print_var(isec_measure);
148  cs_sqr_avg += cs*cs*isec_measure;
149  isec_sum += isec_measure;
150  //DebugOut().fmt("Assembly23: {} {} {} ", ele.idx(), slave_ac_.element_accessor()->id(), isec_measure);
151  pressure_diff(slave_ac_, isec_measure);
152  }
153  if ( ! (slave_ac_.dim() == 2 && master_ac.dim() ==2 ) ) {
154  if( fabs(isec_sum - ele.measure()) > 1E-5) {
155  string out;
156  for(auto & isec : isec_list) {
157  slave_ac_.reinit(isec.second->bulk_ele_idx());
158  out += fmt::format(" {}", slave_ac_.element_accessor().idx());
159  }
160 
161  double diff = (isec_sum - ele.measure())/ele.measure();
162  WarningOut() << print_var(diff)
163  << print_var(ele.idx())
164  << endl
165  << out;
166 
167  }
168  }
169  pressure_diff(master_ac, -isec_sum);
170 
171  //DebugOut().fmt( "cs2: {} d0: {}", cs_sqr_avg, delta_0);
172  master_sigma = master_sigma * (cs_sqr_avg / isec_sum)
173  / isec_sum;
174 
175 
176  add_to_linsys(master_sigma);
177 
178 
179  // 2d-2d
180  //DebugOut() << "2d-2d";
181  if (i < isec_list.size()) {
182  isec_data_list.clear();
183  isec_sum = 0.0;
184  for(; i < isec_list.size(); ++i) {
185  quadrature_.reinit(isec_list[i].second);
187  double isec_measure = quadrature_.measure();
188  isec_sum += isec_measure;
189  //DebugOut().fmt("Assembly22: {} {} {}", ele.idx(), slave_ac_.element_accessor().idx(), isec_measure);
190  pressure_diff(slave_ac_, isec_measure);
191  }
192  pressure_diff(master_ac, -isec_sum);
193 
194  master_sigma = 100 * m_conductivity/ m_crossection / isec_sum;
195 
196  add_to_linsys(master_sigma);
197  }
198 }
199 
201 {
202  // rows
203  for(IsecData &row_ele : isec_data_list) {
204  //columns
205  for(IsecData &col_ele : isec_data_list) {
206 
207 
208  double s = -scale * row_ele.delta * col_ele.delta;
209  //DebugOut().fmt("s: {} {} {} {}", s, scale, row_ele.delta, col_ele.delta);
210  product_ = s * tensor_average(row_ele.dim, col_ele.dim);
211 
212  loc_system_.reset(row_ele.dofs, col_ele.dofs);
213 
214  for(uint i=0; i< row_ele.n_dirichlet; i++)
215  loc_system_.set_solution_row(row_ele.dirichlet_dofs[i], row_ele.dirichlet_sol[i], -1.0);
216  for(uint i=0; i< col_ele.n_dirichlet; i++) loc_system_.set_solution_col(col_ele.dirichlet_dofs[i], col_ele.dirichlet_sol[i]);
217  //ASSERT( arma::norm(product_,2) == 0.0 );
219  // Must have z-coords for every side, can not use averaging vector
220  loc_system_.set_rhs( -s * col_average_[row_ele.dim] * col_ele.ele_z_coord_ );
222 
223  if (fix_velocity_flag) {
224  this->fix_velocity_local(row_ele, col_ele);
225  } else {
226  data_->lin_sys->set_local_system(loc_system_);
227  }
228  }
229  }
230 }
231 
232 
233 void P0_CouplingAssembler::fix_velocity_local(const IsecData &row_ele, const IsecData & col_ele)
234 {
235 
236  uint n_rows = row_ele.vel_dofs.n_rows;
237  uint n_cols = col_ele.dofs.n_rows;
238  arma::vec pressure(n_cols);
239  arma::vec add_velocity(n_rows);
240  double * solution = data_->lin_sys->get_solution_array();
241  for(uint icol=0; icol < n_cols; icol++ ) pressure[icol] = solution[col_ele.dofs[icol]];
242  add_velocity = loc_system_.get_matrix() * pressure - loc_system_.get_rhs();
243  //DebugOut() << "fix_velocity\n" << pressure << add_velocity;
244  for(uint irow=0; irow < n_rows; irow++ ) solution[row_ele.vel_dofs[irow]] += add_velocity[irow] ;
245 }
246 
248 {
249  for(unsigned int i_side=0; i_side < ele_ac.n_sides(); i_side++ ) {
250  dofs[shift+i_side] = ele_ac.edge_row(i_side);
251  Boundary * bcd = ele_ac.element_accessor().side(i_side)->cond();
252 
253  if (bcd) {
254  ElementAccessor<3> b_ele = bcd->element_accessor();
255  auto type = (DarcyMH::EqData::BC_Type)data_->bc_type.value(b_ele.centre(), b_ele);
256  //DebugOut().fmt("bcd id: {} sidx: {} type: {}\n", ele->id(), i_side, type);
257  if (type == DarcyMH::EqData::dirichlet) {
258  //DebugOut().fmt("Dirichlet: {}\n", ele->index());
259  dofs[shift + i_side] = -dofs[shift + i_side];
260  double bc_pressure = data_->bc_pressure.value(b_ele.centre(), b_ele);
261  dirichlet[shift + i_side] = bc_pressure;
262  }
263  }
264  }
265 }
266 
267 
268 /**
269  * P1 connection of different dimensions
270  *
271  * - 20.11. 2014 - very poor convergence, big error in pressure even at internal part of the fracture
272  */
273 
275 /*
276  const IsecList &master_list = master_list_[ele_ac.ele_global_idx()];
277  if (master_list.size() == 0) return; // skip empty masters
278  double master_sigma=data_->sigma.value( ele_ac.element_accessor()->centre(), ele_ac.element_accessor());
279 
280  // set mater sides
281  add_sides(ele_ac, 3, dofs, dirichlet);
282 
283  for(uint i = 0; i < master_list.size(); ++i) {
284  const IntersectionQuadrature &intersec = intersections_[ master_list[i] ];
285  const Element * slave = intersec.slave_iter();
286  ele_ac.reinit(intersec.slave_iter()->index());
287  add_sides(ele_ac, 0, dofs, dirichlet);
288 */
289 
290 /*
291  * Local coordinates on 1D
292  * t0
293  * node 0: 0.0
294  * node 1: 1.0
295  *
296  * base fce points
297  * t0 = 0.0 on side 0 node 0
298  * t0 = 1.0 on side 1 node 1
299  *
300  * Local coordinates on 2D
301  * t0 t1
302  * node 0: 0.0 0.0
303  * node 1: 1.0 0.0
304  * node 2: 0.0 1.0
305  *
306  * base fce points
307  * t0=0.5, t1=0.0 on side 0 nodes (0,1)
308  * t0=0.5, t1=0.5 on side 1 nodes (1,2)
309  * t0=0.0, t1=0.5 on side 2 nodes (2,0)
310  */
311 
312 /*
313 
314  arma::vec point_Y(1);
315  point_Y.fill(1.0);
316  arma::vec point_2D_Y(intersec.map_to_slave(point_Y)); // local coordinates of Y on slave (1, t0, t1)
317  arma::vec point_1D_Y(intersec.map_to_master(point_Y)); // local coordinates of Y on master (1, t0)
318 
319  arma::vec point_X(1);
320  point_X.fill(0.0);
321  arma::vec point_2D_X(intersec.map_to_slave(point_X)); // local coordinates of X on slave (1, t0, t1)
322  arma::vec point_1D_X(intersec.map_to_master(point_X)); // local coordinates of X on master (1, t0)
323 
324  arma::mat base_2D(3, 3);
325  // basis functions are numbered as sides
326  // TODO:
327  // Use RT finite element to evaluate following matrices.
328 
329  // Ravirat - Thomas base functions evaluated in points (0,0), (1,0), (0,1)
330  // 2D RT_i(t0, t1) = a0 + a1*t0 + a2*t1
331  // a0 a1 a2
332  base_2D << 1.0 << 0.0 << -2.0 << arma::endr // RT for side 0
333  << 1.0 << -2.0 << 0.0 << arma::endr // RT for side 1
334  << -1.0 << 2.0 << 2.0 << arma::endr;// RT for side 2
335 
336 
337  arma::mat base_1D(2, 2);
338  // Ravirat - Thomas base functions evaluated in points (0,0), (1,0), (0,1)
339  // 1D RT_i(t0) = a0 + a1 * t0
340  // a0 a1
341  base_1D << 1.0 << -1.0 << arma::endr // RT for side 0,
342  << 0.0 << 1.0 << arma::endr; // RT for side 1,
343 
344 
345 
346  // Consider both 2D and 1D value are defined for the test function
347  // related to the each of 5 DOFs involved in the intersection.
348  // One of these values is always zero.
349  // Compute difference of the 2D and 1D value for every DOF.
350  // Compute value of this difference in both endpoints X,Y of the intersection.
351 
352  arma::vec difference_in_Y(5);
353  arma::vec difference_in_X(5);
354  // slave sides 0,1,2
355  difference_in_Y.subvec(0, 2) = -base_2D * point_2D_Y;
356  difference_in_X.subvec(0, 2) = -base_2D * point_2D_X;
357  // master sides 3,4
358  difference_in_Y.subvec(3, 4) = base_1D * point_1D_Y;
359  difference_in_X.subvec(3, 4) = base_1D * point_1D_X;
360 
361  // applying the Simpson's rule
362  // to the product of two linear functions f, g we get
363  // (b-a)/6 * ( 3*f(a)*g(a) + 3*f(b)*g(b) + 2*f(a)*g(b) + 2*f(b)*g(a) )
364  arma::mat A(5, 5);
365  for (int i = 0; i < 5; ++i) {
366  for (int j = 0; j < 5; ++j) {
367  A(i, j) = -master_sigma * intersec.intersection_true_size() *
368  ( difference_in_Y[i] * difference_in_Y[j]
369  + difference_in_Y[i] * difference_in_X[j]/2
370  + difference_in_X[i] * difference_in_Y[j]/2
371  + difference_in_X[i] * difference_in_X[j]
372  ) * (1.0 / 3.0);
373 
374  }
375  }
376  auto dofs_cp=dofs;
377  data_->lin_sys->set_values( dofs_cp, dofs_cp, A, rhs, dirichlet, dirichlet);
378 
379  }*/
380 }
381 
382 
383 
const arma::vec3 centre() const
arma::vec dirichlet_sol
void set_solution_col(uint loc_col, double solution)
Definition: local_system.cc:91
const arma::mat & get_matrix()
Definition: local_system.hh:65
const arma::vec & get_rhs()
Definition: local_system.hh:66
void add_to_linsys(double scale)
unsigned int uint
void set_matrix(arma::mat matrix)
void reinit(uint loc_ele_idx)
ElementAccessor< 3 > element_accessor()
Boundary * cond() const
Definition: side_impl.hh:71
Wrappers for linear systems based on MPIAIJ and MATIS format.
std::string format(CStringRef format_str, ArgList args)
Definition: format.h:3141
std::shared_ptr< DarcyMH::EqData > AssemblyDataPtr
void add_sides(LocalElementAccessorBase< 3 > ele_ac, unsigned int shift, vector< int > &dofs, vector< double > &dirichlet)
std::vector< std::vector< ILpair > > element_intersections_
bool reinit(const IntersectionLocalBase *isec)
SideIter side(const unsigned int loc_index)
Definition: accessors.hh:137
void pressure_diff(LocalElementAccessorBase< 3 > ele_ac, double delta)
AssemblyDataPtr data_
IntersectionQuadratureP0 quadrature_
void assembly(LocalElementAccessorBase< 3 > ele_ac)
MixedMeshIntersections & mixed_mesh_
void fix_velocity_local(const IsecData &row_ele, const IsecData &col_ele)
void assembly(LocalElementAccessorBase< 3 > ele_ac)
double delta
LocalSystem loc_system_
P0_CouplingAssembler(AssemblyDataPtr data)
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
Definition: accessors.hh:285
LocalSystem::DofVec vel_dofs
void reset()
Resets the matrix, RHS, dofs to zero and clears solution settings.
Definition: local_system.cc:43
double measure() const
Computes the measure of the element.
Definition: accessors.hh:254
LocalElementAccessorBase< 3 > slave_ac_
LocalSystem::DofVec dirichlet_dofs
arma::mat & tensor_average(unsigned int row_dim, unsigned int col_dim)
void set_solution_row(uint loc_row, double solution, double diag=0.0)
Definition: local_system.cc:84
std::vector< arma::vec > col_average_
void set_rhs(arma::vec rhs)
unsigned int n_dirichlet
unsigned int dim
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:246
vector< IsecData > isec_data_list
double ele_z_coord_
#define print_var(var)
Definition: logger.hh:260
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
Definition: accessors.hh:111
void eliminate_solution()
ElementAccessor< 3 > element_accessor()
Definition: boundaries.cc:48
LocalSystem::DofVec dofs