source: CPL/oasis3/trunk/src/lib/psmile/src/prism_enddef_proto.F @ 1677

Last change on this file since 1677 was 1677, checked in by aclsce, 12 years ago

Imported oasis3 (tag ipslcm5a) from cvs server to svn server (igcmg project).

File size: 17.9 KB
Line 
1      SUBROUTINE prism_enddef_proto(kinfo)
2c
3c*    *** Enddef ***   PRISM 1.0
4c
5c     purpose:
6c     --------
7c        beginning of the coupled run
8c
9c     interface:
10c     ----------
11c        kinfo  : output status
12c
13c     lib mp:
14c     -------
15c        mpi-1 or mpi-2
16c
17c     author:
18c     -------
19c        Eric Sevault   - METEO FRANCE
20c        Laurent Terray - CERFACS
21c        Jean Latour    - F.S.E.   (mpi-2)
22c        Arnaud Caubel  - FECIT - removed kmxtag as argument
23c
24c     ----------------------------------------------------------------
25      USE mod_kinds_model
26      USE mod_prism_proto
27      USE mod_comprism_proto
28#if !defined key_noIO
29      USE mod_psmile_io_interfaces
30#endif
31#include <mpif.h>
32c     ----------------------------------------------------------------
33      INTEGER (kind=ip_intwp_p) kinfo, il_err, il_max
34      REAL(kind=ip_double_p) rl_testvar
35      INTEGER (kind=ip_intwp_p) il_ibyt, il_bufsendsize,
36     $    il_bufsizebyt, integer_byte_size, ii, io_size, 
37     $    integer_io_size, il_modbufsizebyt, il_modbufsize,
38     $    il_bufsize, il_errhandler
39c     ----------------------------------------------------------------
40c
41      rl_testvar = 0.0_ip_double_p
42c
43c*    0. First Check
44c     --------------
45c
46      IF (nexit.ne.1) THEN
47          kinfo = CLIM_FastExit
48          WRITE(nulprt,FMT='(A)') 'Enddef - should not be called'
49          GO TO 1010
50      ENDIF
51      kinfo  = CLIM_Ok
52      nexit  = 0
53      il_max = maxval(myport(4,:))
54      ig_maxtype = ig_CLIMmax*kind(rl_testvar)
55      ig_maxtype_field = il_max*kind(rl_testvar)
56c
57      IF (lg_clim_bsend) THEN
58c         Determination of size of dg_bufsend element without kind
59          integer_byte_size = BIT_SIZE(ii)/8
60          INQUIRE (iolength=io_size) ii
61          integer_io_size = io_size
62          INQUIRE (iolength=io_size) rl_testvar
63          il_ibyt = io_size/integer_io_size*integer_byte_size
64c
65c         MPI_BSEND_OVERHEAD/il_ibyt+1 is the number of dg_bufsend
66c         elements needed to store MPI_BSEND_OVERHEAD bytes; there
67c         may be MPI_BSEND_OVERHEAD bytes per message (nports) 
68          il_bufsendsize = (il_max + MPI_BSEND_OVERHEAD/il_ibyt+1)
69     $        * nports
70          WRITE(nulprt,*) 'il_bufsendsize = ', il_bufsendsize
71C 
72C         Test if model has already attached to mpi buffer
73          CALL MPI_Errhandler_set(MPI_COMM_WORLD,
74     $         MPI_ERRORS_RETURN, mpi_err)
75
76          CALL MPI_Buffer_Detach(rl_testvar,il_modbufsizebyt,mpi_err)
77
78          IF (mpi_err .ne. MPI_SUCCESS) THEN
79              il_bufsize = il_bufsendsize
80          ELSE
81              il_modbufsize = il_modbufsizebyt/il_ibyt
82              il_bufsize = il_bufsendsize + il_modbufsize
83              WRITE(nulprt,*) 'il_bufsize = ', il_bufsize
84          ENDIF
85
86          CALL MPI_Errhandler_set(MPI_COMM_WORLD,
87     $        MPI_ERRORS_ARE_FATAL, mpi_err)
88
89          ALLOCATE (dg_bufsend(il_bufsize),  stat = il_err)
90          il_bufsizebyt = il_bufsize * il_ibyt
91          WRITE(nulprt,*) 'il_bufsizebyt = ', il_bufsizebyt
92C 
93C*   Attach a buffer able to contain the number of coupling fields 
94c    (nports) with maximum possible size (il_max) supposed to be declared 
95c    as DOUBLE. 
96          CALL MPI_Buffer_Attach(dg_bufsend, il_bufsizebyt, mpi_err)
97          WRITE(nulprt,*)'Attached buffer of size=',il_bufsizebyt
98c
99      ENDIF
100      ALLOCATE (pkwork_field(il_max), stat = il_err)
101      IF (il_err.ne.CLIM_Ok) WRITE(nulprt,*)
102     $    ' Problem in "pkwork_field" allocation in PRISM_enddef !'
103      pkwork_field(:)=0
104c
105      CALL PRISM_Start_MPI(kinfo)
106c
107c*    1. Allocation of arrays used for exchange of fields and local 
108c*    transformations
109c     -----------------------------------------------------
110      ALLOCATE (rg_field_trans(il_max,nports),stat = il_err)
111      IF (il_err.ne.CLIM_Ok) WRITE(nulprt,*)
112     $    ' Problem in "rg_field_trans" allocation in PRISM_enddef !'
113      rg_field_trans(:,:)=0
114      ALLOCATE (dg_field_trans(il_max,nports),stat = il_err)
115      IF (il_err.ne.CLIM_Ok) WRITE(nulprt,*)
116     $    ' Problem in "dg_field_trans" allocation in PRISM_enddef !'
117      dg_field_trans(:,:)=0
118      ALLOCATE (ig_number(nports),stat = il_err)
119      IF (il_err.ne.CLIM_Ok) WRITE(nulprt,*)
120     $    ' Problem in "ig_number" allocation in PRISM_enddef !'
121      ig_number(:) = 0
122c
123c*    2. Get wall clock differences with other models
124c     -----------------------------------------------
125c
126c       2.1 Initilizations
127c       ------------------
128c
129      IF (lg_oasis_field) THEN
130         delta(0:ncplprocs-1)=0.
131         delte(0:ncplprocs-1)=0.
132      ELSE
133         delta(1:ncplprocs)=0.
134         delte(1:ncplprocs)=0.
135      ENDIF
136c
137c*    correct execution
138c
139      nexit = 1
140#if !defined key_noIO
141      call psmile_def_domains(il_err)
142      call psmile_def_files(il_err)
143      call psmile_def_metadata(il_err)
144#endif
145c
146c     ----------------------------------------------------------------
147c
148 1010 CONTINUE
149      WRITE (nulprt,FMT='(A)') 'Returning from Enddef -- '
150      CALL FLUSH(nulprt)
151      RETURN
152      END
153      SUBROUTINE PRISM_Start_MPI(kinfo)
154c
155c*    *** Start-mpi ***   CLIM 3.0
156c
157c     purpose:
158c     --------
159c        beginning of the coupled run (MPI-2 only)
160c
161c     interface:
162c     ----------
163c        kinfo  : output status
164c
165c     method:
166c     -------
167c        With MPI2, use the new intra-communicator created by the
168c                   MPI_Intercomm_merge in CLIM_Init (models)
169c                   or in inicmc via CLIM_Init (coupler Oasis)
170c
171c     lib mp:
172c     -------
173c        mpi-1 (or mpi-2)
174c
175c     author:
176c     -------
177c        Eric Sevault   - METEO FRANCE
178c        Laurent Terray - CERFACS
179c        Jean Latour    - F.S.E.   (mpi-2)
180c        Sophie Valcke  - CERFACS (08/09/00 -modified)
181c        Arnaud Caubel  - FECIT (08/02 - removed kmxtag as argument)
182c     ----------------------------------------------------------------
183      USE mod_kinds_model
184      USE mod_prism_proto
185      USE mod_comprism_proto
186#include <mpif.h>
187c     ----------------------------------------------------------------
188      INTEGER (kind=ip_intwp_p) kinfo
189c     ----------------------------------------------------------------
190      INTEGER (kind=ip_intwp_p)    ip, info1, info, ilgdt, ipos, is,
191     *             itag1, isdmod, incp, ipotag,nsegid
192      INTEGER (kind=ip_intwp_p)    irempo, il_maxtag, il_start, il_end
193c
194      PARAMETER    (ibuflen=3)
195      INTEGER (kind=ip_intwp_p)    iposbuf, ibuff_mpi(ibuflen) 
196      INTEGER (kind=ip_intwp_p)    istatus(MPI_STATUS_SIZE), imaxbyt
197      LOGICAL      ll_flag
198      INTEGER (kind=ip_intwp_p), allocatable :: ireq(:)
199      INTEGER (kind=ip_intwp_p)              :: ireqrecv
200      INTEGER (kind=ip_intwp_p)              :: irc
201      INTEGER (kind=ip_intwp_p)              :: rpkwork(ig_CLIMmax)
202c     ----------------------------------------------------------------
203c
204c
205c     0. define some variables 
206c     ------------------------
207c
208      WRITE(nulprt,*)'Start - -  '
209c           
210      nexit  = 0
211      istatus(:)=0
212      ibuff_mpi(:)=0
213      call MPI_Attr_get(MPI_COMM_WORLD, MPI_TAG_UB, CLIM_MaxTag, 
214     $    ll_flag,info)
215      if(info.ne.MPI_SUCCESS) then
216        write(nulprt, * )'CLIM_Start_MPI -- MPI_Attr_get failed !',info
217        call MPI_ABORT (mpi_comm, 0, info)
218      endif
219      if ( .not. ll_flag ) then
220         write(nulprt, * ) 'Warning:  MPI_Attr_get did not return '
221         write(nulprt, * ) ' a valid CLIM_MaxTag!'
222         write(nulprt, * ) ' CLIM_MaxTag is set to 32767!'
223         CLIM_MaxTag = 32767
224      else
225         write(nulprt, * ) 'CLIM_MaxTag is now ', CLIM_MaxTag
226!
227!     limit CLIM_MaxTag to 2^^30 - 1 (J. Latour)
228!
229         if ( CLIM_MaxTag .gt. 1073741823 ) CLIM_MaxTag =  1073741823
230         write(nulprt, * ) 'CLIM_MaxTag is now ', CLIM_MaxTag
231!
232      endif
233      il_maxtag = CLIM_MaxTag - 1
234      itag1 = CLIM_MaxTag
235      ilgdt = CLIM_ParSize
236C
237c*    1. broadcast usefull informations
238c     ---------------------------------
239c
240      ibuff_mpi(1) = mynum 
241      ibuff_mpi(2) = modtid(mynum) 
242      ibuff_mpi(3) = nports
243c
244      WRITE(nulprt,FMT='(A,A)') 'Start- send- MODEL ',cmynam
245      WRITE(nulprt,FMT='(A,I9)')'Start- send-  num  :',mynum
246      WRITE(nulprt,FMT='(A,I9)')'Start- send-  tid  :',modtid(mynum)
247      WRITE(nulprt,FMT='(A,I9)')'Start- send-  nport:',nports
248c
249      iposbuf = 0
250      call MPI_Pack ( ibuff_mpi, ibuflen, MPI_INTEGER, pkwork,
251     *                ig_maxtype, iposbuf, mpi_comm, info )
252      call MPI_Pack ( cmynam, CLIM_Clength, MPI_CHARACTER, pkwork,
253     *                ig_maxtype, iposbuf, mpi_comm, info )
254      DO 210 ji=1,nports
255         call MPI_Pack ( cports(ji), CLIM_Clength, MPI_CHARACTER,
256     *           pkwork, ig_maxtype, iposbuf, mpi_comm, info )
257         call MPI_Pack ( myport(1,ji), 5, MPI_INTEGER, pkwork,
258     *                   ig_maxtype, iposbuf, mpi_comm, info )
259         call MPI_Pack ( mydist(1,ji), ilgdt, MPI_INTEGER, pkwork,
260     *                   ig_maxtype, iposbuf, mpi_comm, info )
261 210   CONTINUE
262c
263c*     In the direct case, Oasis process is not a process involved in the 
264c*     coupling
265c
266       IF (lg_oasis_field) THEN
267          il_start = 0
268          il_end = ncplprocs-1
269       ELSE
270          il_start = 1
271          il_end = ncplprocs
272       ENDIF
273       ALLOCATE(ireq(0:ncplprocs))
274       irc = 0
275      DO 220 ji = il_start, il_end
276c
277c       Send to all processors involved in the coupling except itself
278        IF ( ji .NE. mynum ) THEN
279            CALL MPI_Isend ( pkwork, iposbuf, MPI_PACKED, modtid(ji),
280     *                      itag1, mpi_comm, ireq(irc), info )
281            irc = irc+1
282c
283        ENDIF
284 220  CONTINUE
285
286      WRITE (nulprt,FMT='(A,I3,A,I8,A)') 
287     * 'Start - broadcast from mynum = ',mynum,' <MPI ',info,'>'
288c
289c*    3. get these infos from other models and check ports
290c     ----------------------------------------------------
291c
292      imaxbyt=0
293      ireqrecv = 0
294      DO 380 ip=1, ncplprocs-1
295c
296c Blocking receive (to comment if use of non-blocking receive)
297        CALL MPI_IRecv(rpkwork, ig_maxtype, MPI_PACKED, MPI_ANY_SOURCE,
298     *                 itag1, mpi_comm, ireqrecv, info1 )
299        CALL MPI_Wait ( ireqrecv, istatus, info )
300        CALL MPI_Get_count ( istatus, MPI_PACKED, imaxbyt, info )
301        IF ( info .NE. 0 ) THEN
302            kinfo = CLIM_Unpack
303            WRITE(nulprt,FMT='(A,I3,A)')
304     *                 'Import - pb unpack <mpi ',info,'>'
305            GO TO 1010
306        ENDIF
307c END of blocking receive
308c
309c       
310        iposbuf = 0
311        CALL MPI_Unpack ( rpkwork, ig_maxtype, iposbuf, ibuff_mpi,
312     *                      ibuflen, MPI_INTEGER, mpi_comm, info )
313        IF ( info .NE. 0 ) THEN
314            kinfo = CLIM_Unpack
315            WRITE(nulprt,FMT='(A,I3,A)')
316     *                 'Import - pb unpack <mpi ',info,'>'
317            GO TO 1010
318        ENDIF
319c
320        isdmod         = ibuff_mpi(1)
321        modtid(isdmod) = ibuff_mpi(2)
322        irempo         = ibuff_mpi(3)
323c
324        call MPI_Unpack ( rpkwork, ig_maxtype, iposbuf, cnames(isdmod),
325     *                   CLIM_Clength, MPI_CHARACTER, mpi_comm, info )
326        DO 310 ji=1,irempo
327          call MPI_Unpack ( rpkwork, ig_maxtype, iposbuf, clrport(ji),
328     *                   CLIM_Clength, MPI_CHARACTER, mpi_comm, info )
329          call MPI_Unpack( rpkwork, ig_maxtype, iposbuf, irport(1,ji), 
330     *                           5, MPI_INTEGER, mpi_comm, info )
331          call MPI_Unpack( rpkwork, ig_maxtype, iposbuf, irdist(1,ji),
332     *                           ilgdt, MPI_INTEGER, mpi_comm, info )
333 310    CONTINUE
334c
335        WRITE (nulprt,FMT='(A,A)')  'Start - MODEL ',cnames(isdmod)
336        WRITE (nulprt,FMT='(A,I9)') 'Start -   num  :',isdmod
337        WRITE (nulprt,FMT='(A,I9)') 'Start -   tid  :',modtid(isdmod)
338        WRITE (nulprt,FMT='(A,I9)') 'Start -   nport:',irempo
339c
340        ncode(isdmod) = 0
341c
342        DO 350 ji=1,nports
343          DO 340 jj=1,irempo
344            IF (cports(ji).eq.clrport(jj).and.
345     *          myport(1,ji)+irport(1,jj).eq.1) THEN
346                IF ((mydist(CLIM_Strategy,ji).eq.CLIM_Serial.and.
347     *              myport(4,ji).lt.irport(4,jj)).or.
348     *              (irdist(CLIM_Strategy,jj).eq.CLIM_Serial.and.
349     *              irport(4,jj).lt.myport(4,ji))) THEN
350                    kinfo = CLIM_IncSize
351                    WRITE(nulprt,FMT='(A,A,A,I2,A,I2,A,I8,I8)')
352     *                  'Start - WARNING Incompatible sizes - field',
353     *                  cports(ji),'model ',mynum,' and model ',
354     *                  isdmod,': ',myport(4,ji),irport(4,jj)
355                ELSEIF (myport(2,ji).ne.irport(2,jj)) THEN
356                    kinfo = CLIM_BadType
357                    WRITE(nulprt,FMT='(A,A,A,I2,A,I2,A,I2,I2)')
358     *                  'Start - WARNING Incompatible types - field',
359     *                  cports(ji),'model ',mynum,' and model ',
360     *                  isdmod,': ',myport(2,ji),irport(2,jj)
361                ELSE
362                    IF (myport(1,ji).eq.CLIM_Out) THEN
363                        incp   = ig_ntime / ig_frqmin + 1
364                        ipotag = CLIM_MaxTag - (ji-1)*incp - 1
365                    ELSE
366                        incp   = ig_ntime / ig_frqmin + 1
367                        ipotag = CLIM_MaxTag - (jj-1)*incp - 1
368                    ENDIF
369                    IF ((ipotag-incp).lt.il_maxtag) THEN
370                        il_maxtag = ipotag - incp
371                    ENDIF
372                    IF (mydist(CLIM_Strategy,ji).ne.CLIM_Serial.and
373     *               .irdist(CLIM_Strategy,jj).ne.CLIM_Serial) THEN
374                        IF (mydist(CLIM_Strategy,ji).ne.
375     *                      irdist(CLIM_Strategy,jj)) THEN
376                           WRITE (nulprt,FMT=*)'STOP in Clim_Start_MPI'
377                           call flush (nulprt)
378                           call MPI_ABORT (mpi_comm, 0, info)
379                        ELSE
380                            IF (mydist(CLIM_Segments,ji).eq.
381     *                      irdist(CLIM_Segments,jj)) THEN
382                                nsegid = 0                           
383                                DO  is=1,2*mydist(CLIM_Segments,ji)
384                                  IF (mydist(CLIM_Segments+is,ji).eq.
385     *                               irdist(CLIM_Segments+is,jj)) THEN
386                                      nsegid = nsegid + 1
387                                  ENDIF
388                                END DO
389                                IF (nsegid.eq.
390     *                              (2*mydist(CLIM_Segments,ji))) THEN
391                                    nlinks = nlinks + 1
392                                    ipos   = 5
393                                    myport(ipos,ji)=myport(ipos,ji)+ 1
394                                    myport(ipos+myport(ipos,ji),ji)
395     *                                     = nlinks
396                                    mylink(1,nlinks) = isdmod
397                                    mylink(2,nlinks) = modtid(isdmod)
398                                    mylink(3,nlinks) = ipotag
399                                    mylink(4,nlinks) = 1
400                                    mylink(5,nlinks) = 0
401                                    mylink(6,nlinks) = myport(4,ji)
402                                    WRITE(nulprt,FMT='(A,A)')
403     *                                  'Start - LINK ',cports(ji)
404                                    WRITE(nulprt,FMT='(A,I2,A,I1,A,I2,
405     *                                 A,I1,A,I10,A,I10)')
406     *                                  'Start - [model ',
407     *                                  mynum,'/io ',myport(1,ji),'] -
408     *                                    [model ',
409     *                                  isdmod,'/io ',irport(1,jj),']
410     *                                   using tags ',
411     *                                  ipotag,' to ',ipotag-incp+1
412                                    call flush(nulprt)
413                                ENDIF   
414                            ENDIF
415                        ENDIF
416                    ELSE
417                        nlinks = nlinks + 1
418                        ipos   = 5
419                        myport(ipos,ji) = myport(ipos,ji) + 1
420                        myport(ipos+myport(ipos,ji),ji) = nlinks
421                        mylink(1,nlinks) = isdmod
422                        mylink(2,nlinks) = modtid(isdmod)
423                        mylink(3,nlinks) = ipotag
424                        mylink(4,nlinks) = 1
425                        mylink(5,nlinks) = 0
426                        mylink(6,nlinks) = myport(4,ji)
427                        IF (mydist(CLIM_Strategy,ji).eq.CLIM_Serial.and
428     *                   .irdist(CLIM_Strategy,jj).ne.CLIM_Serial) THEN
429                            mylink(4,nlinks) = irdist(CLIM_Segments,jj)
430                            DO 330 is=1,2*irdist(CLIM_Segments,jj)
431                              mylink(4+is,nlinks) = irdist
432     *                                      (CLIM_Segments+is,jj)
433 330                        CONTINUE
434                        ENDIF
435                        WRITE(nulprt,FMT='(A,A)')
436     *                      'Start - LINK ',cports(ji)
437                        WRITE(nulprt,FMT='(A,I2,A,I1,A,I2,A,
438     *                       I1,A,I10,A,I10)')
439     *                      'Start - [model ',
440     *                      mynum,'/io ',myport(1,ji),'] - [model ',
441     *                      isdmod,'/io ',irport(1,jj),'] using tags ',
442     *                      ipotag,' to ',ipotag-incp+1
443                        call flush(nulprt)
444                    ENDIF                   
445                ENDIF
446            ENDIF
447 340      CONTINUE
448 350    CONTINUE
449c
450        WRITE (nulprt,FMT='(A)') 'Start - -'
451 380  CONTINUE
452c
453c     MPI_wait on the above MPI_Isend so to not change the content of pkwork
454c     before the sending is complete.
455      irc = 0
456      DO 420 ji = il_start, il_end
457        IF ( ji .NE. mynum ) THEN
458            CALL MPI_wait ( ireq(irc), istatus, info )
459            WRITE (nulprt, FMT='(A)') 'After MPI_wait'
460            irc = irc+1
461        ENDIF
462 420  CONTINUE
463      DEALLOCATE(ireq)
464c
465c*    correct execution
466c
467      nexit = 1
468c
469c     ----------------------------------------------------------------
470c
471 1010 CONTINUE
472      WRITE (nulprt,FMT='(A)') 'Returning from Start-mpi  -- '
473      CALL FLUSH(nulprt)
474      RETURN
475      END
Note: See TracBrowser for help on using the repository browser.