21 #ifndef _aspect_particle_property_interface_h 22 #define _aspect_particle_property_interface_h 30 #include <deal.II/particles/particle.h> 31 #include <deal.II/particles/particle_handler.h> 32 #include <deal.II/particles/property_pool.h> 33 #include <deal.II/fe/fe_update_flags.h> 43 using namespace ::Particles;
324 initialize_one_particle_property (
const Point<dim> &position,
325 std::vector<double> &particle_properties)
const;
355 update_particle_properties (
const unsigned int data_position,
356 const std::vector<Vector<double>> &solution,
357 const std::vector<std::vector<Tensor<1,dim>>> &gradients,
358 typename ParticleHandler<dim>::particle_iterator_range &particles)
const;
391 update_particle_property (
const unsigned int data_position,
392 const Vector<double> &solution,
393 const std::vector<Tensor<1,dim>> &gradients,
394 typename ParticleHandler<dim>::particle_iterator &particle)
const;
411 need_update ()
const;
423 get_needed_update_flags ()
const;
441 late_initialization_mode ()
const;
457 std::vector<std::pair<std::string, unsigned int>>
458 get_property_information()
const = 0;
480 initialize_one_particle_property (
const Point<dim> &position,
481 std::vector<double> &particle_properties)
const override;
490 std::vector<std::pair<std::string, unsigned int>>
491 get_property_information()
const override;
498 parse_parameters (ParameterHandler &prm)
override;
532 initialize_one_particle (
typename ParticleHandler<dim>::particle_iterator &particle)
const;
541 initialize_late_particle (
const Point<dim> &particle_location,
542 const ParticleHandler<dim> &particle_handler,
544 const typename parallel::distributed::Triangulation<dim>::active_cell_iterator &cell =
typename parallel::distributed::Triangulation<dim>::active_cell_iterator())
const;
558 update_particles (
typename ParticleHandler<dim>::particle_iterator_range &particles,
559 const std::vector<Vector<double>> &solution,
560 const std::vector<std::vector<Tensor<1,dim>>> &gradients)
const;
574 need_update ()
const;
583 get_needed_update_flags ()
const;
590 plugin_name_exists(
const std::string &name)
const;
602 check_plugin_order(
const std::string &first,
const std::string &second)
const;
607 unsigned int get_plugin_index_by_name(
const std::string &name)
const;
623 template <
typename ParticlePropertyType,
624 typename =
typename std::enable_if_t<std::is_base_of<Interface<dim>,ParticlePropertyType>::value>>
627 has_matching_property ()
const;
645 template <
typename ParticlePropertyType,
646 typename =
typename std::enable_if_t<std::is_base_of<Interface<dim>,ParticlePropertyType>::value>>
648 const ParticlePropertyType &
649 get_matching_property ()
const;
659 get_n_property_components ()
const;
670 get_particle_size ()
const;
679 get_data_info()
const;
690 get_property_component_by_name(
const std::string &name)
const;
712 register_particle_property (
const std::string &name,
713 const std::string &description,
714 void (*declare_parameters_function) (ParameterHandler &),
729 parse_parameters (ParameterHandler &prm)
override;
736 void set_particle_world_index(
unsigned int particle_world_index);
768 template <
typename ParticlePropertyType,
typename>
773 return this->
template has_matching_active_plugin<ParticlePropertyType>();
778 template <
typename ParticlePropertyType,
typename>
780 const ParticlePropertyType &
783 return this->
template get_matching_active_plugin<ParticlePropertyType>();
793 #define ASPECT_REGISTER_PARTICLE_PROPERTY(classname,name,description) \ 794 template class classname<2>; \ 795 template class classname<3>; \ 796 namespace ASPECT_REGISTER_PARTICLE_PROPERTY_ ## classname \ 798 aspect::internal::Plugins::RegisterHelper<aspect::Particle::Property::Interface<2>,classname<2>> \ 799 dummy_ ## classname ## _2d (&aspect::Particle::Property::Manager<2>::register_particle_property, \ 800 name, description); \ 801 aspect::internal::Plugins::RegisterHelper<aspect::Particle::Property::Interface<3>,classname<3>> \ 802 dummy_ ## classname ## _3d (&aspect::Particle::Property::Manager<3>::register_particle_property, \ 803 name, description); \ unsigned int n_components() const
std::vector< unsigned int > position_per_field
std::string get_field_name_by_index(const unsigned int field_index) const
StructuredDataLookup< dim > DEAL_II_DEPRECATED
DEAL_II_DEPRECATED const ParticlePropertyType & get_matching_property() const
ParticlePropertyInformation property_information
unsigned int get_components_by_field_index(const unsigned int field_index) const
std::vector< unsigned int > components_per_plugin
void write_plugin_graph(std::ostream &output_stream)
unsigned int get_components_by_field_name(const std::string &name) const
unsigned int get_position_by_field_name(const std::string &name) const
bool fieldname_exists(const std::string &name) const
InitializationModeForLateParticles
std::vector< std::string > field_names
std::vector< unsigned int > components_per_field
unsigned int n_integrator_properties
unsigned int get_field_index_by_name(const std::string &name) const
unsigned int number_of_components
unsigned int number_of_plugins
unsigned int get_position_by_plugin_index(const unsigned int plugin_index) const
unsigned int particle_world_index
std::vector< unsigned int > fields_per_plugin
unsigned int n_fields() const
unsigned int get_components_by_plugin_index(const unsigned int plugin_index) const
unsigned int get_fields_by_plugin_index(const unsigned int plugin_index) const
DEAL_II_DEPRECATED bool has_matching_property() const
unsigned int number_of_fields
std::vector< unsigned int > position_per_plugin
void declare_parameters(ParameterHandler &prm)
unsigned int n_plugins() const
ParticlePropertyInformation()
unsigned int get_position_by_field_index(const unsigned int field_index) const