advecty__slopes.F90 18 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541
  1. !### macro's #####################################################
  2. !
  3. #define TRACEBACK write (gol,'("in ",a," (",a,", line",i5,")")') rname, __FILE__, __LINE__; call goErr
  4. #define IF_NOTOK_RETURN(action) if (status/=0) then; TRACEBACK; action; return; end if
  5. #define IF_ERROR_RETURN(action) if (status> 0) then; TRACEBACK; action; return; end if
  6. !
  7. #include "tm5.inc"
  8. !
  9. !-----------------------------------------------------------------------------
  10. ! TM5 !
  11. !-----------------------------------------------------------------------------
  12. !BOP
  13. !
  14. ! !MODULE: ADVECTY
  15. !
  16. ! !DESCRIPTION:
  17. !\\
  18. !\\
  19. ! !INTERFACE:
  20. !
  21. MODULE ADVECTY
  22. !
  23. ! !USES:
  24. !
  25. USE GO, ONLY : gol, goPr, goErr
  26. USE TM5_DISTGRID, ONLY : DGRID, GET_DISTGRID, UPDATE_HALO, REDUCE
  27. IMPLICIT NONE
  28. PRIVATE
  29. !
  30. ! !PUBLIC MEMBER FUNCTIONS:
  31. !
  32. PUBLIC :: advectyzoom
  33. !
  34. ! !REVISION HISTORY:
  35. ! 2 Feb 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  36. !
  37. ! !REMARKS:
  38. !
  39. !EOP
  40. !------------------------------------------------------------------------
  41. CONTAINS
  42. !--------------------------------------------------------------------------
  43. ! TM5 !
  44. !--------------------------------------------------------------------------
  45. !BOP
  46. !
  47. ! !IROUTINE: ADVECTYZOOM
  48. !
  49. ! !DESCRIPTION: wrapper around DYNAMV
  50. !\\
  51. !\\
  52. ! !INTERFACE:
  53. !
  54. SUBROUTINE ADVECTYZOOM( REGION )
  55. !
  56. ! !USES:
  57. !
  58. use dims, only : lm
  59. use global_data, only : wind_dat, mass_dat
  60. use partools, only : isRoot
  61. #ifdef with_budgets
  62. use budget_global, only : sum_advection
  63. #endif
  64. !
  65. ! !INPUT PARAMETERS:
  66. !
  67. integer, intent(in) :: region
  68. !
  69. ! !REVISION HISTORY:
  70. !
  71. ! written by patrick berkvens and mike botchev, march-june 1999
  72. ! updated and modified by MK, dec 2002
  73. !
  74. ! 2 Feb 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  75. !
  76. !EOP
  77. !------------------------------------------------------------------------
  78. !BOC
  79. integer :: istat
  80. real :: sum_old, sum_new
  81. ! ------ START
  82. #ifdef with_budgets
  83. ! total mass tracer #1
  84. call REDUCE( dgrid(region), mass_dat(region)%rm(:,:,:,1), mass_dat(region)%halo, 'SUM', sum_old, istat)
  85. #endif
  86. CALL DYNAMV(region)
  87. #ifdef with_budgets
  88. call REDUCE( dgrid(region), mass_dat(region)%rm(:,:,:,1), mass_dat(region)%halo, 'SUM', sum_new, istat)
  89. if ( isRoot ) sum_advection(region) = sum_advection(region) + sum_new - sum_old
  90. #endif
  91. END SUBROUTINE ADVECTYZOOM
  92. !EOC
  93. !--------------------------------------------------------------------------
  94. ! TM5 !
  95. !--------------------------------------------------------------------------
  96. !BOP
  97. !
  98. ! !IROUTINE: DYNAMV
  99. !
  100. ! !DESCRIPTION: south-north tracer transport: calculate amount of tracer moved in a south-north
  101. ! advection substep.
  102. !
  103. ! method : slopes scheme
  104. ! reference : Russel and Lerner, 1979
  105. !\\
  106. !\\
  107. ! !INTERFACE:
  108. !
  109. SUBROUTINE DYNAMV( region )!, is, ie, ls, le)
  110. !
  111. ! !USES:
  112. !
  113. use dims, only : im, jm, lm
  114. use dims, only : okdebug, xi, nxi
  115. use dims, only : limits, nloop_max, zero, one
  116. use global_data, only : wind_dat, mass_dat
  117. use meteodata , only : m_dat
  118. use toolbox, only : escape_tm
  119. #ifdef with_budgets
  120. use budget_global, only : budget_flux, jflux1, jflux2, apply_budget_global_flux
  121. #endif
  122. use partools , only : par_reduce
  123. use chem_param, only : ntracet, ra
  124. !
  125. ! !INPUT PARAMETERS:
  126. !
  127. integer,intent(in) :: region
  128. !
  129. ! !REVISION HISTORY:
  130. ! programmed by : mh, mpi HH, feb-jun 1994
  131. ! zoom version written by mike botchev, march-june 1999
  132. ! 2 Feb 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  133. !
  134. !EOP
  135. !------------------------------------------------------------------------
  136. !BOC
  137. real, dimension(:,:,:,:), pointer :: rm, rxm, rym, rzm
  138. real, dimension(:,:,:), pointer :: m, bm, mx
  139. real, dimension(:,:), allocatable :: mnew
  140. real, dimension(:,:), allocatable :: f,pf,fx,fz
  141. integer :: i,j,l,n, is,ie, js,je, ls,le, iee,iss
  142. integer :: i1, i2, j1, j2, jss, jee, status
  143. integer :: imr,imr2,jmr,lmr
  144. real :: sfs,sfzs,sfn,sfzn
  145. real :: my_beta, beta,mxval
  146. integer,dimension(3) :: mxloc
  147. integer :: iloop, nloop(lm(region)), nglob, offsetn
  148. logical :: cfl_ok, SouthPole, NorthPole, SouthBorder, NorthBorder
  149. integer :: lrg, redfact, ixe, ixs
  150. real :: summ
  151. integer :: ns, ne, lss, lee
  152. real :: l_xi, g_xi, l_nxi, g_nxi, g_nl
  153. integer,parameter :: max_nloop = 6
  154. !------ Start
  155. ! Indices
  156. CALL GET_DISTGRID( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2,&
  157. hasSouthPole=SouthPole, hasNorthPole=NorthPole, &
  158. hasSouthBorder=SouthBorder, hasNorthBorder=NorthBorder )
  159. lmr = lm(region)
  160. jmr = jm(region)
  161. is=i1 ; ie=i2 ; js=j1 ; je=j2 ; ls=1 ; le=lmr
  162. ! Work arrays
  163. m => m_dat(region)%data ! halo = 2
  164. rm => mass_dat(region)%rm ! halo = 2
  165. rxm => mass_dat(region)%rxm ! halo = 2
  166. rym => mass_dat(region)%rym ! halo = 2
  167. rzm => mass_dat(region)%rzm ! halo = 2
  168. bm => wind_dat(region)%bm ! halo = 1
  169. CALL UPDATE_HALO( dgrid(region), bm, wind_dat(region)%halo, status)
  170. CALL UPDATE_HALO( dgrid(region), m, m_dat(region)%halo, status)
  171. CALL UPDATE_HALO( dgrid(region), rm, mass_dat(region)%halo, status)
  172. CALL UPDATE_HALO( dgrid(region), rxm, mass_dat(region)%halo, status)
  173. CALL UPDATE_HALO( dgrid(region), rym, mass_dat(region)%halo, status)
  174. CALL UPDATE_HALO( dgrid(region), rzm, mass_dat(region)%halo, status)
  175. ! Reduce work domain if border of zoom region, and take care of poles
  176. if (NorthBorder) then ! j2=jmr by design
  177. if (NorthPole) je=j2-1
  178. endif
  179. if (SouthBorder) then ! j1=1 by design
  180. if (SouthPole) js=j1+1
  181. endif
  182. ! indices for mass update
  183. jss=js
  184. jee=je
  185. if (SouthPole) jss=js-1 ! =j1=1
  186. if (NorthPole) jee=je+1 ! =j2=jmr
  187. ! if requested limit meridional slopes such that no negative
  188. if (limits) then
  189. do n = 1, ntracet
  190. do l = ls, le
  191. do j = js - 1, je + 1
  192. do i = is, ie
  193. rym(i,j,l,n) = max (min (rym(i,j,l,n), rm(i,j,l,n)), -rm(i,j,l,n))
  194. end do
  195. end do
  196. end do
  197. end do
  198. end if
  199. g_nl = 0.0
  200. l_xi = 0.0
  201. l_nxi = 0.0
  202. allocate ( f( i1-1:i2+1, j1-1:j2+1 )) ! halo=1
  203. allocate ( pf( i1-1:i2+1, j1-1:j2+1 )) ! halo=1
  204. allocate ( fx( i1-1:i2+1, j1-1:j2+1 )) ! halo=1
  205. allocate ( fz( i1-1:i2+1, j1-1:j2+1 )) ! halo=1
  206. allocate ( mnew( i1:i2, j1:j2 )) ! halo=0
  207. allocate ( mx(i1-1:i2+1, j1-1:j2+1, ls:le)) ! halo=1
  208. ! ================= Find number of iterations needed (per layer)
  209. LEVELS: do l=ls,le
  210. nloop(l) = 1
  211. beta = 2.0
  212. do while(abs(beta) >= one .and. nloop(l) < max_nloop)
  213. mx(is:ie,js-1:je+1,1) = m(is:ie,js-1:je+1,l) ! use 1st layer of max as a 2D work array
  214. xloop : do iloop = 1, nloop(l)
  215. cfl_ok = .true.
  216. IJ: do j=js,je+1
  217. do i=is,ie
  218. if (bm(i,j,l)>=zero) then
  219. my_beta=bm(i,j,l)/mx(i,j-1,1)
  220. else
  221. my_beta=bm(i,j,l)/mx(i,j,1)
  222. end if
  223. if (abs(my_beta)>=one) exit IJ
  224. end do
  225. end do IJ
  226. my_beta = abs(my_beta)
  227. call Par_Reduce(my_beta, 'MAX', beta, status, all=.true.)
  228. if(status==1) then
  229. write(gol,*) "error par_reduce";call goPr
  230. end if
  231. if (abs(beta)>=one) then
  232. cfl_ok = .false.
  233. exit xloop
  234. end if
  235. ! else...
  236. if (cfl_ok.and.(iloop/=nloop(l))) then
  237. mx(is:ie,jss:jee,1) = mx(is:ie,jss:jee,1) + bm(is:ie,jss:jee,l) - bm(is:ie,jss+1:jee+1,l)
  238. CALL UPDATE_HALO( dgrid(region), mx(:,:,1), 1, status)
  239. end if
  240. end do xloop
  241. if ( .not. cfl_ok ) then
  242. beta = 2.0
  243. bm(:,:,l) = bm(:,:,l)*nloop(l)/(nloop(l)+1) ! reduce mass flux
  244. if ( okdebug) print *,'dynamv: ', i,j,l, &
  245. 'nloop increased to', nloop(l) + 1, beta
  246. nloop(l) = nloop(l) + 1
  247. if(nloop(l) == max_nloop) call escape_tm('nloop too high in dynamv')
  248. end if
  249. end do !while
  250. end do LEVELS
  251. ! advect info
  252. g_nl = max(g_nl, REAL(maxval(nloop)))
  253. ! ================= Loop over tracers and vertical layers to update masses
  254. ! Store init M, so we can get the loop over tracer outside the CFL.
  255. ! Thus we can limit the mpi comm, and switch Trac and levels loops => faster
  256. mx(i1-1:i2+1, j1-1:j2+1, ls:le) = m(i1-1:i2+1, j1-1:j2+1, ls:le)
  257. TRAC: do n=1,ntracet
  258. m(i1-1:i2+1, j1-1:j2+1, ls:le) = mx(i1-1:i2+1, j1-1:j2+1, ls:le) !reset
  259. LEV2: do l=ls,le
  260. CFL: DO iloop = 1,nloop(l)
  261. !-- calculate new air mass distribution
  262. mnew(is:ie,jss:jee) = m(is:ie,jss:jee,l) + &
  263. bm(is:ie,jss:jee,l) - bm(is:ie,jss+1:jee+1,l)
  264. !-- compute all inner fluxes
  265. ! f(*,j,*) is flux entering j-th cell from below?cmk
  266. do j=js+1,je
  267. do i=is,ie
  268. if (bm(i,j,l)>=zero) then
  269. beta=bm(i,j,l)/m(i,j-1,l)
  270. f(i,j)=beta*(rm(i,j-1,l,n)+(one-beta)*rym(i,j-1,l,n))
  271. pf(i,j)=bm(i,j,l)*(beta*beta*rym(i,j-1,l,n)-3.*f(i,j))
  272. fx(i,j)=beta*rxm(i,j-1,l,n)
  273. fz(i,j)=beta*rzm(i,j-1,l,n)
  274. else
  275. beta=bm(i,j,l)/m(i,j,l)
  276. f(i,j)=beta*(rm(i,j,l,n)-(one+beta)*rym(i,j,l,n))
  277. pf(i,j)=bm(i,j,l)*(beta*beta*rym(i,j,l,n)-3.*f(i,j))
  278. fx(i,j)=beta*rxm(i,j,l,n)
  279. fz(i,j)=beta*rzm(i,j,l,n)
  280. end if
  281. l_xi=max(l_xi,abs(beta))
  282. if (abs(beta)>=one) then
  283. l_nxi=l_nxi+1
  284. end if
  285. end do
  286. end do
  287. ! Case of JS
  288. if ( SouthPole ) then ! js=2
  289. do i=is,ie
  290. ! this is the specific for js-1:
  291. fz(i,1)=zero
  292. fx(i,1)=zero
  293. pf(i,1)=zero
  294. f(i,1)=zero
  295. if (bm(i,2,l)>=zero) then ! this is the specific case
  296. beta=bm(i,2,l)/m(i,1,l)
  297. f(i,2)=beta*rm(i,1,l,n)
  298. pf(i,2)=-3.*bm(i,2,l)*f(i,2)
  299. fx(i,2)=zero
  300. fz(i,2)=beta*rzm(i,1,l,n)
  301. else ! like above
  302. beta=bm(i,2,l)/m(i,2,l)
  303. f(i,2)=beta*(rm(i,2,l,n)-(one+beta)*rym(i,2,l,n))
  304. pf(i,2)=bm(i,2,l)*(beta*beta*rym(i,2,l,n)-3.*f(i,2))
  305. fx(i,2)=beta*rxm(i,2,l,n)
  306. fz(i,2)=beta*rzm(i,2,l,n)
  307. end if
  308. l_xi=max(l_xi,abs(beta))
  309. if (abs(beta)>=one) then
  310. l_nxi=l_nxi+1
  311. end if
  312. end do
  313. else ! do like above
  314. j = js
  315. do i=is,ie !no reduced grid allowed
  316. if (bm(i,j,l)>=zero) then
  317. beta=bm(i,j,l)/m(i,j-1,l)
  318. f(i,j)=beta*(rm(i,j-1,l,n)+(one-beta)*rym(i,j-1,l,n))
  319. pf(i,j)=bm(i,j,l)*(beta*beta*rym(i,j-1,l,n)-3.*f(i,j))
  320. fx(i,j)=beta*rxm(i,j-1,l,n)
  321. fz(i,j)=beta*rzm(i,j-1,l,n)
  322. else
  323. beta=bm(i,j,l)/m(i,j,l)
  324. f(i,j)=beta*(rm(i,j,l,n)-(one+beta)*rym(i,j,l,n))
  325. pf(i,j)=bm(i,j,l)*(beta*beta*rym(i,j,l,n)-3.*f(i,j))
  326. fx(i,j)=beta*rxm(i,j,l,n)
  327. fz(i,j)=beta*rzm(i,j,l,n)
  328. end if
  329. l_xi=max(l_xi,abs(beta))
  330. if (abs(beta)>=one) then
  331. l_nxi=l_nxi+1
  332. end if
  333. end do
  334. end if ! compute boundary fluxes south pole...
  335. ! Case of JE+1
  336. if ( NorthPole ) then ! je+1=jmr
  337. do i=is,ie
  338. if (bm(i,jmr,l)>=zero) then !like above
  339. beta=bm(i,jmr,l)/m(i,jmr-1,l)
  340. f(i,jmr)=beta*(rm(i,jmr-1,l,n)+(one-beta)*rym(i,jmr-1,l,n))
  341. pf(i,jmr)=bm(i,jmr,l)*(beta*beta*rym(i,jmr-1,l,n) - &
  342. 3.*f(i,jmr))
  343. fx(i,jmr)=beta*rxm(i,jmr-1,l,n)
  344. fz(i,jmr)=beta*rzm(i,jmr-1,l,n)
  345. else ! this is the specific case
  346. beta=bm(i,jmr,l)/m(i,jmr,l)
  347. ! for solid-body test and polar mixing
  348. f(i,jmr)=beta*rm(i,jmr,l,n)
  349. pf(i,jmr)=-3.*bm(i,jmr,l)*f(i,jmr)
  350. fx(i,jmr)=zero
  351. ! for solid-body test and polar mixing
  352. fz(i,jmr)=beta*rzm(i,jmr,l,n)
  353. end if
  354. l_xi=max(l_xi,abs(beta))
  355. if (abs(beta)>=one) then
  356. l_nxi=l_nxi+1
  357. end if
  358. end do
  359. else !zoom region not touching north pole
  360. j = je+1
  361. do i=is,ie !no reduced grid allowed
  362. if (bm(i,j,l)>=zero) then
  363. beta=bm(i,j,l)/m(i,j-1,l)
  364. f(i,j)=beta*(rm(i,j-1,l,n)+(one-beta)*rym(i,j-1,l,n))
  365. pf(i,j)=bm(i,j,l)*(beta*beta*rym(i,j-1,l,n)-3.*f(i,j))
  366. fx(i,j)=beta*rxm(i,j-1,l,n)
  367. fz(i,j)=beta*rzm(i,j-1,l,n)
  368. else
  369. beta=bm(i,j,l)/m(i,j,l)
  370. f(i,j)=beta*(rm(i,j,l,n)-(one+beta)*rym(i,j,l,n))
  371. pf(i,j)=bm(i,j,l)*(beta*beta*rym(i,j,l,n)-3.*f(i,j))
  372. fx(i,j)=beta*rxm(i,j,l,n)
  373. fz(i,j)=beta*rzm(i,j,l,n)
  374. end if
  375. l_xi=max(l_xi,abs(beta))
  376. if (abs(beta)>=one) then
  377. l_nxi=l_nxi+1
  378. end if
  379. end do
  380. end if ! compute boundary fluxes north pole...
  381. #ifdef with_budgets
  382. !
  383. ! Finished computing fluxes. Now apply
  384. !
  385. !-- first compute INCOMING fluxes from the south...
  386. if (apply_budget_global_flux) then
  387. do i=is,ie
  388. if ((jflux1(region) <= j2 ) .and. (jflux1(region) >= j1 )) then
  389. budget_flux(region)%flux_y1(i,l,n) = budget_flux(region)%flux_y1(i,l,n) + &
  390. f(i,jflux1(region))*1e3/ra(n)
  391. end if
  392. if ((jflux2(region) <= j2 ) .and. (jflux2(region) >= j1 )) then
  393. budget_flux(region)%flux_y2(i,l,n) = budget_flux(region)%flux_y2(i,l,n) + &
  394. f(i,jflux2(region))*1e3/ra(n)
  395. end if
  396. enddo
  397. endif
  398. #endif
  399. !
  400. !-- Calculate new tracer mass, and tracer mass slopes
  401. !
  402. do j=js,je
  403. do i=is,ie
  404. rm(i,j,l,n) = rm(i,j,l,n)+(f(i,j)-f(i,j+1))
  405. rym(i,j,l,n) = rym(i,j,l,n)+(pf(i,j)-pf(i,j+1) &
  406. - (bm(i,j,l)-bm(i,j+1,l))*rym(i,j,l,n) &
  407. + 3.*((bm(i,j,l)+bm(i,j+1,l))*rm(i,j,l,n) &
  408. - (f(i,j)+f(i,j+1))*m(i,j,l)))/mnew(i,j)
  409. if ( limits) rym(i,j,l,n) = max( min(rym(i,j,l,n), &
  410. rm(i,j,l,n)), -rm(i,j,l,n) )
  411. rxm(i,j,l,n)=rxm(i,j,l,n)+(fx(i,j)-fx(i,j+1))
  412. rzm(i,j,l,n)=rzm(i,j,l,n)+(fz(i,j)-fz(i,j+1))
  413. end do
  414. end do !forall
  415. ! case of JS-1 if south pole (js=2)
  416. if ( SouthPole ) then
  417. rm (is:ie, 1,l,n) = rm(is:ie, 1,l,n) - f(is:ie,2)
  418. rzm(is:ie, 1,l,n) = rzm(is:ie, 1,l,n) - fz(is:ie,2)
  419. end if
  420. ! case of JE+1, only if north pole (je+1=jmr)
  421. if ( NorthPole ) then
  422. rm (is:ie,jmr,l,n) = rm(is:ie,jmr,l,n) + f(is:ie,jmr)
  423. rzm(is:ie,jmr,l,n) = rzm(is:ie,jmr,l,n) + fz(is:ie,jmr)
  424. end if
  425. ! store new air mass in m array
  426. m(is:ie,jss:jee,l) = mnew(is:ie,jss:jee)
  427. ! update anything that changes and may be used with +/- index shift,
  428. ! but only if not the last step
  429. if (iloop/=nloop(l)) then
  430. !print*, "PLS-test-halo" (FIXME: j_only is not yet implemented but should reduce comm by a factor of 2 here)
  431. CALL UPDATE_HALO( dgrid(region), m(:,:,l), m_dat(region)%halo, status, j_only=.true.)
  432. CALL UPDATE_HALO( dgrid(region), rm(:,:,l,n), mass_dat(region)%halo, status, j_only=.true.)
  433. CALL UPDATE_HALO( dgrid(region), rxm(:,:,l,n), mass_dat(region)%halo, status, j_only=.true.)
  434. CALL UPDATE_HALO( dgrid(region), rym(:,:,l,n), mass_dat(region)%halo, status, j_only=.true.)
  435. CALL UPDATE_HALO( dgrid(region), rzm(:,:,l,n), mass_dat(region)%halo, status, j_only=.true.)
  436. end if
  437. end do CFL ! iloop
  438. end do LEV2 ! end of l-loop over vertical layers....
  439. end do TRAC ! n-loop
  440. ! restore old bm's
  441. do l=ls,le
  442. bm(:,:,l) = bm(:,:,l)*nloop(l)
  443. end do
  444. ! Gather more alorithm information
  445. CALL PAR_REDUCE( l_xi, 'max', g_xi, status, .true.)
  446. CALL PAR_REDUCE( l_nxi, 'sum', g_nxi, status, .true.)
  447. nloop_max(region,2) = max (nloop_max(region,2), nint (g_nl))
  448. xi(region,2) = max (xi(region,2), g_xi)
  449. nxi(region,2) = nxi(region,2) + nint (g_nxi)
  450. ! Done
  451. deallocate (mx, mnew)
  452. deallocate (f)
  453. deallocate (pf)
  454. deallocate (fx)
  455. deallocate (fz)
  456. nullify(m)
  457. nullify(rm)
  458. nullify(rxm)
  459. nullify(rym)
  460. nullify(rzm)
  461. nullify(bm)
  462. END SUBROUTINE DYNAMV
  463. !EOC
  464. END MODULE ADVECTY