source: LMDZ6/branches/Ocean_skin/libf/phylmd/oasis.F90 @ 3627

Last change on this file since 3627 was 3627, checked in by lguez, 4 years ago

If the ocean skin parameterization is working (passively or actively,
activate_ocean_skin >= 1) and we are coupled to the ocean then
receive bulk salinity of the surface layer of the ocean from the ocean
and feed it to procedure bulk_flux instead of the constant
value 35. If the ocean skin parameterization is working actively
(activate_ocean_skin == 2) and we are coupled to the ocean then send
ocean-air interface temperature to the ocean. We can only send
interface temperature from the previous time-step since communication
with the ocean is before the call to bulk_flux. In module cpl_mod,
define cpl_t_int with rank 1: no dimension for cpl_index because
t_int is only defined over ocean. New dummy argument sss of
procedures cpl_receive_ocean_fields and ocean_cpl_noice. New dummy
argument t_int of cpl_send_ocean_fields. In procedure
surf_ocean, rename local variable s1 to sss and give it the size
klon, which is required by the coupling machinery.

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 20.7 KB
Line 
1!
2MODULE oasis
3!
4! This module contains subroutines for initialization, sending and receiving
5! towards the coupler OASIS3. It also contains some parameters for the coupling.
6!
7! This module should always be compiled. With the coupler OASIS3 available the cpp key
8! CPP_COUPLE should be set and the entier of this file will then be compiled.
9! In a forced mode CPP_COUPLE should not be defined and the compilation ends before
10! the CONTAINS, without compiling the subroutines.
11!
12  USE dimphy
13  USE mod_phys_lmdz_para
14  USE write_field_phy
15
16#ifdef CPP_COUPLE
17! Use of Oasis-MCT coupler
18#if defined CPP_OMCT
19  USE mod_prism
20! Use of Oasis3 coupler
21#else
22  USE mod_prism_proto
23  USE mod_prism_def_partition_proto
24  USE mod_prism_get_proto
25  USE mod_prism_put_proto
26#endif
27#endif
28 
29  IMPLICIT NONE
30 
31  ! Id for fields sent to ocean
32  INTEGER, PARAMETER :: ids_tauxxu = 1
33  INTEGER, PARAMETER :: ids_tauyyu = 2
34  INTEGER, PARAMETER :: ids_tauzzu = 3
35  INTEGER, PARAMETER :: ids_tauxxv = 4
36  INTEGER, PARAMETER :: ids_tauyyv = 5
37  INTEGER, PARAMETER :: ids_tauzzv = 6
38  INTEGER, PARAMETER :: ids_windsp = 7
39  INTEGER, PARAMETER :: ids_shfice = 8
40  INTEGER, PARAMETER :: ids_shfoce = 9
41  INTEGER, PARAMETER :: ids_shftot = 10
42  INTEGER, PARAMETER :: ids_nsfice = 11
43  INTEGER, PARAMETER :: ids_nsfoce = 12
44  INTEGER, PARAMETER :: ids_nsftot = 13
45  INTEGER, PARAMETER :: ids_dflxdt = 14
46  INTEGER, PARAMETER :: ids_totrai = 15
47  INTEGER, PARAMETER :: ids_totsno = 16
48  INTEGER, PARAMETER :: ids_toteva = 17
49  INTEGER, PARAMETER :: ids_icevap = 18
50  INTEGER, PARAMETER :: ids_ocevap = 19
51  INTEGER, PARAMETER :: ids_calvin = 20
52  INTEGER, PARAMETER :: ids_liqrun = 21
53  INTEGER, PARAMETER :: ids_runcoa = 22
54  INTEGER, PARAMETER :: ids_rivflu = 23
55  INTEGER, PARAMETER :: ids_atmco2 = 24
56  INTEGER, PARAMETER :: ids_taumod = 25
57  INTEGER, PARAMETER :: ids_qraioc = 26
58  INTEGER, PARAMETER :: ids_qsnooc = 27
59  INTEGER, PARAMETER :: ids_qraiic = 28
60  INTEGER, PARAMETER :: ids_qsnoic = 29
61  INTEGER, PARAMETER :: ids_t_int = 30
62 
63  INTEGER, PARAMETER :: maxsend    = 30  ! Maximum number of fields to send
64 
65  ! Id for fields received from ocean
66
67  INTEGER, PARAMETER :: idr_sisutw = 1
68  INTEGER, PARAMETER :: idr_icecov = 2
69  INTEGER, PARAMETER :: idr_icealw = 3
70  INTEGER, PARAMETER :: idr_icetem = 4
71  INTEGER, PARAMETER :: idr_curenx = 5
72  INTEGER, PARAMETER :: idr_cureny = 6
73  INTEGER, PARAMETER :: idr_curenz = 7
74  INTEGER, PARAMETER :: idr_oceco2 = 8
75
76  INTEGER, PARAMETER :: idr_sss = 9
77  ! bulk salinity of the surface layer of the ocean, in ppt
78
79  INTEGER, PARAMETER :: maxrecv    = 9  ! Maximum number of fields to receive
80 
81
82  TYPE, PUBLIC ::   FLD_CPL            ! Type for coupling field information
83     CHARACTER(len = 8) ::   name      ! Name of the coupling field   
84     LOGICAL            ::   action    ! To be exchanged or not
85     INTEGER            ::   nid       ! Id of the field
86  END TYPE FLD_CPL
87
88  TYPE(FLD_CPL), DIMENSION(maxsend), SAVE, PUBLIC :: infosend   ! Information for sending coupling fields
89!$OMP THREADPRIVATE(infosend)
90  TYPE(FLD_CPL), DIMENSION(maxrecv), SAVE, PUBLIC :: inforecv   ! Information for receiving coupling fields
91!$OMP THREADPRIVATE(inforecv)
92 
93  LOGICAL,SAVE :: cpl_current
94!$OMP THREADPRIVATE(cpl_current)
95
96#ifdef CPP_COUPLE
97
98CONTAINS
99
100  SUBROUTINE inicma
101!************************************************************************************
102!**** *INICMA*  - Initialize coupled mode communication for atmosphere
103!                 and exchange some initial information with Oasis
104!
105!     Rewrite to take the PRISM/psmile library into account
106!     LF 09/2003
107!
108    USE IOIPSL
109    USE surface_data, ONLY : version_ocean
110    USE carbon_cycle_mod, ONLY : carbon_cycle_cpl
111#ifdef CPP_XIOS
112    USE wxios, ONLY : wxios_context_init
113    USE xios 
114#endif
115    USE print_control_mod, ONLY: lunout
116    USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat, grid_type, unstructured, regular_lonlat
117    USE geometry_mod, ONLY: ind_cell_glo                   
118    USE mod_phys_lmdz_mpi_data, ONLY: klon_mpi_para_nb
119   
120
121
122! Local variables
123!************************************************************************************
124    INTEGER                            :: comp_id
125    INTEGER                            :: ierror, il_commlocal
126    INTEGER                            :: il_part_id
127    INTEGER, ALLOCATABLE               :: ig_paral(:)
128    INTEGER, DIMENSION(2)              :: il_var_nodims
129    INTEGER, DIMENSION(4)              :: il_var_actual_shape
130    INTEGER                            :: il_var_type
131    INTEGER                            :: jf
132    CHARACTER (len = 6)                :: clmodnam
133    CHARACTER (len = 20)               :: modname = 'inicma'
134    CHARACTER (len = 80)               :: abort_message
135    LOGICAL, SAVE                      :: cpl_current_omp
136
137!*    1. Initializations
138!        ---------------
139!************************************************************************************
140    WRITE(lunout,*) ' '
141    WRITE(lunout,*) ' '
142    WRITE(lunout,*) ' ROUTINE INICMA'
143    WRITE(lunout,*) ' **************'
144    WRITE(lunout,*) ' '
145    WRITE(lunout,*) ' '
146
147!
148! Define the model name
149!
150    IF (grid_type==unstructured) THEN
151        clmodnam = 'icosa'                 ! as in $NBMODEL in Cpl/Nam/namcouple.tmp
152    ELSE IF (grid_type==regular_lonlat) THEN
153        clmodnam = 'LMDZ'                  ! as in $NBMODEL in Cpl/Nam/namcouple.tmp
154    ELSE
155        abort_message='Pb : type of grid unknown'
156        CALL abort_physic(modname,abort_message,1)
157    ENDIF
158
159
160!************************************************************************************
161! Define if coupling ocean currents or not
162!************************************************************************************
163!$OMP MASTER
164    cpl_current_omp = .FALSE.
165    CALL getin('cpl_current', cpl_current_omp)
166!$OMP END MASTER
167!$OMP BARRIER
168    cpl_current = cpl_current_omp
169    WRITE(lunout,*) 'Couple ocean currents, cpl_current = ',cpl_current
170
171!************************************************************************************
172! Define coupling variables
173!************************************************************************************
174
175! Atmospheric variables to send
176
177!$OMP MASTER
178    infosend(:)%action = .FALSE.
179
180    infosend(ids_tauxxu)%action = .TRUE. ; infosend(ids_tauxxu)%name = 'COTAUXXU'
181    infosend(ids_tauyyu)%action = .TRUE. ; infosend(ids_tauyyu)%name = 'COTAUYYU'
182    infosend(ids_tauzzu)%action = .TRUE. ; infosend(ids_tauzzu)%name = 'COTAUZZU'
183    infosend(ids_tauxxv)%action = .TRUE. ; infosend(ids_tauxxv)%name = 'COTAUXXV'
184    infosend(ids_tauyyv)%action = .TRUE. ; infosend(ids_tauyyv)%name = 'COTAUYYV'
185    infosend(ids_tauzzv)%action = .TRUE. ; infosend(ids_tauzzv)%name = 'COTAUZZV'
186    infosend(ids_windsp)%action = .TRUE. ; infosend(ids_windsp)%name = 'COWINDSP'
187    infosend(ids_shfice)%action = .TRUE. ; infosend(ids_shfice)%name = 'COSHFICE'
188    infosend(ids_nsfice)%action = .TRUE. ; infosend(ids_nsfice)%name = 'CONSFICE'
189    infosend(ids_dflxdt)%action = .TRUE. ; infosend(ids_dflxdt)%name = 'CODFLXDT'
190    infosend(ids_calvin)%action = .TRUE. ; infosend(ids_calvin)%name = 'COCALVIN'
191   
192    if (activate_ocean_skin == 2) then
193       infosend(ids_t_int)%action = .TRUE.
194       infosend(ids_t_int)%name = 'T_int'
195    end if
196           
197    IF (version_ocean=='nemo') THEN
198        infosend(ids_shftot)%action = .TRUE. ; infosend(ids_shftot)%name = 'COQSRMIX'
199        infosend(ids_nsftot)%action = .TRUE. ; infosend(ids_nsftot)%name = 'COQNSMIX'
200        infosend(ids_totrai)%action = .TRUE. ; infosend(ids_totrai)%name = 'COTOTRAI'
201        infosend(ids_totsno)%action = .TRUE. ; infosend(ids_totsno)%name = 'COTOTSNO'
202        infosend(ids_toteva)%action = .TRUE. ; infosend(ids_toteva)%name = 'COTOTEVA'
203        infosend(ids_icevap)%action = .TRUE. ; infosend(ids_icevap)%name = 'COICEVAP'
204        infosend(ids_liqrun)%action = .TRUE. ; infosend(ids_liqrun)%name = 'COLIQRUN'
205        infosend(ids_taumod)%action = .TRUE. ; infosend(ids_taumod)%name = 'COTAUMOD'
206        IF (carbon_cycle_cpl) THEN
207            infosend(ids_atmco2)%action = .TRUE. ; infosend(ids_atmco2)%name = 'COATMCO2'
208        ENDIF
209        infosend(ids_qraioc)%action = .TRUE. ; infosend(ids_qraioc)%name = 'COQRAIOC'
210        infosend(ids_qsnooc)%action = .TRUE. ; infosend(ids_qsnooc)%name = 'COQSNOOC'
211        infosend(ids_qraiic)%action = .TRUE. ; infosend(ids_qraiic)%name = 'COQRAIIC'
212        infosend(ids_qsnoic)%action = .TRUE. ; infosend(ids_qsnoic)%name = 'COQSNOIC'
213       
214    ELSE IF (version_ocean=='opa8') THEN
215        infosend(ids_shfoce)%action = .TRUE. ; infosend(ids_shfoce)%name = 'COSHFOCE'
216        infosend(ids_nsfoce)%action = .TRUE. ; infosend(ids_nsfoce)%name = 'CONSFOCE'
217        infosend(ids_icevap)%action = .TRUE. ; infosend(ids_icevap)%name = 'COTFSICE'
218        infosend(ids_ocevap)%action = .TRUE. ; infosend(ids_ocevap)%name = 'COTFSOCE'
219        infosend(ids_totrai)%action = .TRUE. ; infosend(ids_totrai)%name = 'COTOLPSU'
220        infosend(ids_totsno)%action = .TRUE. ; infosend(ids_totsno)%name = 'COTOSPSU'
221        infosend(ids_runcoa)%action = .TRUE. ; infosend(ids_runcoa)%name = 'CORUNCOA'
222        infosend(ids_rivflu)%action = .TRUE. ; infosend(ids_rivflu)%name = 'CORIVFLU'
223   ENDIF
224       
225! Oceanic variables to receive
226
227   inforecv(:)%action = .FALSE.
228
229   inforecv(idr_sisutw)%action = .TRUE. ; inforecv(idr_sisutw)%name = 'SISUTESW'
230   inforecv(idr_icecov)%action = .TRUE. ; inforecv(idr_icecov)%name = 'SIICECOV'
231   inforecv(idr_icealw)%action = .TRUE. ; inforecv(idr_icealw)%name = 'SIICEALW'
232   inforecv(idr_icetem)%action = .TRUE. ; inforecv(idr_icetem)%name = 'SIICTEMW'
233
234   if (activate_ocean_skin >= 1) then
235      inforecv(idr_sss)%action = .TRUE.
236      inforecv(idr_sss)%name = 'salinity'
237   end if
238   
239   IF (cpl_current ) THEN
240       inforecv(idr_curenx)%action = .TRUE. ; inforecv(idr_curenx)%name = 'CURRENTX'
241       inforecv(idr_cureny)%action = .TRUE. ; inforecv(idr_cureny)%name = 'CURRENTY'
242       inforecv(idr_curenz)%action = .TRUE. ; inforecv(idr_curenz)%name = 'CURRENTZ'
243   ENDIF
244
245   IF (carbon_cycle_cpl ) THEN
246       inforecv(idr_oceco2)%action = .TRUE. ; inforecv(idr_oceco2)%name = 'SICO2FLX'
247   ENDIF
248
249!************************************************************************************
250! Here we go: psmile initialisation
251!************************************************************************************
252    IF (is_sequential) THEN
253       CALL prism_init_comp_proto (comp_id, clmodnam, ierror)
254       
255       IF (ierror .NE. PRISM_Ok) THEN
256          abort_message=' Probleme init dans prism_init_comp '
257          CALL abort_physic(modname,abort_message,1)
258       ELSE
259          WRITE(lunout,*) 'inicma : init psmile ok '
260       ENDIF
261    ENDIF
262
263    CALL prism_get_localcomm_proto (il_commlocal, ierror)
264!************************************************************************************
265! Domain decomposition
266!************************************************************************************
267    IF (grid_type==unstructured) THEN
268
269      ALLOCATE( ig_paral(klon_mpi_para_nb(mpi_rank) + 2) )
270
271      ig_paral(1) = 4                                      ! points partition for //
272      ig_paral(2) = klon_mpi_para_nb(mpi_rank)             ! nb of local cells
273
274      DO jf=1, klon_mpi_para_nb(mpi_rank)
275        ig_paral(2+jf) = ind_cell_glo(jf)
276      ENDDO
277
278    ELSE IF (grid_type==regular_lonlat) THEN
279
280      ALLOCATE( ig_paral(3) )
281
282      ig_paral(1) = 1                            ! apple partition for //
283      ig_paral(2) = (jj_begin-1)*nbp_lon+ii_begin-1  ! offset
284      ig_paral(3) = (jj_end*nbp_lon+ii_end) - (jj_begin*nbp_lon+ii_begin) + 1
285
286      IF (mpi_rank==mpi_size-1) ig_paral(3)=ig_paral(3)+nbp_lon-1
287    ELSE
288      abort_message='Pb : type of grid unknown'
289      CALL abort_physic(modname,abort_message,1)
290    ENDIF
291
292
293    WRITE(lunout,*) mpi_rank,'ig_paral--->',ig_paral(2),ig_paral(3)
294   
295    ierror=PRISM_Ok
296    CALL prism_def_partition_proto (il_part_id, ig_paral, ierror)
297
298    IF (ierror .NE. PRISM_Ok) THEN
299       abort_message=' Probleme dans prism_def_partition '
300       CALL abort_physic(modname,abort_message,1)
301    ELSE
302       WRITE(lunout,*) 'inicma : decomposition domaine psmile ok '
303    ENDIF
304
305    il_var_nodims(1) = 2                        ! rank of field array (1d or 2d)
306    il_var_nodims(2) = 1                        ! always 1 in current oasis version" doc oasis3mct p18
307
308    il_var_actual_shape(1) = 1                  ! min of 1st dimension (always 1)
309    il_var_actual_shape(2) = nbp_lon            ! max of 1st dimension
310    il_var_actual_shape(3) = 1                  ! min of 2nd dimension (always 1)
311    il_var_actual_shape(4) = nbp_lat            ! max of 2nd dimension
312   
313    il_var_type = PRISM_Real
314
315!************************************************************************************
316! Oceanic Fields to receive
317! Loop over all possible variables
318!************************************************************************************
319    DO jf=1, maxrecv
320       IF (inforecv(jf)%action) THEN
321          CALL prism_def_var_proto(inforecv(jf)%nid, inforecv(jf)%name, il_part_id, &
322               il_var_nodims, PRISM_In, il_var_actual_shape, il_var_type, &
323               ierror)
324          IF (ierror .NE. PRISM_Ok) THEN
325             WRITE(lunout,*) 'inicma : Problem with prism_def_var_proto for field : ',&
326                  inforecv(jf)%name
327             abort_message=' Problem in call to prism_def_var_proto for fields to receive'
328             CALL abort_physic(modname,abort_message,1)
329          ENDIF
330       ENDIF
331    END DO
332   
333!************************************************************************************
334! Atmospheric Fields to send
335! Loop over all possible variables
336!************************************************************************************
337    DO jf=1,maxsend
338       IF (infosend(jf)%action) THEN
339          CALL prism_def_var_proto(infosend(jf)%nid, infosend(jf)%name, il_part_id, &
340               il_var_nodims, PRISM_Out, il_var_actual_shape, il_var_type, &
341               ierror)
342          IF (ierror .NE. PRISM_Ok) THEN
343             WRITE(lunout,*) 'inicma : Problem with prism_def_var_proto for field : ',&
344                  infosend(jf)%name
345             abort_message=' Problem in call to prism_def_var_proto for fields to send'
346             CALL abort_physic(modname,abort_message,1)
347          ENDIF
348       ENDIF
349    END DO
350   
351!************************************************************************************
352! End definition
353!************************************************************************************
354#ifdef CPP_XIOS
355    CALL xios_oasis_enddef()
356#endif
357    CALL prism_enddef_proto(ierror)
358    IF (ierror .NE. PRISM_Ok) THEN
359       abort_message=' Problem in call to prism_endef_proto'
360       CALL abort_physic(modname,abort_message,1)
361    ELSE
362       WRITE(lunout,*) 'inicma : endef psmile ok '
363    ENDIF
364
365#ifdef CPP_XIOS
366!    CALL wxios_context_init()
367#endif
368
369!$OMP END MASTER
370   
371  END SUBROUTINE inicma
372
373!
374!************************************************************************************
375!
376
377  SUBROUTINE fromcpl(ktime, tab_get)
378! ======================================================================
379! L. Fairhead (09/2003) adapted From L.Z.X Li: this subroutine reads the SST
380! and Sea-Ice provided by the coupler. Adaptation to psmile library
381!======================================================================
382!
383    USE print_control_mod, ONLY: lunout
384    USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat
385! Input arguments
386!************************************************************************************
387    INTEGER, INTENT(IN)                               ::  ktime
388
389! Output arguments
390!************************************************************************************
391    REAL, DIMENSION(nbp_lon, jj_nb,maxrecv), INTENT(OUT) :: tab_get
392
393! Local variables
394!************************************************************************************
395    INTEGER                       :: ierror, i
396    INTEGER                       :: istart,iend
397    CHARACTER (len = 20)          :: modname = 'fromcpl'
398    CHARACTER (len = 80)          :: abort_message
399    REAL, DIMENSION(nbp_lon*jj_nb)    :: field
400
401!************************************************************************************
402    WRITE (lunout,*) ' '
403    WRITE (lunout,*) 'Fromcpl: Reading fields from CPL, ktime=',ktime
404    WRITE (lunout,*) ' '
405   
406    istart=ii_begin
407    IF (is_south_pole_dyn) THEN
408       iend=(jj_end-jj_begin)*nbp_lon+nbp_lon
409    ELSE
410       iend=(jj_end-jj_begin)*nbp_lon+ii_end
411    ENDIF
412   
413    DO i = 1, maxrecv
414      IF (inforecv(i)%action .AND. inforecv(i)%nid .NE. -1) THEN
415          field(:) = -99999.
416          CALL prism_get_proto(inforecv(i)%nid, ktime, field(istart:iend), ierror)
417          tab_get(:,:,i) = RESHAPE(field(:),(/nbp_lon,jj_nb/))
418       
419          IF (ierror .NE. PRISM_Ok .AND. ierror.NE.PRISM_Recvd .AND. &
420             ierror.NE.PRISM_FromRest &
421             .AND. ierror.NE.PRISM_Input .AND. ierror.NE.PRISM_RecvOut &
422             .AND. ierror.NE.PRISM_FromRestOut) THEN
423              WRITE (lunout,*)  'Error with receiving filed : ', inforecv(i)%name, ktime   
424              abort_message=' Problem in prism_get_proto '
425              CALL abort_physic(modname,abort_message,1)
426          ENDIF
427      ENDIF
428    END DO
429   
430   
431  END SUBROUTINE fromcpl
432
433!
434!************************************************************************************
435!
436
437  SUBROUTINE intocpl(ktime, last, tab_put)
438! ======================================================================
439! L. Fairhead (09/2003) adapted From L.Z.X Li: this subroutine provides the
440! atmospheric coupling fields to the coupler with the psmile library.
441! IF last time step, writes output fields to binary files.
442! ======================================================================
443!
444!
445    USE print_control_mod, ONLY: lunout
446    USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat
447! Input arguments
448!************************************************************************************
449    INTEGER, INTENT(IN)                              :: ktime
450    LOGICAL, INTENT(IN)                              :: last
451    REAL, DIMENSION(nbp_lon, jj_nb, maxsend), INTENT(IN) :: tab_put
452
453! Local variables
454!************************************************************************************
455    LOGICAL                          :: checkout
456    INTEGER                          :: istart,iend
457    INTEGER                          :: wstart,wend
458    INTEGER                          :: ierror, i
459    REAL, DIMENSION(nbp_lon*jj_nb)       :: field
460    CHARACTER (len = 20),PARAMETER   :: modname = 'intocpl'
461    CHARACTER (len = 80)             :: abort_message
462
463!************************************************************************************
464    checkout=.FALSE.
465
466    WRITE(lunout,*) ' '
467    WRITE(lunout,*) 'Intocpl: sending fields to CPL, ktime= ', ktime
468    WRITE(lunout,*) 'last = ', last
469    WRITE(lunout,*)
470
471
472    istart=ii_begin
473    IF (is_south_pole_dyn) THEN
474       iend=(jj_end-jj_begin)*nbp_lon+nbp_lon
475    ELSE
476       iend=(jj_end-jj_begin)*nbp_lon+ii_end
477    ENDIF
478   
479    IF (checkout) THEN   
480       wstart=istart
481       wend=iend
482       IF (is_north_pole_dyn) wstart=istart+nbp_lon-1
483       IF (is_south_pole_dyn) wend=iend-nbp_lon+1
484       
485       DO i = 1, maxsend
486          IF (infosend(i)%action) THEN
487             field = RESHAPE(tab_put(:,:,i),(/nbp_lon*jj_nb/))
488             CALL writefield_phy(infosend(i)%name,field(wstart:wend),1)
489          END IF
490       END DO
491    END IF
492
493!************************************************************************************
494! PRISM_PUT
495!************************************************************************************
496
497    DO i = 1, maxsend
498      IF (infosend(i)%action .AND. infosend(i)%nid .NE. -1 ) THEN
499          field = RESHAPE(tab_put(:,:,i),(/nbp_lon*jj_nb/))
500          CALL prism_put_proto(infosend(i)%nid, ktime, field(istart:iend), ierror)
501         
502          IF (ierror .NE. PRISM_Ok .AND. ierror.NE.PRISM_Sent .AND. ierror.NE.PRISM_ToRest &
503             .AND. ierror.NE.PRISM_LocTrans .AND. ierror.NE.PRISM_Output .AND. &
504             ierror.NE.PRISM_SentOut .AND. ierror.NE.PRISM_ToRestOut) THEN
505              WRITE (lunout,*) 'Error with sending field :', infosend(i)%name, ktime   
506              abort_message=' Problem in prism_put_proto '
507              CALL abort_physic(modname,abort_message,1)
508          ENDIF
509      ENDIF
510    END DO
511   
512!************************************************************************************
513! Finalize PSMILE for the case is_sequential, if parallel finalization is done
514! from Finalize_parallel in dyn3dpar/parallel.F90
515!************************************************************************************
516
517    IF (last) THEN
518       IF (is_sequential) THEN
519          CALL prism_terminate_proto(ierror)
520          IF (ierror .NE. PRISM_Ok) THEN
521             abort_message=' Problem in prism_terminate_proto '
522             CALL abort_physic(modname,abort_message,1)
523          ENDIF
524       ENDIF
525    ENDIF
526   
527   
528  END SUBROUTINE intocpl
529
530#endif
531 
532END MODULE oasis
Note: See TracBrowser for help on using the repository browser.