|
template<typename ParameterClass > |
| Adapter (ParameterClass const ¶meters, MPI_Comm mpi_comm) |
| Constructor, which sets up the precice Solverinterface.
|
|
void | initialize_precice (VectorType const &dealii_to_precice) |
| Initializes preCICE and passes all relevant data to preCICE.
|
|
void | add_write_surface (dealii::types::boundary_id const surface_id, std::string const &mesh_name, std::vector< std::string > const &write_data_names, WriteDataType const write_data_type, dealii::MatrixFree< dim, double, VectorizedArrayType > const &data, unsigned int const dof_index, unsigned int const write_quad_index) |
|
void | add_read_surface (dealii::MatrixFree< dim, double, VectorizedArrayType > const &data, std::shared_ptr< ContainerInterfaceData< rank, dim, double > > interface_data, std::string const &mesh_name, std::vector< std::string > const &read_data_name) |
|
void | advance (double const computed_timestep_length) |
| Advances preCICE after every timestep.
|
|
void | save_current_state_if_required (std::function< void()> const &save_state) |
| Saves current state of time dependent variables in case of an implicit coupling.
|
|
void | reload_old_state_if_required (std::function< void()> const &reload_old_state) |
| Reloads the previously stored variables in case of an implicit coupling. The current implementation supports subcycling, i.e. previously refers o the last time save_current_state_if_required() has been called.
|
|
void | write_data (std::string const &write_mesh_name, std::string const &write_data_name, VectorType const &write_data, double const computed_timestep_length) |
|
void | read_block_data (std::string const &mesh_name, const std::string &data_name) const |
|
bool | is_coupling_ongoing () const |
| is_coupling_ongoing Calls the preCICE API function isCouplingOnGoing
|
|
bool | is_time_window_complete () const |
| is_time_window_complete Calls the preCICE API function isTimeWindowComplete
|
|
template<int dim, int data_dim, typename VectorType, typename VectorizedArrayType = dealii::VectorizedArray<double>>
class ExaDG::preCICE::Adapter< dim, data_dim, VectorType, VectorizedArrayType >
The Adapter class keeps together with the CouplingInterfaes all functionalities to couple deal.II to other solvers with preCICE i.e. data structures are set up, necessary information is passed to preCICE etc.
template<int dim, int data_dim, typename VectorType , typename VectorizedArrayType >
template<typename ParameterClass >
ExaDG::preCICE::Adapter< dim, data_dim, VectorType, VectorizedArrayType >::Adapter |
( |
ParameterClass const & | parameters, |
|
|
MPI_Comm | mpi_comm ) |
Constructor, which sets up the precice Solverinterface.
- Template Parameters
-
data_dim | Dimension of the coupling data. Equivalent to n_components in the deal.II documentation |
- Parameters
-
[in] | parameters | Parameter class, which hold the data specified in the parameters.prm file |
[in] | dealii_boundary_surface_id | Boundary ID of the triangulation, which is associated with the coupling surface. |
[in] | data | The applied matrix-free object |
[in] | dof_index | Index of the relevant dof_handler in the corresponding MatrixFree object |
[in] | read_quad_index | Index of the quadrature formula in the corresponding MatrixFree object which should be used for data reading |
[in] | is_dirichlet | Boolean to distinguish between Dirichlet type solver (using the DoFs for data reading) and Neumann type solver (using quadrature points for reading) |