Back to home page

MITgcm

 
 

    


File indexing completed on 2019-08-09 05:11:00 UTC

view on githubraw file Latest commit d200a168 on 2019-08-02 06:38:38 UTC
6bfee994e4 Jean*0001 #include "OCN_CPL_OPTIONS.h"
69e21e3ef0 Jean*0002 
79e04e6111 Jean*0003 CBOP 0
                0004 C !ROUTINE: OCN_IMPORT_ATMCONFIG
                0005 
                0006 C !INTERFACE:
69e21e3ef0 Jean*0007       SUBROUTINE OCN_IMPORT_ATMCONFIG( myThid )
                0008 
79e04e6111 Jean*0009 C !DESCRIPTION:
                0010 C     *==========================================================*
                0011 C     | SUBROUTINE OCN_IMPORT_ATMCONFIG
6bfee994e4 Jean*0012 C     | o Routine for importing atmos. config
                0013 C     |   into ocean component.
79e04e6111 Jean*0014 C     *==========================================================*
                0015 C     | This version talks to the MIT Coupler. It uses the
                0016 C     | MIT Coupler "checkpoint 1" library calls.
                0017 C     *==========================================================*
                0018 
                0019 C !USES:
                0020       IMPLICIT NONE
69e21e3ef0 Jean*0021 C     == Global variables ==
                0022 #include "SIZE.h"
                0023 #include "EEPARAMS.h"
                0024 #include "ATMIDS.h"
                0025 #include "OCNCPL.h"
                0026 
79e04e6111 Jean*0027 C !INPUT/OUTPUT PARAMETERS:
69e21e3ef0 Jean*0028 C     == Routine arguments ==
79e04e6111 Jean*0029 C     myThid :: Thread number for this instance of the routine
69e21e3ef0 Jean*0030       INTEGER myThid
79e04e6111 Jean*0031 CEOP
69e21e3ef0 Jean*0032 
6bfee994e4 Jean*0033 C !LOCAL VARIABLES:
                0034 C     == Local variables ==
                0035 C     i,j,bi,bj :: Loop counters
                0036       INTEGER i,j,bi,bj
                0037 
                0038 C-    Initialise land-mask
                0039       DO bj=myByLo(myThid),myByHi(myThid)
d200a16830 Jean*0040        DO bi=myBxLo(myThid),myBxHi(myThid)
6bfee994e4 Jean*0041         DO j=1-OLy,sNy+OLy
                0042          DO i=1-OLx,sNx+OLx
                0043           landMask(i,j,bi,bj) = 0.
                0044          ENDDO
                0045         ENDDO
                0046        ENDDO
                0047       ENDDO
                0048 
                0049 C-    Receive atmos. model configuration info.
79e04e6111 Jean*0050       _BARRIER
                0051       _BEGIN_MASTER( myThid )
6bfee994e4 Jean*0052 
                0053 C     o Import atmosphere model land-mask
79e04e6111 Jean*0054       CALL COMPRECV_R8TILES(
1a1199da4c Jean*0055      I              atmLandName, sNx, OLx, sNy, OLy, 1, nSx, nSy,
                0056      O              landMask )
6bfee994e4 Jean*0057 
79e04e6111 Jean*0058       _END_MASTER( myThid )
                0059       _BARRIER
69e21e3ef0 Jean*0060 
                0061       RETURN
                0062       END