17 #include <deal.II/base/config.h> 19 #ifdef DEAL_II_WITH_P4EST 21 #include <deal.II/lac/vector.h> 22 #include <deal.II/lac/block_vector.h> 23 #include <deal.II/lac/la_parallel_vector.h> 24 #include <deal.II/lac/la_parallel_block_vector.h> 25 #include <deal.II/lac/petsc_parallel_vector.h> 26 #include <deal.II/lac/petsc_parallel_block_vector.h> 27 #include <deal.II/lac/trilinos_vector.h> 28 #include <deal.II/lac/trilinos_parallel_block_vector.h> 30 #include <deal.II/distributed/solution_transfer.h> 31 #include <deal.II/distributed/tria.h> 32 #include <deal.II/dofs/dof_tools.h> 33 #include <deal.II/dofs/dof_accessor.h> 34 #include <deal.II/grid/tria_accessor.h> 35 #include <deal.II/grid/tria_iterator.h> 39 DEAL_II_NAMESPACE_OPEN
46 template <
int dim,
typename VectorType,
typename DoFHandlerType>
49 dof_handler(&dof, typeid(*this).name()),
50 offset (
numbers::invalid_unsigned_int)
54 ExcMessage(
"parallel::distributed::SolutionTransfer requires a parallel::distributed::Triangulation object."));
59 template <
int dim,
typename VectorType,
typename DoFHandlerType>
64 input_vectors = all_in;
65 register_data_attach( get_data_size() * input_vectors.size() );
70 template <
int dim,
typename VectorType,
typename DoFHandlerType>
80 (&dof_handler->get_triangulation())));
86 DoFHandlerType>::pack_callback,
88 std::placeholders::_1,
89 std::placeholders::_2,
90 std::placeholders::_3));
96 template <
int dim,
typename VectorType,
typename DoFHandlerType>
101 std::vector<const VectorType *> all_in(1, &in);
102 prepare_for_coarsening_and_refinement(all_in);
107 template <
int dim,
typename VectorType,
typename DoFHandlerType>
111 std::vector<const VectorType *> all_in(1, &in);
112 prepare_serialization(all_in);
117 template <
int dim,
typename VectorType,
typename DoFHandlerType>
120 (
const std::vector<const VectorType *> &all_in)
122 prepare_for_coarsening_and_refinement (all_in);
127 template <
int dim,
typename VectorType,
typename DoFHandlerType>
131 std::vector<VectorType *> all_in(1, &in);
137 template <
int dim,
typename VectorType,
typename DoFHandlerType>
141 register_data_attach( get_data_size() * all_in.size() );
144 input_vectors.resize(all_in.size());
150 template <
int dim,
typename VectorType,
typename DoFHandlerType>
154 Assert(input_vectors.size()==all_out.size(),
161 (&dof_handler->get_triangulation())));
166 DoFHandlerType>::unpack_callback,
168 std::placeholders::_1,
169 std::placeholders::_2,
170 std::placeholders::_3,
174 for (
typename std::vector<VectorType *>::iterator it=all_out.begin();
179 input_vectors.clear();
184 template <
int dim,
typename VectorType,
typename DoFHandlerType>
188 std::vector<VectorType *> all_out(1, &out);
189 interpolate(all_out);
194 template <
int dim,
typename VectorType,
typename DoFHandlerType>
202 template <
int dim,
typename VectorType,
typename DoFHandlerType>
206 const typename Triangulation<dim,DoFHandlerType::space_dimension>::CellStatus ,
209 typename VectorType::value_type *data_store =
reinterpret_cast<typename VectorType::value_type *
>(data);
211 typename DoFHandlerType::cell_iterator cell(*cell_, dof_handler);
213 const unsigned int dofs_per_cell=cell->get_fe().dofs_per_cell;
215 for (
typename std::vector<const VectorType *>::iterator it=input_vectors.begin();
216 it !=input_vectors.end();
219 cell->get_interpolated_dof_values(*(*it), dofvalues);
220 std::memcpy(data_store, &dofvalues(0),
sizeof(
typename VectorType::value_type)*dofs_per_cell);
221 data_store += dofs_per_cell;
226 template <
int dim,
typename VectorType,
typename DoFHandlerType>
230 const typename Triangulation<dim,DoFHandlerType::space_dimension>::CellStatus ,
232 std::vector<VectorType *> &all_out)
234 typename DoFHandlerType::cell_iterator
235 cell(*cell_, dof_handler);
237 const unsigned int dofs_per_cell=cell->get_fe().dofs_per_cell;
239 const typename VectorType::value_type *data_store =
reinterpret_cast<const typename VectorType::value_type *
>(data);
241 for (
typename std::vector<VectorType *>::iterator it = all_out.begin();
245 std::memcpy(&dofvalues(0), data_store,
sizeof(
typename VectorType::value_type)*dofs_per_cell);
246 cell->set_dof_values_by_interpolation(dofvalues, *(*it));
247 data_store += dofs_per_cell;
257 #include "solution_transfer.inst" 259 DEAL_II_NAMESPACE_CLOSE
void prepare_for_coarsening_and_refinement(const std::vector< const VectorType *> &all_in)
unsigned int get_data_size() const
void deserialize(VectorType &in)
SolutionTransfer(const DoFHandlerType &dof)
void pack_callback(const typename Triangulation< dim, DoFHandlerType::space_dimension >::cell_iterator &cell, const typename Triangulation< dim, DoFHandlerType::space_dimension >::CellStatus status, void *data)
void interpolate(std::vector< VectorType *> &all_out)
unsigned int register_data_attach(const std::size_t size, const std::function< void(const cell_iterator &, const CellStatus, void *)> &pack_callback)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void notify_ready_to_unpack(const unsigned int handle, const std::function< void(const cell_iterator &, const CellStatus, const void *)> &unpack_callback)
SmartPointer< const DoFHandlerType, SolutionTransfer< dim, VectorType, DoFHandlerType > > dof_handler
void prepare_serialization(const VectorType &in)
void unpack_callback(const typename Triangulation< dim, DoFHandlerType::space_dimension >::cell_iterator &cell, const typename Triangulation< dim, DoFHandlerType::space_dimension >::CellStatus status, const void *data, std::vector< VectorType *> &all_out)
::Triangulation< dim, spacedim >::cell_iterator cell_iterator
static ::ExceptionBase & ExcInternalError()