Flow123d
JB_transport-112d700
|
Go to the documentation of this file.
33 #ifndef DARCY_FLOW_LMH_HH
34 #define DARCY_FLOW_LMH_HH
137 <<
"Diverged nonlinear solver. Reason: " << EI_Reason::val
140 <<
"Missing the key 'time', obligatory for the transient problems.");
221 std::shared_ptr<DOFHandlerMultiDim>
dh_;
222 std::shared_ptr<SubDOFHandlerMultiDim>
dh_cr_;
224 std::shared_ptr<SubDOFHandlerMultiDim>
dh_p_;
313 {
eq_fields_->extra_storativity = extra_stor; }
396 {
return *(
eq_data_->lin_sys_schur); }
434 {
return "Flow_Darcy_LMH";}
437 #endif //DARCY_FLOW_LMH_HH
std::vector< bool > save_local_system_
Flag for saving the local system. Currently used only in case of seepage BC.
BCField< 3, FieldValue< 3 >::Scalar > bc_flux
Field< 3, FieldValue< 3 >::Vector > field_ele_velocity
Field< 3, FieldValue< 3 >::Scalar > conductivity
void print_matlab_matrix(string matlab_file)
Print darcy flow matrix in matlab format into a file.
virtual bool zero_time_term(bool time_global=false)
void reset()
Reset data members.
std::shared_ptr< DOFHandlerMultiDim > dh_
full DOF handler represents DOFs of sides, elements and edges
DarcyFlowMHOutput * output_object
Basic time management class.
VectorMPI p_edge_solution
virtual void accept_time_step()
postprocess velocity field (add sources)
std::vector< bool > bc_fluxes_reconstruted
Flag indicating whether the fluxes for seepage BC has been reconstructed already.
void create_linear_system(Input::AbstractRecord rec)
int is_linear
Hack fo BDDC solver.
std::map< LongIdx, LocalSystem > seepage_bc_systems
TYPEDEF_ERR_INFO(EI_Reason, string)
EqFields()
Creation of all fields.
Field< 3, FieldValue< 3 >::Vector > gravity_field
Field< 3, FieldValue< 3 >::Tensor > anisotropy
std::array< unsigned int, 3 > schur_offset_
Index offset in the local system for the Schur complement (of dim = 1,2,3).
std::shared_ptr< LinSys > lin_sys_schur
void solve_time_step(bool output=true)
Solve the problem without moving to next time and without output.
void update_solution() override
static const Input::Type::Record & get_input_type()
virtual void initialize_asm()
Create and initialize assembly objects.
virtual void postprocess()
DECLARE_EXCEPTION(ExcSolverDiverge,<< "Diverged nonlinear solver. Reason: "<< EI_Reason::val)
virtual void assembly_linear_system()
void set_extra_source(const Field< 3, FieldValue< 3 >::Scalar > &extra_src)
Sets external source field (coupling with other equation).
GenericAssemblyBase * mh_matrix_assembly_
Field< 3, FieldValue< 3 >::Scalar > cross_section
BCField< 3, FieldValue< 3 >::Scalar > bc_switch_pressure
static const Input::Type::Record & type_field_descriptor()
void zero_time_step() override
std::vector< char > bc_switch_dirichlet
Idicator of dirichlet or neumann type of switch boundary conditions.
FieldCoords X_
Field holds coordinates for computing of FieldFormulas.
BCField< 3, FieldValue< 3 >::Scalar > bc_robin_sigma
BCField< 3, FieldValue< 3 >::Scalar > bc_switch_piezo_head
std::array< std::vector< unsigned int >, 3 > loc_edge_dofs
BCField< 3, FieldValue< 3 >::Scalar > bc_piezo_head
std::vector< arma::vec > postprocess_solution_
std::shared_ptr< SubDOFHandlerMultiDim > dh_p_
DOF handler represents DOFs of element pressure.
MortarMethod mortar_method_
void solve_nonlinear()
Solve method common to zero_time_step and update solution.
std::vector< LocalSystem > loc_system_
Following data members are stored in vectors, one item for every cell.
std::shared_ptr< Balance > balance_
static const int registrar
Registrar of class to factory.
void set_extra_storativity(const Field< 3, FieldValue< 3 >::Scalar > &extra_stor)
Sets external storarivity field (coupling with other equation).
Basic time management functionality for unsteady (and steady) solvers (class Equation).
static const Input::Type::Selection & get_bc_type_selection()
Return a Selection corresponding to enum BC_Type.
virtual double solved_time() override
Field< 3, FieldValue< 3 >::Scalar > init_piezo_head
Same as previous but used in boundary fields.
BCField< 3, FieldValue< 3 >::Scalar > bc_pressure
BCField< 3, FieldValue< 3 >::Vector > bc_gravity
Holds gravity vector acceptable in FieldModel.
std::shared_ptr< DOFHandlerMultiDim > dh_cr_disc_
DOF handler represents DOFs of sides.
Container for various descendants of FieldCommonBase.
static std::string equation_name()
virtual void initialize_specific()
MortarMethod
Type of experimental Mortar-like method for non-compatible 1d-2d interaction.
std::shared_ptr< EqData > eq_data_
virtual ~DarcyLMH() override
unsigned int nonlinear_iteration_
void initialize() override
Field< 3, FieldValue< 3 >::Scalar > field_ele_pressure
Externally added water source.
std::vector< int > get_component_indices_vec(unsigned int component) const
Get vector of all DOF indices of given component (0..side, 1..element, 2..edge)
void init()
Initialize vectors, ...
Field< 3, FieldValue< 3 >::Scalar > extra_source
Externally added storativity.
LinSys & lin_sys_schur()
Getter for the linear system of the 2. Schur complement.
bool use_steady_assembly_
Field< 3, FieldValue< 3 >::Vector > flux
Field< 3, FieldValue< 3 >::Scalar > sigma
std::array< std::vector< unsigned int >, 3 > loc_side_dofs
static const Input::Type::Selection & get_mh_mortar_selection()
Selection for enum MortarMethod.
Field< 3, FieldValue< 3 >::Scalar > field_ele_piezo_head
void allocate_mh_matrix()
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
bool force_no_neumann_bc
auxiliary flag for switchting Dirichlet like BC
std::shared_ptr< SubDOFHandlerMultiDim > dh_cr_
DOF handler represents DOFs of edges.
Field< 3, FieldValue< 3 >::Scalar > init_pressure
GenericAssemblyBase * reconstruct_schur_assembly_
Field< 3, FieldValue< 3 >::Scalar > water_source_density
FieldCoords & X()
Return coords field.
GenericAssembly< ReadInitCondAssemblyLMH > * read_init_cond_assembly_
general assembly objects, hold assembly objects of appropriate dimension
VectorMPI p_edge_solution_previous
std::vector< LocalConstraint > loc_constraint_
Field< 3, FieldValue< 3 >::Scalar > field_edge_pressure
std::shared_ptr< Balance > balance_
Shared Balance object.
DECLARE_INPUT_EXCEPTION(ExcMissingTimeGovernor,<< "Missing the key 'time', obligatory for the transient problems.")
Field< 3, FieldValue< 3 >::Scalar > storativity
Class template representing a field with values dependent on: point, element, and region.
virtual void read_init_cond_asm()
Call assemble of read_init_cond_assembly_.
BCField< 3, FieldValue< 3 >::Enum > bc_type
std::shared_ptr< EqFields > eq_fields_
Field< 3, FieldValue< 3 >::Scalar > extra_storativity
std::array< unsigned int, 3 > loc_ele_dof
virtual void output_data() override
Write computed fields.
Generic class of assemblation.
DarcyLMH(Mesh &mesh, const Input::Record in_rec, TimeGovernor *tm=nullptr)
CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL.
VectorMPI p_edge_solution_previous_time
unsigned int IntDim
A dimension index type.