SURFEX v7.3
General documentation of Surfex
|
00001 ! ######### 00002 SUBROUTINE WRITESURF_OCEAN_n(HPROGRAM) 00003 ! ######################################## 00004 ! 00005 !!**** *WRITE_OCEAN_n* - writes OCEAN fields 00006 !! 00007 !! PURPOSE 00008 !! ------- 00009 !! 00010 !!** METHOD 00011 !! ------ 00012 !! 00013 !! EXTERNAL 00014 !! -------- 00015 !! 00016 !! 00017 !! IMPLICIT ARGUMENTS 00018 !! ------------------ 00019 !! 00020 !! REFERENCE 00021 !! --------- 00022 !! 00023 !! 00024 !! AUTHOR 00025 !! ------ 00026 !! C. Lebeaupin Brossier *Meteo France* 00027 !! 00028 !! MODIFICATIONS 00029 !! ------------- 00030 !! Original 03/2007 00031 !! Modified 07/2012, P. Le Moigne : CMO1D phasing 00032 !------------------------------------------------------------------------------- 00033 ! 00034 !* 0. DECLARATIONS 00035 ! ------------ 00036 ! 00037 USE MODD_OCEAN_n, ONLY : XSEAT, XSEAS, XSEAU, XSEAV, XSEAE, XSEABATH,& 00038 XSEAHMO, LMERCATOR, XDTFSOL, XDTFNSOL, XKMEL 00039 USE MODD_OCEAN_GRID_n 00040 ! 00041 USE MODI_WRITE_SURF 00042 ! 00043 USE MODD_OCEAN_REL_n, ONLY : XSEAT_REL, XSEAS_REL,XTAU_REL, & 00044 LREL_CUR,LREL_TS,LFLUX_NULL,XQCORR,LFLX_CORR,& 00045 XSEAU_REL,XSEAV_REL,LDIAPYCNAL 00046 ! 00047 USE YOMHOOK ,ONLY : LHOOK, DR_HOOK 00048 USE PARKIND1 ,ONLY : JPRB 00049 ! 00050 IMPLICIT NONE 00051 ! 00052 !* 0.1 Declarations of arguments 00053 ! ------------------------- 00054 ! 00055 CHARACTER(LEN=6), INTENT(IN) :: HPROGRAM ! program calling 00056 00057 ! 00058 !* 0.2 Declarations of local variables 00059 ! ------------------------------- 00060 ! 00061 INTEGER :: IRESP ! IRESP : return-code if a problem appears 00062 CHARACTER(LEN=12) :: YRECFM ! Name of the article to be read 00063 CHARACTER(LEN=4 ) :: YLVL 00064 CHARACTER(LEN=100):: YCOMMENT ! Comment string 00065 CHARACTER(LEN=14) :: YFORM ! Writing format 00066 ! 00067 INTEGER :: JLEVEL ! loop counter on oceanic levels 00068 ! 00069 REAL, DIMENSION(SIZE(XSEAT,1)) :: ZWORK ! 1D array to write data in file 00070 REAL(KIND=JPRB) :: ZHOOK_HANDLE 00071 ! 00072 !------------------------------------------------------------------------------- 00073 IF (LHOOK) CALL DR_HOOK('WRITESURF_OCEAN_N',0,ZHOOK_HANDLE) 00074 ! 00075 !* flag to define if OCEAN is used 00076 ! 00077 YRECFM='SEA_OCEAN' 00078 YCOMMENT='flag to use OCEAN model' 00079 CALL WRITE_SURF(HPROGRAM,YRECFM,LMERCATOR,IRESP,HCOMMENT=YCOMMENT) 00080 ! 00081 IF (.NOT. LMERCATOR .AND. LHOOK) CALL DR_HOOK('WRITESURF_OCEAN_N',1,ZHOOK_HANDLE) 00082 IF (.NOT. LMERCATOR) RETURN 00083 ! 00084 ! Relaxation time 00085 YRECFM='TAU_REL_OC' 00086 YCOMMENT='Relaxation time of ocean model (s)' 00087 CALL WRITE_SURF(HPROGRAM,YRECFM,XTAU_REL,IRESP,HCOMMENT=YCOMMENT) 00088 ! 00089 YRECFM='LREL_CUR_OC' 00090 YCOMMENT='FLAG FOR RELAXATION ON CURRENT' 00091 CALL WRITE_SURF(HPROGRAM,YRECFM,LREL_CUR,IRESP,HCOMMENT=YCOMMENT) 00092 ! 00093 YRECFM='LREL_TS_OC' 00094 YCOMMENT='FLAG FOR RELAXATION ON T,S' 00095 CALL WRITE_SURF(HPROGRAM,YRECFM,LREL_TS,IRESP,HCOMMENT=YCOMMENT) 00096 ! 00097 YRECFM='LFLX_NULL_OC' 00098 YCOMMENT='FLAG FOR ZERO FLUX ' 00099 CALL WRITE_SURF(HPROGRAM,YRECFM,LFLUX_NULL,IRESP,HCOMMENT=YCOMMENT) 00100 00101 YRECFM='LFLX_CORR_OC' 00102 YCOMMENT='FLAG FOR FLUX CORRECTION ' 00103 CALL WRITE_SURF(HPROGRAM,YRECFM,LFLX_CORR,IRESP,HCOMMENT=YCOMMENT) 00104 00105 YRECFM='CORR_FLX_OC' 00106 YCOMMENT='FLUX CORRECTION COEFF (W/M2/K)' 00107 CALL WRITE_SURF(HPROGRAM,YRECFM,XQCORR,IRESP,HCOMMENT=YCOMMENT) 00108 ! 00109 YRECFM='LDIAPYC_OC' 00110 YCOMMENT='FLAG FOR DIAPYCNAL MIXING ' 00111 CALL WRITE_SURF(HPROGRAM,YRECFM,LDIAPYCNAL,IRESP,HCOMMENT=YCOMMENT) 00112 ! 00113 !------------------------------------------------------------------------------- 00114 ! 00115 !* 3. Prognostic fields: 00116 ! ----------------- 00117 ! 00118 ! tendance surface liee au flux non solaire 00119 ! 00120 JLEVEL=1 00121 WRITE(YLVL,'(I4)') JLEVEL 00122 YRECFM='DTFNSOL'//ADJUSTL(YLVL(:LEN_TRIM(YLVL))) 00123 YFORM='(A11,I1.1,A5)' 00124 IF (JLEVEL >= 10) YFORM='(A11,I2.2,A5)' 00125 WRITE(YCOMMENT,FMT=YFORM) 'X_Y_DTFNSOL',JLEVEL,' (°C)' 00126 ZWORK=XDTFNSOL(:) 00127 CALL WRITE_SURF(HPROGRAM,YRECFM,ZWORK,IRESP,HCOMMENT=YCOMMENT) 00128 ! 00129 !* oceanic temperature 00130 ! 00131 DO JLEVEL=NOCKMIN+1,NOCKMAX 00132 WRITE(YLVL,'(I4)') JLEVEL 00133 YRECFM='TEMP_OC'//ADJUSTL(YLVL(:LEN_TRIM(YLVL))) 00134 YFORM='(A11,I1.1,A5)' 00135 IF (JLEVEL >= 10) YFORM='(A11,I2.2,A5)' 00136 WRITE(YCOMMENT,FMT=YFORM) 'X_Y_TEMP_OC',JLEVEL,' (°C)' 00137 ZWORK=XSEAT(:,JLEVEL) 00138 CALL WRITE_SURF(HPROGRAM,YRECFM,ZWORK,IRESP,HCOMMENT=YCOMMENT) 00139 END DO 00140 ! 00141 ! Budget terms lie au flux solaire 00142 DO JLEVEL=NOCKMIN+1,NOCKMAX 00143 WRITE(YLVL,'(I4)') JLEVEL 00144 YRECFM='DTFSOL'//ADJUSTL(YLVL(:LEN_TRIM(YLVL))) 00145 YFORM='(A11,I1.1,A5)' 00146 IF (JLEVEL >= 10) YFORM='(A11,I2.2,A5)' 00147 WRITE(YCOMMENT,FMT=YFORM) 'X_Y_DTFSOL',JLEVEL,' (°C/s)' 00148 ZWORK=XDTFSOL(:,JLEVEL) 00149 CALL WRITE_SURF(HPROGRAM,YRECFM,ZWORK,IRESP,HCOMMENT=YCOMMENT) 00150 END DO 00151 00152 ! 00153 !* relaxation profile for oceanic temperature 00154 DO JLEVEL=NOCKMIN+1,NOCKMAX 00155 WRITE(YLVL,'(I4)') JLEVEL 00156 YRECFM='T_OC_REL'//ADJUSTL(YLVL(:LEN_TRIM(YLVL))) 00157 YFORM='(A11,I1.1,A5)' 00158 IF (JLEVEL >= 10) YFORM='(A11,I2.2,A5)' 00159 WRITE(YCOMMENT,FMT=YFORM) 'X_Y_T_OC_REL',JLEVEL,' (°C)' 00160 ZWORK=XSEAT_REL(:,JLEVEL) 00161 CALL WRITE_SURF(HPROGRAM,YRECFM,ZWORK,IRESP,HCOMMENT=YCOMMENT) 00162 END DO 00163 ! 00164 !* oceanic salinity 00165 ! 00166 DO JLEVEL=NOCKMIN+1,NOCKMAX 00167 WRITE(YLVL,'(I4)') JLEVEL 00168 YRECFM='SALT_OC'//ADJUSTL(YLVL(:LEN_TRIM(YLVL))) 00169 YFORM='(A11,I1.1,A5)' 00170 IF (JLEVEL >= 10) YFORM='(A11,I2.2,A5)' 00171 WRITE(YCOMMENT,FMT=YFORM) 'X_Y_SALT_OC',JLEVEL,'(psu)' 00172 ZWORK=XSEAS(:,JLEVEL) 00173 CALL WRITE_SURF(HPROGRAM,YRECFM,ZWORK,IRESP,HCOMMENT=YCOMMENT) 00174 END DO 00175 ! 00176 !* oceanic salinity profile of relxation 00177 ! 00178 DO JLEVEL=NOCKMIN+1,NOCKMAX 00179 WRITE(YLVL,'(I4)') JLEVEL 00180 YRECFM='S_OC_REL'//ADJUSTL(YLVL(:LEN_TRIM(YLVL))) 00181 YFORM='(A11,I1.1,A5)' 00182 IF (JLEVEL >= 10) YFORM='(A11,I2.2,A5)' 00183 WRITE(YCOMMENT,FMT=YFORM) 'X_Y_S_OC_REL',JLEVEL,'(psu)' 00184 ZWORK=XSEAS_REL(:,JLEVEL) 00185 CALL WRITE_SURF(HPROGRAM,YRECFM,ZWORK,IRESP,HCOMMENT=YCOMMENT) 00186 END DO 00187 ! 00188 !* oceanic zonal current 00189 ! 00190 DO JLEVEL=NOCKMIN+1,NOCKMAX 00191 WRITE(YLVL,'(I4)') JLEVEL 00192 YRECFM='U_OC_REL'//ADJUSTL(YLVL(:LEN_TRIM(YLVL))) 00193 YFORM='(A11,I1.1,A5)' 00194 IF (JLEVEL >= 10) YFORM='(A11,I2.2,A5)' 00195 WRITE(YCOMMENT,FMT=YFORM) 'X_Y_U_OC_REL',JLEVEL,' M/S' 00196 ZWORK=XSEAU_REL(:,JLEVEL) 00197 CALL WRITE_SURF(HPROGRAM,YRECFM,ZWORK,IRESP,HCOMMENT=YCOMMENT) 00198 END DO 00199 ! 00200 DO JLEVEL=NOCKMIN+1,NOCKMAX 00201 WRITE(YLVL,'(I4)') JLEVEL 00202 YRECFM='UCUR_OC'//ADJUSTL(YLVL(:LEN_TRIM(YLVL))) 00203 YFORM='(A11,I1.1,A5)' 00204 IF (JLEVEL >= 10) YFORM='(A11,I2.2,A5)' 00205 WRITE(YCOMMENT,FMT=YFORM) 'X_Y_UCUR_OC',JLEVEL,' (m/s)' 00206 ZWORK=XSEAU(:,JLEVEL) 00207 CALL WRITE_SURF(HPROGRAM,YRECFM,ZWORK,IRESP,HCOMMENT=YCOMMENT) 00208 END DO 00209 ! 00210 !* oceanic meridian current 00211 ! 00212 DO JLEVEL=NOCKMIN+1,NOCKMAX 00213 WRITE(YLVL,'(I4)') JLEVEL 00214 YRECFM='V_OC_REL'//ADJUSTL(YLVL(:LEN_TRIM(YLVL))) 00215 YFORM='(A11,I1.1,A5)' 00216 IF (JLEVEL >= 10) YFORM='(A11,I2.2,A5)' 00217 WRITE(YCOMMENT,FMT=YFORM) 'X_Y_V_OC_REL',JLEVEL,' M/S' 00218 ZWORK=XSEAV_REL(:,JLEVEL) 00219 CALL WRITE_SURF(HPROGRAM,YRECFM,ZWORK,IRESP,HCOMMENT=YCOMMENT) 00220 END DO 00221 ! 00222 DO JLEVEL=NOCKMIN+1,NOCKMAX 00223 WRITE(YLVL,'(I4)') JLEVEL 00224 YRECFM='VCUR_OC'//ADJUSTL(YLVL(:LEN_TRIM(YLVL))) 00225 YFORM='(A11,I1.1,A5)' 00226 IF (JLEVEL >= 10) YFORM='(A11,I2.2,A5)' 00227 WRITE(YCOMMENT,FMT=YFORM) 'X_Y_VCUR_OC',JLEVEL,'(m/s)' 00228 ZWORK=XSEAV(:,JLEVEL) 00229 CALL WRITE_SURF(HPROGRAM,YRECFM,ZWORK,IRESP,HCOMMENT=YCOMMENT) 00230 END DO 00231 ! 00232 !* oceanic turbulent kinetic energy 00233 ! 00234 DO JLEVEL=NOCKMIN+1,NOCKMAX 00235 WRITE(YLVL,'(I4)') JLEVEL 00236 YRECFM='TKE_OC'//ADJUSTL(YLVL(:LEN_TRIM(YLVL))) 00237 YFORM='(A11,I1.1,A5)' 00238 IF (JLEVEL >= 10) YFORM='(A11,I2.2,A5)' 00239 WRITE(YCOMMENT,FMT=YFORM) 'X_Y_TKE_OC',JLEVEL,' (J)' 00240 ZWORK=XSEAE(:,JLEVEL) 00241 CALL WRITE_SURF(HPROGRAM,YRECFM,ZWORK,IRESP,HCOMMENT=YCOMMENT) 00242 END DO 00243 ! 00244 !* oceanic mixing coefficient 00245 ! 00246 DO JLEVEL=NOCKMIN+1,NOCKMAX 00247 WRITE(YLVL,'(I4)') JLEVEL 00248 YRECFM='KMEL_OC'//ADJUSTL(YLVL(:LEN_TRIM(YLVL))) 00249 YFORM='(A11,I1.1,A5)' 00250 IF (JLEVEL >= 10) YFORM='(A11,I2.2,A5)' 00251 WRITE(YCOMMENT,FMT=YFORM) 'X_Y_KMEL_OC',JLEVEL,' (m2/s2)' 00252 ZWORK=XKMEL(:,JLEVEL) 00253 CALL WRITE_SURF(HPROGRAM,YRECFM,ZWORK,IRESP,HCOMMENT=YCOMMENT) 00254 END DO 00255 00256 ! 00257 !------------------------------------------------------------------------------- 00258 ! 00259 !* 4. Semi-prognostic fields: 00260 ! ---------------------- 00261 !* bathymetry indice 00262 ! 00263 DO JLEVEL=NOCKMIN+1,NOCKMAX 00264 WRITE(YLVL,'(I4)') JLEVEL 00265 YRECFM='SEAINDBATH'//ADJUSTL(YLVL(:LEN_TRIM(YLVL))) 00266 YFORM='(A11,I1.1,A5)' 00267 IF (JLEVEL >= 10) YFORM='(A11,I2.2,A5)' 00268 WRITE(YCOMMENT,FMT=YFORM) 'X_Y_SEAINDBATH',JLEVEL,' (J)' 00269 ZWORK=XSEABATH(:,JLEVEL) 00270 CALL WRITE_SURF(HPROGRAM,YRECFM,ZWORK,IRESP,HCOMMENT=YCOMMENT) 00271 END DO 00272 ! 00273 !------------------------------------------------------------------------------- 00274 !Sea Surface Salinity 00275 ! 00276 YRECFM='SSS' 00277 YCOMMENT='SSS' 00278 CALL WRITE_SURF(HPROGRAM,YRECFM,XSEAS(:,NOCKMIN),IRESP,HCOMMENT=YCOMMENT) 00279 !------------------------------------------------------------------------------- 00280 !Sea Surface Salinity 00281 ! 00282 YRECFM='SEA_HMO' 00283 YCOMMENT='X_Y_' 00284 CALL WRITE_SURF(HPROGRAM,YRECFM,XSEAHMO(:),IRESP,HCOMMENT=YCOMMENT) 00285 !------------------------------------------------------------------------------- 00286 IF (LHOOK) CALL DR_HOOK('WRITESURF_OCEAN_N',1,ZHOOK_HANDLE) 00287 ! 00288 END SUBROUTINE WRITESURF_OCEAN_n