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