| Parallel Finite Element - General geometry Ewald-like Method
    Part of Continuum-particle Simulation Suite under MICCOM | 
#include <particle_mesh.h>


| Public Member Functions | |
| ParticleMesh (MeshBase &mesh) | |
| ParticleMesh (MeshBase &mesh, const Real &search_radius_p, const Real &search_radius_e) | |
| ParticleMesh (MeshBase &mesh, PMPeriodicBoundary &pmpb, const Real &search_radius_p, const Real &search_radius_e) | |
| ~ParticleMesh () | |
| void | read_particles_data (const std::string &filename, const std::string &particle_mesh_type) | 
| void | read_particles_data_restart (const std::string &filename, const std::string &particle_mesh_type) | 
| void | read_particles_data_restart (const std::string &filename, const unsigned int &o_step) | 
| void | generate_random_particles (const std::size_t N, const Real bbox_XA, const Real bbox_XB, const Real bbox_YA, const Real bbox_YB, const Real bbox_ZA, const Real bbox_ZB) | 
| void | generate_random_particles (const std::size_t N, const Point &bbox_min, const Point &bbox_max) | 
| std::vector< RigidParticle * > | particles () const | 
| std::size_t | num_particles () const | 
| bool | is_sorted () const | 
| void | reinit () | 
| void | print_particle_info () const | 
| void | set_search_radius (const Real rp, const Real re) | 
| Real | search_radius (const std::string &p_e) const | 
| void | add_periodic_boundary (PMPeriodicBoundary &_periodic_bdry) | 
| PMPeriodicBoundary * | pm_periodic_boundary () | 
| MeshBase & | mesh () | 
| std::size_t | num_mesh_points () const | 
| void | update_mesh (const std::vector< Point > &nodal_pos, const std::vector< Point > &nodal_vel) | 
| void | zero_node_force () | 
| std::vector< Real > | mesh_size () const | 
| void | volume_conservation () | 
| void | initial_particle_center_of_mass (std::vector< Point > ¢er0) const | 
| void | write_particle (const unsigned int &step_id, const unsigned int &o_step, const Real &real_time, const std::vector< std::string > &output_file, unsigned int comm_in_rank) const | 
| void | write_time (const unsigned int &step_id, const unsigned int &o_step, const Real &real_time, unsigned int comm_in_rank) const | 
| void | write_particle_trajectory (const unsigned int &o_step, unsigned int comm_in_rank) const | 
| void | write_surface_node (const unsigned int &o_step, unsigned int comm_in_rank) const | 
| void | write_particle_mesh (const unsigned int &o_step) const | 
| void | write_particle_mesh_restart () const | 
| SerialMesh & | stitched_mesh () | 
| std::map< const std::size_t, std::vector< std::size_t > > | elem_neighbor_list () const | 
| const std::vector< std::size_t > | elem_neighbor_list (const Elem *elem) | 
| const std::vector< std::size_t > | local_elem_neighbor_list (const Elem *elem) | 
| virtual void | clear_kd_tree () | 
| virtual void | build_particle_neighbor_list (const Point &query_pt, const bool is_sorted, std::vector< std::pair< std::size_t, Real > > &IndicesDists) | 
| virtual void | build_particle_neighbor_list () | 
| virtual void | build_particle_neighbor_list_naively () | 
| virtual void | build_elem_neighbor_list (const Elem *elem, const bool is_sorted, std::vector< std::size_t > &n_list) | 
| virtual void | build_elem_neighbor_list () | 
| void | print_elem_neighbor_list (std::ostream &out=libMesh::out) const | 
| Protected Member Functions | |
| virtual void | construct_kd_tree () | 
| libMesh::ParticleMesh< KDDim >::ParticleMesh | ( | MeshBase & | mesh | ) | 
| libMesh::ParticleMesh< KDDim >::ParticleMesh | ( | MeshBase & | mesh, | 
| const Real & | search_radius_p, | ||
| const Real & | search_radius_e | ||
| ) | 
| libMesh::ParticleMesh< KDDim >::ParticleMesh | ( | MeshBase & | mesh, | 
| PMPeriodicBoundary & | pmpb, | ||
| const Real & | search_radius_p, | ||
| const Real & | search_radius_e | ||
| ) | 
| libMesh::ParticleMesh< KDDim >::~ParticleMesh | ( | ) | 
| 
 | inline | 
Add periodic boundary condition
| 
 | virtual | 
Build the neighbor list of partilces for a given elem Output the n_list, which gives the particle ids around elem's centroid within search_radius_e.
Note we don't need the distance between particles and elem centroid here, because this list is only used to evaluate element-wise load vector due to the neighboring particles.

| 
 | virtual | 
FIXME: Should we build this locally and globally(by allgather)???
Build the particle neighbor list for each element within the search radius this will generate the _elem_neighbor_list and _elem_neighbor_list_vector NOTE this takes advantage of KD Tree

| 
 | virtual | 
Build the neighbor list of a particle within the search_radius_p using KD tree. This includes the query particle itself, if the query_pt is also in the particle list. For example, the coords of the query_pt is the same as one of the particles in the list.
IndicesDists output pair<particle_id, distance> Note we need the particle distance to compute their interaction forces.

| 
 | virtual | 
Build the particle-particle neighbor list of particles within the search_radius_p using KD tree(logN complexity). This includes the original particle itself. In fact, the particle neighbor list has been built in "reinit()", so this is NOT USED!

| 
 | virtual | 
NOT ACTUALLY USED NOW, written only for the purpose of comparing with KD-Tree! Build the neighbor list of particles within the search radius directly! (N^2 complexity) This includes the original particle itself. NOTE that the particle list constructed from this function is not sorted!

| 
 | virtual | 
Clear the KD-Tree if need Note, different from the construct_kd_tree, this is a public function.

| 
 | protectedvirtual | 
Build & initialize the KD tree through nanoFlann, if needed.

| 
 | inline | 
Return the information of element neighbor list around the search radius _search_radius_e Mapping between elem id and element neighbor list of particles
| 
 | inline | 
Return the particle neighbor list around a given element around the search radius _search_radius_e
| void libMesh::ParticleMesh< KDDim >::generate_random_particles | ( | const std::size_t | N, | 
| const Real | bbox_XA, | ||
| const Real | bbox_XB, | ||
| const Real | bbox_YA, | ||
| const Real | bbox_YB, | ||
| const Real | bbox_ZA, | ||
| const Real | bbox_ZB | ||
| ) | 
Read chromatin data (cylinder particle) from the local file. mesh_type = "surface_mesh" or "volume_mesh" Generate random particle coordinate data within the bounding box. with r = 1.0, den = 1.0;

| void libMesh::ParticleMesh< KDDim >::generate_random_particles | ( | const std::size_t | N, | 
| const Point & | bbox_min, | ||
| const Point & | bbox_max | ||
| ) | 

| void libMesh::ParticleMesh< KDDim >::initial_particle_center_of_mass | ( | std::vector< Point > & | center0 | ) | const | 

| 
 | inline | 
Check if the neighbor list is constructed sortedly

| 
 | inline | 
Return the particle neighbor list around a given LOCAL element around the search radius _search_radius_e
| 
 | inline | 

| std::vector< Real > libMesh::ParticleMesh< KDDim >::mesh_size | ( | ) | const | 


| std::size_t libMesh::ParticleMesh< KDDim >::num_mesh_points | ( | ) | const | 

| 
 | inline | 
Return the total number of particles

| 
 | inline | 
Return the information of particles

| 
 | inline | 
Retrun the pointer to the periodic boundary for use

| void libMesh::ParticleMesh< KDDim >::print_elem_neighbor_list | ( | std::ostream & | out = libMesh::out | ) | const | 
print out the element neighbor list information

| void libMesh::ParticleMesh< KDDim >::print_particle_info | ( | ) | const | 
print out the particle information

| void libMesh::ParticleMesh< KDDim >::read_particles_data | ( | const std::string & | filename, | 
| const std::string & | particle_mesh_type | ||
| ) | 
Read particle coordinate data from file mesh_type = "surface_mesh" or "volume_mesh" defaulty, we turn the "Electrostatics" off

| void libMesh::ParticleMesh< KDDim >::read_particles_data_restart | ( | const std::string & | filename, | 
| const std::string & | particle_mesh_type | ||
| ) | 
Read particle data from restart mesh files only surface mesh is allowed


| void libMesh::ParticleMesh< KDDim >::read_particles_data_restart | ( | const std::string & | filename, | 
| const unsigned int & | o_step | ||
| ) | 
Read particle data from output_surface_node_*.csv
| void libMesh::ParticleMesh< KDDim >::reinit | ( | ) | 
Should we move this and above function to the ParticleMeshSystem??? Because we need dist-function to calculate the particle-wall force! Reinit the particle-mesh system. This includes: (1) build the particle-particle neighbor list according to search radius; (2) compute the force vectors on particles; (3) build the element-particle neighbor list according to search radius; (4) set the elem_id and proc_id for particles -— not incorporated: particle-wall repulsive force, which will be rebuilt in PMLinearSystem::reinit_particle_mesh()


| Real libMesh::ParticleMesh< KDDim >::search_radius | ( | const std::string & | p_e | ) | const | 
Return the search radius

| 
 | inline | 
Set the values of search radii
| SerialMesh & libMesh::ParticleMesh< KDDim >::stitched_mesh | ( | ) | 


| void libMesh::ParticleMesh< KDDim >::update_mesh | ( | const std::vector< Point > & | nodal_pos, | 
| const std::vector< Point > & | nodal_vel | ||
| ) | 

| void libMesh::ParticleMesh< KDDim >::volume_conservation | ( | ) | 

| void libMesh::ParticleMesh< KDDim >::write_particle | ( | const unsigned int & | step_id, | 
| const unsigned int & | o_step, | ||
| const Real & | real_time, | ||
| const std::vector< std::string > & | output_file, | ||
| unsigned int | comm_in_rank | ||
| ) | const | 


| void libMesh::ParticleMesh< KDDim >::write_particle_mesh | ( | const unsigned int & | o_step | ) | const | 


| void libMesh::ParticleMesh< KDDim >::write_particle_mesh_restart | ( | ) | const | 

| void libMesh::ParticleMesh< KDDim >::write_particle_trajectory | ( | const unsigned int & | o_step, | 
| unsigned int | comm_in_rank | ||
| ) | const | 

| void libMesh::ParticleMesh< KDDim >::write_surface_node | ( | const unsigned int & | o_step, | 
| unsigned int | comm_in_rank | ||
| ) | const | 

| void libMesh::ParticleMesh< KDDim >::write_time | ( | const unsigned int & | step_id, | 
| const unsigned int & | o_step, | ||
| const Real & | real_time, | ||
| unsigned int | comm_in_rank | ||
| ) | const | 

| void libMesh::ParticleMesh< KDDim >::zero_node_force | ( | ) | 

 1.8.11
 1.8.11