← Back to team overview

yade-dev team mailing list archive

[Branch ~yade-pkg/yade/git-trunk] Rev 4162: parallel removal of old interactions

 

------------------------------------------------------------
revno: 4162
committer: Bruno Chareyre <bruno.chareyre@xxxxxxxxxxx>
timestamp: Mon 2014-09-15 12:27:11 +0200
message:
  parallel removal of old interactions
modified:
  core/InteractionContainer.hpp


--
lp:yade
https://code.launchpad.net/~yade-pkg/yade/git-trunk

Your team Yade developers is subscribed to branch lp:yade.
To unsubscribe from this branch go to https://code.launchpad.net/~yade-pkg/yade/git-trunk/+edit-subscription
=== modified file 'core/InteractionContainer.hpp'
--- core/InteractionContainer.hpp	2014-07-03 17:20:40 +0000
+++ core/InteractionContainer.hpp	2014-09-15 10:27:11 +0000
@@ -103,12 +103,35 @@
 		*/
 		template<class T> size_t conditionalyEraseNonReal(const T& t, Scene* rb){
 			// beware iterators here, since erase is invalidating them. We need to iterate carefully, and keep in mind that erasing one interaction is moving the last one to the current position.
-			size_t initSize=currSize;
-		 	for (size_t linPos=0; linPos<currSize;){
-				const shared_ptr<Interaction>& i=linIntrs[linPos];
-				if(!i->isReal() && t.shouldBeErased(i->getId1(),i->getId2(),rb)) erase(i->getId1(),i->getId2(),linPos);
-				else linPos++;}
-			return initSize-currSize;
+			// For the parallel flavor we build the list to be erased in parallel, then it is erased sequentially. Still significant speedup since checking bounds is the most expensive part.
+			#ifdef YADE_OPENMP
+			if (omp_get_max_threads()<=1) {
+			#endif
+				size_t initSize=currSize;
+		 		for (size_t linPos=0; linPos<currSize;){
+					const shared_ptr<Interaction>& i=linIntrs[linPos];
+					if(!i->isReal() && t.shouldBeErased(i->getId1(),i->getId2(),rb)) erase(i->getId1(),i->getId2(),linPos);
+					else linPos++;}
+				return initSize-currSize;
+			#ifdef YADE_OPENMP
+			} else {
+				unsigned nThreads= omp_get_max_threads();
+				assert(nThreads>0);
+				std::vector<std::vector<Vector3i > >toErase;
+				toErase.resize(nThreads,std::vector<Vector3i >());
+				for (unsigned kk=0;  kk<nThreads; kk++) toErase[kk].reserve(1000);//A smarter value than 1000?			
+				size_t initSize=currSize;
+				#pragma omp parallel for schedule(static) num_threads(nThreads)
+				for (size_t linPos=0; linPos<currSize;linPos++){
+					const shared_ptr<Interaction>& i=linIntrs[linPos];
+					if(!i->isReal() && t.shouldBeErased(i->getId1(),i->getId2(),rb)) toErase[omp_get_thread_num()].push_back(Vector3i(i->getId1(),i->getId2(),linPos)) ;
+				}
+				for (int kk=nThreads-1;  kk>=0; kk--)
+					for (int ii(toErase[kk].size()-1); ii>=0;ii--)
+						erase(toErase[kk][ii][0],toErase[kk][ii][1],toErase[kk][ii][2]);
+				return initSize-currSize;
+			}
+		#endif
 		}
 		
 	// we must call Scene's ctor (and from Scene::postLoad), since we depend on the existing BodyContainer at that point.