Flow123d  JS_before_hm-995-g9546b8d
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 "system/index_types.hh"
12 #include "mesh/mesh.h"
14 #include "flow/darcy_flow_mh.hh"
15 #include "la/local_system.hh"
16 #include "fem/dh_cell_accessor.hh"
17 #include <vector>
18 
19 
20 typedef std::shared_ptr<DarcyMH::EqData> AssemblyDataPtr;
21 
22 
24 public:
26 
28  : data_(data),
29  mixed_mesh_(data->mesh->mixed_intersections()),
30  fix_velocity_flag(false)
31  {
32 
33  }
34 
35  virtual ~MortarAssemblyBase() {};
36 
37  // Default assembly is empty to allow dummy implementation for dimensions without coupling.
38  virtual void assembly(FMT_UNUSED const DHCellAccessor& dh_cell) {};
39 
40  void fix_velocity(const DHCellAccessor& dh_cell) {
41  fix_velocity_flag = true;
42  this->assembly(dh_cell);
43  fix_velocity_flag = false;
44  }
45 
46 protected:
51 
52 };
53 
54 
55 struct IsecData {
58  unsigned int dim;
59  double delta;
60  double ele_z_coord_;
63  unsigned int n_dirichlet;
64 };
65 
66 
68 public:
70  void assembly(const DHCellAccessor& dh_cell);
71  void pressure_diff(const DHCellAccessor& dh_cell, double delta);
72  void fix_velocity_local(const IsecData & row_ele, const IsecData &col_ele);
73 private:
74  inline arma::mat & tensor_average(unsigned int row_dim, unsigned int col_dim) {
75  return tensor_average_[4*row_dim + col_dim];
76  }
77 
78  void add_to_linsys(double scale);
79 
81 
82  /// Row matrices to compute element pressure as average of boundary pressures
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 
103  void assembly(const DHCellAccessor& dh_cell);
104  void add_sides(const DHCellAccessor& dh_cell, unsigned int shift, vector<int> &dofs, vector<double> &dirichlet);
105 private:
106 
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
ArmaVec< double, N > vec
Definition: armor.hh:861
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
std::shared_ptr< DarcyMH::EqData > AssemblyDataPtr
LocDofVec dirichlet_dofs
virtual ~MortarAssemblyBase()
Cell accessor allow iterate over DOF handler cells.
vector< double > dirichlet
AssemblyDataPtr data_
IntersectionQuadratureP0 quadrature_
MixedMeshIntersections & mixed_mesh_
double delta
LocalSystem loc_system_
virtual void assembly(FMT_UNUSED const DHCellAccessor &dh_cell)
P1_CouplingAssembler(AssemblyDataPtr data)
LocDofVec vel_dofs
ArmaMat< double, N, M > mat
Definition: armor.hh:864
#define FMT_UNUSED
Definition: posix.h:75
MortarAssemblyBase(AssemblyDataPtr data)
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.
void fix_velocity(const DHCellAccessor &dh_cell)
LocDofVec dofs
Main class for computation of intersection of meshes of combined dimensions.