Salome HOME
d05698790a44b50d78d530221d9296a1f6b4cf7b
[tools/medcoupling.git] / src / ParaMEDMEM / MPIAccessDEC.cxx
1 //  Copyright (C) 2007-2008  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 "MPIAccessDEC.hxx"
21
22 using namespace std;
23
24 namespace ParaMEDMEM
25 {    
26
27   /*!
28     This constructor creates an MPIAccessDEC which has \a source_group as a working side 
29     and  \a target_group as an idle side. 
30     The constructor must be called synchronously on all processors of both processor groups.
31
32     \param source_group working side ProcessorGroup
33     \param target_group lazy side ProcessorGroup
34     \param Asynchronous Communication mode (default asynchronous)
35     \param nStepBefore Number of Time step needed for the interpolation before current time
36     \param nStepAfter Number of Time step needed for the interpolation after current time
37
38   */
39
40   MPIAccessDEC::MPIAccessDEC( const ProcessorGroup& source_group,
41                               const ProcessorGroup& target_group,
42                               bool Asynchronous )
43   {
44
45     ProcessorGroup * union_group = source_group.fuse(target_group) ;  
46     int i ;
47     std::set<int> procs;
48     for ( i = 0 ; i < union_group->size() ; i++ )
49       {
50         procs.insert(i) ;
51       }
52     MPIProcessorGroup *mpilg = (MPIProcessorGroup *)&source_group;
53     _MPI_union_group = new ParaMEDMEM::MPIProcessorGroup( union_group->getCommInterface(),procs,mpilg->getWorldComm());
54     delete union_group ;
55     _my_rank = _MPI_union_group->myRank() ;
56     _group_size = _MPI_union_group->size() ;
57     _MPI_access = new MPIAccess( _MPI_union_group ) ;
58     _asynchronous = Asynchronous ;
59     _time_messages = new vector< vector< TimeMessage > > ;
60     _time_messages->resize( _group_size ) ;
61     _out_of_time = new vector< bool > ;
62     _out_of_time->resize( _group_size ) ;
63     _data_messages_recv_count = new vector< int > ;
64     _data_messages_recv_count->resize( _group_size ) ;
65     for ( i = 0 ; i < _group_size ; i++ )
66       {
67         (*_out_of_time)[i] = false ;
68         (*_data_messages_recv_count)[i] = 0 ;
69       }
70     _data_messages_type = new vector< MPI_Datatype > ;
71     _data_messages_type->resize( _group_size ) ;
72     _data_messages = new vector< vector< void * > > ;
73     _data_messages->resize( _group_size ) ;
74     _time_interpolator = NULL ;
75     _map_of_send_buffers = new map< int , SendBuffStruct * > ;
76     cout << "MPIAccessDEC" << _my_rank << " Asynchronous " << _asynchronous << endl ;
77   }
78
79   MPIAccessDEC::~MPIAccessDEC()
80   {
81     checkFinalSent() ;
82     checkFinalRecv() ;
83     delete _MPI_union_group ;
84     delete _MPI_access ;
85     if ( _time_interpolator )
86       delete _time_interpolator ;
87     if ( _time_messages )
88       delete _time_messages ;
89     if ( _out_of_time )
90       delete _out_of_time ;
91     if ( _data_messages_recv_count )
92       delete _data_messages_recv_count ;
93     if ( _data_messages_type )
94       delete _data_messages_type ;
95     if ( _data_messages )
96       delete _data_messages ;
97     if ( _map_of_send_buffers )
98       delete _map_of_send_buffers ;
99   } 
100
101   void MPIAccessDEC::setTimeInterpolator( TimeInterpolationMethod aTimeInterp ,
102                                           double InterpPrecision, int nStepBefore,
103                                           int nStepAfter )
104   {
105     cout << "MPIAccessDEC::SetTimeInterpolator" << _my_rank << " Asynchronous "
106          << _asynchronous << " TimeInterpolationMethod " << aTimeInterp
107          << " InterpPrecision " << InterpPrecision << " nStepBefore " << nStepBefore
108          << " nStepAfter " << nStepAfter << endl ;
109     if ( _time_interpolator )
110       delete _time_interpolator ;
111     switch ( aTimeInterp )
112       {
113       case WithoutTimeInterp :
114         _time_interpolator = NULL ;
115         _n_step_before = 0 ;
116         _n_step_after = 0 ;
117         break ;
118       case LinearTimeInterp :
119         _time_interpolator = new LinearTimeInterpolator( InterpPrecision , nStepBefore ,
120                                                          nStepAfter ) ;
121         _n_step_before = nStepBefore ;
122         _n_step_after = nStepAfter ;
123         int i ;
124         for ( i = 0 ; i < _group_size ; i++ )
125           {
126             (*_time_messages)[ i ].resize( _n_step_before + _n_step_after ) ;
127             (*_data_messages)[ i ].resize( _n_step_before + _n_step_after ) ;
128             int j ;
129             for ( j = 0 ; j < _n_step_before + _n_step_after ; j++ )
130               {
131                 (*_time_messages)[ i ][ j ].time = -1 ;
132                 (*_time_messages)[ i ][ j ].deltatime = -1 ;
133                 (*_data_messages)[ i ][ j ] = NULL ;
134               }
135           }
136         break ;
137       }
138   }
139
140   /*!
141     Send sendcount datas from sendbuf[offset] with type sendtype to target of IntraCommunicator
142     (Internal Protected method)
143
144     Returns the request identifier SendRequestId
145
146   */
147   int MPIAccessDEC::send( void* sendbuf, int sendcount , int offset ,
148                           MPI_Datatype sendtype , int target , int &SendRequestId )
149   {
150     int sts ;
151     if ( _asynchronous )
152       {
153         if ( sendtype == MPI_INT )
154           {
155             sts = _MPI_access->ISend( &((int *) sendbuf)[offset] , sendcount , sendtype ,
156                                       target , SendRequestId ) ;
157           }
158         else
159           {
160             sts = _MPI_access->ISend( &((double *) sendbuf)[offset] , sendcount , sendtype ,
161                                       target , SendRequestId ) ;
162           }
163       }
164     else
165       {
166         if ( sendtype == MPI_INT )
167           {
168             sts = _MPI_access->send( &((int *) sendbuf)[offset] , sendcount , sendtype ,
169                                      target , SendRequestId ) ;
170           }
171         else
172           {
173             sts = _MPI_access->send( &((double *) sendbuf)[offset] , sendcount , sendtype ,
174                                      target , SendRequestId ) ;
175           }
176       }
177     return sts ;
178   }
179
180   /*!
181     Receive recvcount datas to recvbuf[offset] with type recvtype from target of IntraCommunicator
182     (Internal Protected method)
183
184     Returns the request identifier RecvRequestId
185
186   */
187   int MPIAccessDEC::recv( void* recvbuf, int recvcount , int offset ,
188                           MPI_Datatype recvtype , int target , int &RecvRequestId )
189   {
190     int sts ;
191     if ( _asynchronous )
192       {
193         if ( recvtype == MPI_INT )
194           {
195             sts = _MPI_access->IRecv( &((int *) recvbuf)[offset] , recvcount , recvtype ,
196                                       target , RecvRequestId ) ;
197           }
198         else
199           {
200             sts = _MPI_access->IRecv( &((double *) recvbuf)[offset] , recvcount , recvtype ,
201                                       target , RecvRequestId ) ;
202           }
203       }
204     else
205       {
206         if ( recvtype == MPI_INT )
207           {
208             sts = _MPI_access->recv( &((int *) recvbuf)[offset] , recvcount , recvtype ,
209                                      target , RecvRequestId ) ;
210           }
211         else
212           {
213             sts = _MPI_access->recv( &((double *) recvbuf)[offset] , recvcount , recvtype ,
214                                      target , RecvRequestId ) ;
215           }
216       }
217     return sts ;
218   }
219
220   /*!
221     Send sendcount datas from sendbuf[offset] with type sendtype to target of IntraCommunicator
222     Receive recvcount datas to recvbuf[offset] with type recvtype from target of IntraCommunicator
223     (Internal Protected method)
224
225     Returns the request identifier SendRequestId
226     Returns the request identifier RecvRequestId
227
228   */
229   int MPIAccessDEC::sendRecv( void* sendbuf, int sendcount , int sendoffset ,
230                               MPI_Datatype sendtype ,
231                               void* recvbuf, int recvcount , int recvoffset ,
232                               MPI_Datatype recvtype , int target ,
233                               int &SendRequestId , int &RecvRequestId )
234   {
235     int sts ;
236     if ( _asynchronous )
237       {
238         if ( sendtype == MPI_INT )
239           {
240             if ( recvtype == MPI_INT )
241               {
242                 sts = _MPI_access->ISendRecv( &((int *) sendbuf)[sendoffset] , sendcount ,
243                                               sendtype , target , SendRequestId ,
244                                               &((int *) recvbuf)[recvoffset] , recvcount ,
245                                               recvtype , target , RecvRequestId ) ;
246               }
247             else
248               {
249                 sts = _MPI_access->ISendRecv( &((int *) sendbuf)[sendoffset] , sendcount ,
250                                               sendtype , target , SendRequestId ,
251                                               &((double *) recvbuf)[recvoffset] ,
252                                               recvcount , recvtype , target , RecvRequestId ) ;
253               }
254           }
255         else
256           {
257             if ( recvtype == MPI_INT )
258               {
259                 sts = _MPI_access->ISendRecv( &((double *) sendbuf)[sendoffset] , sendcount ,
260                                               sendtype , target , SendRequestId ,
261                                               &((int *) recvbuf)[recvoffset] ,
262                                               recvcount , recvtype , target , RecvRequestId ) ;
263               }
264             else
265               {
266                 sts = _MPI_access->ISendRecv( &((double *) sendbuf)[sendoffset] , sendcount ,
267                                               sendtype , target , SendRequestId ,
268                                               &((double *) recvbuf)[recvoffset] ,
269                                               recvcount , recvtype , target , RecvRequestId ) ;
270               }
271           }
272       }
273     else
274       {
275         if ( sendtype == MPI_INT )
276           {
277             if ( recvtype == MPI_INT )
278               {
279                 sts = _MPI_access->sendRecv( &((int *) sendbuf)[sendoffset] , sendcount ,
280                                              sendtype , target , SendRequestId ,
281                                              &((int *) recvbuf)[recvoffset] , recvcount ,
282                                              recvtype , target , RecvRequestId ) ;
283               }
284             else
285               {
286                 sts = _MPI_access->sendRecv( &((int *) sendbuf)[sendoffset] , sendcount ,
287                                              sendtype , target , SendRequestId ,
288                                              &((double *) recvbuf)[recvoffset] ,
289                                              recvcount , recvtype , target , RecvRequestId ) ;
290               }
291           }
292         else
293           {
294             if ( recvtype == MPI_INT )
295               {
296                 sts = _MPI_access->sendRecv( &((double *) sendbuf)[sendoffset] , sendcount ,
297                                              sendtype , target , SendRequestId ,
298                                              &((int *) recvbuf)[recvoffset] ,
299                                              recvcount , recvtype , target , RecvRequestId ) ;
300               }
301             else
302               {
303                 cout << "SendRecv" << _my_rank << " target " << target << " sendbuf "
304                      << &((double *) sendbuf)[sendoffset] << " sendcount " << sendcount
305                      << " recvbuf " << &((double *) recvbuf)[recvoffset] << " recvcount "
306                      << recvcount << endl ;
307                 sts = _MPI_access->sendRecv( &((double *) sendbuf)[sendoffset] , sendcount ,
308                                              sendtype , target , SendRequestId ,
309                                              &((double *) recvbuf)[recvoffset] ,
310                                              recvcount , recvtype , target , RecvRequestId ) ;
311               }
312           }
313       }
314     return sts ;
315   }
316
317   /*!
318     Send sendcount datas from sendbuf[offset] with type sendtype to all targets of IntraCommunicator
319     Receive recvcount datas to recvbuf[offset] with type recvtype from all targets of IntraCommunicator
320
321   */
322   int MPIAccessDEC::allToAll( void* sendbuf, int sendcount, MPI_Datatype sendtype ,
323                               void* recvbuf, int recvcount, MPI_Datatype recvtype )
324   {
325     if ( _time_interpolator )
326       {
327         return allToAllTime( sendbuf, sendcount, sendtype , recvbuf, recvcount, recvtype ) ;
328       }
329     int sts ;
330     int target ;
331     int sendoffset = 0 ;
332     int recvoffset = 0 ;
333     int SendRequestId ;
334     int RecvRequestId ;
335
336     //Free of SendBuffers 
337     if ( _asynchronous )
338       checkSent() ;
339
340     //DoSend + DoRecv : SendRecv
341     SendBuffStruct * aSendDataStruct = NULL ;
342     if ( _asynchronous && sendbuf )
343       {
344         aSendDataStruct = new SendBuffStruct ;
345         aSendDataStruct->SendBuffer = sendbuf ;
346         aSendDataStruct->Counter = 0 ;
347         aSendDataStruct->DataType = sendtype ;
348       }
349     for ( target = 0 ; target < _group_size ; target++ )
350       {
351         sts = sendRecv( sendbuf , sendcount , sendoffset , sendtype ,
352                         recvbuf , recvcount , recvoffset , recvtype ,
353                         target , SendRequestId , RecvRequestId ) ;
354         if ( _asynchronous && sendbuf && sendcount )
355           {
356             aSendDataStruct->Counter += 1 ;
357             (*_map_of_send_buffers)[ SendRequestId ] = aSendDataStruct ;
358           }
359         sendoffset += sendcount ;
360         recvoffset += recvcount ;
361       }
362     if ( !_asynchronous && sendbuf )
363       {
364         if ( sendtype == MPI_INT )
365           {
366             delete [] (int *) sendbuf ;
367           }
368         else
369           {
370             delete [] (double *) sendbuf ;
371           }
372       }
373     return sts ;
374   }
375
376   /*!
377     Send sendcounts[target] datas from sendbuf[sdispls[target]] with type sendtype to all targets of IntraCommunicator
378     Receive recvcounts[target] datas to recvbuf[rdispls[target]] with type recvtype from all targets of IntraCommunicator
379
380   */
381   int MPIAccessDEC::allToAllv( void* sendbuf, int* sendcounts, int* sdispls,
382                                MPI_Datatype sendtype ,
383                                void* recvbuf, int* recvcounts, int* rdispls,
384                                MPI_Datatype recvtype )
385   {
386     if ( _time_interpolator )
387       {
388         return allToAllvTime( sendbuf, sendcounts, sdispls, sendtype ,
389                               recvbuf, recvcounts, rdispls, recvtype ) ;
390       }
391     int sts ;
392     int target ;
393     int SendRequestId ;
394     int RecvRequestId ;
395
396     //Free of SendBuffers 
397     if ( _asynchronous )
398       {
399         checkSent() ;
400       }
401
402     //DoSend + DoRecv : SendRecv
403     SendBuffStruct * aSendDataStruct = NULL ;
404     if ( _asynchronous && sendbuf )
405       {
406         aSendDataStruct = new SendBuffStruct ;
407         aSendDataStruct->SendBuffer = sendbuf ;
408         aSendDataStruct->Counter = 0 ;
409         aSendDataStruct->DataType = sendtype ;
410       }
411     for ( target = 0 ; target < _group_size ; target++ )
412       {
413         if ( sendcounts[target] || recvcounts[target] )
414           {
415             sts = sendRecv( sendbuf , sendcounts[target] , sdispls[target] , sendtype ,
416                             recvbuf , recvcounts[target] , rdispls[target] , recvtype ,
417                             target , SendRequestId , RecvRequestId ) ;
418             if ( _asynchronous && sendbuf && sendcounts[target])
419               {
420                 aSendDataStruct->Counter += 1 ;
421                 (*_map_of_send_buffers)[ SendRequestId ] = aSendDataStruct ;
422               }
423           }
424       }
425     if ( !_asynchronous && sendbuf )
426       {
427         if ( sendtype == MPI_INT )
428           {
429             delete [] (int *) sendbuf ;
430           }
431         else
432           {
433             delete [] (double *) sendbuf ;
434           }
435       }
436     return sts ;
437   }
438
439   /*
440     MPIAccessDEC and the management of SendBuffers :
441     =================================================
442
443     . In the collective communications collectives we send only parts of
444     the same buffer to each "target". So in asynchronous mode it is
445     necessary that all parts are free before to delete/free the
446     buffer.
447
448     . We assume that buffers are allocated with a new double[]. so a
449     delete [] is done.
450
451     . The structure SendBuffStruct permit to keep the adress of the buffer
452     and to manage a reference counter of that buffer. It contains
453     also MPI_Datatype for the delete [] (double *) ... when the counter
454     is null.
455
456     . The map _MapOfSendBuffers etablish the correspondance between each
457     RequestId given by a MPI_Access->ISend(...) and a SendBuffStruct
458     for each "target" of a part of the buffer.
459
460     . All that concerns only asynchronous Send. In synchronous mode,
461     we delete senbuf just after the Send.
462   */
463
464   /*
465     MPIAccessDEC and the management of RecvBuffers :
466     =================================================
467
468     If there is no interpolation, no special action is done.
469
470     With interpolation for each target :
471     ------------------------------------
472     . We have _time_messages[target] which is a vector of TimesMessages.
473     We have 2 TimesMessages in our case with a linear interpolation.
474     They contain the previous time(t0)/deltatime and the last
475     time(t1)/deltatime.
476
477     . We have _data_messages[target] which is a vector of DatasMessages.
478     We have 2 DatasMessages in our case with a linear interpolation.
479     They contain the previous datas at time(t0)/deltatime and at last
480     time(t1)/deltatime.
481
482     . At time _t(t*) of current processus we do the interpolation of
483     the values of the 2 DatasMessages which are returned in the part of
484     recvbuf corresponding to the target with t0 < t* <= t1.
485
486     . Because of the difference of "deltatimes" between processes, we
487     may have t0 < t1 < t* and there is an extrapolation.
488
489     . The vectors _out_of_time, _DataMessagesRecvCount and _DataMessagesType
490     contain for each target true if t* > last t1, recvcount and
491     MPI_Datatype for the finalize of messages at the end.
492   */
493
494   /*!
495     Send a TimeMessage to all targets of IntraCommunicator
496     Receive the TimeMessages from targets of IntraCommunicator if necessary.
497
498     Send sendcount datas from sendbuf[offset] with type sendtype to all targets of IntraCommunicator
499     Returns recvcount datas to recvbuf[offset] with type recvtype after an interpolation
500     with datas received from all targets of IntraCommunicator.
501
502   */
503   int MPIAccessDEC::allToAllTime( void* sendbuf, int sendcount , MPI_Datatype sendtype ,
504                                   void* recvbuf, int recvcount , MPI_Datatype recvtype )
505   {
506     int sts ;
507     int target ;
508     int sendoffset = 0 ;
509     int SendTimeRequestId ;
510     int SendDataRequestId ;
511
512     if ( _time_interpolator == NULL )
513       {
514         return MPI_ERR_OTHER ;
515       }
516
517     //Free of SendBuffers 
518     if ( _asynchronous )
519       {
520         checkSent() ;
521       }
522
523     //DoSend : Time + SendBuff
524     SendBuffStruct * aSendTimeStruct = NULL ;
525     SendBuffStruct * aSendDataStruct = NULL ;
526     if ( sendbuf && sendcount )
527       {
528         TimeMessage * aSendTimeMessage = new TimeMessage ;
529         if ( _asynchronous )
530           {
531             aSendTimeStruct = new SendBuffStruct ;
532             aSendTimeStruct->SendBuffer = aSendTimeMessage ;
533             aSendTimeStruct->Counter = 0 ;
534             aSendTimeStruct->DataType = _MPI_access->timeType() ;
535             aSendDataStruct = new SendBuffStruct ;
536             aSendDataStruct->SendBuffer = sendbuf ;
537             aSendDataStruct->Counter = 0 ;
538             aSendDataStruct->DataType = sendtype ;
539           }
540         aSendTimeMessage->time = _t ;
541         aSendTimeMessage->deltatime = _dt ;
542         for ( target = 0 ; target < _group_size ; target++ )
543           {
544             sts = send( aSendTimeMessage , 1 , 0 , _MPI_access->timeType() , target ,
545                         SendTimeRequestId ) ;
546             sts = send( sendbuf , sendcount , sendoffset , sendtype , target , SendDataRequestId ) ;
547             if ( _asynchronous )
548               {
549                 aSendTimeStruct->Counter += 1 ;
550                 (*_map_of_send_buffers)[ SendTimeRequestId ] = aSendTimeStruct ;
551                 aSendDataStruct->Counter += 1 ;
552                 (*_map_of_send_buffers)[ SendDataRequestId ] = aSendDataStruct ;
553               }
554             sendoffset += sendcount ;
555           }
556         if ( !_asynchronous )
557           {
558             delete aSendTimeMessage ;
559             if ( sendtype == MPI_INT )
560               {
561                 delete [] (int *) sendbuf ;
562               }
563             else
564               {
565                 delete [] (double *) sendbuf ;
566               }
567           }
568       }
569
570     //CheckTime + DoRecv + DoInterp
571     if ( recvbuf && recvcount )
572       {
573         for ( target = 0 ; target < _group_size ; target++ )
574           {
575             int recvsize = recvcount*_MPI_access->extent( recvtype ) ;
576             checkTime( recvcount , recvtype , target , false ) ;
577             //===========================================================================
578             //TODO : it is assumed actually that we have only 1 timestep before nad after
579             //===========================================================================
580             if ( _time_interpolator && (*_time_messages)[target][0].time != -1 )
581               {
582                 if ( (*_out_of_time)[target] )
583                   {
584                     cout << " =====================================================" << endl
585                          << "Recv" << _my_rank << " <-- target " << target << " t0 "
586                          << (*_time_messages)[target][0].time << " < t1 "
587                          << (*_time_messages)[target][1].time << " < t* " << _t << endl
588                          << " =====================================================" << endl ;
589                   }
590                 if ( recvtype == MPI_INT )
591                   {
592                     _time_interpolator->doInterp( (*_time_messages)[target][0].time,
593                                                   (*_time_messages)[target][1].time, _t, recvcount ,
594                                                   _n_step_before, _n_step_after,
595                                                   (int **) &(*_data_messages)[target][0],
596                                                   (int **) &(*_data_messages)[target][1],
597                                                   &((int *)recvbuf)[target*recvcount] ) ;
598                   }
599                 else
600                   {
601                     _time_interpolator->doInterp( (*_time_messages)[target][0].time,
602                                                   (*_time_messages)[target][1].time, _t, recvcount ,
603                                                   _n_step_before, _n_step_after,
604                                                   (double **) &(*_data_messages)[target][0],
605                                                   (double **) &(*_data_messages)[target][1],
606                                                   &((double *)recvbuf)[target*recvcount] ) ;
607                   }
608               }
609             else
610               {
611                 char * buffdest = (char *) recvbuf ;
612                 char * buffsrc = (char *) (*_data_messages)[target][1] ;
613                 memcpy( &buffdest[target*recvsize] , buffsrc , recvsize ) ;
614               }
615           }
616       }
617
618     return sts ;
619   }
620
621   int MPIAccessDEC::allToAllvTime( void* sendbuf, int* sendcounts, int* sdispls,
622                                    MPI_Datatype sendtype ,
623                                    void* recvbuf, int* recvcounts, int* rdispls,
624                                    MPI_Datatype recvtype )
625   {
626     int sts ;
627     int target ;
628     int SendTimeRequestId ;
629     int SendDataRequestId ;
630
631     if ( _time_interpolator == NULL )
632       {
633         return MPI_ERR_OTHER ;
634       }
635
636     //Free of SendBuffers 
637     if ( _asynchronous )
638       {
639         checkSent() ;
640       }
641
642     /*
643       . DoSend :
644       + We create a TimeMessage (look at that structure in MPI_Access).
645       + If we are in asynchronous mode, we create two structures SendBuffStruct
646       aSendTimeStruct and aSendDataStruct that we fill.
647       + We fill the structure aSendTimeMessage with time/deltatime of
648       the current process. "deltatime" must be nul if it is the last step of
649       Time.
650       + After that for each "target", we Send the TimeMessage and the part
651       of sendbuf corresponding to that target.
652       + If we are in asynchronous mode, we increment the counter and we add
653       aSendTimeStruct and aSendDataStruct to _MapOfSendBuffers with the
654       identifiers SendTimeRequestId and SendDataRequestId returned by
655       MPI_Access->Send(...).
656       + And if we are in synchronous mode we delete the SendMessages.
657     */
658     //DoSend : Time + SendBuff
659     SendBuffStruct * aSendTimeStruct = NULL ;
660     SendBuffStruct * aSendDataStruct = NULL ;
661     if ( sendbuf )
662       {
663         TimeMessage * aSendTimeMessage = new TimeMessage ;
664         if ( _asynchronous )
665           {
666             aSendTimeStruct = new SendBuffStruct ;
667             aSendTimeStruct->SendBuffer = aSendTimeMessage ;
668             aSendTimeStruct->Counter = 0 ;
669             aSendTimeStruct->DataType = _MPI_access->timeType() ;
670             aSendDataStruct = new SendBuffStruct ;
671             aSendDataStruct->SendBuffer = sendbuf ;
672             aSendDataStruct->Counter = 0 ;
673             aSendDataStruct->DataType = sendtype ;
674           }
675         aSendTimeMessage->time = _t ;
676         aSendTimeMessage->deltatime = _dt ;
677         for ( target = 0 ; target < _group_size ; target++ )
678           {
679             if ( sendcounts[target] )
680               {
681                 sts = send( aSendTimeMessage , 1 , 0 , _MPI_access->timeType() , target ,
682                             SendTimeRequestId ) ;
683                 sts = send( sendbuf , sendcounts[target] , sdispls[target] , sendtype , target ,
684                             SendDataRequestId ) ;
685                 if ( _asynchronous )
686                   {
687                     aSendTimeStruct->Counter += 1 ;
688                     (*_map_of_send_buffers)[ SendTimeRequestId ] = aSendTimeStruct ;
689                     aSendDataStruct->Counter += 1 ;
690                     (*_map_of_send_buffers)[ SendDataRequestId ] = aSendDataStruct ;
691                   }
692               }
693           }
694         if ( !_asynchronous )
695           {
696             delete aSendTimeMessage ;
697             if ( sendtype == MPI_INT )
698               {
699                 delete [] (int *) sendbuf ;
700               }
701             else
702               {
703                 delete [] (double *) sendbuf ;
704               }
705           }
706       }
707
708     /*
709       . CheckTime + DoRecv + DoInterp
710       + For each target we call CheckTime
711       + If there is a TimeInterpolator and if the TimeMessage of the target
712       is not the first, we call the interpolator which return its
713       results in the part of the recv buffer corresponding to the "target".
714       + If not, there is a copy of received datas for that first step of time
715       in the part of the recv buffer corresponding to the "target".
716     */
717     //CheckTime + DoRecv + DoInterp
718     if ( recvbuf )
719       {
720         for ( target = 0 ; target < _group_size ; target++ )
721           {
722             if ( recvcounts[target] )
723               {
724                 int recvsize = recvcounts[target]*_MPI_access->extent( recvtype ) ;
725                 checkTime( recvcounts[target] , recvtype , target , false ) ;
726                 //===========================================================================
727                 //TODO : it is assumed actually that we have only 1 timestep before nad after
728                 //===========================================================================
729                 if ( _time_interpolator && (*_time_messages)[target][0].time != -1 )
730                   {
731                     if ( (*_out_of_time)[target] )
732                       {
733                         cout << " =====================================================" << endl
734                              << "Recv" << _my_rank << " <-- target " << target << " t0 "
735                              << (*_time_messages)[target][0].time << " < t1 "
736                              << (*_time_messages)[target][1].time << " < t* " << _t << endl
737                              << " =====================================================" << endl ;
738                       }
739                     if ( recvtype == MPI_INT )
740                       {
741                         _time_interpolator->doInterp( (*_time_messages)[target][0].time,
742                                                       (*_time_messages)[target][1].time, _t,
743                                                       recvcounts[target] , _n_step_before, _n_step_after,
744                                                       (int **) &(*_data_messages)[target][0],
745                                                       (int **) &(*_data_messages)[target][1],
746                                                       &((int *)recvbuf)[rdispls[target]] ) ;
747                       }
748                     else
749                       {
750                         _time_interpolator->doInterp( (*_time_messages)[target][0].time,
751                                                       (*_time_messages)[target][1].time, _t,
752                                                       recvcounts[target] , _n_step_before, _n_step_after,
753                                                       (double **) &(*_data_messages)[target][0],
754                                                       (double **) &(*_data_messages)[target][1],
755                                                       &((double *)recvbuf)[rdispls[target]] ) ;
756                       }
757                   }
758                 else
759                   {
760                     char * buffdest = (char *) recvbuf ;
761                     char * buffsrc = (char *) (*_data_messages)[target][1] ;
762                     memcpy( &buffdest[rdispls[target]*_MPI_access->extent( recvtype )] , buffsrc ,
763                             recvsize ) ;
764                   }
765               }
766           }
767       }
768
769     return sts ;
770   }
771
772   /*
773     . CheckTime(recvcount , recvtype , target , UntilEnd)
774     + At the beginning, we read the first TimeMessage in
775     &(*_TimeMessages)[target][1] and the first DataMessage
776     in the allocated buffer (*_DataMessages)[target][1].
777     + deltatime of TimesMessages must be nul if it is the last one.
778     + While : _t(t*) is the current time of the processus.
779     "while _t(t*) is greater than the time of the "target"
780     (*_TimeMessages)[target][1].time and
781     (*_TimeMessages)[target][1].deltatime is not nul",
782     So at the end of the while we have :
783     _t(t*) <= (*_TimeMessages)[target][1].time with
784     _t(t*) > (*_TimeMessages)[target][0].time
785     or we have the last TimeMessage of the "target".
786     + If it is the finalization of the recv of TimeMessages and
787     DataMessages (UntilEnd value is true), we execute the while
788     until (*_TimeMessages)[target][1].deltatime is nul.
789     + In the while :
790     We copy the last TimeMessage in the previoud TimeMessage and
791     we read a new TimeMessage
792     We delete the previous DataMessage.
793     We copy the last DataMessage pointer in the previous one.
794     We allocate a new last DataMessage buffer
795     (*_DataMessages)[target][1] and we read the corresponding
796     datas in that buffe.
797     + If the current time of the current process is greater than the
798     last time (*_TimeMessages)[target][1].time du target, we give
799     a true value to (*_OutOfTime)[target].
800     (*_TimeMessages)[target][1].deltatime is nul.
801   */
802   int MPIAccessDEC::checkTime( int recvcount , MPI_Datatype recvtype , int target ,
803                                bool UntilEnd )
804   {
805     int sts = MPI_SUCCESS ;
806     int RecvTimeRequestId ;
807     int RecvDataRequestId ;
808     //Pour l'instant on cherche _time_messages[target][0] < _t <= _time_messages[target][1]
809     //===========================================================================
810     //TODO : it is assumed actually that we have only 1 timestep before and after
811     //       instead of _n_step_before and _n_step_after ...
812     //===========================================================================
813     (*_data_messages_recv_count)[target] = recvcount ;
814     (*_data_messages_type)[target] = recvtype ;
815     if ( (*_time_messages)[target][1].time == -1 )
816       {
817         (*_time_messages)[target][0] = (*_time_messages)[target][1] ;
818         sts = recv( &(*_time_messages)[target][1] , 1 , _MPI_access->timeType() ,
819                     target , RecvTimeRequestId ) ;
820         (*_data_messages)[target][0] = (*_data_messages)[target][1] ;
821         if ( recvtype == MPI_INT )
822           {
823             (*_data_messages)[target][1] = new int[recvcount] ;
824           }
825         else
826           {
827             (*_data_messages)[target][1] = new double[recvcount] ;
828           }
829         sts = recv( (*_data_messages)[target][1] , recvcount , recvtype , target ,
830                     RecvDataRequestId ) ;
831       }
832     else
833       {
834         while ( ( _t > (*_time_messages)[target][1].time || UntilEnd ) &&
835                 (*_time_messages)[target][1].deltatime != 0 )
836           {
837             (*_time_messages)[target][0] = (*_time_messages)[target][1] ;
838             sts = recv( &(*_time_messages)[target][1] , 1 , _MPI_access->timeType() ,
839                         target , RecvTimeRequestId ) ;
840             if ( UntilEnd )
841               {
842                 cout << "CheckTime" << _my_rank << " TimeMessage target " << target
843                      << " RecvTimeRequestId " << RecvTimeRequestId << " MPITag "
844                      << _MPI_access->recvMPITag(target) << endl ;
845               }
846             if ( recvtype == MPI_INT )
847               {
848                 delete [] (int *) (*_data_messages)[target][0] ;
849               }
850             else
851               {
852                 delete [] (double *) (*_data_messages)[target][0] ;
853               }
854             (*_data_messages)[target][0] = (*_data_messages)[target][1] ;
855             if ( recvtype == MPI_INT )
856               {
857                 (*_data_messages)[target][1] = new int[recvcount] ;
858               }
859             else
860               {
861                 (*_data_messages)[target][1] = new double[recvcount] ;
862               }
863             sts = recv( (*_data_messages)[target][1] , recvcount , recvtype , target ,
864                         RecvDataRequestId ) ;
865             if ( UntilEnd )
866               {
867                 cout << "CheckTime" << _my_rank << " DataMessage target " << target
868                      << " RecvDataRequestId " << RecvDataRequestId << " MPITag "
869                      << _MPI_access->recvMPITag(target) << endl ;
870               }
871           }
872
873         if ( _t > (*_time_messages)[target][0].time &&
874              _t <= (*_time_messages)[target][1].time )
875           {
876           }
877         else
878           {
879             (*_out_of_time)[target] = true ;
880           }
881       }
882     return sts ;
883   }
884
885   /*
886     . CheckSent() :
887     + call  SendRequestIds of MPI_Access in order to get all
888     RequestIds of SendMessages of all "targets".
889     + For each RequestId, CheckSent call "Test" of MPI_Access in order
890     to know if the buffer is "free" (flag = true). If it is the
891     FinalCheckSent (WithWait = true), we call Wait instead of Test.
892     + If the buffer is "free", the counter of the structure SendBuffStruct
893     (from _MapOfSendBuffers) is decremented.
894     + If that counter is nul we delete the TimeMessage or the
895     SendBuffer according to the DataType.
896     + And we delete the structure SendBuffStruct before the suppression
897     (erase) of that item of _MapOfSendBuffers
898   */
899   int MPIAccessDEC::checkSent(bool WithWait)
900   {
901     int sts = MPI_SUCCESS ;
902     int flag = WithWait ;
903     int size = _MPI_access->sendRequestIdsSize() ;
904     int * ArrayOfSendRequests = new int[ size ] ;
905     int nSendRequest = _MPI_access->sendRequestIds( size , ArrayOfSendRequests ) ;
906     bool SendTrace = false ;
907     int i ;
908     for ( i = 0 ; i < nSendRequest ; i++ )
909       {
910         if ( WithWait )
911           {
912             cout << "CheckSent" << _my_rank << " " << i << "./" << nSendRequest
913                  << " SendRequestId " << ArrayOfSendRequests[i] << " MPITarget "
914                  << _MPI_access->MPITarget(ArrayOfSendRequests[i]) << " MPITag "
915                  << _MPI_access->MPITag(ArrayOfSendRequests[i]) << " Wait :" << endl ;
916             sts = _MPI_access->wait( ArrayOfSendRequests[i] ) ;
917           }
918         else
919           {
920             sts = _MPI_access->test( ArrayOfSendRequests[i] , flag ) ;
921           }
922         if ( flag )
923           {
924             _MPI_access->deleteRequest( ArrayOfSendRequests[i] ) ;
925             if ( SendTrace )
926               {
927                 cout << "CheckSent" << _my_rank << " " << i << "./" << nSendRequest
928                      << " SendRequestId " << ArrayOfSendRequests[i]
929                      << " flag " << flag
930                      << " Counter " << (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ]->Counter
931                      << " DataType " << (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ]->DataType
932                      << endl ;
933               }
934             (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ]->Counter -= 1 ;
935             if ( SendTrace )
936               {
937                 if ( (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ]->DataType == 
938                      _MPI_access->timeType() )
939                   {
940                     cout << "CheckTimeSent" << _my_rank << " Request " ;
941                   }
942                 else
943                   {
944                     cout << "CheckDataSent" << _my_rank << " Request " ;
945                   }
946                 cout << ArrayOfSendRequests[i]
947                      << " _map_of_send_buffers->SendBuffer "
948                      << (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ]->SendBuffer
949                      << " Counter " << (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ]->Counter
950                      << endl ;
951               }
952             if ( (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ]->Counter  == 0 )
953               {
954                 if ( SendTrace )
955                   {
956                     cout << "CheckSent" << _my_rank << " SendRequestId " << ArrayOfSendRequests[i]
957                          << " Counter " << (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ]->Counter
958                          << " flag " << flag << " SendBuffer "
959                          << (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ]->SendBuffer
960                          << " deleted. Erase in _map_of_send_buffers :" << endl ;
961                   }
962                 if ( (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ]->DataType ==
963                      _MPI_access->timeType() )
964                   {
965                     delete (TimeMessage * ) (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ]->SendBuffer ;
966                   }
967                 else
968                   {
969                     if ( (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ]->DataType == MPI_INT )
970                       {
971                         delete [] (int *) (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ]->SendBuffer ;
972                       }
973                     else
974                       {
975                         delete [] (double *) (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ]->SendBuffer ;
976                       }
977                   }
978                 delete (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ] ;
979               }
980             if ( SendTrace )
981               {
982                 cout << "CheckSent" << _my_rank << " Erase in _map_of_send_buffers SendRequestId "
983                      << ArrayOfSendRequests[i] << endl ;
984               }
985             (*_map_of_send_buffers).erase( ArrayOfSendRequests[i] ) ;
986           }
987         else if ( SendTrace )
988           {
989             cout << "CheckSent" << _my_rank << " " << i << "./" << nSendRequest
990                  << " SendRequestId " << ArrayOfSendRequests[i]
991                  << " flag " << flag
992                  << " Counter " << (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ]->Counter
993                  << " DataType " << (*_map_of_send_buffers)[ ArrayOfSendRequests[i] ]->DataType
994                  << endl ;
995           }
996       }
997     if ( SendTrace )
998       {
999         _MPI_access->check() ;
1000       }
1001     delete [] ArrayOfSendRequests ;
1002     return sts ;
1003   }
1004
1005   int MPIAccessDEC::checkFinalRecv()
1006   {
1007     int sts = MPI_SUCCESS ;
1008     if ( _time_interpolator )
1009       {
1010         int target ;
1011         for ( target = 0 ; target < _group_size ; target++ )
1012           {
1013             if ( (*_data_messages)[target][0] != NULL )
1014               {
1015                 sts = checkTime( (*_data_messages_recv_count)[target] , (*_data_messages_type)[target] ,
1016                                  target , true ) ;
1017                 if ( (*_data_messages_type)[target] == MPI_INT )
1018                   {
1019                     delete [] (int *) (*_data_messages)[target][0] ;
1020                   }
1021                 else
1022                   {
1023                     delete [] (double *) (*_data_messages)[target][0] ;
1024                   }
1025                 (*_data_messages)[target][0] = NULL ;
1026                 if ( (*_data_messages)[target][1] != NULL )
1027                   {
1028                     if ( (*_data_messages_type)[target] == MPI_INT )
1029                       {
1030                         delete [] (int *) (*_data_messages)[target][1] ;
1031                       }
1032                     else
1033                       {
1034                         delete [] (double *) (*_data_messages)[target][1] ;
1035                       }
1036                     (*_data_messages)[target][1] = NULL ;
1037                   }
1038               }
1039           }
1040       }
1041     return sts ;
1042   }
1043
1044   ostream & operator<< (ostream & f ,const TimeInterpolationMethod & interpolationmethod )
1045   {
1046     switch (interpolationmethod)
1047       {
1048       case WithoutTimeInterp :
1049         f << " WithoutTimeInterpolation ";
1050         break;
1051       case LinearTimeInterp :
1052         f << " LinearTimeInterpolation ";
1053         break;
1054       default :
1055         f << " UnknownTimeInterpolation ";
1056         break;
1057       }
1058
1059     return f;
1060   }
1061 }