33 #ifndef DARCY_FLOW_MH_HH
34 #define DARCY_FLOW_MH_HH
64 template<
unsigned int dim,
unsigned int spacedim>
class FE_RT0;
65 template<
unsigned int degree,
unsigned int dim,
unsigned int spacedim>
class FE_P_disc;
66 template<
unsigned int dim,
unsigned int spacedim>
class MappingP1;
67 template<
unsigned int dim,
unsigned int spacedim>
class FEValues;
68 template<
unsigned int dim,
unsigned int spacedim>
class FESideValues;
69 template<
unsigned int dim>
class QGauss;
145 <<
"Diverged nonlinear solver. Reason: " << EI_Reason::val
394 arma::mat master_map(1,2), slave_map(1,3);
395 master_map.fill(1.0 / 2);
396 slave_map.fill(1.0 / 3);
407 unsigned int &ele_type,
409 arma::vec &dirichlet);
482 static const int registrar;
488 #endif //DARCY_FLOW_MH_HH
void get_solution_vector(double *&vec, unsigned int &vec_size) override
void make_serial_scatter()
Abstract base class for equation clasess.
Container for various descendants of FieldCommonBase.
void assembly_mh_matrix(MultidimAssembler ma)
Field< 3, FieldValue< 3 >::Scalar > init_pressure
virtual void initialize_specific()
void get_parallel_solution_vector(Vec &vector) override
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
std::vector< unsigned int > dirichlet_edge
virtual void assembly_source_term()
Source term is implemented differently in LMH version.
vector< vector< arma::mat > > tensor_average
Row matrices to compute element pressure as average of boundary pressures.
virtual void output_data() override
Write computed fields.
Class template representing a field with values dependent on: point, element, and region...
void assembly(LinSys &ls)
static const int registrar
Registrar of class to factory.
Field< 3, FieldValue< 3 >::Scalar > sigma
static const Input::Type::Selection & get_mh_mortar_selection()
Selection for enum MortarMethod.
void initialize() override
void assembly(LinSys &ls)
std::vector< char > bc_switch_dirichlet
Idicator of dirichlet or neumann type of switch boundary conditions.
static const Input::Type::Record & type_field_descriptor()
vector< IsecList >::const_iterator ml_it_
P1_CouplingAssembler(const DarcyMH &darcy)
void pressure_diff(int i_ele, vector< int > &dofs, unsigned int &ele_type, double &delta, arma::vec &dirichlet)
vector< double > dirichlet
P0_CouplingAssembler(const DarcyMH &darcy)
bool solution_changed_for_scatter
void solve_nonlinear()
Solve method common to zero_time_step and update solution.
std::shared_ptr< EqData > data_
void mat_count_off_proc_values(Mat m, Vec v)
class to manage local indices on sub-domain to global indices on domain
MortarMethod mortar_method_
vector< unsigned int > IsecList
Symmetric Gauss-Legendre quadrature formulae on simplices.
unsigned int nonlinear_iteration_
void set_mesh_data_for_bddc(LinSys_BDDC *bddc_ls)
BCField< 3, FieldValue< 3 >::Scalar > bc_flux
void update_solution() override
Discontinuous Lagrangean finite element on dim dimensional simplex.
std::shared_ptr< arma::mat > local_matrix
BCField< 3, FieldValue< 3 >::Scalar > bc_pressure
virtual void postprocess()
const MH_DofHandler & get_mh_dofhandler()
Field< 3, FieldValue< 3 >::Scalar > water_source_density
BCField< 3, FieldValue< 3 >::Scalar > bc_switch_pressure
Field< 3, FieldValue< 3 >::Scalar > storativity
std::shared_ptr< Balance > balance_
Virtual class for construction and partitioning of a distributed sparse graph.
Raviart-Thomas element of order 0.
const vector< IsecList > & master_list_
bool use_steady_assembly_
Affine mapping between reference and actual cell.
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
static const Input::Type::Selection & get_bc_type_selection()
Return a Selection corresponding to enum BC_Type.
const vector< Intersection > & intersections_
virtual void assembly_linear_system()
const vector< Intersection > & intersections_
DarcyMH(Mesh &mesh, const Input::Record in_rec)
CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL.
double delta_0
measure of master element, should be sum of intersection measures
void create_linear_system(Input::AbstractRecord rec)
Initialize global_row_4_sub_row.
MortarMethod
Type of experimental Mortar-like method for non-compatible 1d-2d interaction.
DECLARE_EXCEPTION(ExcSolverDiverge,<< "Diverged nonlinear solver. Reason: "<< EI_Reason::val)
std::shared_ptr< Balance > balance
Vector classes to support both Iterator, index and Id access and creating co-located vectors...
Field< 3, FieldValue< 3 >::TensorFixed > anisotropy
unsigned int water_balance_idx_
index of water balance within the Balance object.
Abstract linear system class.
Calculates finite element data on the actual cell.
TYPEDEF_ERR_INFO(EI_Reason, string)
std::vector< std::shared_ptr< AssemblyBase > > MultidimAssembler
EqData()
Creation of all fields.
void zero_time_step() override
void add_sides(const Element *ele, unsigned int shift, vector< int > &dofs, vector< double > &dirichlet)
virtual void read_initial_condition()
Field< 3, FieldValue< 3 >::Scalar > cross_section
BCField< 3, FieldValue< 3 >::Enum > bc_type
DarcyFlowMHOutput * output_object
Field< 3, FieldValue< 3 >::Scalar > conductivity
void set_solution(double time, double *solution, double precision)
virtual void prepare_new_time_step()
postprocess velocity field (add sources)
virtual double solution_precision() const
BCField< 3, FieldValue< 3 >::Scalar > bc_robin_sigma
virtual bool zero_time_term()
virtual void setup_time_term()
static const Input::Type::Record & get_input_type()
Solver based on Multilevel BDDC - using corresponding class of OpenFTL package.
Calculates finite element data on a side.
Mixed-hybrid of steady Darcy flow with sources and variable density.