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