NimbleSM
NimbleSM is a solid mechanics simulation code for dynamic systems
Loading...
Searching...
No Matches
nimble::BvhContactManager Class Reference

#include <bvh_contact_manager.h>

Inheritance diagram for nimble::BvhContactManager:

Public Member Functions

 BvhContactManager (std::shared_ptr< ContactInterface > interface, nimble::DataManager &data_manager, std::size_t _overdecomposition, std::string splitting_alg)
 
 BvhContactManager (const BvhContactManager &)=delete
 
 BvhContactManager (BvhContactManager &&)=delete
 
BvhContactManageroperator= (const BvhContactManager &)=delete
 
BvhContactManageroperator= (BvhContactManager &&)=delete
 
 ~BvhContactManager ()=default
 

Protected Member Functions

void ComputeParallelContactForce (int step, bool debug_output, nimble::Viewify< 2 > contact_force) override
 
void ComputeBoundingVolumes ()
 

Protected Attributes

bvh::collision_world m_world
 
bvh::collision_object * m_nodes
 
bvh::collision_object * m_faces
 
bvh::split_algorithm m_split_alg = bvh::split_algorithm::geom_axis
 
std::vector< NarrowphaseResultm_last_results
 

Constructor & Destructor Documentation

◆ BvhContactManager() [1/3]

nimble::BvhContactManager::BvhContactManager ( std::shared_ptr< ContactInterface > interface,
nimble::DataManager & data_manager,
std::size_t _overdecomposition,
std::string splitting_alg )
218 : ParallelContactManager(interface, data_manager),
219 m_world{_overdecomposition},
220 m_nodes{&m_world.create_collision_object()},
221 m_faces{&m_world.create_collision_object()}
222{
223 m_world.set_narrowphase_functor<ContactEntity>(NarrowphaseFunc{this});
224 if (splitting_alg == "geom_axis")
225 m_split_alg = bvh::split_algorithm::geom_axis;
226 else if (splitting_alg == "ml_geom_axis")
227 m_split_alg = bvh::split_algorithm::ml_geom_axis;
228 else if (splitting_alg == "clustering")
229 m_split_alg = bvh::split_algorithm::clustering;
230 else
231 throw std::runtime_error("invalid splitting algorithm");
232}
bvh::collision_object * m_nodes
Definition bvh_contact_manager.h:97
bvh::collision_world m_world
Definition bvh_contact_manager.h:96
bvh::split_algorithm m_split_alg
Definition bvh_contact_manager.h:100
bvh::collision_object * m_faces
Definition bvh_contact_manager.h:98

◆ BvhContactManager() [2/3]

nimble::BvhContactManager::BvhContactManager ( const BvhContactManager & )
delete

◆ BvhContactManager() [3/3]

nimble::BvhContactManager::BvhContactManager ( BvhContactManager && )
delete

◆ ~BvhContactManager()

nimble::BvhContactManager::~BvhContactManager ( )
default

Member Function Documentation

◆ ComputeBoundingVolumes()

void nimble::BvhContactManager::ComputeBoundingVolumes ( )
protected
236{
237 Kokkos::parallel_for( contact_nodes_d_.extent( 0 ), KOKKOS_LAMBDA( int _i ) {
238 auto &node = contact_nodes_d_( _i );
239 const double inflation_length = node.inflation_factor * node.char_len_;
240 ContactEntity::vertex* v = reinterpret_cast<ContactEntity::vertex*>(&node.coord_1_x_);
241 node.kdop_ = bvh::bphase_kdop::from_sphere(*v, inflation_length);
242 } );
243 Kokkos::parallel_for( contact_faces_d_.extent( 0 ), KOKKOS_LAMBDA( int _i ) {
244 auto &face = contact_faces_d_( _i );
245 const double inflation_length = face.inflation_factor * face.char_len_;
246 ContactEntity::vertex* v = reinterpret_cast<ContactEntity::vertex*>(&face.coord_1_x_);
247 face.kdop_ = bvh::bphase_kdop::from_vertices(v, v + 3, inflation_length);
248 } );
249}

◆ ComputeParallelContactForce()

void nimble::BvhContactManager::ComputeParallelContactForce ( int step,
bool debug_output,
nimble::Viewify< 2 > contact_force )
overrideprotected
253{
254 auto model_ptr = this->data_manager_.GetModelData();
255 auto model_data = &dynamic_cast<nimble_kokkos::ModelData &>(*model_ptr);
256
257 auto field_ids = this->data_manager_.GetFieldIDs();
258 auto displacement_d = model_data->GetDeviceVectorNodeData(field_ids.displacement);
259
260 auto contact_force_h = model_data->GetHostVectorNodeData(field_ids.contact_force);
261 auto contact_force_d = model_data->GetDeviceVectorNodeData(field_ids.contact_force);
262 Kokkos::deep_copy(contact_force_d, 0.0);
263
264 this->ApplyDisplacements(displacement_d);
265
266 total_search_time.Start();
267 m_world.start_iteration();
268
269 // Update collision objects, this will build the trees
271
272 this->startTimer("BVH::SetEntity::Nodes");
273 m_nodes->set_entity_data(contact_nodes_d_, m_split_alg);
274 this->stopTimer("BVH::SetEntity::Nodes");
275
276 this->startTimer("BVH::SetEntity::Faces");
277 m_faces->set_entity_data(contact_faces_d_, m_split_alg);
278 this->stopTimer("BVH::SetEntity::Faces");
279
280 m_nodes->init_broadphase();
281 m_faces->init_broadphase();
282
283 m_faces->broadphase(*m_nodes);
284
285 m_last_results.clear();
286 m_faces->for_each_result<NarrowphaseResult>(
287 [this](const NarrowphaseResult& _res) { m_last_results.emplace_back(_res); });
288
289 m_world.finish_iteration();
290 total_search_time.Stop();
291
292 total_enforcement_time.Start();
293 // Update contact entities
294 for (auto&& r : m_last_results) {
295 if (r.node) {
296 if (r.local_index >= contact_nodes_d_.extent( 0 ))
297 std::cerr << "contact node index " << r.local_index << " is out of bounds (" << contact_nodes_d_.extent( 0 ) << ")\n";
298 auto& node = contact_nodes_d_(r.local_index);
299 node.set_contact_status(true);
300 node.SetNodalContactForces(r.contact_force);
301 node.ScatterForceToContactManagerForceVector(force_d_);
302 } else {
303 if (r.local_index >= contact_faces_d_.extent( 0 ))
304 std::cerr << "contact face index " << r.local_index << " is out of bounds (" << contact_faces_d_.extent( 0 ) << ")\n";
305 auto& face = contact_faces_d_(r.local_index);
306 face.set_contact_status(true);
307 face.SetNodalContactForces(r.contact_force, r.bary);
308 face.ScatterForceToContactManagerForceVector(force_d_);
309 }
310 }
311 total_num_contacts += m_last_results.size();
312 total_enforcement_time.Stop();
313
314 this->GetForces(contact_force_d);
315 Kokkos::deep_copy(contact_force_h, contact_force_d);
316
317 constexpr int vector_dim = 3;
318 auto myVectorCommunicator = this->data_manager_.GetVectorCommunicator();
319 myVectorCommunicator->VectorReduction(vector_dim, contact_force_h);
320}
void ComputeBoundingVolumes()
Definition bvh_contact_manager.cc:235
std::vector< NarrowphaseResult > m_last_results
Definition bvh_contact_manager.h:102

◆ operator=() [1/2]

BvhContactManager & nimble::BvhContactManager::operator= ( BvhContactManager && )
delete

◆ operator=() [2/2]

BvhContactManager & nimble::BvhContactManager::operator= ( const BvhContactManager & )
delete

Member Data Documentation

◆ m_faces

bvh::collision_object* nimble::BvhContactManager::m_faces
protected

◆ m_last_results

std::vector<NarrowphaseResult> nimble::BvhContactManager::m_last_results
protected

◆ m_nodes

bvh::collision_object* nimble::BvhContactManager::m_nodes
protected

◆ m_split_alg

bvh::split_algorithm nimble::BvhContactManager::m_split_alg = bvh::split_algorithm::geom_axis
protected

◆ m_world

bvh::collision_world nimble::BvhContactManager::m_world
protected

The documentation for this class was generated from the following files: