Flow123d  release_3.0.0-863-g23f23ed
mortar_assembly.hh
Go to the documentation of this file.
1 /*
2  * mortar_assembly.hh
3  *
4  * Created on: Feb 22, 2017
5  * Author: jb
6  */
7 
8 #ifndef SRC_FLOW_MORTAR_ASSEMBLY_HH_
9 #define SRC_FLOW_MORTAR_ASSEMBLY_HH_
10 
11 #include "mesh/side_impl.hh"
12 #include "mesh/mesh.h"
14 #include "flow/darcy_flow_mh.hh"
15 #include "la/local_system.hh"
16 #include <vector>
17 
18 
19 typedef std::shared_ptr<DarcyMH::EqData> AssemblyDataPtr;
20 
21 
23 public:
25 
27  : data_(data),
28  mixed_mesh_(data->mesh->mixed_intersections()),
29  fix_velocity_flag(false)
30  {
31 
32  }
33 
34  virtual ~MortarAssemblyBase() {};
35 
36  // Default assembly is empty to allow dummy implementation for dimensions without coupling.
37  virtual void assembly(LocalElementAccessorBase<3> ele_ac) {};
38 
40  fix_velocity_flag = true;
41  this->assembly(ele_ac);
42  fix_velocity_flag = false;
43  }
44 
45 protected:
50 
51 };
52 
53 
54 struct IsecData {
57  unsigned int dim;
58  double delta;
59  double ele_z_coord_;
61  arma::vec dirichlet_sol;
62  unsigned int n_dirichlet;
63 };
64 
65 
67 public:
70  void pressure_diff(LocalElementAccessorBase<3> ele_ac, double delta);
71  void fix_velocity_local(const IsecData & row_ele, const IsecData &col_ele);
72 private:
73  inline arma::mat & tensor_average(unsigned int row_dim, unsigned int col_dim) {
74  return tensor_average_[4*row_dim + col_dim];
75  }
76 
77  void add_to_linsys(double scale);
78 
80 
81  /// Row matrices to compute element pressure as average of boundary pressures
85  arma::mat product_;
87 
88 };
89 
90 
91 
93 public:
95 : MortarAssemblyBase(data),
96  rhs(5),
97  dofs(5),
98  dirichlet(5)
99  {
100  rhs.zeros();
101  }
102 
104  void add_sides(LocalElementAccessorBase<3> ele_ac, unsigned int shift, vector<int> &dofs, vector<double> &dirichlet);
105 private:
106 
107  arma::vec rhs;
110 };
111 
112 
113 
114 #endif /* SRC_FLOW_MORTAR_ASSEMBLY_HH_ */
std::vector< arma::mat > tensor_average_
Row matrices to compute element pressure as average of boundary pressures.
arma::vec dirichlet_sol
virtual void assembly(LocalElementAccessorBase< 3 > ele_ac)
std::shared_ptr< DarcyMH::EqData > AssemblyDataPtr
virtual ~MortarAssemblyBase()
vector< double > dirichlet
AssemblyDataPtr data_
arma::Col< LongIdx > DofVec
Definition: local_system.hh:30
IntersectionQuadratureP0 quadrature_
MixedMeshIntersections & mixed_mesh_
double delta
LocalSystem loc_system_
P1_CouplingAssembler(AssemblyDataPtr data)
LocalSystem::DofVec vel_dofs
void fix_velocity(LocalElementAccessorBase< 3 > ele_ac)
MortarAssemblyBase(AssemblyDataPtr data)
LocalElementAccessorBase< 3 > slave_ac_
LocalSystem::DofVec dirichlet_dofs
arma::mat & tensor_average(unsigned int row_dim, unsigned int col_dim)
vector< unsigned int > IsecList
std::vector< arma::vec > col_average_
unsigned int n_dirichlet
unsigned int dim
vector< IsecData > isec_data_list
double ele_z_coord_
mixed-hybrid model of linear Darcy flow, possibly unsteady.
Main class for computation of intersection of meshes of combined dimensions.
LocalSystem::DofVec dofs