convection.F90 24 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811
  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: CONVECTION
  15. !
  16. ! !DESCRIPTION: methods to init and apply convection
  17. !\\
  18. !\\
  19. ! !INTERFACE:
  20. !
  21. MODULE CONVECTION
  22. !
  23. ! !USES:
  24. !
  25. use GO, only : gol, goPr, goErr
  26. IMPLICIT NONE
  27. PRIVATE
  28. !
  29. ! !PUBLIC MEMBER FUNCTIONS:
  30. !
  31. public :: Convection_Init, Convection_Done
  32. public :: Convec
  33. !
  34. ! !PRIVATE DATA MEMBERS:
  35. !
  36. character(len=*), parameter :: mname = 'convection'
  37. !
  38. ! !REVISION HISTORY:
  39. ! 30 Mar 2010 - P. Le Sager - re-instate test on wet dep flag in convec
  40. ! 17 Nov 2011 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  41. !
  42. ! !REMARKS:
  43. !
  44. !EOP
  45. !---------------------------------------------------------------------------
  46. CONTAINS
  47. !--------------------------------------------------------------------------
  48. ! TM5 !
  49. !--------------------------------------------------------------------------
  50. !BOP
  51. !
  52. ! !IROUTINE: CONVECTION_INIT
  53. !
  54. ! !DESCRIPTION: set %used flag of releveant meteo to True
  55. !\\
  56. !\\
  57. ! !INTERFACE:
  58. !
  59. SUBROUTINE CONVECTION_INIT( status )
  60. !
  61. ! !USES:
  62. !
  63. use dims , only : nregions
  64. use Meteo , only : Set
  65. use Meteodata , only : entu_dat, entd_dat, detu_dat, detd_dat
  66. use Meteodata , only : mfw_dat, omega_dat
  67. #ifndef without_wet_deposition
  68. use Meteodata , only : cp_dat
  69. #endif
  70. !
  71. ! !OUTPUT PARAMETERS:
  72. !
  73. integer, intent(out) :: status
  74. !
  75. ! !REVISION HISTORY:
  76. !
  77. !EOP
  78. !------------------------------------------------------------------------
  79. !BOC
  80. character(len=*), parameter :: rname = mname//'/Convection_Init'
  81. integer :: region
  82. ! loop over all (zoom) regions:
  83. do region = 1, nregions
  84. ! entrements/detrements
  85. call Set( entu_dat(region), status, used=.true. )
  86. call Set( entd_dat(region), status, used=.true. )
  87. call Set( detu_dat(region), status, used=.true. )
  88. call Set( detd_dat(region), status, used=.true. )
  89. #ifndef without_wet_deposition
  90. ! convective precip:
  91. call Set( cp_dat(region), status, used=.true. )
  92. #endif
  93. ! omega is used when creating updr/downdr from raw data;
  94. ! computed from mfw :
  95. call Set( mfw_dat(region), status, used=.true. )
  96. call Set( omega_dat(region), status, used=.true. )
  97. end do
  98. ! ok
  99. status = 0
  100. END SUBROUTINE CONVECTION_INIT
  101. !EOC
  102. !--------------------------------------------------------------------------
  103. ! TM5 !
  104. !--------------------------------------------------------------------------
  105. !BOP
  106. !
  107. ! !IROUTINE: CONVECTION_DONE
  108. !
  109. ! !DESCRIPTION:
  110. !\\
  111. !\\
  112. ! !INTERFACE:
  113. !
  114. SUBROUTINE CONVECTION_DONE( status )
  115. !
  116. ! !OUTPUT PARAMETERS:
  117. !
  118. integer, intent(out) :: status
  119. !
  120. ! !REVISION HISTORY:
  121. !
  122. !EOP
  123. !------------------------------------------------------------------------
  124. !BOC
  125. character(len=*), parameter :: rname = mname//'/Convection_Done'
  126. ! nothing to do
  127. status = 0
  128. END SUBROUTINE CONVECTION_DONE
  129. !EOC
  130. !--------------------------------------------------------------------------
  131. ! TM5 !
  132. !--------------------------------------------------------------------------
  133. !BOP
  134. !
  135. ! !IROUTINE: CONVEC
  136. !
  137. ! !DESCRIPTION:
  138. !
  139. !**** convec - mix tracers by convection v 8.5
  140. !
  141. ! programmed by mh mpi HH 1-oct-1991
  142. ! modified by mh mpi HH 11-jun-1994
  143. !
  144. ! purpose
  145. ! -------
  146. ! mix tracers vertically by convection
  147. !
  148. ! interface
  149. ! ---------
  150. ! call convec
  151. !
  152. ! method
  153. ! ------
  154. ! the matrix conv is applied simultaneously to rm, rxm and rym.
  155. ! The diagonal elements of conv are applied to rzm.
  156. !
  157. ! externals
  158. ! ---------
  159. ! subroutines: tstamp
  160. !
  161. ! reference
  162. ! ---------
  163. ! see manual
  164. !-----------------------------------------------------------------------
  165. ! Following the method by Walter Guelle (1997) convective removal of
  166. ! tracers is added with some little changes:
  167. ! a vector cvsfac is added containing scavenging efficiencies (0-1)
  168. ! per tracer.
  169. ! aj, knmi may 1998
  170. ! Maarten Krol, july 2000, implemented zoom version
  171. ! lmax_conv is now set as the maximum layer to which convection
  172. ! is taken into account.
  173. !
  174. ! Second moments have been added: Bram Bregman, August 2004
  175. !\\
  176. !\\
  177. ! !INTERFACE:
  178. !
  179. SUBROUTINE CONVEC( region, status )
  180. !
  181. ! !USES:
  182. !
  183. use dims, only : tref !, dx, dy, xref, yref
  184. use dims, only : itau, okdebug, kdebug, nconv, lmax_conv
  185. use dims, only : nregions
  186. use dims, only : revert
  187. use global_data, only : region_dat, conv_dat, mass_dat
  188. use meteodata, only : entu_dat, entd_dat, detu_dat, detd_dat
  189. use meteodata, only : cp_dat
  190. use meteodata, only : m_dat
  191. use chem_param, only : ra, ntracet
  192. #ifdef with_budgets
  193. #ifndef without_wet_deposition
  194. use budget_global, only : nzon_vg
  195. use wet_deposition,only : sum_wetdep, buddep_dat
  196. #endif
  197. #endif
  198. use toolbox, only : lvlpress
  199. #ifndef without_wet_deposition
  200. use wet_deposition,only : cvsfac, cp_scale
  201. #endif
  202. use datetime, only : tstamp
  203. #ifdef with_tendencies
  204. use tracer_data, only : PLC_AddValue, plc_ipr_lwdep, plc_ipr_tcnvd, plc_kg_from_tm
  205. #endif
  206. use tm5_distgrid, only : dgrid, Get_DistGrid, reduce
  207. use partools, only : isRoot
  208. !
  209. ! !INPUT PARAMETERS:
  210. !
  211. integer,intent(in) :: region
  212. !
  213. ! !OUTPUT PARAMETERS:
  214. !
  215. integer,intent(out) :: status
  216. !
  217. ! !REVISION HISTORY:
  218. ! 17 Nov 2011 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  219. !
  220. !EOP
  221. !------------------------------------------------------------------------
  222. !BOC
  223. character(len=*), parameter :: rname = mname//'/CONVEC'
  224. real,dimension(:,:,:,:), pointer :: rm
  225. #ifdef slopes
  226. real,dimension(:,:,:,:), pointer :: rxm,rym,rzm
  227. #ifdef secmom
  228. real,dimension(:,:,:,:), pointer :: rxxm,rxym,rxzm,ryym,ryzm,rzzm
  229. #endif
  230. #endif
  231. real,dimension(:,:,:) , pointer :: m
  232. real,dimension(:,:,:) , pointer :: entu, entd, detu, detd
  233. #ifndef without_diffusion
  234. real,dimension(:,:,:) , pointer :: dkg
  235. #endif
  236. integer,dimension(:,:) , pointer :: cloud_base, cloud_top, cloud_lfs, zoomed
  237. integer,dimension(3) :: mxloc
  238. ! arrays for convection matrix
  239. ! conv1: convection matrix, lbdcv: wet removal matrix
  240. real,dimension(lmax_conv,lmax_conv) :: conv1
  241. real,dimension(lmax_conv,lmax_conv) :: lbdcv
  242. real,dimension(0:lmax_conv,lmax_conv) :: f
  243. real,dimension(0:lmax_conv,lmax_conv) :: fu
  244. real,dimension(0:lmax_conv) :: amu,amd
  245. real,dimension(lmax_conv) :: new_mass
  246. real,dimension(lmax_conv) :: zwetrm, zcnvd
  247. real,dimension(lmax_conv) :: zrm
  248. real,dimension(lmax_conv) :: srm
  249. #ifdef slopes
  250. real,dimension(lmax_conv) :: zrxm,zrym
  251. real,dimension(lmax_conv) :: srxm,srym,srzm
  252. #ifdef secmom
  253. real,dimension(lmax_conv) :: zrxxm,zrxym,zrxzm,zryym,zryzm
  254. real,dimension(lmax_conv) :: srxxm,srxym,srxzm,sryym,sryzm,srzzm
  255. #endif
  256. #endif
  257. integer :: n,l,k,j,i, i1, i2, j1, j2
  258. integer :: info,ld,li,kk, lshallowtop, nzv
  259. real :: scavf,dt,zxi,wetpct
  260. real :: minvalue, mxval,mnval,mxvald
  261. real :: minvalue_all, mxval_all,mnval_all,mxvald_all
  262. real :: sceffdeep, cp_value
  263. ! for global wetdep budget of tracer #1
  264. real :: g_sum_wet
  265. real, allocatable :: l_sum_wet(:,:)
  266. !====================== start ==================================
  267. sceffdeep = 1.0
  268. ! rain scale on this resolution: e.g. on 6x4 this read 1x1/(6x4)
  269. ! thus 1 mm/hr becomes 1mm/day
  270. ! cp_scaler = cp_scale*xref(region)*yref(region)/(dx*dy)
  271. !CMK removed after test study: hardly resolution dependent;;
  272. !-------- local indices
  273. call Get_DistGrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  274. #ifdef with_budgets
  275. #ifndef without_wet_deposition
  276. allocate(l_sum_wet(i1:i2,j1:j2))
  277. l_sum_wet = 0.0
  278. #endif
  279. #endif
  280. !-------- shortcuts
  281. m => m_dat(region)%data
  282. rm => mass_dat(region)%rm
  283. #ifdef slopes
  284. rxm => mass_dat(region)%rxm
  285. rym => mass_dat(region)%rym
  286. rzm => mass_dat(region)%rzm
  287. #ifdef secmom
  288. rxxm => mass_dat(region)%rxxm
  289. rxym => mass_dat(region)%rxym
  290. rxzm => mass_dat(region)%rxzm
  291. ryym => mass_dat(region)%ryym
  292. ryzm => mass_dat(region)%ryzm
  293. rzzm => mass_dat(region)%rzzm
  294. #endif
  295. #endif
  296. zoomed => region_dat(region)%zoomed
  297. entu => entu_dat(region)%data
  298. detu => detu_dat(region)%data
  299. entd => entd_dat(region)%data
  300. detd => detd_dat(region)%data
  301. cloud_top => conv_dat(region)%cloud_top
  302. cloud_base => conv_dat(region)%cloud_base
  303. cloud_lfs => conv_dat(region)%cloud_lfs
  304. #ifndef without_diffusion
  305. dkg => conv_dat(region)%dkg
  306. #endif
  307. if ( okdebug ) call tstamp(kdebug,itau,'convec ')
  308. if ( okdebug ) print *,'convec: Convection called (region)', region
  309. !----------- timestep and position in packed arrays
  310. dt = nconv/(2.*tref(region))
  311. if ( okdebug ) then
  312. minvalue = 2.0
  313. mxvald = -1.0
  314. mxval = 0.0
  315. mnval = 2.0
  316. end if
  317. !------------ start main loop over surface cells
  318. !
  319. ! PLS : re-incorporate the test on wet dep, from older cy3base/src.
  320. ! It hides cvsfac and cp_scale when wet dep is off.
  321. !$OMP PARALLEL &
  322. !$OMP default (none) &
  323. #ifdef with_budgets
  324. !$OMP shared (nzon_vg) &
  325. #endif
  326. !$OMP shared (region, ntracetloc, sceffdeep, dt, &
  327. !$OMP g_sum_wet, cloud_base, cloud_top, cloud_lfs, cp_dat, &
  328. !$OMP zoomed, rm, rxm, rym, rzm, m, entu, entd, detu, detd, &
  329. !$OMP dkg, revert, isr, ier, jsr, jer) &
  330. #ifndef without_wet_deposition
  331. !$OMP shared ( cvsfac,cp_scale ) &
  332. #ifdef with_budgets
  333. !$OMP shared ( buddep_dat,ra ) &
  334. #endif
  335. #endif
  336. !$OMP private (i, j, k,kk,l,n, li, ld, n, nzv, zxi, scavf, &
  337. !$OMP conv1, lbdcv, f, fu,amu, amd, zrm, zrxm, &
  338. !$OMP zrym, zwetrm, zcnvd, srm, srxm, srym, srzm, cp_value)
  339. do j = j1, j2
  340. do i = i1, i2
  341. #ifdef with_zoom
  342. if (zoomed(i,j) /= region) cycle
  343. #endif
  344. ! calculate convection matrix and directly apply to rm(i,j,*,1:ntracet)
  345. ! Set all amu's, amd's entu' etc to 0 for gridcells without clouds
  346. !
  347. ! updraft - amu is positive upward
  348. !
  349. f(:,:)=0.
  350. fu(:,:)=0.
  351. amu(:)=0.
  352. li = cloud_top(i,j)
  353. do k=1,li
  354. amu(k)=amu(k-1)+entu(i,j,k)-detu(i,j,k)
  355. if ( amu(k) > 0. ) then
  356. !
  357. ! we limit zxi from below at 0. in case
  358. ! there are inconsistent en/detrainment
  359. ! rates. The case zxi>1 should not occur
  360. ! since eu and du are >=0.
  361. !
  362. zxi=max(0.,1.-detu(i,j,k)/(amu(k-1)+entu(i,j,k)))
  363. else
  364. ! set massflux at upper boundary to strictly 0.
  365. amu(k)=0.
  366. zxi=0.
  367. end if
  368. do kk=1,k-1
  369. fu(k,kk)=fu(k-1,kk)*zxi !all previous levels
  370. end do
  371. !cmk: note: the fu(k-1,k) values are zero (see manual Heimann)
  372. fu(k,k)=entu(i,j,k)/m(i,j,k)*zxi
  373. end do !k
  374. !
  375. ! downdraft - amd is negative downward
  376. !
  377. amd(:)=0.
  378. ld = cloud_lfs(i,j)
  379. do k=ld,2,-1
  380. amd(k-1)=amd(k)-entd(i,j,k)+detd(i,j,k)
  381. if ( amd(k-1) < 0. ) then
  382. zxi=max(0.,1.+detd(i,j,k)/(amd(k)-entd(i,j,k)))
  383. else
  384. amd(k-1)=0.
  385. zxi=0.
  386. end if
  387. do kk=k+1,ld !
  388. f(k-1,kk)=f(k,kk)*zxi !f is here only fd downdraftmatrix
  389. end do !kk
  390. f(k-1,k)=-entd(i,j,k)/m(i,j,k)*zxi !note the negatives for j
  391. end do !k
  392. !
  393. ! add coefficients from up and downdraft together
  394. !
  395. do k=1,lmax_conv-1
  396. do kk=1,lmax_conv
  397. f(k,kk)=fu(k,kk)+f(k,kk)
  398. end do !kk
  399. ! add coefficient from subsidence
  400. !CMK SH f(k,k+1)=f(k,k+1)-(amu(k)+amd(k))/m(i,j,k+1)
  401. ! reported by Sander Houweling.........
  402. f(k,k+1)=f(k,k+1)-amu(k)/m(i,j,k+1)
  403. f(k,k )=f(k,k )-amd(k)/m(i,j,k )
  404. !
  405. #ifndef without_diffusion
  406. ! diffusion terms
  407. ! no diffusion in the stratosphere k>lmax_conv?
  408. f(k,k ) = f(k,k ) + dkg(i,j,k)/m(i,j,k )
  409. f(k,k+1) = f(k,k+1) - dkg(i,j,k)/m(i,j,k+1)
  410. #endif
  411. end do
  412. ! >>> wet removal >>>
  413. ! initialize wet-removal
  414. lbdcv(:,:) = 0.0
  415. #ifndef without_wet_deposition
  416. ! fill scaveging stuff if necessary
  417. !
  418. ! WARNING: For matrices f(k,kk), fu(k,kk), and fd(k,kk), index k
  419. ! corresponds to the upper boundary of level k
  420. ! and index j to level j.
  421. ! For matrices g(i,j,k,kk) and lbdcv(i,j,k,kk),
  422. ! both indices k and j
  423. ! correspond to respective levels and not to boundaries.
  424. !
  425. cp_value = cp_dat(region)%data(i,j,1) ! rain in m/s in region
  426. do k=2,lmax_conv
  427. do kk=1,k-1
  428. !cmk lbdcv(k,kk)=sceffdeep*cp_eff(region)%surf(i,j)* &
  429. !cmk dt*(fu(k-1,kk) - fu(k,kk))
  430. ! scale factor between 0 - 1, with 1 for cp > cp_scaler
  431. !CMK OLD lbdcv(k,kk)=sceffdeep*min(cp_value, cp_scaler)/cp_scaler * &
  432. lbdcv(k,kk)=sceffdeep*(1.0 - exp(-cp_value/cp_scale))* &
  433. dt*(fu(k-1,kk) - fu(k,kk))
  434. enddo !kk
  435. enddo !k
  436. ! CMK end
  437. !
  438. ! which ensures that the rate of scavenging is equal to dn/dt=-s Mu n
  439. !
  440. ! TvN (8 July 2004)
  441. #endif
  442. ! <<<
  443. !---- WG end 21/04/98 -----
  444. ! generate forward matrix
  445. do k=1,lmax_conv
  446. do kk=1,lmax_conv
  447. conv1(k,kk)=-dt*(f(k-1,kk)-f(k,kk))
  448. end do !kk
  449. conv1(k,k) = conv1(k,k) + 1.0
  450. end do !k
  451. call fastminv(conv1,lmax_conv)
  452. ! some checks...can be removed later!
  453. !if(okdebug) then
  454. ! minvalue = min(minvalue,minval(conv1))
  455. ! do l=1,lmax_conv
  456. ! mxval = max(mxval,sum(conv1(:,l)))
  457. ! mnval = min(mnval,sum(conv1(:,l)))
  458. ! end do
  459. !
  460. ! ! test deviations for a uniform mixing ratio of 1.....ppb
  461. !
  462. ! do l=1,lmax_conv
  463. ! new_mass(l) = dot_product(conv1(l,:),m(i,j,1:lmax_conv)*1e-9)
  464. ! end do
  465. ! mxvald = max(mxvald, maxval(abs((m(i,j,1:lmax_conv)*1e-9-new_mass)
  466. ! /(m(i,j,1:lmax_conv)*1e-9))))
  467. !end if !okdebug
  468. !
  469. ! directly apply the matrix:
  470. ! copy i,j
  471. if(revert==-1) conv1 = transpose(conv1)
  472. ! loop over local tracers
  473. TRACER : do n = 1, ntracet
  474. ! copy tracer concentrations for this column in srm/srxm/etc:
  475. do l=1,lmax_conv
  476. srm(l) = rm(i,j,l,n)
  477. #ifdef slopes
  478. srxm(l) = rxm(i,j,l,n)
  479. srym(l) = rym(i,j,l,n)
  480. srzm(l) = rzm(i,j,l,n)
  481. #ifdef secmom
  482. srxxm(l) = rxxm(i,j,l,n)
  483. srxym(l) = rxym(i,j,l,n)
  484. srxzm(l) = rxzm(i,j,l,n)
  485. sryym(l) = ryym(i,j,l,n)
  486. sryzm(l) = ryzm(i,j,l,n)
  487. srzzm(l) = rzzm(i,j,l,n)
  488. #endif
  489. #endif
  490. end do
  491. ! init budgets to zero:
  492. do l=1,lmax_conv
  493. zwetrm(l) = 0.0
  494. zcnvd (l) = 0.0
  495. end do
  496. ! init column of concentrations to zero:
  497. do l=1,lmax_conv
  498. zrm(l)=0.
  499. #ifdef slopes
  500. zrxm(l) = 0.0
  501. zrym(l) = 0.0
  502. #ifdef secmom
  503. zrxxm(l) = 0.0
  504. zrxym(l) = 0.0
  505. zrxzm(l) = 0.0
  506. zryym(l) = 0.0
  507. zryzm(l) = 0.0
  508. #endif
  509. #endif
  510. end do
  511. ! apply convection matrix:
  512. do l=1,lmax_conv
  513. #ifdef slopes
  514. srzm(l)=conv1(l,l)*srzm(l) !rzm is directly calculated
  515. srxm(l)=conv1(l,l)*srxm(l) !no convection -> no multip. !WP-Emma van der Veen
  516. srym(l)=conv1(l,l)*srym(l) !no convection -> no multip. !WP-Emma van der Veen
  517. #ifdef secmom
  518. srzzm(l)=conv1(l,l)*srzzm(l) !rzzm is directly calculated
  519. #endif
  520. #endif
  521. do k=1,lmax_conv
  522. ! species dependent removal....
  523. #ifndef without_wet_deposition
  524. scavf = min( lbdcv(l,k)*cvsfac(n), conv1(l,k) )
  525. #else
  526. scavf = 0.0
  527. #endif
  528. ! for budgets:
  529. zcnvd (l) = zcnvd (l) + conv1(l,k) * srm(k)
  530. zwetrm(l) = zwetrm(l) + scavf * srm(k) ! loss is positive
  531. ! store new concentrations in zrm/zrxm/etc
  532. ! using old concentrations in srm/srxm/etc
  533. zrm(l) = zrm(l) + (conv1(l,k)-scavf) * srm(k)
  534. #ifdef slopes
  535. !WP-vdV! zrxm (l) = zrxm (l) + (conv1(l,k)-scavf) * srxm(k)
  536. !WP-vdV! zrym (l) = zrym (l) + (conv1(l,k)-scavf) * srym(k)
  537. #ifdef secmom
  538. zrxxm(l) = zrxxm(l) + (conv1(l,k)-scavf) * srxxm(k)
  539. zrxym(l) = zrxym(l) + (conv1(l,k)-scavf) * srxym(k)
  540. zrxzm(l) = zrxzm(l) + (conv1(l,k)-scavf) * srxzm(k)
  541. zryym(l) = zryym(l) + (conv1(l,k)-scavf) * sryym(k)
  542. zryzm(l) = zryzm(l) + (conv1(l,k)-scavf) * sryzm(k)
  543. #endif
  544. #endif
  545. end do !k
  546. ! substract original tracer mass
  547. ! VH, TvN (2007-11-16)
  548. zcnvd(l)=zcnvd(l)-srm(l)
  549. end do !l
  550. ! restore concentrations:
  551. do l=1,lmax_conv
  552. rm(i,j,l,n)=zrm(l)
  553. #ifdef slopes
  554. !WP-vdV! rxm(i,j,l,n)=zrxm(l)
  555. !WP-vdV! rym(i,j,l,n)=zrym(l)
  556. rxm(i,j,l,n)=srxm(l) !VDV
  557. rym(i,j,l,n)=srym(l) !VDv
  558. rzm(i,j,l,n)=srzm(l)
  559. #ifdef secmom
  560. rxxm(i,j,l,n)=zrxxm(l)
  561. rxym(i,j,l,n)=zrxym(l)
  562. rxzm(i,j,l,n)=srxzm(l)
  563. ryym(i,j,l,n)=zryym(l)
  564. ryzm(i,j,l,n)=zryzm(l)
  565. rzzm(i,j,l,n)=srzzm(l)
  566. #endif
  567. #endif
  568. end do
  569. #ifdef with_budgets
  570. #ifndef without_wet_deposition
  571. if (cvsfac(n).gt.0.) then
  572. do l=1,lmax_conv
  573. nzv=nzon_vg(l)
  574. buddep_dat(region)%cp(i,j,nzv,n) = buddep_dat(region)%cp(i,j,nzv,n) + &
  575. zwetrm(l)/ra(n)*1.e3 !kg=> moles
  576. if ( n == 1 ) l_sum_wet(i,j) = l_sum_wet(i,j) + zwetrm(l) !in kg!
  577. end do !l
  578. end if
  579. #endif
  580. #endif
  581. #ifdef with_tendencies
  582. ! Add conv/diff and wet deposition budgets in kg to tendencies:
  583. do l = 1, lmax_conv
  584. call PLC_AddValue( region, plc_ipr_tcnvd, i, j, l, n, &
  585. zcnvd(l)*plc_kg_from_tm(n), status )
  586. IF_NOTOK_RETURN(status=1)
  587. call PLC_AddValue( region, plc_ipr_lwdep, i, j, l, n, &
  588. zwetrm(l)*plc_kg_from_tm(n), status )
  589. IF_NOTOK_RETURN(status=1)
  590. end do
  591. #endif
  592. end do TRACER !ntracet
  593. end do !i
  594. end do !j
  595. !$OMP END PARALLEL
  596. #ifndef without_wet_deposition
  597. #ifdef with_budgets
  598. call REDUCE( dgrid(region), l_sum_wet, 0, 'SUM', g_sum_wet, status)
  599. IF_NOTOK_RETURN(status=1)
  600. if ( isRoot ) sum_wetdep(region) = sum_wetdep(region) + g_sum_wet
  601. deallocate( l_sum_wet )
  602. #endif
  603. #endif
  604. !--- Cleanup
  605. nullify(m)
  606. nullify(rm)
  607. #ifdef slopes
  608. nullify(rxm)
  609. nullify(rym)
  610. nullify(rzm)
  611. #ifdef secmom
  612. nullify(rxxm)
  613. nullify(rxym)
  614. nullify(rxzm)
  615. nullify(ryym)
  616. nullify(ryzm)
  617. nullify(rzzm)
  618. #endif
  619. #endif
  620. nullify(zoomed)
  621. nullify(entu)
  622. nullify(detu)
  623. nullify(entd)
  624. nullify(detd)
  625. nullify(cloud_top)
  626. nullify(cloud_base)
  627. nullify(cloud_lfs)
  628. #ifndef without_diffusion
  629. nullify(dkg)
  630. #endif
  631. ! ok
  632. status = 0
  633. END SUBROUTINE CONVEC
  634. !EOC
  635. #ifdef with_lapack
  636. !-----------------------------------------------------------------------
  637. ! compute inverse matrix using LU decompositon and forward-backward
  638. ! substitution without pivoting (for easy vectorization)
  639. !
  640. ! LAPACK version, WP, Jan-2006
  641. !-----------------------------------------------------------------------
  642. subroutine fastminv(a,lm)
  643. use toolbox, only: escape_tm
  644. ! input/output
  645. integer,intent(in) :: lm
  646. real,dimension(lm,lm),intent(inout) :: a
  647. ! local
  648. integer :: istat
  649. real,dimension(lm) :: work
  650. real,dimension(lm) :: ipiv
  651. call DGETRF( lm, lm, a, lm, ipiv, istat )
  652. if(istat.ne.0) call escape_tm('fastminv: Failed to make LU decomposition ')
  653. call DGETRI( lm, a, lm, ipiv, work, lm, istat )
  654. if(istat.ne.0) call escape_tm('fastminv: Failed to make inverse matrix')
  655. end subroutine fastminv
  656. #else
  657. !-----------------------------------------------------------------------
  658. ! compute inverse matrix using LU decompositon and forward-backward
  659. ! substitution without pivoting (for easy vectorization)
  660. !
  661. ! Source: EMR and NR, sec. edition (p. 35ff).
  662. !
  663. ! mh, {Sat, Jun 4, 1994} 16:27:15
  664. !-----------------------------------------------------------------------
  665. subroutine fastminv(a,lm)
  666. ! input/output
  667. integer,intent(in) :: lm
  668. real,dimension(lm,lm),intent(inout) :: a
  669. ! local
  670. real,dimension(lm) :: y
  671. real,dimension(lm) :: b
  672. real,dimension(lm,lm) :: bb
  673. integer :: i,j,k
  674. do j=1,lm
  675. do i=1,j
  676. bb(i,j)=a(i,j)
  677. do k=1,i-1
  678. bb(i,j)=bb(i,j)-bb(i,k)*bb(k,j)
  679. end do
  680. end do
  681. do i=j+1,lm
  682. bb(i,j)=a(i,j)
  683. do k=1,j-1
  684. bb(i,j)=bb(i,j)-bb(i,k)*bb(k,j)
  685. end do
  686. bb(i,j)=bb(i,j)/bb(j,j)
  687. end do
  688. end do
  689. do k=1,lm
  690. do i=1,lm
  691. b(i)=0.
  692. end do
  693. b(k)=1.
  694. do i=1,lm
  695. y(i)=b(i)
  696. do j=1,i-1
  697. y(i)=y(i)-bb(i,j)*y(j)
  698. end do
  699. end do
  700. do i=lm,1,-1
  701. a(i,k)=y(i)
  702. do j=i+1,lm
  703. a(i,k)=a(i,k)-bb(i,j)*a(j,k)
  704. end do
  705. a(i,k)=a(i,k)/bb(i,i)
  706. end do
  707. end do
  708. end subroutine fastminv
  709. #endif
  710. end module convection