Salome HOME
Merge from V6_main_20120808 08Aug12
[modules/med.git] / src / ParaMEDMEM / OverlapElementLocator.cxx
1 // Copyright (C) 2007-2012  CEA/DEN, EDF R&D
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 #include "OverlapElementLocator.hxx"
21
22 #include "CommInterface.hxx"
23 #include "Topology.hxx"
24 #include "BlockTopology.hxx"
25 #include "ParaFIELD.hxx"
26 #include "ParaMESH.hxx"
27 #include "ProcessorGroup.hxx"
28 #include "MPIProcessorGroup.hxx"
29 #include "OverlapInterpolationMatrix.hxx"
30 #include "MEDCouplingFieldDouble.hxx"
31 #include "MEDCouplingFieldDiscretization.hxx"
32 #include "DirectedBoundingBox.hxx"
33 #include "InterpKernelAutoPtr.hxx"
34
35 #include <limits>
36
37 using namespace std;
38
39 namespace ParaMEDMEM 
40
41   OverlapElementLocator::OverlapElementLocator(const ParaFIELD *sourceField, const ParaFIELD *targetField, const ProcessorGroup& group)
42     : _local_source_field(sourceField),
43       _local_target_field(targetField),
44       _local_source_mesh(0),
45       _local_target_mesh(0),
46       _domain_bounding_boxes(0),
47       _group(group)
48   { 
49     if(_local_source_field)
50       _local_source_mesh=_local_source_field->getSupport()->getCellMesh();
51     if(_local_target_field)
52       _local_target_mesh=_local_target_field->getSupport()->getCellMesh();
53     _comm=getCommunicator();
54     computeBoundingBoxes();
55   }
56
57   OverlapElementLocator::~OverlapElementLocator()
58   {
59     delete [] _domain_bounding_boxes;
60   }
61
62   const MPI_Comm *OverlapElementLocator::getCommunicator() const
63   {
64     const MPIProcessorGroup* group=static_cast<const MPIProcessorGroup*>(&_group);
65     return group->getComm();
66   }
67
68   void OverlapElementLocator::computeBoundingBoxes()
69   {
70     CommInterface comm_interface=_group.getCommInterface();
71     const MPIProcessorGroup* group=static_cast<const MPIProcessorGroup*> (&_group);
72     _local_space_dim=0;
73     if(_local_source_mesh)
74       _local_space_dim=_local_source_mesh->getSpaceDimension();
75     else
76       _local_space_dim=_local_target_mesh->getSpaceDimension();
77     //
78     const MPI_Comm* comm = group->getComm();
79     int bbSize=2*2*_local_space_dim;//2 (for source/target) 2 (min/max)
80     _domain_bounding_boxes=new double[bbSize*_group.size()];
81     INTERP_KERNEL::AutoPtr<double> minmax=new double[bbSize];
82     //Format minmax : Xmin_src,Xmax_src,Ymin_src,Ymax_src,Zmin_src,Zmax_src,Xmin_trg,Xmax_trg,Ymin_trg,Ymax_trg,Zmin_trg,Zmax_trg
83     if(_local_source_mesh)
84       _local_source_mesh->getBoundingBox(minmax);
85     else
86       {
87         for(int i=0;i<_local_space_dim;i++)
88           {
89             minmax[i*2]=std::numeric_limits<double>::max();
90             minmax[i*2+1]=-std::numeric_limits<double>::max();
91           }
92       }
93     if(_local_target_mesh)
94       _local_target_mesh->getBoundingBox(minmax+2*_local_space_dim);
95     else
96       {
97         for(int i=0;i<_local_space_dim;i++)
98           {
99             minmax[i*2+2*_local_space_dim]=std::numeric_limits<double>::max();
100             minmax[i*2+1+2*_local_space_dim]=-std::numeric_limits<double>::max();
101           }
102       }
103     comm_interface.allGather(minmax, bbSize, MPI_DOUBLE,
104                              _domain_bounding_boxes,bbSize, MPI_DOUBLE, 
105                              *comm);
106   
107     // Computation of all pairs needing an interpolation pairs are duplicated now !
108     
109     _proc_pairs.clear();//first is source second is target
110     _proc_pairs.resize(_group.size());
111     for(int i=0;i<_group.size();i++)
112       for(int j=0;j<_group.size();j++)
113         {
114           if(intersectsBoundingBox(i,j))
115             _proc_pairs[i].push_back(j);
116         }
117
118     // OK now let's assigning as balanced as possible, job to each proc of group
119     std::vector< std::vector< std::pair<int,int> > > pairsToBeDonePerProc(_group.size());
120     int i=0;
121     for(std::vector< std::vector< int > >::const_iterator it1=_proc_pairs.begin();it1!=_proc_pairs.end();it1++,i++)
122       for(std::vector< int >::const_iterator it2=(*it1).begin();it2!=(*it1).end();it2++)
123         {
124           if(pairsToBeDonePerProc[i].size()<=pairsToBeDonePerProc[*it2].size())//it includes the fact that i==*it2
125             pairsToBeDonePerProc[i].push_back(std::pair<int,int>(i,*it2));
126           else
127             pairsToBeDonePerProc[*it2].push_back(std::pair<int,int>(i,*it2));
128         }
129     //Keeping todo list of current proc. _to_do_list contains a set of pair where at least _group.myRank() appears once.
130     //This proc will be in charge to perform interpolation of any of element of '_to_do_list'
131     //If _group.myRank()==myPair.first, current proc should fetch target mesh of myPair.second (if different from _group.myRank()).
132     //If _group.myRank()==myPair.second, current proc should fetch source mesh of myPair.second.
133     
134     int myProcId=_group.myRank();
135     _to_do_list=pairsToBeDonePerProc[myProcId];
136
137     //Feeding now '_procs_to_send'. A same id can appears twice. The second parameter in pair means what to send true=source, false=target
138     _procs_to_send.clear();
139     for(int i=_group.size()-1;i>=0;i--)
140       if(i!=myProcId)
141         {
142           const std::vector< std::pair<int,int> >& anRemoteProcToDoList=pairsToBeDonePerProc[i];
143           for(std::vector< std::pair<int,int> >::const_iterator it=anRemoteProcToDoList.begin();it!=anRemoteProcToDoList.end();it++)
144             {
145               if((*it).first==myProcId)
146                 _procs_to_send.push_back(std::pair<int,bool>(i,true));
147               if((*it).second==myProcId)
148                 _procs_to_send.push_back(std::pair<int,bool>(i,false));
149             }
150         }
151   }
152
153   /*!
154    * The aim of this method is to perform the communication to get data corresponding to '_to_do_list' attribute.
155    * The principle is the following : if proc n1 and n2 need to perform a cross sending with n1<n2, then n1 will send first and receive then.
156    */
157   void OverlapElementLocator::exchangeMeshes(OverlapInterpolationMatrix& matrix)
158   {
159     int myProcId=_group.myRank();
160     //starting to receive every procs whose id is lower than myProcId.
161     std::vector< std::pair<int,int> > toDoListForFetchRemaining;
162     for(std::vector< std::pair<int,int> >::const_iterator it=_to_do_list.begin();it!=_to_do_list.end();it++)
163       {
164         if((*it).first!=(*it).second)
165           {
166             if((*it).first==myProcId)
167               {
168                 if((*it).second<myProcId)
169                   receiveRemoteMesh((*it).second,false);
170                 else
171                   toDoListForFetchRemaining.push_back(std::pair<int,int>((*it).first,(*it).second));
172               }
173             else
174               {//(*it).second==myProcId
175                 if((*it).first<myProcId)
176                   receiveRemoteMesh((*it).first,true);
177                 else
178                   toDoListForFetchRemaining.push_back(std::pair<int,int>((*it).first,(*it).second));
179               }
180           }
181       }
182     //sending source or target mesh to remote procs
183     for(std::vector< std::pair<int,bool> >::const_iterator it2=_procs_to_send.begin();it2!=_procs_to_send.end();it2++)
184       sendLocalMeshTo((*it2).first,(*it2).second,matrix);
185     //fetching remaining meshes
186     for(std::vector< std::pair<int,int> >::const_iterator it=toDoListForFetchRemaining.begin();it!=toDoListForFetchRemaining.end();it++)
187       {
188         if((*it).first!=(*it).second)
189           {
190             if((*it).first==myProcId)
191               receiveRemoteMesh((*it).second,false);
192             else//(*it).second==myProcId
193               receiveRemoteMesh((*it).first,true);
194           }
195       }
196   }
197   
198   std::string OverlapElementLocator::getSourceMethod() const
199   {
200     return _local_source_field->getField()->getDiscretization()->getStringRepr();
201   }
202
203   std::string OverlapElementLocator::getTargetMethod() const
204   {
205     return _local_target_field->getField()->getDiscretization()->getStringRepr();
206   }
207
208   const MEDCouplingPointSet *OverlapElementLocator::getSourceMesh(int procId) const
209   {
210     int myProcId=_group.myRank();
211     if(myProcId==procId)
212       return _local_source_mesh;
213     std::pair<int,bool> p(procId,true);
214     std::map<std::pair<int,bool>, MEDCouplingAutoRefCountObjectPtr< MEDCouplingPointSet > >::const_iterator it=_remote_meshes.find(p);
215     return (*it).second;
216   }
217
218   const DataArrayInt *OverlapElementLocator::getSourceIds(int procId) const
219   {
220     int myProcId=_group.myRank();
221     if(myProcId==procId)
222       return 0;
223     std::pair<int,bool> p(procId,true);
224     std::map<std::pair<int,bool>, MEDCouplingAutoRefCountObjectPtr< DataArrayInt > >::const_iterator it=_remote_elems.find(p);
225     return (*it).second;
226   }
227
228   const MEDCouplingPointSet *OverlapElementLocator::getTargetMesh(int procId) const
229   {
230     int myProcId=_group.myRank();
231     if(myProcId==procId)
232       return _local_target_mesh;
233     std::pair<int,bool> p(procId,false);
234     std::map<std::pair<int,bool>, MEDCouplingAutoRefCountObjectPtr< MEDCouplingPointSet > >::const_iterator it=_remote_meshes.find(p);
235     return (*it).second;
236   }
237
238   const DataArrayInt *OverlapElementLocator::getTargetIds(int procId) const
239   {
240     int myProcId=_group.myRank();
241     if(myProcId==procId)
242       return 0;
243     std::pair<int,bool> p(procId,false);
244     std::map<std::pair<int,bool>, MEDCouplingAutoRefCountObjectPtr< DataArrayInt > >::const_iterator it=_remote_elems.find(p);
245     return (*it).second;
246   }
247
248   bool OverlapElementLocator::intersectsBoundingBox(int isource, int itarget) const
249   {
250     const double *source_bb=_domain_bounding_boxes+isource*2*2*_local_space_dim;
251     const double *target_bb=_domain_bounding_boxes+itarget*2*2*_local_space_dim+2*_local_space_dim;
252
253     for (int idim=0; idim < _local_space_dim; idim++)
254       {
255         const double eps = -1e-12;//tony to change
256         bool intersects = (target_bb[idim*2]<source_bb[idim*2+1]+eps)
257           && (source_bb[idim*2]<target_bb[idim*2+1]+eps);
258         if (!intersects)
259           return false; 
260       }
261     return true;
262   }
263
264   /*!
265    * This methods sends local source if 'sourceOrTarget'==True to proc 'procId'.
266    * This methods sends local target if 'sourceOrTarget'==False to proc 'procId'.
267    *
268    * This method prepares the matrix too, for matrix assembling and future matrix-vector computation.
269    */
270   void OverlapElementLocator::sendLocalMeshTo(int procId, bool sourceOrTarget, OverlapInterpolationMatrix& matrix) const
271   {
272    vector<int> elems;
273    //int myProcId=_group.myRank();
274    const double *distant_bb=0;
275    MEDCouplingPointSet *local_mesh=0;
276    const ParaFIELD *field=0;
277    if(sourceOrTarget)//source for local but target for distant
278      {
279        distant_bb=_domain_bounding_boxes+procId*2*2*_local_space_dim+2*_local_space_dim;
280        local_mesh=_local_source_mesh;
281        field=_local_source_field;
282      }
283    else//target for local but source for distant
284      {
285        distant_bb=_domain_bounding_boxes+procId*2*2*_local_space_dim;
286        local_mesh=_local_target_mesh;
287        field=_local_target_field;
288      }
289    local_mesh->getCellsInBoundingBox(distant_bb,getBoundingBoxAdjustment(),elems);
290    DataArrayInt *idsToSend;
291    MEDCouplingPointSet *send_mesh=static_cast<MEDCouplingPointSet *>(field->getField()->buildSubMeshData(&elems[0],&elems[elems.size()],idsToSend));
292    if(sourceOrTarget)
293      matrix.keepTracksOfSourceIds(procId,idsToSend);//Case#1 in Step2 of main algorithm.
294    else
295      matrix.keepTracksOfTargetIds(procId,idsToSend);//Case#0 in Step2 of main algorithm.
296    sendMesh(procId,send_mesh,idsToSend);
297    send_mesh->decrRef();
298    idsToSend->decrRef();
299   }
300
301   /*!
302    * This method recieves source remote mesh on proc 'procId' if sourceOrTarget==True
303    * This method recieves target remote mesh on proc 'procId' if sourceOrTarget==False
304    */
305   void OverlapElementLocator::receiveRemoteMesh(int procId, bool sourceOrTarget)
306   {
307     DataArrayInt *da=0;
308     MEDCouplingPointSet *m=0;
309     receiveMesh(procId,m,da);
310     std::pair<int,bool> p(procId,sourceOrTarget);
311     _remote_meshes[p]=m;
312     _remote_elems[p]=da;
313   }
314
315   void OverlapElementLocator::sendMesh(int procId, const MEDCouplingPointSet *mesh, const DataArrayInt *idsToSend) const
316   {
317     CommInterface comInterface=_group.getCommInterface();
318     // First stage : exchanging sizes
319     vector<double> tinyInfoLocalD;//tinyInfoLocalD not used for the moment
320     vector<int> tinyInfoLocal;
321     vector<string> tinyInfoLocalS;
322     mesh->getTinySerializationInformation(tinyInfoLocalD,tinyInfoLocal,tinyInfoLocalS);
323     const MPI_Comm *comm=getCommunicator();
324     //
325     int lgth[2];
326     lgth[0]=tinyInfoLocal.size();
327     lgth[1]=idsToSend->getNbOfElems();
328     comInterface.send(&lgth,2,MPI_INT,procId,1140,*_comm);
329     comInterface.send(&tinyInfoLocal[0],tinyInfoLocal.size(),MPI_INT,procId,1141,*comm);
330     //
331     DataArrayInt *v1Local=0;
332     DataArrayDouble *v2Local=0;
333     mesh->serialize(v1Local,v2Local);
334     comInterface.send(v1Local->getPointer(),v1Local->getNbOfElems(),MPI_INT,procId,1142,*comm);
335     comInterface.send(v2Local->getPointer(),v2Local->getNbOfElems(),MPI_DOUBLE,procId,1143,*comm);
336     //finished for mesh, ids now
337     comInterface.send((int *)idsToSend->getConstPointer(),lgth[1],MPI_INT,procId,1144,*comm);
338     //
339     v1Local->decrRef();
340     v2Local->decrRef();
341   }
342
343   void OverlapElementLocator::receiveMesh(int procId, MEDCouplingPointSet* &mesh, DataArrayInt *&ids) const
344   {
345     int lgth[2];
346     MPI_Status status;
347     const MPI_Comm *comm=getCommunicator();
348     CommInterface comInterface=_group.getCommInterface();
349     comInterface.recv(lgth,2,MPI_INT,procId,1140,*_comm,&status);
350     std::vector<int> tinyInfoDistant(lgth[0]);
351     ids=DataArrayInt::New();
352     ids->alloc(lgth[1],1);
353     comInterface.recv(&tinyInfoDistant[0],lgth[0],MPI_INT,procId,1141,*comm,&status);
354     mesh=MEDCouplingPointSet::BuildInstanceFromMeshType((MEDCouplingMeshType)tinyInfoDistant[0]);
355     std::vector<std::string> unusedTinyDistantSts;
356     vector<double> tinyInfoDistantD(1);//tinyInfoDistantD not used for the moment
357     DataArrayInt *v1Distant=DataArrayInt::New();
358     DataArrayDouble *v2Distant=DataArrayDouble::New();
359     mesh->resizeForUnserialization(tinyInfoDistant,v1Distant,v2Distant,unusedTinyDistantSts);
360     comInterface.recv(v1Distant->getPointer(),v1Distant->getNbOfElems(),MPI_INT,procId,1142,*comm,&status);
361     comInterface.recv(v2Distant->getPointer(),v2Distant->getNbOfElems(),MPI_DOUBLE,procId,1143,*comm,&status);
362     mesh->unserialization(tinyInfoDistantD,tinyInfoDistant,v1Distant,v2Distant,unusedTinyDistantSts);
363     //finished for mesh, ids now
364     comInterface.recv(ids->getPointer(),lgth[1],MPI_INT,procId,1144,*comm,&status);
365     //
366     v1Distant->decrRef();
367     v2Distant->decrRef();
368   }
369 }