plasim_dummy.f90 25 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817
  1. module pumamod
  2. use resmod
  3. !
  4. ! Parameter
  5. !
  6. parameter(NPRO = NPRO_ATM) ! Number of processes (resmod)
  7. parameter(NLEV=NLEV_ATM)
  8. parameter(NLAT = NLAT_ATM) ! Number of latitudes (resmod)
  9. parameter(NLON = NLAT + NLAT) ! Number of longitudes
  10. parameter(NLPP = NLAT / NPRO) ! Latitudes per process
  11. parameter(NHOR = NLON * NLPP) ! Horizontal part
  12. parameter(NUGP = NLON * NLAT) ! Number of gridpoints
  13. parameter(NTRU = (NLON-1) / 3) ! Triangular truncation
  14. parameter(NTP1 = NTRU + 1) ! Truncation + 1
  15. parameter(NRSP =(NTRU+1)*(NTRU+2)) ! No of real global modes
  16. parameter(NCSP = NRSP / 2) ! No of complex global modes
  17. parameter(NSPP = (NRSP+NPRO-1)/NPRO) ! Modes per process
  18. parameter(NESP = NSPP * NPRO) ! Dim of spectral fields
  19. parameter(NROOT = 0) ! Master node
  20. parameter(SOLAR_DAY=86400.)
  21. parameter(TFREEZE=271.25)
  22. parameter(CRHOI=920.)
  23. !
  24. real (kind=8) :: gwd(NHOR)
  25. !
  26. real :: ccls(NHOR)
  27. real :: cctaux(NHOR,0:13) = 0.
  28. real :: cctauy(NHOR,0:13) = 0.
  29. real :: ccprl(NHOR,0:13) = 0.
  30. real :: ccprc(NHOR,0:13) = 0.
  31. real :: ccevap(NHOR,0:13) = 0.
  32. real :: ccroff(NHOR,0:13) = 0.
  33. real :: cciced(NHOR,0:13) = 0.
  34. real :: ccsnow(NHOR,0:13) = 0.
  35. real :: ccsst(NHOR,0:13) = 0.
  36. real :: ccswr(NHOR,0:13) = 0.
  37. real :: cclwr(NHOR,0:13) = 0.
  38. real :: ccshf(NHOR,0:13) = 0.
  39. real :: cclhf(NHOR,0:13) = 0.
  40. !
  41. real :: cltaux(NHOR) = 0.
  42. real :: cltauy(NHOR) = 0.
  43. real :: clprl(NHOR) = 0.
  44. real :: clprc(NHOR) = 0.
  45. real :: clevap(NHOR) = 0.
  46. real :: clroff(NHOR) = 0.
  47. real :: cliced(NHOR) = 0.
  48. real :: clsnow(NHOR) = 0.
  49. real :: clsst(NHOR) = 0.
  50. real :: clswr(NHOR) = 0.
  51. real :: cllwr(NHOR) = 0.
  52. real :: clshf(NHOR) = 0.
  53. real :: cllhf(NHOR) = 0.
  54. !
  55. real :: yls(NHOR) = 0.
  56. real :: ytaux(NHOR) = 0.
  57. real :: ytauy(NHOR) = 0.
  58. real :: yprl(NHOR) = 0.
  59. real :: yprc(NHOR) = 0.
  60. real :: yevap(NHOR) = 0.
  61. real :: ysst(NHOR) = 0.
  62. real :: ypme(NHOR) = 0.
  63. real :: yroff(NHOR) = 0.
  64. real :: yheat(NHOR) = 0.
  65. real :: yicesnow(NHOR)= 0.
  66. real :: yfldo(NHOR) = 0.
  67. real :: ymld(NHOR) = 50.
  68. !
  69. real :: ols(NHOR) = 0.
  70. real :: otaux(NHOR) = 0.
  71. real :: otauy(NHOR) = 0.
  72. real :: oprl(NHOR) = 0.
  73. real :: oprc(NHOR) = 0.
  74. real :: oevap(NHOR) = 0.
  75. real :: osst(NHOR) = 0.
  76. real :: opme(NHOR) = 0.
  77. real :: oheat(NHOR) = 0.
  78. real :: oroff(NHOR) = 0.
  79. real :: oofl(NHOR) = 0.
  80. real :: oicesnow(NHOR)= 0.
  81. real :: oice(NHOR) = 0.
  82. real :: osnow(NHOR) = 0.
  83. !
  84. ! global integer
  85. !
  86. integer :: ndatim(7) = -1
  87. integer :: nstep = 0 ! time step
  88. integer :: naccuout = 0 ! counter for accumulated output
  89. integer :: nrestart = 0 ! switch for restart
  90. integer :: n_run_years = 0
  91. integer :: n_run_months = 0
  92. integer :: n_run_days = -1
  93. integer :: n_start_step = 0
  94. integer :: ntspd = 32
  95. integer :: n_days_per_month = 30
  96. integer :: n_days_per_year = 360
  97. integer :: n_start_year = 1
  98. integer :: n_start_month = 1
  99. integer :: nafter = 32
  100. integer :: naomod = 320
  101. integer :: nlsg = 1
  102. integer :: ngui = 0
  103. !
  104. ! Parallel Stuff
  105. !
  106. integer :: mpinfo = 0
  107. integer :: mypid = 0
  108. integer :: myworld = 0
  109. integer :: nproc = NPRO
  110. !
  111. ! needed to compile surfmod
  112. !
  113. real :: so(NSPP),dls(NHOR)
  114. !
  115. end module pumamod
  116. !
  117. program climate
  118. use pumamod
  119. !
  120. call clini
  121. !
  122. mocd = n_run_months ! month countdown
  123. nscd = n_run_days * ntspd ! step countdown
  124. if (nrestart == 0) nstep = n_start_step ! timestep since 01-01-0001
  125. !
  126. call updatim(nstep) ! set date & time array ndatim
  127. !
  128. do while (mocd > 0 .or. nscd > 0) ! main loop
  129. call clstep
  130. if(mod(nstep,nafter) == 0) then
  131. call clout
  132. endif
  133. iyea = ndatim(1) ! current year
  134. imon = ndatim(2) ! current month
  135. nstep = nstep + 1
  136. call updatim(nstep) ! set date & time array ndatim
  137. if (imon /= ndatim(2)) then
  138. mocd = mocd - 1 ! next month
  139. write(nud,"('Completed month ',I2.2,'-',I4.4)") imon,iyea
  140. endif
  141. if (nscd > 0) nscd = nscd - 1
  142. if (nscd == 0) exit
  143. enddo ! month countdown
  144. !
  145. call clstop
  146. !
  147. stop
  148. end
  149. !
  150. ! ================
  151. ! SUBROUTINE CLINI
  152. ! ================
  153. !
  154. subroutine clini
  155. use pumamod
  156. !
  157. ! version identifier (date)
  158. !
  159. character(len=80) :: version = '25.07.2008 by Larry'
  160. !
  161. integer :: ih(8)
  162. !
  163. real :: zin(NLON*NLAT)
  164. real :: zcltaux(NLON*NLAT,0:13) = 0.
  165. real :: zcltauy(NLON*NLAT,0:13) = 0.
  166. real :: zclprl(NLON*NLAT,0:13) = 0.
  167. real :: zclprc(NLON*NLAT,0:13) = 0.
  168. real :: zclevap(NLON*NLAT,0:13) = 0.
  169. real :: zclroff(NLON*NLAT,0:13) = 0.
  170. real :: zcliced(NLON*NLAT,0:13) = 0.
  171. real :: zclsnow(NLON*NLAT,0:13) = 0.
  172. real :: zclsst(NLON*NLAT,0:13) = 0.
  173. real :: zclswr(NLON*NLAT,0:13) = 0.
  174. real :: zcllwr(NLON*NLAT,0:13) = 0.
  175. real :: zclshf(NLON*NLAT,0:13) = 0.
  176. real :: zcllhf(NLON*NLAT,0:13) = 0.
  177. real :: zlsm(NLON*NLAT) = 0.
  178. !
  179. logical :: lrestart = .false.
  180. !
  181. character(len=80) :: climatefile='climate.txt'
  182. character(len=80) :: surface='surface.txt'
  183. namelist/clpar/ntspd,nafter,n_run_days,n_run_months,n_run_years &
  184. & ,n_start_year,n_start_month,nlsg,climatefile &
  185. & ,surface
  186. !
  187. ! get process id
  188. !
  189. call get_mpi_info(nproc,mypid)
  190. !
  191. ! read and print namelist and distribute it
  192. !
  193. if (mypid == NROOT) then
  194. open(12,file='cl_namelist',form='formatted')
  195. read(12,clpar)
  196. write(nud,'(/," *******************************************")')
  197. write(nud,'(" * CLIMATE ",a30," *")') trim(version)
  198. write(nud,'(" *******************************************")')
  199. write(nud,'(" * Namelist CLPAR from <cl_namelist> *")')
  200. write(nud,'(" *******************************************")')
  201. write(nud,clpar)
  202. close(12)
  203. n_start_step = ntspd * (n_start_year * n_days_per_year &
  204. + (n_start_month-1) * n_days_per_month)
  205. if(n_run_years > 0) n_run_months=n_run_years*12
  206. endif
  207. !
  208. call mpbci(ntspd)
  209. call mpbci(nlsg)
  210. call mpbci(nafter)
  211. call mpbci(n_days_per_month)
  212. call mpbci(n_days_per_year)
  213. call mpbci(n_run_days)
  214. call mpbci(n_run_months)
  215. call mpbci(n_run_years)
  216. call mpbci(n_start_step)
  217. call mpbci(n_start_year)
  218. call mpbci(n_start_month)
  219. !
  220. call calini(n_days_per_month,n_days_per_year,n_start_step,ntspd &
  221. & ,solar_day,mypid)
  222. !
  223. if(mypid == NROOT) then
  224. call restart_ini(lrestart,'climate_restart')
  225. if (lrestart) then
  226. nrestart = 1
  227. endif
  228. endif
  229. call mpbci(nrestart)
  230. !
  231. if (nrestart == 0) then
  232. !
  233. ! new start (read start file)
  234. !
  235. if (mypid == NROOT) then
  236. open(20,FILE=trim(climatefile),FORM='formatted')
  237. do
  238. read(20,'(8I10)',IOSTAT=iostat) ih(:)
  239. if (iostat /= 0) exit
  240. read(20,'(8E12.6)',IOSTAT=iostat) zin(:)
  241. if (iostat /= 0) exit
  242. if (ih(1) == 169) then ! Code 169 for SST
  243. jmon=ih(3)
  244. if (jmon > 12) jmon = mod(jmon/100,100)
  245. write(nud,*) 'SST month ',jmon,' found in ',trim(climatefile)
  246. zclsst(:,jmon) = zin(:)
  247. elseif (ih(1) == 180) then ! Code 180 for taux
  248. jmon=ih(3)
  249. if (jmon > 12) jmon = mod(jmon/100,100)
  250. write(nud,*) 'TAUX month ',jmon,' found in ',trim(climatefile)
  251. zcltaux(:,jmon) = zin(:)
  252. elseif (ih(1) == 181) then ! Code 181 for tauy
  253. jmon=ih(3)
  254. if (jmon > 12) jmon = mod(jmon/100,100)
  255. write(nud,*) 'TAUY month ',jmon,' found in ',trim(climatefile)
  256. zcltauy(:,jmon) = zin(:)
  257. elseif (ih(1) == 182) then ! Code 182 for evap
  258. jmon=ih(3)
  259. if (jmon > 12) jmon = mod(jmon/100,100)
  260. write(nud,*) 'EVAP month ',jmon,' found in ',trim(climatefile)
  261. zclevap(:,jmon) = zin(:)
  262. elseif (ih(1) == 142) then ! Code 142 for prl
  263. jmon=ih(3)
  264. if (jmon > 12) jmon = mod(jmon/100,100)
  265. write(nud,*) 'PRL month ',jmon,' found in ',trim(climatefile)
  266. zclprl(:,jmon) = zin(:)
  267. elseif (ih(1) == 143) then ! Code 143 for prc
  268. jmon=ih(3)
  269. if (jmon > 12) jmon = mod(jmon/100,100)
  270. write(nud,*) 'PRC month ',jmon,' found in ',trim(climatefile)
  271. zclprc(:,jmon) = zin(:)
  272. elseif (ih(1) == 160) then ! Code 160 for roff
  273. jmon=ih(3)
  274. if (jmon > 12) jmon = mod(jmon/100,100)
  275. write(nud,*) 'ROFF month ',jmon,' found in ',trim(climatefile)
  276. zclroff(:,jmon) = zin(:)
  277. elseif (ih(1) == 211) then ! Code 211 for iced
  278. jmon=ih(3)
  279. if (jmon > 12) jmon = mod(jmon/100,100)
  280. write(nud,*) 'ICED month ',jmon,' found in ',trim(climatefile)
  281. zcliced(:,jmon) = zin(:)
  282. elseif (ih(1) == 141) then ! Code 141 for snow
  283. jmon=ih(3)
  284. if (jmon > 12) jmon = mod(jmon/100,100)
  285. write(nud,*) 'SNOW month ',jmon,' found in ',trim(climatefile)
  286. zclsnow(:,jmon) = zin(:)
  287. elseif (ih(1) == 176) then ! Code 176 for swr
  288. jmon=ih(3)
  289. if (jmon > 12) jmon = mod(jmon/100,100)
  290. write(nud,*) 'SWR month ',jmon,' found in ',trim(climatefile)
  291. zclswr(:,jmon) = zin(:)
  292. elseif (ih(1) == 177) then ! Code 177 for lwr
  293. jmon=ih(3)
  294. if (jmon > 12) jmon = mod(jmon/100,100)
  295. write(nud,*) 'LWR month ',jmon,' found in ',trim(climatefile)
  296. zcllwr(:,jmon) = zin(:)
  297. elseif (ih(1) == 146) then ! Code 146 for shf
  298. jmon=ih(3)
  299. if (jmon > 12) jmon = mod(jmon/100,100)
  300. write(nud,*) 'SHF month ',jmon,' found in ',trim(climatefile)
  301. zclshf(:,jmon) = zin(:)
  302. elseif (ih(1) == 147) then ! Code 147 for lhf
  303. jmon=ih(3)
  304. if (jmon > 12) jmon = mod(jmon/100,100)
  305. write(nud,*) 'LHF month ',jmon,' found in ',trim(climatefile)
  306. zcllhf(:,jmon) = zin(:)
  307. elseif (ih(1) == 172) then ! Code 172 for Land-Sea-Mask
  308. write(nud,*) 'land-sea-mask found in ',trim(climatefile)
  309. zlsm(:) = zin(:)
  310. endif
  311. enddo
  312. close(20)
  313. !
  314. ! read sst and ice from surface.txt (if exist) and replace
  315. !
  316. open(20,FILE=trim(surface),FORM='formatted')
  317. do
  318. read(20,'(8I10)',IOSTAT=iostat) ih(:)
  319. if (iostat /= 0) exit
  320. read(20,'(8E12.6)',IOSTAT=iostat) zin(:)
  321. if (iostat /= 0) exit
  322. if (ih(1) == 169) then ! Code 169 for SST
  323. jmon=ih(3)
  324. if (jmon > 12) jmon = mod(jmon/100,100)
  325. write(nud,*) 'SST month ',jmon,' found in ',trim(surface)
  326. zclsst(:,jmon) = zin(:)
  327. elseif (ih(1) == 211) then ! Code 211 for iced
  328. jmon=ih(3)
  329. if (jmon > 12) jmon = mod(jmon/100,100)
  330. write(nud,*) 'ICED month ',jmon,' found in ',trim(surface)
  331. zcliced(:,jmon) = zin(:)
  332. elseif (ih(1) == 172) then ! Code 172 for Land-Sea-Mask
  333. write(nud,*) 'land-sea-mask found in ',trim(surface)
  334. zlsm(:) = zin(:)
  335. endif
  336. enddo
  337. close(20)
  338. endif ! mypid
  339. zclsst(:,0)=zclsst(:,12)
  340. zclsst(:,13)=zclsst(:,1)
  341. zcltaux(:,0)=zcltaux(:,12)
  342. zcltaux(:,13)=zcltaux(:,1)
  343. zcltauy(:,0)=zcltauy(:,12)
  344. zcltauy(:,13)=zcltauy(:,1)
  345. zclprl(:,0)=zclprl(:,12)
  346. zclprl(:,13)=zclprl(:,1)
  347. zclprc(:,0)=zclprc(:,12)
  348. zclprc(:,13)=zclprc(:,1)
  349. zclevap(:,0)=zclevap(:,12)
  350. zclevap(:,13)=zclevap(:,1)
  351. zclroff(:,0)=zclroff(:,12)
  352. zclroff(:,13)=zclroff(:,1)
  353. zcliced(:,0)=zcliced(:,12)
  354. zcliced(:,13)=zcliced(:,1)
  355. zclsnow(:,0)=zclsnow(:,12)
  356. zclsnow(:,13)=zclsnow(:,1)
  357. zclswr(:,0)=zclswr(:,12)
  358. zclswr(:,13)=zclswr(:,1)
  359. zcllwr(:,0)=zcllwr(:,12)
  360. zcllwr(:,13)=zcllwr(:,1)
  361. zclshf(:,0)=zclshf(:,12)
  362. zclshf(:,13)=zclshf(:,1)
  363. zcllhf(:,0)=zcllhf(:,12)
  364. zcllhf(:,13)=zcllhf(:,1)
  365. !
  366. call mpscgp(zclsst,ccsst,14)
  367. call mpscgp(zcltaux,cctaux,14)
  368. call mpscgp(zcltauy,cctauy,14)
  369. call mpscgp(zclevap,ccevap,14)
  370. call mpscgp(zclprl,ccprl,14)
  371. call mpscgp(zclprc,ccprc,14)
  372. call mpscgp(zclroff,ccroff,14)
  373. call mpscgp(zcliced,cciced,14)
  374. call mpscgp(zclsnow,ccsnow,14)
  375. call mpscgp(zclswr,ccswr,14)
  376. call mpscgp(zcllwr,cclwr,14)
  377. call mpscgp(zclshf,ccshf,14)
  378. call mpscgp(zcllhf,cclhf,14)
  379. call mpscgp(zls,ccls,1)
  380. !
  381. ! make clsst >= tfreeze
  382. !
  383. do jm=0,13
  384. ccsst(:,jm)=AMAX1(ccsst(:,jm),TFREEZE)
  385. enddo
  386. !
  387. ! initialize
  388. !
  389. naccuout=0
  390. !
  391. call clget
  392. !
  393. yls(:)=real(nint(ccls(:)))
  394. ytaux(:)=cltaux(:)
  395. ytauy(:)=cltauy(:)
  396. yprl(:)=clprl(:)
  397. yprc(:)=clprc(:)
  398. yevap(:)=clevap(:)
  399. ysst(:)=clsst(:)
  400. ypme(:)=clprl(:)+clprc(:)+clevap(:)
  401. yroff(:)=clroff(:)
  402. yicesnow(:)=cliced(:)*CRHOI/1000.+clsnow(:)
  403. yheat(:)=clswr(:)+cllwr(:)+clshf(:)+cllhf(:)
  404. !
  405. else
  406. !
  407. ! restart from restart file
  408. !
  409. if (mypid == NROOT) then
  410. call get_restart_integer('nstep' ,nstep)
  411. call get_restart_integer('naccu' ,naccuout)
  412. endif
  413. !
  414. call mpgetgp('ccls',ccls,NHOR,1)
  415. call mpgetgp('ccsst',ccsst,NHOR,14)
  416. call mpgetgp('cctaux',cctaux,NHOR,14)
  417. call mpgetgp('cctauy',cctauy,NHOR,14)
  418. call mpgetgp('ccprl',ccprl,NHOR,14)
  419. call mpgetgp('ccprc',ccprc,NHOR,14)
  420. call mpgetgp('ccevap',ccevap,NHOR,14)
  421. call mpgetgp('ccroff',ccroff,NHOR,14)
  422. call mpgetgp('cciced',cciced,NHOR,14)
  423. call mpgetgp('ccsnow',ccsnow,NHOR,14)
  424. call mpgetgp('ccswr',ccswr,NHOR,14)
  425. call mpgetgp('cclwr',cclwr,NHOR,14)
  426. call mpgetgp('ccshf',ccshf,NHOR,14)
  427. call mpgetgp('cclhf',cclhf,NHOR,14)
  428. call mpgetgp('ols',ols,NHOR,1)
  429. call mpgetgp('otaux',otaux,NHOR,1)
  430. call mpgetgp('otauy',otauy,NHOR,1)
  431. call mpgetgp('oprl',oprl,NHOR,1)
  432. call mpgetgp('oprc',oprc,NHOR,1)
  433. call mpgetgp('oevap',oevap,NHOR,1)
  434. call mpgetgp('osst',osst,NHOR,1)
  435. call mpgetgp('opme',opme,NHOR,1)
  436. call mpgetgp('oicesnow',oicesnow,NHOR,1)
  437. call mpgetgp('oice',oice,NHOR,1)
  438. call mpgetgp('osnow',osnow,NHOR,1)
  439. call mpgetgp('oheat',oheat,NHOR,1)
  440. call mpgetgp('oofl',oofl,NHOR,1)
  441. call mpgetgp('oroff',oroff,NHOR,1)
  442. yls(:)=ccls(:)
  443. !
  444. endif
  445. !
  446. call mpbci(nstep)
  447. call mpbci(naccuout)
  448. !
  449. ! open output file
  450. !
  451. if (mypid == NROOT) then
  452. open(31,file='cl_output',form='unformatted')
  453. endif
  454. !
  455. if (nlsg > 0) then
  456. call mpgagp(zlsm,yls,1)
  457. call ntomin(nstep,ndatim(5),ndatim(4),ndatim(3),ndatim(2) &
  458. & ,ndatim(1))
  459. if(mypid == nroot) then
  460. !
  461. ! check if day_per_year are right (lsg coupling needs 360)
  462. !
  463. kdpy=n_days_per_year
  464. if(kdpy .ne. 360.) then
  465. write(nud,*) '!ERROR! for LSG coupling you need to set '
  466. write(nud,*) ' n_days_per_year in plasim namelist INP '
  467. write(nud,*) ' to 360 !'
  468. write(nud,*) ' at the moment, n_days_per_year= ',kdpy
  469. write(nud,*) 'Model stoped!'
  470. stop 'ERROR'
  471. endif
  472. call clsgini(ndatim,ntspd,naomod,nlsg,ngui,zlsm,nlon,nlat)
  473. endif
  474. endif
  475. !
  476. if (mypid == NROOT) then
  477. open(31,file='cl_output',form='unformatted')
  478. endif
  479. !
  480. return
  481. end subroutine clini
  482. ! ===================================================================
  483. ! SUBROUTINE clstep
  484. ! ===================================================================
  485. subroutine clstep
  486. use pumamod
  487. !
  488. real :: zsst(nlon*nlat)
  489. real :: ztaux(nlon*nlat)
  490. real :: ztauy(nlon*nlat)
  491. real :: zpme(nlon*nlat)
  492. real :: zroff(nlon*nlat)
  493. real :: zice(nlon*nlat)
  494. real :: zheat(nlon*nlat)
  495. real :: zfldo(nlon*nlat)
  496. real :: zmld(nlon*nlat)
  497. !
  498. ! dbug arrays
  499. !
  500. real,allocatable :: zprf1(:),zprf2(:),zprf3(:),zprf4(:),zprf5(:)
  501. !
  502. ! get climate
  503. !
  504. call clget
  505. !
  506. ! set climate (exept sst which is set later according to oceanmod)
  507. !
  508. yls(:)=ccls(:)
  509. ytaux(:)=cltaux(:)
  510. ytauy(:)=cltauy(:)
  511. yprl(:)=clprl(:)
  512. yprc(:)=clprc(:)
  513. yevap(:)=clevap(:)
  514. ypme(:)=clprl(:)+clprc(:)+clevap(:)
  515. yroff(:)=clroff(:)
  516. yicesnow(:)=cliced(:)*CRHOI/1000.+clsnow(:)
  517. yheat(:)=clswr(:)+cllwr(:)+clshf(:)+cllhf(:)
  518. !
  519. ols(:)=ccls(:)+ols(:)
  520. otaux(:)=cltaux(:)+otaux(:)
  521. otauy(:)=cltauy(:)+otauy(:)
  522. oprl(:)=clprl(:)+oprl(:)
  523. oprc(:)=clprc(:)+oprc(:)
  524. oevap(:)=clevap(:)+oevap(:)
  525. osst(:)=clsst(:)+osst(:)
  526. opme(:)=clprl(:)+clprc(:)+clevap(:)+opme(:)
  527. oroff(:)=clroff(:)+oroff(:)
  528. oicesnow(:)=cliced(:)*CRHOI/1000.+clsnow(:)+oicesnow(:)
  529. oice(:)=cliced(:)+oice(:)
  530. osnow(:)=clsnow(:)+osnow(:)
  531. oheat(:)=clswr(:)+cllwr(:)+clshf(:)+cllhf(:)+oheat(:)
  532. naccuout=naccuout+1
  533. !
  534. ! do the lsg coupling
  535. !
  536. if (nlsg > 0) then
  537. if(nlsg < 2) call mpgagp(zsst,ysst,1)
  538. call mpgagp(ztaux,ytaux,1)
  539. call mpgagp(ztauy,ytauy,1)
  540. call mpgagp(zpme,ypme,1)
  541. call mpgagp(zroff,yroff,1)
  542. call mpgagp(zice,yicesnow,1)
  543. call mpgagp(zheat,yheat,1)
  544. !! call mpgagp(zmld,ymld,1)
  545. call ntomin(nstep,ndatim(5),ndatim(4),ndatim(3),ndatim(2) &
  546. & ,ndatim(1))
  547. if(mypid == nroot) then
  548. call clsgstep(ndatim,nstep,zsst,ztaux,ztauy,zpme,zroff,zice &
  549. & ,zheat,zfldo)
  550. endif
  551. if(nlsg < 2) then
  552. if(mod(nstep,naomod) == naomod-1) call mpscgp(zfldo,yfldo,1)
  553. else
  554. yfldo(:)=0.
  555. endif
  556. oofl(:)=oofl(:)+yfldo(:)
  557. endif
  558. !
  559. ! if lsg coupling due to ahfl and osst, replace ssts
  560. !
  561. if(nlsg > 1) then
  562. if(mod(nstep,naomod) == naomod-1) call mpscgp(zsst,ysst,1)
  563. endif
  564. !
  565. ! set sst according to oceanmod
  566. !
  567. ysst(:)=clsst(:)
  568. !
  569. return
  570. end subroutine clstep
  571. ! ===================================
  572. ! SUBROUTINE CLSTOP
  573. ! ===================================
  574. subroutine clstop
  575. use pumamod
  576. !
  577. ! write restart file
  578. !
  579. if (mypid == NROOT) then
  580. call restart_prepare('climate_status')
  581. call put_restart_integer('nstep',nstep)
  582. call put_restart_integer('naccu',naccuout)
  583. endif
  584. !
  585. call mpputgp('ccls',ccls,NHOR,1)
  586. call mpputgp('ccsst',ccsst,NHOR,14)
  587. call mpputgp('cctaux',cctaux,NHOR,14)
  588. call mpputgp('cctauy',cctauy,NHOR,14)
  589. call mpputgp('ccprl',ccprl,NHOR,14)
  590. call mpputgp('ccprc',ccprc,NHOR,14)
  591. call mpputgp('ccevap',ccevap,NHOR,14)
  592. call mpputgp('ccroff',ccroff,NHOR,14)
  593. call mpputgp('cciced',cciced,NHOR,14)
  594. call mpputgp('ccsnow',ccsnow,NHOR,14)
  595. call mpputgp('ccswr',ccswr,NHOR,14)
  596. call mpputgp('cclwr',cclwr,NHOR,14)
  597. call mpputgp('ccshf',ccshf,NHOR,14)
  598. call mpputgp('cclhf',cclhf,NHOR,14)
  599. call mpputgp('ols',ols,NHOR,1)
  600. call mpputgp('otaux',otaux,NHOR,1)
  601. call mpputgp('otauy',otauy,NHOR,1)
  602. call mpputgp('oprl',oprl,NHOR,1)
  603. call mpputgp('oprc',oprc,NHOR,1)
  604. call mpputgp('oevap',oevap,NHOR,1)
  605. call mpputgp('osst',osst,NHOR,1)
  606. call mpputgp('opme',opme,NHOR,1)
  607. call mpputgp('oicesnow',oicesnow,NHOR,1)
  608. call mpputgp('oice',oice,NHOR,1)
  609. call mpputgp('osnow',osnow,NHOR,1)
  610. call mpputgp('oheat',oheat,NHOR,1)
  611. call mpputgp('oofl',oofl,NHOR,1)
  612. call mpputgp('oroff',oroff,NHOR,1)
  613. !
  614. ! close output file
  615. !
  616. if (mypid == NROOT) then
  617. close(31)
  618. endif
  619. !
  620. ! stop lsg coupling
  621. !
  622. if(nlsg > 0) then
  623. if(mypid == NROOT) then
  624. call clsgstop
  625. endif
  626. endif
  627. !
  628. return
  629. end subroutine clstop
  630. ! ===================================
  631. ! SUBROUTINE CLOUT
  632. ! ===================================
  633. subroutine clout
  634. use pumamod
  635. !
  636. integer ih(8)
  637. !
  638. zn=1./real(naccuout)
  639. ols(:)=ols(:)*zn
  640. osst(:)=osst(:)*zn
  641. otaux(:)=otaux(:)*zn
  642. otauy(:)=otauy(:)*zn
  643. oprl(:)=oprl(:)*zn
  644. oprc(:)=oprc(:)*zn
  645. oevap(:)=oevap(:)*zn
  646. oroff(:)=oroff(:)*zn
  647. oice(:)=oice(:)*zn
  648. oofl(:)=oofl(:)*zn
  649. oicesnow(:)=oicesnow(:)*zn
  650. opme(:)=opme(:)*zn
  651. osnow(:)=osnow(:)*zn
  652. oheat(:)=oheat(:)*zn
  653. !
  654. call writegp(31,oprl,142,1)
  655. call writegp(31,oprc,143,1)
  656. call writegp(31,oroff,160,1)
  657. call writegp(31,osst,169,1)
  658. call writegp(31,ols,172,1)
  659. call writegp(31,otaux,180,1)
  660. call writegp(31,otauy,181,1)
  661. call writegp(31,oevap,182,1)
  662. call writegp(31,oice,211,1)
  663. call writegp(31,oofl,222,1)
  664. call writegp(31,oicesnow,69,1)
  665. call writegp(31,oheat,263,1)
  666. call writegp(31,opme,260,1)
  667. call writegp(31,osnow,141,1)
  668. !
  669. naccuout=0
  670. ols(:)=0.
  671. osst(:)=0.
  672. otaux(:)=0.
  673. otauy(:)=0.
  674. oprl(:)=0.
  675. oprc(:)=0.
  676. oevap(:)=0.
  677. oroff(:)=0.
  678. oice(:)=0.
  679. oofl(:)=0.
  680. oicesnow(:)=0.
  681. opme(:)=0.
  682. osnow(:)=0.
  683. oheat(:)=0.
  684. !
  685. return
  686. end subroutine clout
  687. !====================================================================
  688. ! SUBROUTINE CLGET
  689. !====================================================================
  690. subroutine clget
  691. use pumamod
  692. ! *********************
  693. ! * get annual cycle *
  694. ! *********************
  695. ! modified for 14 months version of SST
  696. ! jm2: 0=Dec, 1-12:Months, 13=Jan
  697. jstep=nstep+1 !advance time step (ts = ts(t)+ dtsdt)
  698. call ntomin(jstep,nmin,nhour,nday,nmonth,nyear)
  699. !
  700. jm1 = nmonth
  701. if (nday > 15) then ! Does not work for other planets
  702. jm2=jm1+1
  703. else
  704. jm2=jm1-1
  705. endif
  706. ! ****************************************
  707. !* interpolate
  708. ! ****************************************
  709. zgw2 = abs((nday-1)*1440+nhour*60+nmin-21600)/43200.
  710. zgw1 = 1.0 - zgw2
  711. clsst(:) = zgw1*ccsst(:,jm1)+zgw2*ccsst(:,jm2)
  712. cltaux(:) = zgw1*cctaux(:,jm1)+zgw2*cctaux(:,jm2)
  713. cltauy(:) = zgw1*cctauy(:,jm1)+zgw2*cctauy(:,jm2)
  714. clprl(:) = zgw1*ccprl(:,jm1)+zgw2*ccprl(:,jm2)
  715. clprc(:) = zgw1*ccprc(:,jm1)+zgw2*ccprc(:,jm2)
  716. clevap(:) = zgw1*ccevap(:,jm1)+zgw2*ccevap(:,jm2)
  717. clroff(:) = zgw1*ccroff(:,jm1)+zgw2*ccroff(:,jm2)
  718. cliced(:) = zgw1*cciced(:,jm1)+zgw2*cciced(:,jm2)
  719. clsnow(:) = zgw1*ccsnow(:,jm1)+zgw2*ccsnow(:,jm2)
  720. clswr(:) = zgw1*ccswr(:,jm1)+zgw2*ccswr(:,jm2)
  721. cllwr(:) = zgw1*cclwr(:,jm1)+zgw2*cclwr(:,jm2)
  722. clshf(:) = zgw1*ccshf(:,jm1)+zgw2*ccshf(:,jm2)
  723. cllhf(:) = zgw1*cclhf(:,jm1)+zgw2*cclhf(:,jm2)
  724. clsnow(:)=max(clsnow(:),0.)
  725. cliced(:)=max(cliced(:),0.)
  726. clprl(:)=max(clprl(:),0.)
  727. clprc(:)=max(clprc(:),0.)
  728. clroff(:)=max(clroff(:),0.)
  729. clevap(:)=min(clevap(:),0.)
  730. where(cliced(:)==0.) clsnow=0.
  731. return
  732. end subroutine clget
  733. ! ==================
  734. ! SUBROUTINE UPDATIM
  735. ! ==================
  736. subroutine updatim(kstep)
  737. use pumamod
  738. if (n_days_per_year == 365) then
  739. call step2cal(kstep,ntspd,ndatim)
  740. else
  741. call step2cal30(kstep,ntspd,ndatim)
  742. endif
  743. return
  744. end
  745. subroutine writegp(kunit,pf,kcode,klev)
  746. use pumamod
  747. real :: pf(NHOR)
  748. real :: zf(NUGP)
  749. integer(kind=4) :: ihead(8)
  750. real(kind=4) :: zzf(NUGP)
  751. istep = nstep
  752. call ntomin(istep,nmin,nhour,nday,nmonth,nyear)
  753. ihead(1) = kcode
  754. ihead(2) = klev
  755. ihead(3) = nday + 100 * nmonth + 10000 * nyear
  756. ihead(4) = nmin + 100 * nhour
  757. ihead(5) = NLON
  758. ihead(6) = NLAT
  759. ihead(7) = 1
  760. ihead(8) = n_days_per_year
  761. call mpgagp(zf,pf,1)
  762. if (mypid == NROOT) then
  763. write(kunit) ihead
  764. zzf(:)=zf(:)
  765. write(kunit) zzf
  766. endif
  767. return
  768. end
  769. !
  770. ! dummies
  771. !
  772. subroutine get_surf_array(yn,pa,k1,k2,kread)
  773. character (len=*) :: yn ! array name
  774. real :: pa(k1,k2) ! array to receive values
  775. integer :: k1 ! 1. dim of array
  776. integer :: k2 ! 2. dim of array
  777. integer :: kread ! flag to indicate success
  778. return
  779. end
  780. subroutine get_3d_year(yn,pa,kdim,klev,kread)
  781. character (len=*) :: yn
  782. real :: pa(kdim,klev,0:13)
  783. return
  784. end
  785. subroutine get_surf_year(yn,pa,kdim,kread)
  786. character (len=*) :: yn
  787. integer :: ih(8)
  788. real :: pa(kdim,0:13)
  789. return
  790. end