outmea.f 9.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395
  1. ! ------------------------------------------------------------------
  2. !
  3. subroutine outmea(pt,ps,putot,pvtot,pzeta,psice,pflukhea,pflukwat,&
  4. & pub,pvb,pw,pconvad,ptaux,ptauy,ptbound,pfluhea, &
  5. & ppsi,pfluwat,ksum,knt)
  6. use lsgvar
  7. implicit none
  8. integer :: ksum,knt
  9. ! ------------------------------------------------------------------
  10. !
  11. !**** *outmea*
  12. !
  13. ! by U. Mikolajewicz 12/87.
  14. !
  15. ! Purpose.
  16. ! --------
  17. ! *outmea writes output data on tapeunit *nopost*.
  18. ! This file is further processed in the postprocessor.
  19. ! The file is defined during the run.
  20. ! The time for the next output is read from file and the
  21. ! name of the next output file is generated.
  22. !
  23. ! Input.
  24. ! ------
  25. ! common blocks /lsgfie/ and /lsgsur/
  26. !
  27. ! Output.
  28. ! -------
  29. ! file with output data.
  30. !
  31. ! Interface.
  32. ! ----------
  33. ! *call* *outmea*
  34. !
  35. ! Calls subroutines *datfrnt*.
  36. ! ------------------------------------------------------------------
  37. !
  38. !
  39. ! Dimension of parameters.
  40. ! -----------------------
  41. real (kind=8) :: pt(ien,jen,ken)
  42. real (kind=8) :: ps(ien,jen,ken)
  43. real (kind=8) :: putot(ien,jen,ken)
  44. real (kind=8) :: pvtot(ien,jen,ken)
  45. real (kind=8) :: pw(ien,jen,ken)
  46. real (kind=8) :: pub(ien,jen)
  47. real (kind=8) :: pvb(ien,jen)
  48. real (kind=8) :: pzeta(ien,jen)
  49. real (kind=8) :: psice(ien,jen)
  50. real (kind=8) :: pflukhea(ien,jen)
  51. real (kind=8) :: pflukwat(ien,jen)
  52. real (kind=8) :: pconvad(ien,jen,ken)
  53. real (kind=8) :: ptaux(ien,jen)
  54. real (kind=8) :: ptauy(ien,jen)
  55. real (kind=8) :: ptbound(ien,jen)
  56. real (kind=8) :: pfluhea(ien,jen)
  57. real (kind=8) :: ppsi(ien,jen)
  58. real (kind=8) :: pfluwat(ien,jen)
  59. !
  60. ! Declaration of local variables.
  61. ! -----------------------------
  62. character(len=8) filemea
  63. integer :: i,j,k,ji,iyy,idd,iyear,khelp,i3dfields,i2dfields
  64. integer :: jr,jlev
  65. real (kind=8) :: zero,zdum,zcode
  66. !
  67. ! 32bit variables for file <filemea>
  68. !
  69. integer (kind=4) :: jddr(512)
  70. real (kind=4) :: yddr(512)
  71. real (kind=4) :: ylen,yprel,ycode,ylev,ydum
  72. real (kind=4) :: yhelp(ien,jen)
  73. !
  74. ! arrays of default kind
  75. !
  76. integer :: iddr(512)
  77. real (kind=8) :: zddr(512)
  78. !
  79. !* 1. Actualize *iddr*.
  80. ! -----------------
  81. zero=0.
  82. zdum=zero
  83. !
  84. ! *iddr* and *zddr* are header fields for output file.
  85. !
  86. do ji=1,512
  87. iddr(ji)=nddr(ji)
  88. zddr(ji)=oddr(ji)
  89. end do
  90. iddr(404)=0
  91. iddr(22)=6*ken+8
  92. iddr(16)=6*ken+8
  93. iddr(11)=-9999
  94. iddr(12)=0
  95. iddr(508)=ksum
  96. call datfrnt(nt,iyy,idd)
  97. call datfrnt(nt-ksum+1,iddr(510),iddr(509))
  98. iddr(512)=iyy
  99. iddr(511)=idd
  100. ! iyear=iyy+1 ! wrong year
  101. iyear=iyy
  102. khelp=mod(iyear,100000)
  103. !
  104. ! Coupled with PlaSim: name without year to be edited by coupled script
  105. write (filemea,'(a)') 'LSG_outm'
  106. !
  107. ! V2.2 2005/09/14
  108. ! number of fields stored in iddr(16):
  109. ! 6 fields more than by Uwe M.,
  110. ! where 8 horzontal fields were stored (iddr(16)=6*ken+8):
  111. ! (2 more + 4 forcing fields):
  112. ! code -27: ppsi
  113. ! code -65: pfluwat
  114. ! (evtl. tracers here)
  115. ! code 52: ptaux
  116. ! code 53: ptauy
  117. ! code 92: ptbound (ts or ta, depending on nsmix)
  118. ! code 18: pfluhea (=pflukheat, depending on mix)
  119. !
  120. i3dfields=6
  121. i2dfields=8+6
  122. iddr(16)=i3dfields*ken+i2dfields
  123. iddr(22)=i3dfields*ken+i2dfields
  124. !
  125. 9880 format ( ' OUTMEA: output average ',a,' on file ',a, &
  126. & ' YYYYY -MMDD :',i5.2,i6.4,' +'/ &
  127. & ' OUTMEA: output average of',i6,' timesteps')
  128. write (no6,9880) '(real*4)',filemea,iyy,idd,ksum
  129. !
  130. ! Pot. temperature in Kelvin.
  131. !
  132. do jr=1,ken
  133. iddr(22+2*jr-1)=-2
  134. end do
  135. do jr=1,ken
  136. !
  137. ! Salinity.
  138. !
  139. iddr(22+2*ken+2*jr-1)=-5
  140. iddr(22+2*ken+2*jr)=nint(du(jr))
  141. !
  142. ! Velocities.
  143. !
  144. iddr(22+4*ken+2*jr-1)=-3
  145. iddr(22+4*ken+2*jr)=nint(du(jr))
  146. iddr(22+6*ken+2*jr-1)=-4
  147. iddr(22+6*ken+2*jr)=nint(du(jr))
  148. iddr(22+8*ken+2*jr-1)=-7
  149. iddr(22+8*ken+2*jr)=nint(dw(jr))
  150. end do
  151. !
  152. ! Barotropic velocities.
  153. !
  154. iddr(22+10*ken+1)=-37
  155. iddr(22+10*ken+2)=-100
  156. iddr(22+10*ken+3)=-38
  157. iddr(22+10*ken+4)=-100
  158. !
  159. ! Surface elevation.
  160. !
  161. iddr(22+10*ken+5)=-1
  162. iddr(22+10*ken+6)=-100
  163. !
  164. ! Ice thickness.
  165. !
  166. iddr(22+10*ken+7)=-13
  167. iddr(22+10*ken+8)=-100
  168. !
  169. ! Topography in vector-points.
  170. !
  171. iddr(22+10*ken+9)=-99
  172. iddr(22+10*ken+10)=-100
  173. !
  174. ! Topography in scalar-points.
  175. !
  176. iddr(22+10*ken+11)=-98
  177. iddr(22+10*ken+12)=-100
  178. !
  179. ! Fresh water fluxes due to newtonian coupling.
  180. !
  181. iddr(22+10*ken+13)=-67
  182. iddr(22+10*ken+14)=-100
  183. !
  184. ! Heat fluxes due to newtonian coupling.
  185. !
  186. iddr(22+10*ken+15)=-68
  187. iddr(22+10*ken+16)=-100
  188. !
  189. ! Convective adjustment.
  190. !
  191. iddr(22+10*ken+17)=-69
  192. iddr(22+10*ken+18)=-100
  193. !
  194. ! Horizontal stream function.
  195. !
  196. iddr(22+10*ken+19)=-27
  197. iddr(22+10*ken+20)=-100
  198. !
  199. ! Freshwater flux.
  200. !
  201. iddr(22+10*ken+21)=-65
  202. iddr(22+10*ken+22)=-100
  203. !
  204. do k=2,ken
  205. if (iddr(22+10*ken+22+(k-1)*2)>=300) cycle
  206. iddr(22+10*ken+21+(k-1)*2)=-69
  207. iddr(22+10*ken+22+(k-1)*2)=nint(dw(k-1))
  208. end do
  209. !
  210. !* 2. Write output on file *filemea*.
  211. ! -------------------------------
  212. !
  213. ! Write *iddr* and *zddr* on file.
  214. !
  215. ! write (no6,*) " open mean output file ",filemea
  216. open (nopost,file=filemea,access="sequential",form="unformatted")
  217. rewind nopost
  218. jddr(:)=iddr(:)
  219. yddr(:)=zddr(:)
  220. write (nopost) jddr
  221. write (nopost) yddr
  222. ylen=ien*jen
  223. yprel=6.
  224. !
  225. ! Write temperature.
  226. !
  227. zcode=-2
  228. do jlev=1,ken
  229. ycode=-2.
  230. write (nopost) yprel,ylen,ycode,yddr(10+jlev-1),ydum,ydum
  231. yhelp(:,:)=pt(:,:,jlev)+tkelvin
  232. write (nopost) yhelp(:,:)
  233. end do
  234. !
  235. ! Salinity.
  236. !
  237. zcode=-5
  238. do jlev=1,ken
  239. ycode=-5.
  240. write (nopost) yprel,ylen,ycode,yddr(10+jlev-1),ydum,ydum
  241. yhelp(:,:)=ps(:,:,jlev)
  242. write (nopost) yhelp(:,:)
  243. end do
  244. !
  245. ! Velocities.
  246. !
  247. zcode=-3
  248. do jlev=1,ken
  249. ycode=-3.
  250. write (nopost) yprel,ylen,ycode,yddr(10+jlev-1),ydum,ydum
  251. yhelp(:,:)=putot(:,:,jlev)
  252. write (nopost) yhelp(:,:)
  253. end do
  254. !
  255. zcode=-4
  256. do jlev=1,ken
  257. ycode=-4.
  258. write (nopost) yprel,ylen,ycode,yddr(10+jlev-1),ydum,ydum
  259. yhelp(:,:)=pvtot(:,:,jlev)
  260. write (nopost) yhelp(:,:)
  261. end do
  262. !
  263. zcode=-7
  264. do jlev=1,ken
  265. ycode=-7.
  266. write (nopost) yprel,ylen,ycode,yddr(10+jlev-1),ydum,ydum
  267. yhelp(:,:)=pw(:,:,jlev)
  268. write (nopost) yhelp
  269. end do
  270. !
  271. ! Now 2-d fields:
  272. !
  273. ! 2-d fields output in real*4
  274. ylev=-100.
  275. ycode=-37.
  276. write (nopost) yprel,ylen,ycode,ylev,ydum,ydum
  277. yhelp(:,:)=pub(:,:)
  278. write (nopost) yhelp
  279. ycode=-38.
  280. write (nopost) yprel,ylen,ycode,ylev,ydum,ydum
  281. yhelp(:,:)=pvb(:,:)
  282. write (nopost) yhelp
  283. !
  284. ! Surface elevation and ice thickness.
  285. !
  286. ycode=-1.
  287. write (nopost) yprel,ylen,ycode,ylev,ydum,ydum
  288. yhelp(:,:)=pzeta(:,:)
  289. write (nopost) yhelp
  290. !
  291. ycode=-13.
  292. write (nopost) yprel,ylen,ycode,ylev,ydum,ydum
  293. yhelp(:,:)=psice(:,:)
  294. write (nopost) yhelp
  295. !
  296. ! Horizontal stream function.
  297. !
  298. ycode=-27.
  299. write (nopost) yprel,ylen,ycode,ylev,ydum,ydum
  300. yhelp(:,:)=ppsi(:,:)
  301. write (nopost) yhelp
  302. !
  303. ! Store topography in vector-points.
  304. !
  305. ycode=-99.
  306. write (nopost) yprel,ylen,ycode,ylev,ydum,ydum
  307. yhelp(:,:)=depth(:,:)
  308. write (nopost) yhelp
  309. !
  310. ! Store topography in scalar-points.
  311. !
  312. ycode=-98.
  313. write (nopost) yprel,ylen,ycode,ylev,ydum,ydum
  314. yhelp(:,:)=depp(:,:)
  315. write (nopost) yhelp
  316. !
  317. ! Store *pfluwat*.
  318. !
  319. ycode=-65.
  320. write (nopost) yprel,ylen,ycode,ylev,ydum,ydum
  321. yhelp(:,:)=pfluwat(:,:)
  322. write (nopost) yhelp
  323. !
  324. ! Store *pflukwat*.
  325. !
  326. ycode=-67.
  327. write (nopost) yprel,ylen,ycode,ylev,ydum,ydum
  328. yhelp(:,:)=pflukwat(:,:)
  329. write (nopost) yhelp
  330. !
  331. ! Store *pflukhea*.
  332. !
  333. ycode=-68.
  334. write (nopost) yprel,ylen,ycode,ylev,ydum,ydum
  335. yhelp(:,:)=pflukhea(:,:)
  336. write (nopost) yhelp
  337. !
  338. ! Store convective adjustments.
  339. !
  340. ycode=-66.
  341. ylev=-100.
  342. write (nopost) yprel,ylen,ycode,ylev,ydum,ydum
  343. yhelp(:,:)=pconvad(:,:,1)
  344. write (nopost) yhelp
  345. ycode=-69.
  346. do k=2,ken
  347. ylev=dw(k-1)
  348. write (nopost) yprel,ylen,ycode,ylev,ydum,ydum
  349. yhelp(:,:)=pconvad(:,:,k)
  350. write (nopost) yhelp
  351. end do
  352. !
  353. ! Store *ptaux*
  354. !
  355. ycode=52.
  356. ylev=-100.
  357. write (nopost) yprel,ylen,ycode,ylev,ydum,ydum
  358. yhelp(:,:)=ptaux(:,:)*rhonul
  359. write (nopost) yhelp
  360. !
  361. ! Store *ptauy*
  362. !
  363. ycode=53.
  364. ylev=-100.
  365. write (nopost) yprel,ylen,ycode,ylev,ydum,ydum
  366. yhelp(:,:)=ptauy(:,:)*rhonul
  367. write (nopost) yhelp
  368. !
  369. ! Store *ptbound*
  370. !
  371. ycode=-92.
  372. ylev=-100.
  373. write (nopost) yprel,ylen,ycode,ylev,ydum,ydum
  374. yhelp(:,:)=ptbound(:,:)+tkelvin
  375. write (nopost) yhelp
  376. !
  377. ! Store *pfluhea*.
  378. !
  379. ycode=-18.
  380. write (nopost) yprel,ylen,ycode,ylev,ydum,ydum
  381. yhelp(:,:)=pfluhea(:,:)
  382. write (nopost) yhelp
  383. !
  384. !* 3. Store output on file *filemea*.
  385. ! -------------------------------
  386. close (nopost)
  387. return
  388. end subroutine outmea