SURFEX v7.3
General documentation of Surfex
 All Classes Files Functions Variables Typedefs
/home/dasprezs/EXPORT_v7_3/src/SURFEX/read_oceann.F90
Go to the documentation of this file.
00001 !     #########
00002       SUBROUTINE READ_OCEAN_n(HPROGRAM)
00003 !     #########################################
00004 !
00005 !!****  *READ_OCEAN_n* - read oceanic variables
00006 !!
00007 !!
00008 !!    PURPOSE
00009 !!    -------
00010 !!
00011 !!**  METHOD
00012 !!    ------
00013 !!
00014 !!    EXTERNAL
00015 !!    --------
00016 !!
00017 !!
00018 !!    IMPLICIT ARGUMENTS
00019 !!    ------------------
00020 !!
00021 !!    REFERENCE
00022 !!    ---------
00023 !!
00024 !!
00025 !!    AUTHOR
00026 !!    ------
00027 !!      C. Lebeaupin Brossier   *Meteo France*  
00028 !!
00029 !!    MODIFICATIONS
00030 !!    -------------
00031 !!      Original    04/2007 
00032 !!      Modofied    07/2012, P. Le Moigne : CMO1D phasing
00033 !-------------------------------------------------------------------------------
00034 !
00035 !*       0.    DECLARATIONS
00036 !              ------------
00037 !
00038 USE MODD_SURF_PAR, ONLY : XUNDEF
00039 USE MODD_OCEAN_n
00040 USE MODD_DIAG_OCEAN_n
00041 USE MODD_OCEAN_GRID_n, ONLY : NOCKMIN,NOCKMAX
00042 !
00043 !
00044 USE MODI_READ_SURF
00045 USE MODI_OCEAN_MERCATORVERGRID
00046 !
00047 USE MODD_OCEAN_REL_n
00048 !
00049 USE YOMHOOK   ,ONLY : LHOOK,   DR_HOOK
00050 USE PARKIND1  ,ONLY : JPRB
00051 !
00052 USE MODI_GET_TYPE_DIM_n
00053 !
00054 IMPLICIT NONE
00055 !
00056 !*       0.1   Declarations of arguments
00057 !              -------------------------
00058 !
00059  CHARACTER(LEN=6),  INTENT(IN)  :: HPROGRAM ! calling program
00060 !
00061 !*       0.2   Declarations of local variables
00062 !              -------------------------------
00063 !
00064 INTEGER           :: ILU          ! 1D physical dimension
00065 !
00066 INTEGER           :: IRESP          ! Error code after redding
00067 !
00068  CHARACTER(LEN=4)  :: YLVL
00069 !
00070  CHARACTER(LEN=12) :: YRECFM         ! Name of the article to be read
00071  CHARACTER(LEN=14) :: YFORM          ! Writing format
00072 REAL, DIMENSION(:),ALLOCATABLE  :: ZWORK      ! 1D array to write data in file
00073 !
00074 INTEGER :: JLEVEL ! loop counter on oceanic levels
00075 INTEGER :: J      ! loop counter on sea grid points
00076 !
00077 INTEGER           :: IVERSION       ! surface version
00078 REAL(KIND=JPRB) :: ZHOOK_HANDLE
00079 !
00080 !-------------------------------------------------------------------------------
00081 IF (LHOOK) CALL DR_HOOK('READ_OCEAN_N',0,ZHOOK_HANDLE)
00082 NOCTCOUNT=0
00083 !
00084 YRECFM='VERSION'
00085  CALL READ_SURF(HPROGRAM,YRECFM,IVERSION,IRESP)
00086 !
00087 !* flag to use or not Ocean model
00088 !
00089 IF (IVERSION<=3) THEN
00090    LMERCATOR=.FALSE.
00091 ELSE
00092    YRECFM='SEA_OCEAN'
00093    CALL READ_SURF(HPROGRAM,YRECFM,LMERCATOR,IRESP)
00094 ENDIF
00095 !
00096 IF (.NOT. LMERCATOR) THEN
00097   ALLOCATE(XSEAT(0,0))
00098   ALLOCATE(XSEAT_REL(0,0))
00099   ALLOCATE(XSEAS(0,0))
00100   ALLOCATE(XSEAS_REL(0,0))
00101   ALLOCATE(XSEAU_REL(0,0))
00102   ALLOCATE(XSEAV_REL(0,0))
00103   ALLOCATE(XSEAU(0,0))
00104   ALLOCATE(XSEAV(0,0))
00105   ALLOCATE(XSEAE(0,0))
00106   ALLOCATE(XSEABATH(0,0))
00107   ALLOCATE(XSEAHMO(0))
00108   ALLOCATE(XLE        (0,0))
00109   ALLOCATE(XLK        (0,0))
00110   ALLOCATE(XKMEL      (0,0))
00111   ALLOCATE(XKMELM     (0,0))
00112   ALLOCATE(XSEATEND   (0))
00113   ALLOCATE(XDTFSOL(0,0))
00114   ALLOCATE(XDTFNSOL(0))
00115   IF (LHOOK) CALL DR_HOOK('READ_OCEAN_N',1,ZHOOK_HANDLE)
00116   RETURN
00117 ENDIF
00118 !
00119 !-------------------------------------------------------------------------------
00120  CALL OCEAN_MERCATORVERGRID
00121 ! Relaxation time and logical
00122 YRECFM='TAU_REL_OC'
00123  CALL READ_SURF(HPROGRAM,YRECFM,XTAU_REL,IRESP)
00124 !
00125 YRECFM='LREL_CUR_OC'
00126  CALL READ_SURF(HPROGRAM,YRECFM,LREL_CUR,IRESP)
00127 
00128 YRECFM='LREL_TS_OC'
00129  CALL READ_SURF(HPROGRAM,YRECFM,LREL_TS,IRESP)
00130 YRECFM='LFLX_NULL_OC'
00131  CALL READ_SURF(HPROGRAM,YRECFM,LFLUX_NULL,IRESP)
00132 YRECFM='LFLX_CORR_OC'
00133  CALL READ_SURF(HPROGRAM,YRECFM,LFLX_CORR,IRESP)
00134 YRECFM='CORR_FLX_OC'
00135  CALL READ_SURF(HPROGRAM,YRECFM,XQCORR,IRESP)
00136 YRECFM='LDIAPYC_OC'
00137  CALL READ_SURF(HPROGRAM,YRECFM,LDIAPYCNAL,IRESP)
00138 !
00139 !* 1D physical dimension
00140 !
00141 YRECFM='SIZE_SEA'
00142  CALL GET_TYPE_DIM_n('SEA   ',ILU)
00143 !
00144 !*       2.     Prognostic fields:
00145 !               -----------------
00146 !
00147 ALLOCATE(ZWORK(ILU))
00148 !* oceanic temperature
00149 !
00150 ALLOCATE(XSEAT(ILU,NOCKMIN:NOCKMAX))
00151 !
00152 DO JLEVEL=NOCKMIN+1,NOCKMAX
00153   WRITE(YLVL,'(I4)') JLEVEL
00154   YRECFM='TEMP_OC'//ADJUSTL(YLVL(:LEN_TRIM(YLVL)))
00155   CALL READ_SURF(HPROGRAM,YRECFM,ZWORK(:),IRESP)
00156   XSEAT(:,JLEVEL)=ZWORK(:)
00157 END DO
00158 XSEAT(:,NOCKMIN)=XSEAT(:,NOCKMIN+1)
00159 !
00160 !* relaxation profile for oceanic temperature
00161 !
00162 ALLOCATE(XSEAT_REL(ILU,NOCKMIN:NOCKMAX))
00163 !
00164 DO JLEVEL=NOCKMIN+1,NOCKMAX
00165   WRITE(YLVL,'(I4)') JLEVEL
00166   YRECFM='T_OC_REL'//ADJUSTL(YLVL(:LEN_TRIM(YLVL)))
00167   CALL READ_SURF(HPROGRAM,YRECFM,ZWORK(:),IRESP)
00168   XSEAT_REL(:,JLEVEL)=ZWORK(:)
00169 END DO
00170 XSEAT_REL(:,NOCKMIN)=XSEAT_REL(:,NOCKMIN+1)
00171 !
00172 !* oceanic salinity
00173 !
00174 ALLOCATE(XSEAS(ILU,NOCKMIN:NOCKMAX))
00175 !
00176 DO JLEVEL=NOCKMIN+1,NOCKMAX
00177   WRITE(YLVL,'(I4)') JLEVEL
00178   YRECFM='SALT_OC'//ADJUSTL(YLVL(:LEN_TRIM(YLVL)))
00179   CALL READ_SURF(HPROGRAM,YRECFM,ZWORK(:),IRESP)
00180   XSEAS(:,JLEVEL)=ZWORK(:)
00181 END DO
00182 XSEAS(:,NOCKMIN)=XSEAS(:,NOCKMIN+1)
00183 !
00184 !* oceanic salinity profile of relaxation
00185 !
00186 ALLOCATE(XSEAS_REL(ILU,NOCKMIN:NOCKMAX))
00187 !
00188 DO JLEVEL=NOCKMIN+1,NOCKMAX
00189   WRITE(YLVL,'(I4)') JLEVEL
00190   YRECFM='S_OC_REL'//ADJUSTL(YLVL(:LEN_TRIM(YLVL)))
00191   CALL READ_SURF(HPROGRAM,YRECFM,ZWORK(:),IRESP)
00192   XSEAS_REL(:,JLEVEL)=ZWORK(:)
00193 END DO
00194 XSEAS_REL(:,NOCKMIN)=XSEAS_REL(:,NOCKMIN+1)
00195 !
00196 !* oceanic current
00197 !
00198 ALLOCATE(XSEAU_REL(ILU,NOCKMIN:NOCKMAX))
00199 ALLOCATE(XSEAV_REL(ILU,NOCKMIN:NOCKMAX))
00200 !
00201 DO JLEVEL=NOCKMIN+1,NOCKMAX
00202   WRITE(YLVL,'(I4)') JLEVEL
00203   YRECFM='U_OC_REL'//ADJUSTL(YLVL(:LEN_TRIM(YLVL)))
00204   CALL READ_SURF(HPROGRAM,YRECFM,ZWORK(:),IRESP)
00205   XSEAU_REL(:,JLEVEL)=ZWORK(:)
00206 END DO
00207 XSEAU_REL(:,NOCKMIN)=XSEAU_REL(:,NOCKMIN+1)
00208 !
00209 DO JLEVEL=NOCKMIN+1,NOCKMAX
00210   WRITE(YLVL,'(I4)') JLEVEL
00211   YRECFM='V_OC_REL'//ADJUSTL(YLVL(:LEN_TRIM(YLVL)))
00212   CALL READ_SURF(HPROGRAM,YRECFM,ZWORK(:),IRESP)
00213   XSEAV_REL(:,JLEVEL)=ZWORK(:)
00214 END DO
00215 XSEAV_REL(:,NOCKMIN)=XSEAV_REL(:,NOCKMIN+1)
00216 !
00217 ALLOCATE(XSEAU(ILU,NOCKMIN:NOCKMAX))
00218 ALLOCATE(XSEAV(ILU,NOCKMIN:NOCKMAX))
00219 !
00220 DO JLEVEL=NOCKMIN+1,NOCKMAX
00221   WRITE(YLVL,'(I4)') JLEVEL
00222   YRECFM='UCUR_OC'//ADJUSTL(YLVL(:LEN_TRIM(YLVL)))
00223   CALL READ_SURF(HPROGRAM,YRECFM,ZWORK(:),IRESP)
00224   XSEAU(:,JLEVEL)=ZWORK(:)
00225 END DO
00226 DO JLEVEL=NOCKMIN+1,NOCKMAX
00227   WRITE(YLVL,'(I4)') JLEVEL
00228   YRECFM='VCUR_OC'//ADJUSTL(YLVL(:LEN_TRIM(YLVL)))
00229   CALL READ_SURF(HPROGRAM,YRECFM,ZWORK(:),IRESP)
00230   XSEAV(:,JLEVEL)=ZWORK(:)
00231 END DO
00232 XSEAU(:,NOCKMIN)=XSEAU(:,NOCKMIN+1)
00233 XSEAV(:,NOCKMIN)=XSEAV(:,NOCKMIN+1)
00234 !
00235 !* oceanic turbulent kinetic energy
00236 !
00237 ALLOCATE(XSEAE(ILU,NOCKMIN:NOCKMAX))
00238 !
00239 DO JLEVEL=NOCKMIN+1,NOCKMAX
00240   WRITE(YLVL,'(I4)') JLEVEL
00241   YRECFM='TKE_OC'//ADJUSTL(YLVL(:LEN_TRIM(YLVL)))
00242   CALL READ_SURF(HPROGRAM,YRECFM,ZWORK(:),IRESP)
00243   XSEAE(:,JLEVEL)=ZWORK(:)
00244 END DO
00245 XSEAE(:,NOCKMIN)=XSEAE(:,NOCKMIN+1)
00246 !
00247 !
00248 !-------------------------------------------------------------------------------
00249 !
00250 !*       4.     Semi-prognostic fields:
00251 !               ----------------------
00252 !
00253 !* bathymetry indice
00254 !
00255 ALLOCATE(XSEABATH(ILU,NOCKMIN:NOCKMAX))
00256 !
00257 DO JLEVEL=NOCKMIN+1,NOCKMAX
00258   WRITE(YLVL,'(I4)') JLEVEL
00259   YRECFM='SEAINDBATH'//ADJUSTL(YLVL(:LEN_TRIM(YLVL)))
00260   CALL READ_SURF(HPROGRAM,YRECFM,ZWORK(:),IRESP)
00261   XSEABATH(:,JLEVEL)=ZWORK(:)
00262 END DO
00263 XSEABATH(:,NOCKMIN)=1.
00264 !
00265 !-------------------------------------------------------------------------------
00266 !Complete undefined values for the oceanic 1D model convergence
00267 DO J=1,ILU
00268   DO JLEVEL=NOCKMIN+2,NOCKMAX
00269     IF (XSEABATH(J,JLEVEL)==0.) THEN
00270       XSEAT(J,JLEVEL)=XSEAT(J,JLEVEL-1)
00271       XSEAS(J,JLEVEL)=XSEAS(J,JLEVEL-1)
00272       XSEAU(J,JLEVEL)=XSEAU(J,JLEVEL-1)
00273       XSEAV(J,JLEVEL)=XSEAV(J,JLEVEL-1)
00274       XSEAE(J,JLEVEL)=XSEAE(J,JLEVEL-1)
00275       !
00276       XSEAT_REL(J,JLEVEL)=XSEAT_REL(J,JLEVEL-1)
00277       XSEAS_REL(J,JLEVEL)=XSEAS_REL(J,JLEVEL-1)
00278       XSEAU_REL(J,JLEVEL)=XSEAU_REL(J,JLEVEL-1)
00279       XSEAV_REL(J,JLEVEL)=XSEAV_REL(J,JLEVEL-1)      
00280     ENDIF
00281   ENDDO
00282 ENDDO
00283 !
00284 DEALLOCATE(ZWORK)
00285 !-------------------------------------------------------------------------------
00286 ALLOCATE(XSEAHMO(ILU))
00287 YRECFM='SEA_HMO'
00288  CALL READ_SURF(HPROGRAM,YRECFM,XSEAHMO(:),IRESP)
00289 !
00290 !-------------------------------------------------------------------------------
00291 ALLOCATE(XLE        (SIZE(XSEAT,1),NOCKMIN:NOCKMAX))
00292 ALLOCATE(XLK        (SIZE(XSEAT,1),NOCKMIN:NOCKMAX))
00293 ALLOCATE(XKMEL      (SIZE(XSEAT,1),NOCKMIN:NOCKMAX))
00294 ALLOCATE(XKMELM     (SIZE(XSEAT,1),NOCKMIN:NOCKMAX))
00295 XLE(:,:)    =XUNDEF
00296 XLK(:,:)    =XUNDEF
00297 XKMEL(:,:)  =XUNDEF
00298 XKMELM(:,:) =XUNDEF
00299 !
00300 ALLOCATE(XSEATEND   (SIZE(XSEAT,1)))
00301 XSEATEND(:) =XUNDEF
00302 !
00303 ALLOCATE(XDTFSOL(ILU,NOCKMIN:NOCKMAX))
00304 ALLOCATE(XDTFNSOL(ILU))
00305 !
00306 XDTFSOL(:,:) = XUNDEF 
00307 XDTFNSOL(:) = XUNDEF 
00308 !
00309 IF (LHOOK) CALL DR_HOOK('READ_OCEAN_N',1,ZHOOK_HANDLE)
00310 !
00311 !------------------------------------------------------------------------------
00312 END SUBROUTINE READ_OCEAN_n