MODULE cpl_oasis3 !!====================================================================== !! *** MODULE cpl_oasis *** !! Coupled O/A : coupled ocean-atmosphere case using OASIS3-MCT !!===================================================================== !! History : !! 9.0 ! 04-06 (R. Redler, NEC Laboratories Europe, Germany) Original code !! " " ! 04-11 (R. Redler, NEC Laboratories Europe; N. Keenlyside, W. Park, IFM-GEOMAR, Germany) revision !! " " ! 04-11 (V. Gayler, MPI M&D) Grid writing !! " " ! 05-08 (R. Redler, W. Park) frld initialization, paral(2) revision !! " " ! 05-09 (R. Redler) extended to allow for communication over root only !! " " ! 06-01 (W. Park) modification of physical part !! " " ! 06-02 (R. Redler, W. Park) buffer array fix for root exchange !! 3.4 ! 11-11 (C. Harris) Changes to allow mutiple category fields !!---------------------------------------------------------------------- !!---------------------------------------------------------------------- !! 'key_oasis3' coupled Ocean/Atmosphere via OASIS3 !!---------------------------------------------------------------------- !! cpl_init : initialization of coupled mode communication !! cpl_define : definition of grid and fields !! cpl_snd : snd out fields in coupled mode !! cpl_rcv : receive fields in coupled mode !! cpl_finalize : finalize the coupled mode communication !!---------------------------------------------------------------------- #if defined key_oasis3 || defined key_oasis3mct USE mod_oasis ! OASIS3-MCT module #endif USE par_oce ! ocean parameters USE dom_oce ! ocean space and time domain USE in_out_manager ! I/O manager USE lbclnk ! ocean lateral boundary conditions (or mpp link) #if defined key_cpl_rootexchg USE lib_mpp, only : mppsync USE lib_mpp, only : mppscatter,mppgather #endif IMPLICIT NONE PRIVATE PUBLIC cpl_init PUBLIC cpl_define PUBLIC cpl_snd PUBLIC cpl_rcv PUBLIC cpl_freq PUBLIC cpl_finalize #if defined key_mpp_mpi INCLUDE 'mpif.h' #endif INTEGER, PARAMETER :: localRoot = 0 LOGICAL :: commRank ! true for ranks doing OASIS communication #if defined key_cpl_rootexchg LOGICAL :: rootexchg =.true. ! logical switch #else LOGICAL :: rootexchg =.false. ! logical switch #endif INTEGER, PUBLIC :: OASIS_Rcv = 1 !: return code if received field INTEGER, PUBLIC :: OASIS_idle = 0 !: return code if nothing done by oasis INTEGER :: ncomp_id ! id returned by oasis_init_comp INTEGER :: nerror ! return error code #if ! defined key_oasis3 ! OASIS Variables not used. defined only for compilation purpose INTEGER :: OASIS_Out = -1 INTEGER :: OASIS_REAL = -1 INTEGER :: OASIS_Ok = -1 INTEGER :: OASIS_In = -1 INTEGER :: OASIS_Sent = -1 INTEGER :: OASIS_SentOut = -1 INTEGER :: OASIS_ToRest = -1 INTEGER :: OASIS_ToRestOut = -1 INTEGER :: OASIS_Recvd = -1 INTEGER :: OASIS_RecvOut = -1 INTEGER :: OASIS_FromRest = -1 INTEGER :: OASIS_FromRestOut = -1 #endif INTEGER, PUBLIC, PARAMETER :: nmaxfld=40 ! Maximum number of coupling fields INTEGER, PUBLIC, PARAMETER :: nmaxcat=5 ! Maximum number of coupling fields INTEGER, PUBLIC, PARAMETER :: nmaxcpl=5 ! Maximum number of coupling fields TYPE, PUBLIC :: FLD_CPL !: Type for coupling field information LOGICAL :: laction ! To be coupled or not CHARACTER(len = 8) :: clname ! Name of the coupling field CHARACTER(len = 1) :: clgrid ! Grid type REAL(wp) :: nsgn ! Control of the sign change INTEGER, DIMENSION(nmaxcat,nmaxcpl) :: nid ! Id of the field (no more than 9 categories and 9 extrena models) INTEGER :: nct ! Number of categories in field INTEGER :: ncplmodel ! Maximum number of models to/from which this variable may be sent/received END TYPE FLD_CPL TYPE(FLD_CPL), DIMENSION(nmaxfld), PUBLIC :: srcv, ssnd !: Coupling fields REAL(wp), DIMENSION(:,:), ALLOCATABLE :: exfld ! Temporary buffer for receiving REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: tbuf ! Temporary buffer for sending / receiving INTEGER, PUBLIC :: localComm !!---------------------------------------------------------------------- !! NEMO/OPA 3.3 , NEMO Consortium (2010) !! $Id$ !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) !!---------------------------------------------------------------------- CONTAINS SUBROUTINE cpl_init( kl_comm ) !!------------------------------------------------------------------- !! *** ROUTINE cpl_init *** !! !! ** Purpose : Initialize coupled mode communication for ocean !! exchange between AGCM, OGCM and COUPLER. (OASIS3 software) !! !! ** Method : OASIS3 MPI communication !!-------------------------------------------------------------------- INTEGER, INTENT(out) :: kl_comm ! local communicator of the model !!-------------------------------------------------------------------- ! WARNING: No write in numout in this routine !============================================ !------------------------------------------------------------------ ! 1st Initialize the OASIS system for the application !------------------------------------------------------------------ CALL oasis_init_comp ( ncomp_id, 'toyoce', nerror ) IF ( nerror /= OASIS_Ok ) & CALL oasis_abort (ncomp_id, 'cpl_init', 'Failure in oasis_init_comp') !------------------------------------------------------------------ ! 3rd Get an MPI communicator for OPA local communication !------------------------------------------------------------------ CALL oasis_get_localcomm ( kl_comm, nerror ) IF ( nerror /= OASIS_Ok ) & CALL oasis_abort (ncomp_id, 'cpl_init','Failure in oasis_get_localcomm' ) localComm = kl_comm ! END SUBROUTINE cpl_init SUBROUTINE cpl_define( krcv, ksnd, kcplmodel ) !!------------------------------------------------------------------- !! *** ROUTINE cpl_define *** !! !! ** Purpose : Define grid and field information for ocean !! exchange between AGCM, OGCM and COUPLER. (OASIS3 software) !! !! ** Method : OASIS3 MPI communication !!-------------------------------------------------------------------- INTEGER, INTENT(in) :: krcv, ksnd ! Number of received and sent coupling fields INTEGER, INTENT(in) :: kcplmodel ! Maximum number of models to/from which NEMO is potentialy sending/receiving data ! INTEGER :: id_part INTEGER :: paral(5) ! OASIS3 box partition INTEGER :: ishape(2,2) ! shape of arrays passed to PSMILe INTEGER :: ji,jc,jm ! local loop indicees CHARACTER(LEN=64) :: zclname CHARACTER(LEN=2) :: cli2 !!-------------------------------------------------------------------- IF(lwp) WRITE(numout,*) IF(lwp) WRITE(numout,*) 'cpl_define : initialization in coupled ocean/atmosphere case' IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~~~~' IF(lwp) WRITE(numout,*) commRank = .false. IF ( rootexchg ) THEN IF ( nproc == localRoot ) commRank = .true. ELSE commRank = .true. ENDIF IF( kcplmodel > nmaxcpl ) THEN CALL oasis_abort ( ncomp_id, 'cpl_define', 'kcplmodel is larger than nmaxcpl, increase nmaxcpl') ; RETURN ENDIF ! ! ... Define the shape for the area that excludes the halo ! For serial configuration (key_mpp_mpi not being active) ! nl* is set to the global values 1 and jp*glo. ! ishape(:,1) = (/ 1, nlei-nldi+1 /) ishape(:,2) = (/ 1, nlej-nldj+1 /) ! ! ! ----------------------------------------------------------------- ! ... Define the partition ! ----------------------------------------------------------------- IF ( rootexchg ) THEN paral(1) = 2 ! box partitioning paral(2) = 0 ! NEMO lower left corner global offset paral(3) = jpiglo ! local extent in i paral(4) = jpjglo ! local extent in j paral(5) = jpiglo ! global extent in x ELSE paral(1) = 2 ! box partitioning paral(2) = jpiglo * (nldj-1+njmpp-1) + (nldi-1+nimpp-1) ! NEMO lower left corner global offset paral(3) = nlei-nldi+1 ! local extent in i paral(4) = nlej-nldj+1 ! local extent in j paral(5) = jpiglo ! global extent in x IF( ln_ctl ) THEN WRITE(numout,*) ' multiexchg: paral (1:5)', paral WRITE(numout,*) ' multiexchg: jpi, jpj =', jpi, jpj WRITE(numout,*) ' multiexchg: nldi, nlei, nimpp =', nldi, nlei, nimpp WRITE(numout,*) ' multiexchg: nldj, nlej, njmpp =', nldj, nlej, njmpp ENDIF ENDIF IF ( commRank ) CALL oasis_def_partition ( id_part, paral, nerror ) ! ... Allocate memory for data exchange ! ALLOCATE(exfld(paral(3), paral(4)), stat = nerror) IF( nerror > 0 ) THEN CALL oasis_abort ( ncomp_id, 'cpl_define', 'Failure in allocating exfld') ; RETURN ENDIF IF ( rootexchg ) THEN ! Should possibly use one of the work arrays for tbuf really ALLOCATE(tbuf(jpi, jpj, jpnij), stat = nerror) IF( nerror > 0 ) THEN CALL oasis_abort ( ncomp_id, 'cpl_define', 'Failure in allocating tbuf') ; RETURN ENDIF ENDIF ! IF (commRank ) THEN ! ! ... Announce send variables. ! ssnd(:)%ncplmodel = kcplmodel ! DO ji = 1, ksnd IF ( ssnd(ji)%laction ) THEN IF( ssnd(ji)%nct > nmaxcat ) THEN CALL oasis_abort ( ncomp_id, 'cpl_define', 'Number of categories of '// & & TRIM(ssnd(ji)%clname)//' is larger than nmaxcat, increase nmaxcat' ) RETURN ENDIF DO jc = 1, ssnd(ji)%nct DO jm = 1, kcplmodel IF ( ssnd(ji)%nct .GT. 1 ) THEN WRITE(cli2,'(i2.2)') jc zclname = TRIM(ssnd(ji)%clname)//'_cat'//cli2 ELSE zclname = ssnd(ji)%clname ENDIF IF ( kcplmodel > 1 ) THEN WRITE(cli2,'(i2.2)') jm zclname = 'model'//cli2//'_'//TRIM(zclname) ENDIF #if defined key_agrif IF( agrif_fixed() /= 0 ) THEN zclname=TRIM(Agrif_CFixed())//'_'//TRIM(zclname) END IF #endif IF( ln_ctl ) WRITE(numout,*) "Define", ji, jc, jm, " "//TRIM(zclname), " for ", OASIS_Out CALL oasis_def_var (ssnd(ji)%nid(jc,jm), zclname, id_part , (/ 2, 0 /), & & OASIS_Out , ishape , OASIS_REAL, nerror ) IF ( nerror /= OASIS_Ok ) THEN WRITE(numout,*) 'Failed to define transient ', ji, jc, jm, " "//TRIM(zclname) CALL oasis_abort ( ssnd(ji)%nid(jc,jm), 'cpl_define', 'Failure in oasis_def_var' ) ENDIF IF( ln_ctl .AND. ssnd(ji)%nid(jc,jm) /= -1 ) WRITE(numout,*) "variable defined in the namcouple" IF( ln_ctl .AND. ssnd(ji)%nid(jc,jm) == -1 ) WRITE(numout,*) "variable NOT defined in the namcouple" END DO END DO ENDIF END DO ! ! ... Announce received variables. ! srcv(:)%ncplmodel = kcplmodel ! DO ji = 1, krcv IF ( srcv(ji)%laction ) THEN IF( srcv(ji)%nct > nmaxcat ) THEN CALL oasis_abort ( ncomp_id, 'cpl_define', 'Number of categories of '// & & TRIM(srcv(ji)%clname)//' is larger than nmaxcat, increase nmaxcat' ) RETURN ENDIF DO jc = 1, srcv(ji)%nct DO jm = 1, kcplmodel IF ( srcv(ji)%nct .GT. 1 ) THEN WRITE(cli2,'(i2.2)') jc zclname = TRIM(srcv(ji)%clname)//'_cat'//cli2 ELSE zclname = srcv(ji)%clname ENDIF IF ( kcplmodel > 1 ) THEN WRITE(cli2,'(i2.2)') jm zclname = 'model'//cli2//'_'//TRIM(zclname) ENDIF #if defined key_agrif IF( agrif_fixed() /= 0 ) THEN zclname=TRIM(Agrif_CFixed())//'_'//TRIM(zclname) END IF #endif IF( ln_ctl ) WRITE(numout,*) "Define", ji, jc, jm, " "//TRIM(zclname), " for ", OASIS_In CALL oasis_def_var (srcv(ji)%nid(jc,jm), zclname, id_part , (/ 2, 0 /), & & OASIS_In , ishape , OASIS_REAL, nerror ) IF ( nerror /= OASIS_Ok ) THEN WRITE(numout,*) 'Failed to define transient ', ji, jc, jm, " "//TRIM(zclname) CALL oasis_abort ( srcv(ji)%nid(jc,jm), 'cpl_define', 'Failure in oasis_def_var' ) ENDIF IF( ln_ctl .AND. srcv(ji)%nid(jc,jm) /= -1 ) WRITE(numout,*) "variable defined in the namcouple" IF( ln_ctl .AND. srcv(ji)%nid(jc,jm) == -1 ) WRITE(numout,*) "variable NOT defined in the namcouple" END DO END DO ENDIF END DO ! ENDIF ! commRank=true !------------------------------------------------------------------ ! End of definition phase !------------------------------------------------------------------ IF ( commRank ) THEN CALL oasis_enddef(nerror) IF( nerror /= OASIS_Ok ) CALL oasis_abort ( ncomp_id, 'cpl_define', 'Failure in oasis_enddef') ENDIF ! END SUBROUTINE cpl_define SUBROUTINE cpl_snd( kid, kstep, pdata, kinfo ) !!--------------------------------------------------------------------- !! *** ROUTINE cpl_snd *** !! !! ** Purpose : - At each coupling time-step,this routine sends fields !! like sst or ice cover to the coupler or remote application. !!---------------------------------------------------------------------- INTEGER , INTENT(in ) :: kid ! variable index in the array INTEGER , INTENT( out) :: kinfo ! OASIS3 info argument INTEGER , INTENT(in ) :: kstep ! ocean time-step in seconds REAL(wp), DIMENSION(:,:,:), INTENT(in ) :: pdata !! INTEGER :: jn,jc,jm ! local loop index !!-------------------------------------------------------------------- ! ! snd data to OASIS3 ! DO jc = 1, ssnd(kid)%nct DO jm = 1, ssnd(kid)%ncplmodel IF( ssnd(kid)%nid(jc,jm) /= -1 ) THEN IF ( rootexchg ) THEN ! ! collect data on the local root process ! CALL mppgather (pdata(:,:,jc),localRoot,tbuf) CALL mppsync IF ( nproc == localRoot ) THEN DO jn = 1, jpnij exfld(nimppt(jn)-1+nldit(jn):nimppt(jn)+nleit(jn)-1,njmppt(jn)-1+nldjt(jn):njmppt(jn)+nlejt(jn)-1)= & tbuf(nldit(jn):nleit(jn),nldjt(jn):nlejt(jn),jn) ENDDO ! snd data to OASIS3 CALL oasis_put ( ssnd(kid)%nid(jc,jm), kstep, exfld, kinfo ) ENDIF ELSE ! snd data to OASIS3 CALL oasis_put ( ssnd(kid)%nid(jc,jm), kstep, pdata(nldi:nlei, nldj:nlej,jc), kinfo ) ENDIF IF ( ln_ctl ) THEN IF ( kinfo == OASIS_Sent .OR. kinfo == OASIS_ToRest .OR. & & kinfo == OASIS_SentOut .OR. kinfo == OASIS_ToRestOut ) THEN WRITE(numout,*) '****************' WRITE(numout,*) 'oasis_put: Outgoing ', ssnd(kid)%clname WRITE(numout,*) 'oasis_put: ivarid ', ssnd(kid)%nid(jc,jm) WRITE(numout,*) 'oasis_put: kstep ', kstep WRITE(numout,*) 'oasis_put: info ', kinfo WRITE(numout,*) ' - Minimum value is ', MINVAL(pdata(:,:,jc)) WRITE(numout,*) ' - Maximum value is ', MAXVAL(pdata(:,:,jc)) WRITE(numout,*) ' - Sum value is ', SUM(pdata(:,:,jc)) WRITE(numout,*) '****************' ENDIF ENDIF ENDIF ENDDO ENDDO ! END SUBROUTINE cpl_snd SUBROUTINE cpl_rcv( kid, kstep, pdata, pmask, kinfo ) !!--------------------------------------------------------------------- !! *** ROUTINE cpl_rcv *** !! !! ** Purpose : - At each coupling time-step,this routine receives fields !! like stresses and fluxes from the coupler or remote application. !!---------------------------------------------------------------------- INTEGER , INTENT(in ) :: kid ! variable index in the array INTEGER , INTENT(in ) :: kstep ! ocean time-step in seconds REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pdata ! IN to keep the value if nothing is done REAL(wp), DIMENSION(:,:,:), INTENT(in ) :: pmask ! coupling mask INTEGER , INTENT( out) :: kinfo ! OASIS3 info argument !! INTEGER :: jn,jc,jm ! local loop index LOGICAL :: llaction, llfisrt !!-------------------------------------------------------------------- ! ! receive local data from OASIS3 on every process ! kinfo = OASIS_idle ! DO jc = 1, srcv(kid)%nct llfisrt = .TRUE. DO jm = 1, srcv(kid)%ncplmodel IF( srcv(kid)%nid(jc,jm) /= -1 ) THEN ! ! receive data from OASIS3 ! IF ( commRank ) CALL oasis_get ( srcv(kid)%nid(jc,jm), kstep, exfld, kinfo ) IF ( rootexchg ) CALL MPI_BCAST ( kinfo, 1, MPI_INTEGER, localRoot, localComm, nerror ) llaction = kinfo == OASIS_Recvd .OR. kinfo == OASIS_FromRest .OR. & & kinfo == OASIS_RecvOut .OR. kinfo == OASIS_FromRestOut IF ( ln_ctl ) WRITE(numout,*) "llaction, kinfo, kstep, ivarid: " , llaction, kinfo, kstep, srcv(kid)%nid(jc,jm) IF ( llaction ) THEN kinfo = OASIS_Rcv IF( llfisrt ) THEN IF ( rootexchg ) THEN ! distribute data to processes ! IF ( nproc == localRoot ) THEN DO jn = 1, jpnij tbuf(nldit(jn):nleit(jn),nldjt(jn):nlejt(jn),jn)= & exfld(nimppt(jn)-1+nldit(jn):nimppt(jn)+nleit(jn)-1,njmppt(jn)-1+nldjt(jn):njmppt(jn)+nlejt(jn)-1) ! NOTE: we are missing combining this with pmask (see else below) ENDDO ENDIF CALL mppscatter (tbuf,localRoot,pdata(:,:,jc)) CALL mppsync ELSE pdata(nldi:nlei, nldj:nlej, jc) = exfld(:,:) * pmask(nldi:nlei,nldj:nlej,jm) ENDIF llfisrt = .FALSE. ELSE pdata(nldi:nlei,nldj:nlej,jc) = pdata(nldi:nlei,nldj:nlej,jc) + exfld(:,:) * pmask(nldi:nlei,nldj:nlej,jm) ENDIF IF ( ln_ctl ) THEN WRITE(numout,*) '****************' WRITE(numout,*) 'oasis_get: Incoming ', srcv(kid)%clname WRITE(numout,*) 'oasis_get: ivarid ' , srcv(kid)%nid(jc,jm) WRITE(numout,*) 'oasis_get: kstep', kstep WRITE(numout,*) 'oasis_get: info ', kinfo WRITE(numout,*) ' - Minimum value is ', MINVAL(pdata(:,:,jc)) WRITE(numout,*) ' - Maximum value is ', MAXVAL(pdata(:,:,jc)) WRITE(numout,*) ' - Sum value is ', SUM(pdata(:,:,jc)) WRITE(numout,*) '****************' ENDIF ENDIF ENDIF ENDDO !--- Fill the overlap areas and extra hallows (mpp) !--- check periodicity conditions (all cases) IF( .not. llfisrt ) CALL lbc_lnk( pdata(:,:,jc), srcv(kid)%clgrid, srcv(kid)%nsgn ) ENDDO ! END SUBROUTINE cpl_rcv INTEGER FUNCTION cpl_freq( kid ) !!--------------------------------------------------------------------- !! *** ROUTINE cpl_freq *** !! !! ** Purpose : - send back the coupling frequency for a particular field !!---------------------------------------------------------------------- INTEGER,INTENT(in) :: kid ! variable index !! INTEGER :: info INTEGER, DIMENSION(1) :: itmp !!---------------------------------------------------------------------- #if defined key_oasis3 CALL oasis_get_freqs(kid, 1, itmp, info) cpl_freq = itmp(1) #endif #if defined key_oasis3mct cpl_freq = namflddti( kid ) #endif ! END FUNCTION cpl_freq SUBROUTINE cpl_finalize !!--------------------------------------------------------------------- !! *** ROUTINE cpl_finalize *** !! !! ** Purpose : - Finalizes the coupling. If MPI_init has not been !! called explicitly before cpl_init it will also close !! MPI communication. !!---------------------------------------------------------------------- ! DEALLOCATE( exfld ) IF ( rootexchg ) DEALLOCATE ( tbuf ) IF (nstop == 0) THEN CALL oasis_terminate( nerror ) ELSE CALL oasis_abort( ncomp_id, "cpl_finalize", "NEMO ABORT STOP" ) ENDIF ! END SUBROUTINE cpl_finalize #if ! defined key_oasis3 !!---------------------------------------------------------------------- !! No OASIS Library OASIS3 Dummy module... !!---------------------------------------------------------------------- SUBROUTINE oasis_init_comp(k1,cd1,k2) CHARACTER(*), INTENT(in ) :: cd1 INTEGER , INTENT( out) :: k1,k2 k1 = -1 ; k2 = -1 WRITE(numout,*) 'oasis_init_comp: Error you sould not be there...', cd1 END SUBROUTINE oasis_init_comp SUBROUTINE oasis_abort(k1,cd1,cd2) INTEGER , INTENT(in ) :: k1 CHARACTER(*), INTENT(in ) :: cd1,cd2 WRITE(numout,*) 'oasis_abort: Error you sould not be there...', cd1, cd2 END SUBROUTINE oasis_abort SUBROUTINE oasis_get_localcomm(k1,k2) INTEGER , INTENT( out) :: k1,k2 k1 = -1 ; k2 = -1 WRITE(numout,*) 'oasis_get_localcomm: Error you sould not be there...' END SUBROUTINE oasis_get_localcomm SUBROUTINE oasis_def_partition(k1,k2,k3) INTEGER , INTENT( out) :: k1,k3 INTEGER , INTENT(in ) :: k2(5) k1 = k2(1) ; k3 = k2(5) WRITE(numout,*) 'oasis_def_partition: Error you sould not be there...' END SUBROUTINE oasis_def_partition SUBROUTINE oasis_def_var(k1,cd1,k2,k3,k4,k5,k6,k7) CHARACTER(*), INTENT(in ) :: cd1 INTEGER , INTENT(in ) :: k2,k3(2),k4,k5(2,2),k6 INTEGER , INTENT( out) :: k1,k7 k1 = -1 ; k7 = -1 WRITE(numout,*) 'oasis_def_var: Error you sould not be there...', cd1 END SUBROUTINE oasis_def_var SUBROUTINE oasis_enddef(k1) INTEGER , INTENT( out) :: k1 k1 = -1 WRITE(numout,*) 'oasis_enddef: Error you sould not be there...' END SUBROUTINE oasis_enddef SUBROUTINE oasis_put(k1,k2,p1,k3) REAL(wp), DIMENSION(:,:), INTENT(in ) :: p1 INTEGER , INTENT(in ) :: k1,k2 INTEGER , INTENT( out) :: k3 k3 = -1 WRITE(numout,*) 'oasis_put: Error you sould not be there...' END SUBROUTINE oasis_put SUBROUTINE oasis_get(k1,k2,p1,k3) REAL(wp), DIMENSION(:,:), INTENT( out) :: p1 INTEGER , INTENT(in ) :: k1,k2 INTEGER , INTENT( out) :: k3 p1(1,1) = -1. ; k3 = -1 WRITE(numout,*) 'oasis_get: Error you sould not be there...' END SUBROUTINE oasis_get SUBROUTINE oasis_get_freqs(k1,k2,k3,k4) INTEGER , INTENT(in ) :: k1,k2 INTEGER, DIMENSION(1), INTENT( out) :: k3 INTEGER , INTENT( out) :: k4 k3(1) = k1 ; k4 = k2 WRITE(numout,*) 'oasis_get_freqs: Error you sould not be there...' END SUBROUTINE oasis_get_freqs SUBROUTINE oasis_terminate(k1) INTEGER , INTENT( out) :: k1 k1 = -1 WRITE(numout,*) 'oasis_terminate: Error you sould not be there...' END SUBROUTINE oasis_terminate #endif !!===================================================================== END MODULE cpl_oasis3