m_fixhycom_eco_metno.F90 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424
  1. module m_fixhycom_eco_metno
  2. !Ehouarn: March 2011
  3. !fixanalysis: remapping of tracers after physical analysis!
  4. !use of remapping subroutines embedded in hycom to interpolate!
  5. !biogeochemical tracer on the analysis grid (dp)!
  6. !Remapping is realized after correction of negative anlaysis dp!
  7. contains
  8. subroutine hybgen_weno_coefs(s,dp,lc,ci,kk,ks,thin)
  9. implicit none
  10. integer kk,ks
  11. logical lc(kk)
  12. real s(kk,ks),dp(kk),ci(kk,ks,2),thin
  13. !
  14. !-----------------------------------------------------------------------
  15. ! 1) coefficents for remaping from one set of vertical cells to another.
  16. ! method: monotonic WENO-like alternative to PPM across each input cell
  17. ! a second order polynomial approximation of the profiles
  18. ! using a WENO reconciliation of the slopes to compute the
  19. ! interfacial values
  20. !
  21. ! REFERENCE?
  22. !
  23. ! 2) input arguments:
  24. ! s - initial scalar fields in pi-layer space
  25. ! dp - initial layer thicknesses (>=thin)
  26. ! lc - use PCM for selected layers
  27. ! kk - number of layers
  28. ! ks - number of fields
  29. ! thin - layer thickness (>0) that can be ignored
  30. !
  31. ! 3) output arguments:
  32. ! ci - coefficents for hybgen_weno_remap
  33. ! ci.1 is value at interface above
  34. ! ci.2 is value at interface below
  35. !
  36. ! 4) Laurent Debreu, Grenoble.
  37. ! Alan J. Wallcraft, Naval Research Laboratory, July 2008.
  38. !-----------------------------------------------------------------------
  39. !
  40. real, parameter :: dsmll=1.0e-8
  41. !
  42. integer j,i
  43. real q,q01,q02,q001,q002
  44. real qdpjm(kk),qdpjmjp(kk),dpjm2jp(kk)
  45. real zw(kk+1,3)
  46. !compute grid metrics
  47. do j=2,kk-1
  48. qdpjm( j) = 1.0/(dp(j-1) + dp(j))
  49. qdpjmjp(j) = 1.0/(dp(j-1) + dp(j) + dp(j+1))
  50. dpjm2jp(j) = dp(j-1) + 2.0*dp(j) + dp(j+1)
  51. enddo !j
  52. j=kk
  53. qdpjm( j) = 1.0/(dp(j-1) + dp(j))
  54. !
  55. do i= 1,ks
  56. do j=2,kk
  57. zw(j,3) = qdpjm(j)*(s(j,i)-s(j-1,i))
  58. enddo !j
  59. j = 1 !PCM first layer
  60. ci(j,i,1) = s(j,i)
  61. ci(j,i,2) = s(j,i)
  62. zw(j, 1) = 0.0
  63. zw(j, 2) = 0.0
  64. do j=2,kk-1
  65. if (lc(j) .or. dp(j).le.thin) then !use PCM
  66. ci(j,i,1) = s(j,i)
  67. ci(j,i,2) = s(j,i)
  68. zw(j, 1) = 0.0
  69. zw(j, 2) = 0.0
  70. else
  71. q001 = dp(j)*zw(j+1,3)
  72. q002 = dp(j)*zw(j, 3)
  73. if (q001*q002 < 0.0) then
  74. q001 = 0.0
  75. q002 = 0.0
  76. endif
  77. q01 = dpjm2jp(j)*zw(j+1,3)
  78. q02 = dpjm2jp(j)*zw(j, 3)
  79. if (abs(q001) > abs(q02)) then
  80. q001 = q02
  81. endif
  82. if (abs(q002) > abs(q01)) then
  83. q002 = q01
  84. endif
  85. q = (q001-q002)*qdpjmjp(j)
  86. q001 = q001-q*dp(j+1)
  87. q002 = q002+q*dp(j-1)
  88. ci(j,i,2) = s(j,i)+q001
  89. ci(j,i,1) = s(j,i)-q002
  90. zw( j,1) = (2.0*q001-q002)**2
  91. zw( j,2) = (2.0*q002-q001)**2
  92. endif !PCM:WEND
  93. enddo !j
  94. j = kk !PCM last layer
  95. ci(j,i,1) = s(j,i)
  96. ci(j,i,2) = s(j,i)
  97. zw(j, 1) = 0.0
  98. zw(j, 2) = 0.0
  99. do j=2,kk
  100. q002 = max(zw(j-1,2),dsmll)
  101. q001 = max(zw(j, 1),dsmll)
  102. zw(j,3) = (q001*ci(j-1,i,2)+q002*ci(j,i,1))/(q001+q002)
  103. enddo !j
  104. zw( 1,3) = 2.0*s( 1,i)-zw( 2,3) !not used?
  105. zw(kk+1,3) = 2.0*s(kk,i)-zw(kk,3) !not used?
  106. do j=2,kk-1
  107. if (.not.(lc(j) .or. dp(j).le.thin)) then !don't use PCM
  108. q01 = zw(j+1,3)-s(j,i)
  109. q02 = s(j,i)-zw(j,3)
  110. q001 = 2.0*q01
  111. q002 = 2.0*q02
  112. if (q01*q02 < 0.0) then
  113. q01 = 0.0
  114. q02 = 0.0
  115. elseif (abs(q01) > abs(q002)) then
  116. q01 = q002
  117. elseif (abs(q02) > abs(q001)) then
  118. q02 = q001
  119. endif
  120. ci(j,i,1) = s(j,i)-q02
  121. ci(j,i,2) = s(j,i)+q01
  122. endif !PCM:WEND
  123. enddo !j
  124. enddo !i
  125. return
  126. end subroutine hybgen_weno_coefs
  127. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  128. subroutine hybgen_weno_remap(si,pi,dpi,ci,&
  129. so,po,ki,ko,ks,thin)
  130. implicit none
  131. !
  132. integer ki,ko,ks
  133. real si(ki,ks),pi(ki+1),dpi(ki),ci(ki,ks,2),&
  134. so(ko,ks),po(ko+1),thin
  135. !
  136. !-----------------------------------------------------------------------
  137. ! 1) remap from one set of vertical cells to another.
  138. ! method: monotonic WENO-like alternative to PPM across each input cell
  139. ! a second order polynomial approximation of the profiles
  140. ! using a WENO reconciliation of the slopes to compute the
  141. ! interfacial values
  142. ! the output is the average of the interpolation
  143. ! profile across each output cell.
  144. !
  145. ! REFERENCE?
  146. !
  147. ! 2) input arguments:
  148. ! si - initial scalar fields in pi-layer space
  149. ! pi - initial layer interface depths (non-negative)
  150. ! pi( 1) is the surface
  151. ! pi(ki+1) is the bathymetry
  152. ! pi(k+1) >= pi(k)
  153. ! dpi - initial layer thicknesses (dpi(k)=pi(k+1)-pi(k))
  154. ! ci - coefficents from hybgen_weno_coefs
  155. ! ci.1 is value at interface above
  156. ! ci.2 is value at interface below
  157. ! ki - number of input layers
  158. ! ko - number of output layers
  159. ! ks - number of fields
  160. ! po - target interface depths (non-negative)
  161. ! po( 1) is the surface
  162. ! po(ko+1) is the bathymetry (== pi(ki+1))
  163. ! po(k+1) >= po(k)
  164. ! thin - layer thickness (>0) that can be ignored
  165. !
  166. ! 3) output arguments:
  167. ! so - scalar fields in po-layer space
  168. !
  169. ! 4) Laurent Debreu, Grenoble.
  170. ! Alan J. Wallcraft, Naval Research Laboratory, Aug. 2007.
  171. !-----------------------------------------------------------------------
  172. !
  173. integer i,k,l,lb,lt
  174. real dpb,dpt,qb0,qb1,qb2,qt0,qt1,qt2,xb,xt,zb,zt,zx,o
  175. real*8 sz
  176. !
  177. zx=pi(ki+1) !maximum depth
  178. zb=max(po(1),pi(1))
  179. lb=1
  180. do while (pi(lb+1).lt.zb .and. lb.lt.ki)
  181. lb=lb+1
  182. enddo
  183. do k= 1,ko !output layers
  184. zt = zb
  185. zb = min(po(k+1),zx)
  186. ! write(lp,*) 'k,zt,zb = ',k,zt,zb
  187. lt=lb !top will always correspond to bottom of previous
  188. !find input layer containing bottom output interface
  189. do while (pi(lb+1).lt.zb .and. lb.lt.ki)
  190. lb=lb+1
  191. enddo
  192. if (zb-zt.le.thin .or. zt.ge.zx) then
  193. if (k.ne.1) then
  194. !
  195. ! --- thin or bottomed layer, values taken from layer above
  196. !
  197. do i= 1,ks
  198. so(k,i) = so(k-1,i)
  199. enddo !i
  200. else !thin surface layer
  201. do i= 1,ks
  202. so(k,i) = si(k,i)
  203. enddo !i
  204. endif
  205. else
  206. ! form layer averages.
  207. !
  208. ! if (pi(lb).gt.zt) then
  209. ! write(lp,*) 'bad lb = ',lb
  210. ! stop
  211. ! endif
  212. xt=(zt-pi(lt))/max(dpi(lt),thin)
  213. xb=(zb-pi(lb))/max(dpi(lb),thin)
  214. if (lt.ne.lb) then !multiple layers
  215. dpt = pi(lt+1)-zt
  216. dpb = zb-pi(lb)
  217. qt1 = xt*(xt-1.0)
  218. qt2 = qt1+xt
  219. qt0 = 1.0-qt1-qt2
  220. qb1 = (xb-1.0)**2
  221. qb2 = qb1-1.0+xb
  222. qb0 = 1.0-qb1-qb2
  223. do i= 1,ks
  224. o = si((lt+lb)/2,i) !offset to reduce round-off
  225. sz = dpt*(qt0*(si(lt,i) -o) + &
  226. qt1*(ci(lt,i,1)-o) + &
  227. qt2*(ci(lt,i,2)-o) )
  228. do l=lt+1,lb-1
  229. sz = sz+dpi(l)*(si(l,i) - o)
  230. enddo !l
  231. sz = sz + dpb*(qb0*(si(lb,i) -o) + &
  232. qb1*(ci(lb,i,1)-o) + &
  233. qb2*(ci(lb,i,2)-o) )
  234. so(k,i) = o + sz/(zb-zt) !zb-zt>=thin
  235. enddo !i
  236. else !single layer
  237. qt1 = xb**2 + xt**2 + xb*xt + 1.0 - 2.0*(xb+xt)
  238. qt2 = qt1 - 1.0 + (xb+xt)
  239. qt0 = 1.0 - qt1 - qt2
  240. do i= 1,ks
  241. sz=qt0*(si(lt,i) -o) + &
  242. qt1*(ci(lt,i,1)-o) + &
  243. qt2*(ci(lt,i,2)-o)
  244. so(k,i) = o + sz
  245. enddo !i
  246. endif !layers
  247. endif !thin:std layer
  248. enddo !k
  249. return
  250. end subroutine hybgen_weno_remap
  251. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!1
  252. integer function tracr_get_incr(char2)
  253. !function returns the number of the tracer.
  254. !char => integer
  255. implicit none
  256. character(len=2) :: char2
  257. tracr_get_incr=-1
  258. select case (char2)
  259. case ('01')
  260. tracr_get_incr=1
  261. case ('02')
  262. tracr_get_incr=2
  263. case ('03')
  264. tracr_get_incr=3
  265. case ('04')
  266. tracr_get_incr=4
  267. case ('05')
  268. tracr_get_incr=5
  269. case ('06')
  270. tracr_get_incr=6
  271. case ('07')
  272. tracr_get_incr=7
  273. case ('08')
  274. tracr_get_incr=8
  275. case ('09')
  276. tracr_get_incr=9
  277. case ('10')
  278. tracr_get_incr=10
  279. case ('11')
  280. tracr_get_incr=11
  281. case default
  282. print *,'tracer unknown',char2
  283. end select
  284. return
  285. end function
  286. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!1
  287. integer function compute_kisop(temp,sal,nz)
  288. !function defines which layers are isopycnal
  289. implicit none
  290. integer::nz
  291. real,dimension(1:nz)::temp,sal
  292. real::eps,tmp
  293. integer::k
  294. eps=0.1
  295. do k=2,nz-1
  296. if (sig(temp(k),sal(k))-sig_ref(k).lt.0.)then
  297. tmp=(sig(temp(k),sal(k))-sig_ref(k))/(sig_ref(k)-sig_ref(k-1))
  298. if (tmp.gt.eps)then
  299. compute_kisop=k
  300. return
  301. else
  302. tmp=(sig(temp(k),sal(k))-sig_ref(k))/(sig_ref(k+1)-sig_ref(k))
  303. if(tmp.gt.eps)then
  304. compute_kisop=k
  305. return
  306. endif
  307. endif
  308. endif
  309. enddo
  310. compute_kisop=nz
  311. return
  312. end function
  313. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!1
  314. real function sig(t,s)
  315. !function returns the value of sigma_0
  316. !according to T and S
  317. implicit none
  318. real::t,s
  319. real :: c1,c2,c3,c4,c5,c6,c7
  320. c1=-1.36471E-01
  321. c2= 4.68181E-02
  322. c3= 8.07004E-01
  323. c4=-7.45353E-03
  324. c5=-2.94418E-03
  325. c6= 3.43570E-05
  326. c7= 3.48658E-05
  327. sig=(c1+c3*s+t*(c2+c5*s+t*(c4+c7*s+c6*t)))
  328. return
  329. end function
  330. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  331. real function sig_ref(k)
  332. !function return the value of the target density for a given layer
  333. implicit none
  334. integer::k
  335. select case (k)
  336. case (1)
  337. sig_ref=0.1
  338. case (2)
  339. sig_ref=0.2
  340. case (3)
  341. sig_ref=0.3
  342. case (4)
  343. sig_ref=0.4
  344. case (5)
  345. sig_ref=0.5
  346. case (6)
  347. sig_ref=24.05
  348. case (7)
  349. sig_ref=24.96
  350. case (8)
  351. sig_ref=25.68
  352. case (9)
  353. sig_ref=26.05
  354. case (10)
  355. sig_ref=26.30
  356. case (11)
  357. sig_ref=26.60
  358. case (12)
  359. sig_ref=26.83
  360. case (13)
  361. sig_ref=27.03
  362. case (14)
  363. sig_ref=27.20
  364. case (15)
  365. sig_ref=27.33
  366. case (16)
  367. sig_ref=27.46
  368. case (17)
  369. sig_ref=27.55
  370. case (18)
  371. sig_ref=27.66
  372. case (19)
  373. sig_ref=27.74
  374. case (20)
  375. sig_ref=27.82
  376. case (21)
  377. sig_ref=27.90
  378. case (22)
  379. sig_ref=27.97
  380. case (23)
  381. sig_ref=28.01
  382. case (24)
  383. sig_ref=28.04
  384. case (25)
  385. sig_ref=28.07
  386. case (26)
  387. sig_ref=28.09
  388. case (27)
  389. sig_ref=28.11
  390. case (28)
  391. sig_ref=28.13
  392. end select
  393. return
  394. end function
  395. end module