tm5_prism.F90 62 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514151515161517151815191520152115221523152415251526152715281529153015311532153315341535153615371538153915401541154215431544154515461547154815491550155115521553155415551556155715581559156015611562156315641565156615671568156915701571157215731574157515761577157815791580158115821583158415851586158715881589159015911592159315941595159615971598159916001601160216031604160516061607160816091610161116121613161416151616161716181619162016211622162316241625162616271628162916301631163216331634163516361637163816391640164116421643164416451646164716481649165016511652165316541655165616571658165916601661166216631664166516661667166816691670167116721673167416751676167716781679168016811682168316841685168616871688168916901691169216931694169516961697169816991700170117021703170417051706170717081709171017111712171317141715171617171718171917201721172217231724172517261727172817291730173117321733173417351736173717381739174017411742174317441745174617471748174917501751175217531754175517561757175817591760176117621763176417651766176717681769177017711772177317741775177617771778177917801781178217831784178517861787178817891790179117921793179417951796179717981799180018011802180318041805180618071808180918101811181218131814181518161817181818191820182118221823182418251826182718281829183018311832183318341835183618371838183918401841184218431844184518461847184818491850185118521853185418551856185718581859186018611862186318641865186618671868186918701871
  1. !
  2. #define TRACEBACK write (gol,'("in ",a," (",a,", line",i5,")")') rname, __FILE__, __LINE__; call goErr
  3. !
  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. #define PRISM_ERR write (gol,'("from call to PRISM routine")'); call goErr
  8. #define IF_PRISM_NOTOK_RETURN(action) if (status/=OASIS_OK) then; PRISM_ERR; TRACEBACK; action; return; end if
  9. !
  10. #include "tm5.inc"
  11. !
  12. !-----------------------------------------------------------------------------
  13. ! TM5 !
  14. !-----------------------------------------------------------------------------
  15. !BOP
  16. !
  17. ! !MODULE: TM5_PRISM
  18. !
  19. ! !DESCRIPTION:
  20. !\\
  21. !\\
  22. ! !INTERFACE:
  23. !
  24. MODULE TM5_PRISM
  25. !
  26. ! !USES:
  27. !
  28. USE MOD_OASIS
  29. USE GO, ONLY : gol, goPr, goErr, TDate
  30. IMPLICIT NONE
  31. PRIVATE
  32. !
  33. ! !PUBLIC MEMBER FUNCTIONS:
  34. !
  35. public :: TM5_Prism_Init ! read control parameters from rc file
  36. public :: TM5_Prism_Init2 ! prism init: grids, partition, coupled variables
  37. public :: TM5_Prism_Done ! deallocate
  38. public :: SetPrismTime ! current time/date into prism format
  39. public :: InqCplVar ! check if a coupled variable exists
  40. public :: Init, Done, Setup, Remap ! methods for TshRemap object
  41. !
  42. ! !PUBLIC DATA MEMBERS:
  43. !
  44. character(len=6), public, parameter :: comp_name = 'ctm5mp'
  45. integer, public :: PRISM_start_date(6) ! prism reference start time
  46. integer, public :: comp_id ! tm5 id for prism
  47. integer, public :: exchange_period ! coupling exchange period in hours
  48. integer, public :: ifs_cpl_nlev ! number of levels for coupling with IFS
  49. integer, public :: ifs_cpl_nlev_cutoff ! reduced number of levels for coupling with IFS,
  50. ! applied for aerosols
  51. !
  52. logical, public :: refine_levels ! do we send/receive all IFS levels? (then we need to "refine them" here)
  53. logical, public :: coupled_to_lpj ! is TM5 coupled to LPJ-Guess?
  54. !
  55. ! !PRIVATE DATA MEMBERS:
  56. !
  57. character(len=*), parameter :: mname = 'TM5_Prism'
  58. integer, parameter :: len_grid_name = 4 ! oasis3: fixed length for gridname
  59. !
  60. integer :: ifs_npes ! ifs number of procs
  61. integer :: ifs_shT ! ifs spectral fields : resolution
  62. integer :: ifs_shn ! ifs spectral fields : number of coeff
  63. integer :: ifs_shn_recv ! ifs spectral fields : number of coeff received
  64. integer, parameter :: wp = SELECTED_REAL_KIND(12,307) ! working precision = double
  65. character(len=256) :: error_message
  66. !
  67. ! !PUBLIC TYPES:
  68. !
  69. TYPE, PUBLIC :: TshRemap
  70. !
  71. ! remapping: for each element in received grid, identify target indices:
  72. !
  73. ! remap(:,:,1) : 1 = real part, 2 = imag part
  74. ! remap(:,:,2) : index in triangle : 1=(0,0) 2=(0,1) 3=(0,2) ... np=(m,n)
  75. ! remap(:,:,3) : level
  76. !
  77. ! Info array has values : m*1000 + n + lev/100 and negative for imaginary part
  78. !
  79. ! receive spectral info every timestep to avoid memory growth
  80. type(TDate) :: t ! store time of current info
  81. integer, pointer :: remap(:,:,:) ! remapping indices
  82. integer :: shT ! truncation
  83. END TYPE TshRemap
  84. TYPE, PUBLIC :: TCplVar ! --- TM Coupled Variable ---
  85. character(len=128) :: name ! tm5 internal name
  86. character(len=128) :: cpl_name ! short name used in rcfile
  87. integer, pointer :: var_id(:) ! list of id return by oasis3 (one per level)
  88. character(len=128), pointer :: var_name(:)! list of names used in the namcouple (one per level)
  89. logical :: serial ! serial transfer
  90. character(len=3) :: intent ! in or out
  91. character(len=2) :: grid_type ! spectral or gridpoint
  92. integer :: region !
  93. integer :: nlev !
  94. real :: west_deg, dlon_deg, south_deg, dlat_deg ! lon/lat grids
  95. integer :: nlon, nlat
  96. integer :: shT, shn, shn_recv ! spectral info
  97. integer :: itr ! tracer index if any
  98. logical :: cache ! cache
  99. type(TDate) :: cache_tmid
  100. real, pointer :: cache_data(:,:,:) => null()
  101. END TYPE TCplVar
  102. ! --- var -----------------------------------
  103. integer, parameter :: maxcplvar = 165 ! max nb of coupled fields
  104. type(TCplVar), public, save :: CplVar(maxcplvar) ! array of coupled fields
  105. integer, public :: ncplvar ! actual nb of coupled fields
  106. integer :: region_glb, region_sfc
  107. character(len=1024) :: prism_get_list ! comma seperated lists of coupled fields
  108. character(len=1024) :: prism_put_list
  109. integer, allocatable :: part_id(:)
  110. integer, allocatable :: var_shape(:,:)
  111. integer :: sp_part_id
  112. integer, allocatable :: sp_var_shape(:)
  113. !
  114. ! !INTERFACE:
  115. !
  116. interface Init
  117. module procedure shremap_Init
  118. end interface
  119. interface Done
  120. module procedure shremap_Done
  121. end interface
  122. interface Setup
  123. module procedure shremap_Setup
  124. end interface
  125. interface Remap
  126. module procedure shremap_Remap
  127. end interface
  128. interface SetPrismTime
  129. module procedure SetPrismTime_date
  130. end interface
  131. !
  132. ! !REVISION HISTORY:
  133. ! 10 Sep 2013 - Ph. Le Sager - cleanup, document, remove oasis4 stuff (can
  134. ! always be retrieved with svn if needed)
  135. ! 8 Oct 2013 - Ph. Le Sager - dummy CO2 exchange with LPJ-Guess
  136. !
  137. ! !REMARKS:
  138. !
  139. !EOP
  140. !------------------------------------------------------------------------
  141. CONTAINS
  142. !--------------------------------------------------------------------------
  143. ! TM5 !
  144. !--------------------------------------------------------------------------
  145. !BOP
  146. !
  147. ! !IROUTINE: TM5_PRISM_INIT
  148. !
  149. ! !DESCRIPTION: read coupling information from rc file.
  150. !\\
  151. !\\
  152. ! !INTERFACE:
  153. !
  154. SUBROUTINE TM5_Prism_Init( rcfile, status )
  155. !
  156. ! !USES:
  157. !
  158. use GO, only : TrcFile, Init, Done, ReadRc
  159. !
  160. ! !INPUT PARAMETERS:
  161. !
  162. character(len=*), intent(in) :: rcfile
  163. !
  164. ! !OUTPUT PARAMETERS:
  165. !
  166. integer, intent(out) :: status
  167. !
  168. ! !REMARKS:
  169. !
  170. !EOP
  171. !------------------------------------------------------------------------
  172. !BOC
  173. character(len=*), parameter :: rname = mname//'/TM5_Prism_Init'
  174. type(TrcFile) :: rcF
  175. ! --- begin -----------------------------------
  176. ! * settings from rcfile
  177. call Init( rcF, rcfile, status )
  178. IF_NOTOK_RETURN(status=1)
  179. call ReadRc( rcF, 'ifs.cpl.nlev', ifs_cpl_nlev, status )
  180. IF_NOTOK_RETURN(status=1)
  181. call ReadRc( rcf, 'ifs.cpl.nlev.cutoff', ifs_cpl_nlev_cutoff, status )
  182. IF_NOTOK_RETURN(status=1)
  183. call ReadRc( rcF, 'ifs.nproc', ifs_npes, status )
  184. IF_NOTOK_RETURN(status=1)
  185. call ReadRc( rcF, 'ifs.shresol', ifs_shT, status )
  186. IF_NOTOK_RETURN(status=1)
  187. call ReadRc( rcF, 'cpl.exchange.period', exchange_period, status )
  188. IF_NOTOK_RETURN(status=1)
  189. call ReadRc( rcF, 'prism.get', prism_get_list, status )
  190. IF_NOTOK_RETURN(status=1)
  191. call ReadRc( rcF, 'prism.put', prism_put_list, status )
  192. IF_NOTOK_RETURN(status=1)
  193. call ReadRc( rcF, 'coupled_to_lpjguess', coupled_to_lpj, status, default=.false. )
  194. IF_ERROR_RETURN(status=1)
  195. select case (ifs_cpl_nlev)
  196. case (91,62)
  197. refine_levels=.true.
  198. case(34,31)
  199. refine_levels=.false.
  200. case default
  201. write(gol,*) " Wrong (sub)set of levels is exchanged between IFS and TM5 " ; call goErr
  202. write(gol,*) " Either (34 out of) 91, or (31 out of) 62" ; call goErr
  203. status=1
  204. TRACEBACK; return
  205. end select
  206. call Done( rcF, status )
  207. IF_NOTOK_RETURN(status=1)
  208. ! * spectral grids
  209. ifs_shn = (ifs_shT+1)*(ifs_shT+2)/2 ! number of coeff
  210. !PRIOR ! because of partioning at ifs, more spectral coeff are received;
  211. !PRIOR ! assume that the maximum number of coeff on an ifs process is
  212. !PRIOR ! 2 more than an equal distribution:
  213. !PRIOR ifs_shn_recv = ( ceiling(real(ifs_shn)/real(ifs_npes)) + 2 ) * ifs_npes
  214. ifs_shn_recv = ifs_shn !PLS: with MCT, the exact number of coeff is needed
  215. ! check that compilation was perform as expected with respect to optical feedback
  216. #ifndef with_ecearth_optics
  217. if (index(prism_put_list,'AOD_01') /=0) then
  218. write(gol,*) "Feedback of aerosols optical properties requires 'with_ecearth_optics' cpp"; call goErr
  219. write(gol,*) "You must recompile TM5-MP with cpp defs 'with_ecearth_optics'"; call goErr
  220. write(gol,*) "This can be done in the config-build.xml file with the TM5_MDEFS_FFLAGS key."; call goErr
  221. status=1
  222. TRACEBACK; return
  223. endif
  224. #endif
  225. status = 0
  226. END SUBROUTINE TM5_PRISM_INIT
  227. !EOC
  228. !--------------------------------------------------------------------------
  229. ! TM5 !
  230. !--------------------------------------------------------------------------
  231. !BOP
  232. !
  233. ! !IROUTINE: TM5_PRISM_INIT2
  234. !
  235. ! !DESCRIPTION: prism grid writing, prism partition, init coupled variables
  236. !\\
  237. !\\
  238. ! !INTERFACE:
  239. !
  240. subroutine TM5_Prism_Init2( nreg_all, nreg, ireg_sfc, lli, levi, status )
  241. !
  242. ! !USES:
  243. !
  244. use Grid, only : TllGridInfo, TLevelInfo
  245. use Grid, only : TshGridInfo, Init, Done
  246. #ifdef parallel_cplng
  247. use tm5_distgrid, only : dgrid, Get_DistGrid
  248. #endif
  249. use GO, only : NewDate, goReadFromLine
  250. use dims, only : lm
  251. use chem_param, only : names, io3, ich4, ico2
  252. use partools, only : isRoot
  253. #ifdef with_m7
  254. use chem_param, only : inus_n, iso4nus, iais_n, iso4ais, ibcais, ipomais
  255. use chem_param, only : iacs_n, iso4acs, ibcacs, ipomacs, issacs, iduacs
  256. use chem_param, only : icos_n, iso4cos, ibccos, ipomcos, isscos, iducos
  257. use chem_param, only : iaii_n, ibcaii, ipomaii, iaci_n, iduaci, icoi_n, iducoi
  258. use chem_param, only : ino3_a, imsa
  259. #endif
  260. !
  261. ! !INPUT PARAMETERS:
  262. !
  263. integer, intent(in) :: nreg_all
  264. integer, intent(in) :: nreg
  265. integer, intent(in) :: ireg_sfc
  266. type(TllGridInfo), intent(in) :: lli(1:nreg_all)
  267. type(TLevelInfo), intent(in) :: levi
  268. !
  269. ! !OUTPUT PARAMETERS:
  270. !
  271. integer, intent(out) :: status
  272. !
  273. ! !REVISION HISTORY:
  274. ! 8 Oct 2013 - Ph. Le Sager - coupling with LPJG
  275. !
  276. ! !REMARKS:
  277. !
  278. !EOP
  279. !------------------------------------------------------------------------
  280. !BOC
  281. character(len=*), parameter :: rname = mname//'/TM5_Prism_Init2'
  282. integer :: ireg, region, id, ivar, ilev, i1, j1, i2, j2
  283. integer, parameter :: nc = 4 ! corners
  284. integer :: nx, ny
  285. integer :: i, j, k, m, n, z
  286. character(len=len_grid_name) :: point_name
  287. real(ip_realwp_p), allocatable :: lons(:,:), lats(:,:)
  288. real(ip_realwp_p), allocatable :: clons(:,:,:), clats(:,:,:)
  289. real(ip_realwp_p), allocatable :: area(:,:)
  290. integer, allocatable :: mask(:,:)
  291. character(len=len_grid_name) :: sp_point_name
  292. real(ip_realwp_p), allocatable :: sp_lons(:,:), sp_lats(:,:)
  293. integer, allocatable :: sp_mask(:,:)
  294. real(ip_realwp_p) :: sp_dlon, sp_dlat
  295. #ifdef parallel_cplng
  296. integer :: part_val(1:5)
  297. #else
  298. integer :: part_val(1:3)
  299. #endif
  300. integer :: sp_part_val(1:3)
  301. character(len=128) :: cpl_name
  302. integer :: var_nodims(2)
  303. integer :: var_type
  304. integer(kind=ip_intwp_p) :: var_intent
  305. logical :: write_grid
  306. #ifdef parallel_cplng
  307. type(TllGridInfo) :: local_lli ! local Lat-Lon grid info
  308. #endif
  309. ! --- begin -----------------------------------
  310. write (gol,'("initialize prism coupling:")'); call goPr
  311. write (gol,'(" component : ",a)') trim(comp_name); call goPr
  312. ! store in module variables
  313. region_glb = 1
  314. region_sfc = ireg_sfc
  315. ! storage for variable shape:
  316. allocate( part_id(nreg_all) )
  317. allocate( var_shape(4,nreg_all) )
  318. allocate( sp_var_shape(4) )
  319. ! init to zero on all pe's
  320. part_id = 0
  321. var_shape = 0
  322. sp_part_id = 0
  323. sp_var_shape = 0
  324. ! ============ oasis3 grid writing =================
  325. write_grid=.false. !! HARDCODED !!
  326. ! Define the grids by master proc
  327. if ( isroot .and. write_grid ) then
  328. call oasis_start_grids_writing( status )
  329. ! **** lon/lat grid
  330. do region = 1, nreg_all
  331. ! name of grid points
  332. if ( region == region_glb ) then
  333. ! global region
  334. point_name = 'CTM3'
  335. else if ( region == region_sfc ) then
  336. ! global surface grid:
  337. point_name = 'CTM1'
  338. else
  339. ! global grids only ...
  340. cycle
  341. end if
  342. write (gol,'(" define points ",a," ...")') point_name; call goPr
  343. ! grid size:
  344. nx = lli(region)%nlon
  345. ny = lli(region)%nlat
  346. write (gol,'(" region : ",i6)') region; call goPr
  347. write (gol,'(" size : ",2i6)') nx, ny; call goPr
  348. allocate( lons(nx,ny) )
  349. allocate( lats(nx,ny) )
  350. allocate( clons(nx,ny,nc) )
  351. allocate( clats(nx,ny,nc) )
  352. allocate( area(nx,ny) )
  353. allocate( mask(nx,ny) )
  354. ! set lon/lat grid (grid cell centers):
  355. do i = 1, nx
  356. lons(i,:) = lli(region)%lon_deg(i)
  357. end do
  358. do j = 1, ny
  359. lats(:,j) = lli(region)%lat_deg(j)
  360. end do
  361. call oasis_write_grid( point_name, nx, ny, lons, lats )
  362. ! set corner lon/lat:
  363. ! 3 o o 2
  364. ! 4 o o 1
  365. do i = 1, nx
  366. clons(i,:,1) = lli(region)%blon_deg(i )
  367. clons(i,:,2) = lli(region)%blon_deg(i )
  368. clons(i,:,3) = lli(region)%blon_deg(i-1)
  369. clons(i,:,4) = lli(region)%blon_deg(i-1)
  370. end do
  371. do j = 1, ny
  372. clats(:,j,1) = lli(region)%blat_deg(j-1)
  373. clats(:,j,2) = lli(region)%blat_deg(j )
  374. clats(:,j,3) = lli(region)%blat_deg(j )
  375. clats(:,j,4) = lli(region)%blat_deg(j-1)
  376. end do
  377. write (gol,'(" write corners ...")'); call goPr
  378. call oasis_write_corner( point_name, nx, ny, nc, clons, clats )
  379. ! land/sea mask
  380. mask = 0 ! not masked; gives warnings about 'sea-world' cplout ...
  381. write (gol,'(" write mask ...")'); call goPr
  382. call oasis_write_mask( point_name, nx, ny, mask )
  383. do j = 1, ny
  384. area(:,j) = lli(region)%area_m2(j)
  385. end do
  386. write (gol,'(" write area ...")'); call goPr
  387. call oasis_write_area( point_name, nx, ny, area )
  388. deallocate( lons )
  389. deallocate( lats )
  390. deallocate( clons )
  391. deallocate( clats )
  392. deallocate( area )
  393. deallocate( mask )
  394. end do ! regions
  395. ! *** SPECTRAL GRID
  396. write(sp_point_name, '("C",i3.3)') ifs_shT
  397. write (gol,'(" define points ",a," ...")') trim(sp_point_name); call goPr
  398. allocate( sp_lons(2*ifs_shn,1) )
  399. allocate( sp_lats(2*ifs_shn,1) )
  400. allocate( sp_mask(2*ifs_shn,1) )
  401. ! Triangular storage:
  402. !
  403. ! NSMAX * * .. *
  404. ! :
  405. !
  406. ! 1 * *
  407. ! n 0 *
  408. !
  409. ! 0 1 .. NSMAX
  410. ! m "wavenumber"
  411. !
  412. ! dummy locations : (n*2+z+0.5)*dlon, (m+0.5)*dlat
  413. ! where z=0 is real part and z=1 is imaginary part
  414. ! dummy spacing:
  415. sp_dlon = 0.1 ! degree
  416. sp_dlat = 0.1 ! degree
  417. ! index in coeff array:
  418. k = 0
  419. ! loop over global wavenumbers:
  420. do m = 0, ifs_shT
  421. ! loop over complex coeff:
  422. do n = m, ifs_shT
  423. ! loop over real/complex
  424. do z = 0, 1
  425. ! next coeff:
  426. k = k + 1
  427. ! cell centers:
  428. sp_lons(k,1) = -180.0 + (n*2+z+0.5) * sp_dlon
  429. sp_lats(k,1) = -90.0 + (m +0.5) * sp_dlat
  430. ! no mask:
  431. sp_mask(k,1) = 0 ! not masked
  432. end do ! re,im
  433. end do ! n
  434. end do ! m
  435. call oasis_write_grid( sp_point_name, 2*ifs_shn, 1, sp_lons, sp_lats )
  436. call oasis_write_mask( sp_point_name, 2*ifs_shn, 1, sp_mask )
  437. deallocate( sp_lons )
  438. deallocate( sp_lats )
  439. deallocate( sp_mask )
  440. call oasis_terminate_grids_writing()
  441. else
  442. write (gol,'(" not necessary to write grids (not root) ...")'); call goPr
  443. end if ! root and write_grid
  444. ! ============ Partition =================
  445. write (gol,'(" define partitions ...")'); call goPr
  446. ! *** LAT/LON
  447. reg: do region = 1, nreg_all
  448. if ( (region /= region_glb) .and. (region /= region_sfc) ) cycle ! global grids only
  449. nx = lli(region)%nlon
  450. ny = lli(region)%nlat
  451. #ifdef parallel_cplng
  452. call Get_DistGrid( dgrid(region), I_STRT=i1, I_STOP=i2 , J_STRT=j1, J_STOP=j2 )
  453. ! store shape:
  454. var_shape(1:2,region) = (/1,i2-i1+1/)
  455. var_shape(3:4,region) = (/1,j2-j1+1/)
  456. part_val(1) = 2
  457. part_val(2) = i1-1+(j1-1)*nx
  458. part_val(3) = i2-i1+1
  459. part_val(4) = j2-j1+1
  460. part_val(5) = nx
  461. #else
  462. ! store shape:
  463. var_shape(1:2,region) = (/1,nx/)
  464. var_shape(3:4,region) = (/1,ny/)
  465. part_val(1) = 0 ! serial partition
  466. part_val(2) = 0
  467. part_val(3) = 0
  468. if ( isroot ) part_val(3) = nx*ny ! total grid size
  469. #endif
  470. status = OASIS_OK ! <-- status was not set by PRISM_Def_Partition_Proto (is it still true in OASIS3-MCT?)
  471. call oasis_def_partition( part_id(region), part_val, status )
  472. if (status/=OASIS_OK) then
  473. write (error_message,'("from OASIS_DEF_PARTITION : ",i6)') status
  474. TRACEBACK
  475. call oasis_abort( comp_id, rname, error_message )
  476. end if
  477. #ifdef parallel_cplng
  478. write (gol,'(" partition : ",i6," ; ",5i6)') part_id(region), part_val; call goPr
  479. #else
  480. write (gol,'(" partition : ",i6," ; ",3i6)') part_id(region), part_val; call goPr
  481. #endif
  482. end do reg
  483. ! *** SPECTRAL
  484. sp_part_val(1) = 0
  485. sp_part_val(2) = 0
  486. sp_part_val(3) = 0
  487. if ( isroot ) sp_part_val(3) = ifs_shn*2 ! total grid size
  488. !status = OASIS_OK ! <-- status was not set by PRISM_Def_Partition_Proto (is it still true in OASIS3-MCT?)
  489. call oasis_def_partition( sp_part_id, sp_part_val, status )
  490. if (status/=OASIS_OK) then
  491. write (error_message,'("from OASIS_DEF_PARTITION : ",i6)') status
  492. TRACEBACK
  493. call oasis_abort( comp_id, rname, error_message )
  494. end if
  495. write (gol,'(" partition : ",i6," ; ",3i6)') sp_part_id, sp_part_val; call goPr
  496. ! -------------------------------------------------------------------------
  497. ! * CONFIGURE COUPLING FIELDS
  498. ! -------------------------------------------------------------------------
  499. !
  500. ! (0) DEFAULT
  501. !
  502. write (gol,'(" init cplvar list ...")'); call goPr
  503. var_nodims(1) = 2 ! rank of coupling field
  504. var_nodims(2) = 1 ! number of bundles in coupling field (always 1)
  505. var_type = OASIS_Real
  506. do ivar = 1, size(CplVar)
  507. CplVar(ivar)%cpl_name = ''
  508. nullify( CplVar(ivar)%var_id )
  509. nullify( CplVar(ivar)%var_name )
  510. CplVar(ivar)%grid_type = 'll' ! lon/lat
  511. CplVar(ivar)%region = -1
  512. CplVar(ivar)%nlev = -1
  513. #ifdef parallel_cplng
  514. CplVar(ivar)%serial = .false.
  515. #else
  516. CplVar(ivar)%serial = .true.
  517. #endif
  518. CplVar(ivar)%intent = 'xxx'
  519. CplVar(ivar)%cache = .false. ! at init, flag means "will use cache"
  520. CplVar(ivar)%cache_tmid = NewDate(year=9999)
  521. nullify( CplVar(ivar)%cache_data )
  522. end do
  523. ! Above init should be same as:
  524. ! CplVar = TCplVar('','',null(),null(),.true.,'xxx','ll',-1,-1,0.,0.,0.,0.,0,0,0,0,0,.false.,NewDate(year=9999),null())
  525. ! We could also just set pointers to => null() in declaration l.117
  526. ncplvar = 0 ! no fields defined yet
  527. !
  528. ! (1) METEO VARIABLES
  529. !
  530. write (gol,'(" set meteo cplvars ...")'); call goPr
  531. write (gol,'(" list : ",a)') trim(prism_get_list); call goPr
  532. ivar = 0
  533. GET: DO
  534. if ( len_trim(prism_get_list) == 0 ) exit ! leave if empty
  535. ! extract next name from list
  536. call goReadFromLine( prism_get_list, cpl_name, status, sep=',' )
  537. IF_NOTOK_RETURN(status=1)
  538. write (gol,'(" extracted ",a," ...")') trim(cpl_name); call goPr
  539. ivar = ivar + 1
  540. if ( ivar > maxcplvar ) then
  541. write (gol,'("ivar exceeds maxcplvar ",i6)') maxcplvar; call goErr
  542. TRACEBACK; status=1; return
  543. end if
  544. CplVar(ivar)%cpl_name = trim(cpl_name)
  545. CplVar(ivar)%intent = 'in'
  546. select case ( trim(cpl_name) )
  547. case ( 'SPINF' )
  548. CplVar(ivar)%name = 'spinf2d'
  549. CplVar(ivar)%grid_type = 'sh'
  550. CplVar(ivar)%serial = .true.
  551. CplVar(ivar)%nlev = 1
  552. case ( 'LNSP' )
  553. CplVar(ivar)%name = 'LNSP'
  554. CplVar(ivar)%grid_type = 'sh' ! spherical harmonics
  555. CplVar(ivar)%serial = .true.
  556. CplVar(ivar)%nlev = 1
  557. CplVar(ivar)%cache = .true.
  558. case ( 'VOR' )
  559. CplVar(ivar)%name = 'VO'
  560. CplVar(ivar)%grid_type = 'sh' ! spherical harmonics
  561. CplVar(ivar)%serial = .true.
  562. CplVar(ivar)%nlev = ifs_cpl_nlev
  563. CplVar(ivar)%cache = .true.
  564. case ( 'DIV' )
  565. CplVar(ivar)%name = 'D'
  566. CplVar(ivar)%grid_type = 'sh' ! spherical harmonics
  567. CplVar(ivar)%serial = .true.
  568. CplVar(ivar)%nlev = ifs_cpl_nlev
  569. CplVar(ivar)%cache = .true.
  570. case ( 'SPRES' )
  571. CplVar(ivar)%name = 'sp'
  572. CplVar(ivar)%region = region_glb
  573. CplVar(ivar)%nlev = 1
  574. CplVar(ivar)%cache = .true.
  575. case ( 'TMP' )
  576. CplVar(ivar)%name = 'T'
  577. CplVar(ivar)%region = region_glb
  578. CplVar(ivar)%nlev = ifs_cpl_nlev
  579. case ( 'HUM' )
  580. CplVar(ivar)%name = 'Q'
  581. CplVar(ivar)%region = region_glb
  582. CplVar(ivar)%nlev = ifs_cpl_nlev
  583. case ( 'OROG' )
  584. CplVar(ivar)%name = 'oro'
  585. CplVar(ivar)%region = region_sfc
  586. CplVar(ivar)%nlev = 1
  587. case ( 'CLW' )
  588. CplVar(ivar)%name = 'CLWC'
  589. CplVar(ivar)%region = region_glb
  590. CplVar(ivar)%nlev = ifs_cpl_nlev
  591. case ( 'CIW')
  592. CplVar(ivar)%name = 'CIWC'
  593. CplVar(ivar)%region = region_glb
  594. CplVar(ivar)%nlev = ifs_cpl_nlev
  595. case ( 'CC' )
  596. CplVar(ivar)%name = 'CC'
  597. CplVar(ivar)%region = region_glb
  598. CplVar(ivar)%nlev = ifs_cpl_nlev
  599. case ( 'CCO' )
  600. CplVar(ivar)%name = 'CCO'
  601. CplVar(ivar)%region = region_glb
  602. CplVar(ivar)%nlev = ifs_cpl_nlev
  603. case ( 'CCU' )
  604. CplVar(ivar)%name = 'CCU'
  605. CplVar(ivar)%region = region_glb
  606. CplVar(ivar)%nlev = ifs_cpl_nlev
  607. case ( 'UMF' )
  608. CplVar(ivar)%name = 'UDMF'
  609. CplVar(ivar)%region = region_glb
  610. CplVar(ivar)%nlev = ifs_cpl_nlev
  611. case ( 'DMF' )
  612. CplVar(ivar)%name = 'DDMF'
  613. CplVar(ivar)%region = region_glb
  614. CplVar(ivar)%nlev = ifs_cpl_nlev
  615. case ( 'UDR' )
  616. CplVar(ivar)%name = 'UDDR'
  617. CplVar(ivar)%region = region_glb
  618. CplVar(ivar)%nlev = ifs_cpl_nlev
  619. case ( 'DDR' )
  620. CplVar(ivar)%name = 'DDDR'
  621. CplVar(ivar)%region = region_glb
  622. CplVar(ivar)%nlev = ifs_cpl_nlev
  623. case ( 'LSMSK' )
  624. CplVar(ivar)%name = 'lsm'
  625. CplVar(ivar)%region = region_sfc
  626. CplVar(ivar)%nlev = 1
  627. case ( 'ALB' )
  628. CplVar(ivar)%name = 'albedo'
  629. CplVar(ivar)%region = region_sfc
  630. CplVar(ivar)%nlev = 1
  631. case ( 'SR' )
  632. CplVar(ivar)%name = 'sr'
  633. CplVar(ivar)%region = region_sfc
  634. CplVar(ivar)%nlev = 1
  635. case ( 'CI' )
  636. CplVar(ivar)%name = 'ci'
  637. CplVar(ivar)%region = region_sfc
  638. CplVar(ivar)%nlev = 1
  639. case ( 'SST__' )
  640. CplVar(ivar)%name = 'sst'
  641. CplVar(ivar)%region = region_sfc
  642. CplVar(ivar)%nlev = 1
  643. case ( 'U10M' )
  644. CplVar(ivar)%name = 'u10m'
  645. CplVar(ivar)%region = region_sfc
  646. CplVar(ivar)%nlev = 1
  647. case ( 'V10M' )
  648. CplVar(ivar)%name = 'v10m'
  649. CplVar(ivar)%region = region_sfc
  650. CplVar(ivar)%nlev = 1
  651. case ( 'SRC' )
  652. CplVar(ivar)%name = 'src'
  653. CplVar(ivar)%region = region_sfc
  654. CplVar(ivar)%nlev = 1
  655. case ( 'D2M' )
  656. CplVar(ivar)%name = 'd2m'
  657. CplVar(ivar)%region = region_sfc
  658. CplVar(ivar)%nlev = 1
  659. case ( 'T2M' )
  660. CplVar(ivar)%name = 't2m'
  661. CplVar(ivar)%region = region_sfc
  662. CplVar(ivar)%nlev = 1
  663. case ( 'SSHF' )
  664. CplVar(ivar)%name = 'sshf'
  665. CplVar(ivar)%region = region_sfc
  666. CplVar(ivar)%nlev = 1
  667. case ( 'SLHF' )
  668. CplVar(ivar)%name = 'slhf'
  669. CplVar(ivar)%region = region_sfc
  670. CplVar(ivar)%nlev = 1
  671. case ( 'EWSS' )
  672. CplVar(ivar)%name = 'ewss'
  673. CplVar(ivar)%region = region_sfc
  674. CplVar(ivar)%nlev = 1
  675. case ( 'NSSS' )
  676. CplVar(ivar)%name = 'nsss'
  677. CplVar(ivar)%region = region_sfc
  678. CplVar(ivar)%nlev = 1
  679. case ( 'CP' )
  680. CplVar(ivar)%name = 'cp'
  681. CplVar(ivar)%region = region_sfc
  682. CplVar(ivar)%nlev = 1
  683. case ( 'LSP' )
  684. CplVar(ivar)%name = 'lsp'
  685. CplVar(ivar)%region = region_sfc
  686. CplVar(ivar)%nlev = 1
  687. case ( 'SSR' )
  688. CplVar(ivar)%name = 'ssr'
  689. CplVar(ivar)%region = region_sfc
  690. CplVar(ivar)%nlev = 1
  691. case ( 'SKT__')
  692. CplVar(ivar)%name = 'skt'
  693. CplVar(ivar)%region = region_sfc
  694. CplVar(ivar)%nlev = 1
  695. case ( 'SF___' )
  696. CplVar(ivar)%name = 'sf'
  697. CplVar(ivar)%region = region_sfc
  698. CplVar(ivar)%nlev = 1
  699. case ( 'SD' )
  700. CplVar(ivar)%name = 'sd'
  701. CplVar(ivar)%region = region_sfc
  702. CplVar(ivar)%nlev = 1
  703. case ( 'SWVL1' )
  704. CplVar(ivar)%name = 'swvl1'
  705. CplVar(ivar)%region = region_sfc
  706. CplVar(ivar)%nlev = 1
  707. case ( 'TV01' )
  708. CplVar(ivar)%name = 'tv01'
  709. CplVar(ivar)%region = region_sfc
  710. CplVar(ivar)%nlev = 1
  711. case ( 'TV02' )
  712. CplVar(ivar)%name = 'tv02'
  713. CplVar(ivar)%region = region_sfc
  714. CplVar(ivar)%nlev = 1
  715. case ( 'TV03' )
  716. CplVar(ivar)%name = 'tv03'
  717. CplVar(ivar)%region = region_sfc
  718. CplVar(ivar)%nlev = 1
  719. case ( 'TV04' )
  720. CplVar(ivar)%name = 'tv04'
  721. CplVar(ivar)%region = region_sfc
  722. CplVar(ivar)%nlev = 1
  723. case ( 'TV05' )
  724. CplVar(ivar)%name = 'tv05'
  725. CplVar(ivar)%region = region_sfc
  726. CplVar(ivar)%nlev = 1
  727. case ( 'TV06' )
  728. CplVar(ivar)%name = 'tv06'
  729. CplVar(ivar)%region = region_sfc
  730. CplVar(ivar)%nlev = 1
  731. case ( 'TV07' )
  732. CplVar(ivar)%name = 'tv07'
  733. CplVar(ivar)%region = region_sfc
  734. CplVar(ivar)%nlev = 1
  735. case ( 'TV09' )
  736. CplVar(ivar)%name = 'tv09'
  737. CplVar(ivar)%region = region_sfc
  738. CplVar(ivar)%nlev = 1
  739. case ( 'TV10' )
  740. CplVar(ivar)%name = 'tv10'
  741. CplVar(ivar)%region = region_sfc
  742. CplVar(ivar)%nlev = 1
  743. case ( 'TV11' )
  744. CplVar(ivar)%name = 'tv11'
  745. CplVar(ivar)%region = region_sfc
  746. CplVar(ivar)%nlev = 1
  747. case ( 'TV13' )
  748. CplVar(ivar)%name = 'tv13'
  749. CplVar(ivar)%region = region_sfc
  750. CplVar(ivar)%nlev = 1
  751. case ( 'TV16' )
  752. CplVar(ivar)%name = 'tv16'
  753. CplVar(ivar)%region = region_sfc
  754. CplVar(ivar)%nlev = 1
  755. case ( 'TV17' )
  756. CplVar(ivar)%name = 'tv17'
  757. CplVar(ivar)%region = region_sfc
  758. CplVar(ivar)%nlev = 1
  759. case ( 'TV18' )
  760. CplVar(ivar)%name = 'tv18'
  761. CplVar(ivar)%region = region_sfc
  762. CplVar(ivar)%nlev = 1
  763. case ( 'TV19' )
  764. CplVar(ivar)%name = 'tv19'
  765. CplVar(ivar)%region = region_sfc
  766. CplVar(ivar)%nlev = 1
  767. case ( 'CVL' )
  768. CplVar(ivar)%name = 'cvl'
  769. CplVar(ivar)%region = region_sfc
  770. CplVar(ivar)%nlev = 1
  771. case ( 'CVH' )
  772. CplVar(ivar)%name = 'cvh'
  773. CplVar(ivar)%region = region_sfc
  774. CplVar(ivar)%nlev = 1
  775. case default
  776. write (gol,'("unsupported cpl_name : ",a)') trim(cpl_name); call goErr
  777. TRACEBACK; status=1; return
  778. end select
  779. ! check
  780. if ( CplVar(ivar)%nlev < 1 ) then
  781. write (gol,'("found nlev ",i6," in cplvar ",i4," (",a,")")') CplVar(ivar)%nlev, ivar, trim(cpl_name); call goErr
  782. TRACEBACK; status=1; return
  783. end if
  784. ! storage per level
  785. allocate( CplVar(ivar)%var_id (CplVar(ivar)%nlev) )
  786. allocate( CplVar(ivar)%var_name(CplVar(ivar)%nlev) )
  787. ! name used in namcouple
  788. if ( CplVar(ivar)%nlev == 1 ) then
  789. ilev = 1
  790. write (CplVar(ivar)%var_name(ilev),'("C_",a)') trim(CplVar(ivar)%cpl_name)
  791. write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev, trim(CplVar(ivar)%var_name(ilev)); call goPr
  792. else
  793. do ilev = 1, CplVar(ivar)%nlev
  794. write (CplVar(ivar)%var_name(ilev),'("C_",a,".L",i3.3)') trim(CplVar(ivar)%cpl_name), ilev
  795. write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev, trim(CplVar(ivar)%var_name(ilev)); call goPr
  796. end do
  797. end if
  798. ! store latest number:
  799. ncplvar = ivar
  800. END DO GET ! list coupled var to get
  801. !
  802. ! (2) VARIABLES SENT TO OASIS
  803. !
  804. write (gol,'(" set concentration/optics cplvars ...")'); call goPr
  805. write (gol,'(" list : ",a)') trim(prism_put_list); call goPr
  806. ivar = ncplvar
  807. PUT: DO
  808. if ( len_trim(prism_put_list) == 0 ) exit
  809. ! extract next name from list
  810. call goReadFromLine( prism_put_list, cpl_name, status, sep=',' )
  811. IF_NOTOK_RETURN(status=1)
  812. write (gol,'(" extracted ",a," ...")') trim(cpl_name); call goPr
  813. ivar = ivar + 1
  814. if ( ivar > maxcplvar ) then
  815. write (gol,'("ivar exceeds maxcplvar ",i6)') maxcplvar; call goErr
  816. TRACEBACK; status=1; return
  817. end if
  818. CplVar(ivar)%cpl_name = trim(cpl_name)
  819. #ifdef parallel_cplng
  820. CplVar(ivar)%serial = .false.
  821. #else
  822. CplVar(ivar)%serial = .true.
  823. #endif
  824. CplVar(ivar)%intent = 'out'
  825. CplVar(ivar)%region = region_glb
  826. select case ( trim(cpl_name) )
  827. case ( 'O3', 'CH4' )
  828. ! send whole atmosphere for ozone and methane
  829. if (.not.refine_levels) then ! in both cases this should be ifs_cpl_nlev!
  830. CplVar(ivar)%nlev = lm(region_glb)
  831. else
  832. CplVar(ivar)%nlev = ifs_cpl_nlev
  833. endif
  834. case default
  835. ! for aerosols, do not send all levels in stratosphere
  836. ! this works, for refine_levels true or false
  837. CplVar(ivar)%nlev = ifs_cpl_nlev_cutoff
  838. end select
  839. select case ( trim(cpl_name) )
  840. case ( 'O3' )
  841. CplVar(ivar)%itr = io3
  842. case ( 'CH4' )
  843. CplVar(ivar)%itr = ich4
  844. #ifdef with_m7
  845. case ( 'N1' )
  846. CplVar(ivar)%itr = inus_n
  847. case ( 'SU1' )
  848. CplVar(ivar)%itr = iso4nus
  849. case ( 'N2' )
  850. CplVar(ivar)%itr = iais_n
  851. case ( 'SU2' )
  852. CplVar(ivar)%itr = iso4ais
  853. case ( 'BC2' )
  854. CplVar(ivar)%itr = ibcais
  855. case ( 'OC2' )
  856. CplVar(ivar)%itr = ipomais
  857. case ( 'N3' )
  858. CplVar(ivar)%itr = iacs_n
  859. case ( 'SU3' )
  860. CplVar(ivar)%itr = iso4acs
  861. case ( 'BC3' )
  862. CplVar(ivar)%itr = ibcacs
  863. case ( 'OC3' )
  864. CplVar(ivar)%itr = ipomacs
  865. case ( 'SS3' )
  866. CplVar(ivar)%itr = issacs
  867. case ( 'DU3' )
  868. CplVar(ivar)%itr = iduacs
  869. case ( 'N4' )
  870. CplVar(ivar)%itr = icos_n
  871. case ( 'SU4' )
  872. CplVar(ivar)%itr = iso4cos
  873. case ( 'BC4' )
  874. CplVar(ivar)%itr = ibccos
  875. case ( 'OC4' )
  876. CplVar(ivar)%itr = ipomcos
  877. case ( 'SS4' )
  878. CplVar(ivar)%itr = isscos
  879. case ( 'DU4' )
  880. CplVar(ivar)%itr = iducos
  881. case ( 'N5' )
  882. CplVar(ivar)%itr = iaii_n
  883. case ( 'BC5' )
  884. CplVar(ivar)%itr = ibcaii
  885. case ( 'OC5' )
  886. CplVar(ivar)%itr = ipomaii
  887. case ( 'N6' )
  888. CplVar(ivar)%itr = iaci_n
  889. case ( 'DU6' )
  890. CplVar(ivar)%itr = iduaci
  891. case ( 'N7' )
  892. CplVar(ivar)%itr = icoi_n
  893. case ( 'DU7' )
  894. CplVar(ivar)%itr = iducoi
  895. case ( 'NO3' )
  896. CplVar(ivar)%itr = ino3_a
  897. case ( 'MSA' )
  898. CplVar(ivar)%itr = imsa
  899. #endif
  900. #ifdef with_ecearth_optics
  901. case ( 'AOD_01', 'AOD_02', 'AOD_03', 'AOD_04', 'AOD_05', 'AOD_06', 'AOD_07', &
  902. 'AOD_08', 'AOD_09', 'AOD_10', 'AOD_11', 'AOD_12', 'AOD_13', 'AOD_14', &
  903. 'SSA_01', 'SSA_02', 'SSA_03', 'SSA_04', 'SSA_05', 'SSA_06', 'SSA_07', &
  904. 'SSA_08', 'SSA_09', 'SSA_10', 'SSA_11', 'SSA_12', 'SSA_13', 'SSA_14', &
  905. 'ASF_01', 'ASF_02', 'ASF_03', 'ASF_04', 'ASF_05', 'ASF_06', 'ASF_07', &
  906. 'ASF_08', 'ASF_09', 'ASF_10', 'ASF_11', 'ASF_12', 'ASF_13', 'ASF_14' )
  907. #endif
  908. case default
  909. write (gol,'("unsupported cpl_name : ",a)') trim(cpl_name); call goErr
  910. TRACEBACK; status=1; return
  911. end select
  912. ! set name:
  913. select case ( trim(cpl_name) )
  914. #ifdef with_ecearth_optics
  915. case ( 'AOD_01', 'AOD_02', 'AOD_03', 'AOD_04', 'AOD_05', 'AOD_06', 'AOD_07', &
  916. 'AOD_08', 'AOD_09', 'AOD_10', 'AOD_11', 'AOD_12', 'AOD_13', 'AOD_14', &
  917. 'SSA_01', 'SSA_02', 'SSA_03', 'SSA_04', 'SSA_05', 'SSA_06', 'SSA_07', &
  918. 'SSA_08', 'SSA_09', 'SSA_10', 'SSA_11', 'SSA_12', 'SSA_13', 'SSA_14', &
  919. 'ASF_01', 'ASF_02', 'ASF_03', 'ASF_04', 'ASF_05', 'ASF_06', 'ASF_07', &
  920. 'ASF_08', 'ASF_09', 'ASF_10', 'ASF_11', 'ASF_12', 'ASF_13', 'ASF_14' )
  921. write(CplVar(ivar)%name,'(a)') trim(cpl_name)
  922. #endif
  923. case default
  924. CplVar(ivar)%name = names(CplVar(ivar)%itr)
  925. end select
  926. ! check ..
  927. if ( CplVar(ivar)%nlev < 1 ) then
  928. write (gol,'("found nlev ",i6," in cplvar ",i4," (",a,")")') CplVar(ivar)%nlev, ivar, trim(cpl_name); call goErr
  929. TRACEBACK; status=1; return
  930. end if
  931. ! storage per level
  932. allocate( CplVar(ivar)%var_id (CplVar(ivar)%nlev) )
  933. allocate( CplVar(ivar)%var_name(CplVar(ivar)%nlev) )
  934. ! name used in namcouple
  935. if ( CplVar(ivar)%nlev == 1 ) then
  936. ilev = 1
  937. write (CplVar(ivar)%var_name(ilev),'(a,"TM5")') trim(CplVar(ivar)%cpl_name)
  938. write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev, trim(CplVar(ivar)%var_name(ilev)); call goPr
  939. else
  940. do ilev = 1, CplVar(ivar)%nlev
  941. write (CplVar(ivar)%var_name(ilev),'("C_",a,".L",i3.3)') trim(CplVar(ivar)%cpl_name), ilev
  942. end do
  943. write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev-1, trim(CplVar(ivar)%var_name(ilev-1)); call goPr
  944. end if
  945. ncplvar = ivar
  946. END DO PUT ! list with coupled names of var to send
  947. !
  948. ! (3) COUPLING WITH LPJG
  949. !
  950. if (coupled_to_lpj) then
  951. ! Sent to LPJ-Guess
  952. ivar = ncplvar + 1
  953. if ( ivar > maxcplvar ) then
  954. write (gol,'("ivar exceeds maxcplvar ",i6)') maxcplvar; call goErr
  955. TRACEBACK; status=1; return
  956. end if
  957. CplVar(ivar)%cpl_name = 'co2'
  958. CplVar(ivar)%itr = ico2
  959. CplVar(ivar)%name = names(CplVar(ivar)%itr)
  960. #ifdef parallel_cplng
  961. CplVar(ivar)%serial = .false.
  962. #else
  963. CplVar(ivar)%serial = .true.
  964. #endif
  965. CplVar(ivar)%intent = 'out'
  966. CplVar(ivar)%region = region_glb
  967. CplVar(ivar)%nlev = 1
  968. allocate( CplVar(ivar)%var_id (1) )
  969. allocate( CplVar(ivar)%var_name(1) )
  970. CplVar(ivar)%var_name(1) = "CO2_TM5"
  971. write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev, trim(CplVar(ivar)%var_name(1)); call goPr
  972. ! Received from LPJ-Guess
  973. ivar = ivar + 1
  974. if ( ivar > maxcplvar ) then
  975. write (gol,'("ivar exceeds maxcplvar ",i6)') maxcplvar; call goErr
  976. TRACEBACK; status=1; return
  977. end if
  978. CplVar(ivar)%cpl_name = 'c_flux'
  979. CplVar(ivar)%name = 'c_flux'
  980. #ifdef parallel_cplng
  981. CplVar(ivar)%serial = .false.
  982. #else
  983. CplVar(ivar)%serial = .true.
  984. #endif
  985. CplVar(ivar)%intent = 'in'
  986. CplVar(ivar)%region = region_glb
  987. CplVar(ivar)%nlev = 1
  988. CplVar(ivar)%itr = 999
  989. allocate( CplVar(ivar)%var_id (1) )
  990. allocate( CplVar(ivar)%var_name(1) )
  991. CplVar(ivar)%var_name(1) = "TM5_CFLX"
  992. write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev, trim(CplVar(ivar)%var_name(1)); call goPr
  993. ivar = ivar + 1
  994. if ( ivar > maxcplvar ) then
  995. write (gol,'("ivar exceeds maxcplvar ",i6)') maxcplvar; call goErr
  996. TRACEBACK; status=1; return
  997. end if
  998. CplVar(ivar)%cpl_name = 'c_npp'
  999. CplVar(ivar)%name = 'c_npp'
  1000. #ifdef parallel_cplng
  1001. CplVar(ivar)%serial = .false.
  1002. #else
  1003. CplVar(ivar)%serial = .true.
  1004. #endif
  1005. CplVar(ivar)%intent = 'in'
  1006. CplVar(ivar)%region = region_glb
  1007. CplVar(ivar)%nlev = 1
  1008. CplVar(ivar)%itr = 999
  1009. allocate( CplVar(ivar)%var_id (1) )
  1010. allocate( CplVar(ivar)%var_name(1) )
  1011. CplVar(ivar)%var_name(1) = "TM5_CNPP"
  1012. write (gol,'(" cplvar ",2i4," ",a)') ivar, ilev, trim(CplVar(ivar)%var_name(1)); call goPr
  1013. ncplvar = ivar
  1014. end if
  1015. !
  1016. ! (4) DEFINE OASIS VARIABLES
  1017. !
  1018. write (gol,'(" define oasis variables ...")'); call goPr
  1019. do ivar = 1, ncplvar
  1020. ireg = CplVar(ivar)%region
  1021. select case ( CplVar(ivar)%intent )
  1022. case ( 'in' )
  1023. var_intent = OASIS_In
  1024. case ( 'out' )
  1025. var_intent = OASIS_Out
  1026. case default
  1027. write (gol,'("unsupported intent : ",a)') trim(CplVar(ivar)%intent); call goErr
  1028. TRACEBACK; status=1; return
  1029. end select
  1030. do ilev = 1, CplVar(ivar)%nlev
  1031. select case ( CplVar(ivar)%grid_type )
  1032. case ( 'sh' )
  1033. !DBG write (gol,'(" ",i4," (",i3,") spectral variable ",a," ...")') ivar, ilev, trim(CplVar(ivar)%var_name(ilev)); call goPr
  1034. !DBG write (gol,'(" region : ",i6)' ) ireg ; call goPr
  1035. !DBG write (gol,'(" partition : ",i6)' ) sp_part_id ; call goPr
  1036. !DBG write (gol,'(" nodims : ",2i6)' ) var_nodims ; call goPr
  1037. !DBG write (gol,'(" shape : ",4i6)' ) sp_var_shape(1:4) ; call goPr
  1038. call oasis_def_var( &
  1039. CplVar(ivar)%var_id(ilev), &
  1040. trim(CplVar(ivar)%var_name(ilev)), &
  1041. sp_part_id, var_nodims, var_intent, &
  1042. sp_var_shape(1:4), var_type, status )
  1043. IF_PRISM_NOTOK_RETURN(status=1)
  1044. case ( 'll' )
  1045. !DBG write (gol,'(" ",i4," (",i3,") gridded variable ",a," ...")') ivar, ilev, trim(CplVar(ivar)%var_name(ilev)); call goPr
  1046. !DBG write (gol,'(" region : ",i6)' ) ireg ; call goPr
  1047. !DBG write (gol,'(" partition : ",i6)' ) part_id(ireg) ; call goPr
  1048. !DBG write (gol,'(" nodims : ",2i6)' ) var_nodims ; call goPr
  1049. !DBG write (gol,'(" shape : ",4i6)' ) var_shape(1:4,ireg) ; call goPr
  1050. call oasis_def_var( &
  1051. CplVar(ivar)%var_id(ilev), &
  1052. trim(CplVar(ivar)%var_name(ilev)), &
  1053. part_id(ireg), var_nodims, var_intent, &
  1054. var_shape(1:4,ireg), var_type, status )
  1055. IF_PRISM_NOTOK_RETURN(status=1)
  1056. case default
  1057. write (gol,'("unsupported grid_type : ",a)') trim(CplVar(ivar)%grid_type); call goErr
  1058. TRACEBACK; status=1; return
  1059. end select
  1060. end do ! levels
  1061. end do ! var
  1062. !
  1063. ! (4) STORE GRID INFO
  1064. !
  1065. write (gol,'("add grid info to cplvars ...")'); call goPr
  1066. do ivar = 1, ncplvar
  1067. do ilev = 1, CplVar(ivar)%nlev
  1068. select case ( CplVar(ivar)%grid_type )
  1069. case ( 'sh' )
  1070. CplVar(ivar)%shT = ifs_shT
  1071. CplVar(ivar)%shn = ifs_shn
  1072. CplVar(ivar)%shn_recv = ifs_shn_recv
  1073. if ( CplVar(ivar)%cache ) allocate( CplVar(ivar)%cache_data(2,CplVar(ivar)%shn,CplVar(ivar)%nlev) )
  1074. case ( 'll' )
  1075. #ifdef parallel_cplng
  1076. call Get_DistGrid( dgrid(CplVar(ivar)%region), lli=local_lli)
  1077. CplVar(ivar)%west_deg = local_lli%lon_deg(1)
  1078. CplVar(ivar)%south_deg = local_lli%lat_deg(1)
  1079. CplVar(ivar)%dlon_deg = local_lli%dlon_deg
  1080. CplVar(ivar)%dlat_deg = local_lli%dlat_deg
  1081. CplVar(ivar)%nlon = local_lli%nlon
  1082. CplVar(ivar)%nlat = local_lli%nlat
  1083. if ( CplVar(ivar)%cache ) &
  1084. allocate( CplVar(ivar)%cache_data(CplVar(ivar)%nlon, CplVar(ivar)%nlat, CplVar(ivar)%nlev) )
  1085. #else
  1086. CplVar(ivar)%west_deg = lli(CplVar(ivar)%region)%lon_deg(1)
  1087. CplVar(ivar)%south_deg = lli(CplVar(ivar)%region)%lat_deg(1)
  1088. CplVar(ivar)%dlon_deg = lli(CplVar(ivar)%region)%dlon_deg
  1089. CplVar(ivar)%dlat_deg = lli(CplVar(ivar)%region)%dlat_deg
  1090. CplVar(ivar)%nlon = lli(CplVar(ivar)%region)%nlon
  1091. CplVar(ivar)%nlat = lli(CplVar(ivar)%region)%nlat
  1092. if ( CplVar(ivar)%cache ) &
  1093. allocate( CplVar(ivar)%cache_data(CplVar(ivar)%nlon,CplVar(ivar)%nlat,CplVar(ivar)%nlev) )
  1094. #endif
  1095. case default
  1096. write (gol,'("unsupported grid_type : ",a)') trim(CplVar(ivar)%grid_type); call goErr
  1097. TRACEBACK; status=1; return
  1098. end select
  1099. end do ! levels
  1100. end do ! var
  1101. !
  1102. ! (5) FINALISE
  1103. !
  1104. call oasis_enddef( status )
  1105. if (status/=OASIS_OK) then
  1106. write (error_message,'("from OASIS_ENDDEF : ",i6)') status
  1107. TRACEBACK
  1108. call oasis_abort( comp_id, rname, error_message )
  1109. end if
  1110. if (isRoot) then
  1111. write (gol,'(" ")' ) ; call goPr
  1112. write (gol,'("initialized oasis coupling:")' ) ; call goPr
  1113. write (gol,'(" component : ",a)' ) trim(comp_name) ; call goPr
  1114. write (gol,'(" real kind : ",i4)' ) wp ; call goPr
  1115. write (gol,'(" ")' ) ; call goPr
  1116. end if
  1117. status = 0
  1118. END SUBROUTINE TM5_PRISM_INIT2
  1119. !EOC
  1120. !--------------------------------------------------------------------------
  1121. ! TM5 !
  1122. !--------------------------------------------------------------------------
  1123. !BOP
  1124. !
  1125. ! !IROUTINE: TM5_PRISM_DONE
  1126. !
  1127. ! !DESCRIPTION: cleanup (ie deallocate).
  1128. !\\
  1129. !\\
  1130. ! !INTERFACE:
  1131. !
  1132. SUBROUTINE TM5_Prism_Done( status )
  1133. !
  1134. ! !OUTPUT PARAMETERS:
  1135. !
  1136. integer, intent(out) :: status
  1137. !
  1138. ! !REVISION HISTORY:
  1139. ! 10 Sep 2013 - Ph. Le Sager -
  1140. !
  1141. ! !REMARKS:
  1142. !
  1143. !EOP
  1144. !------------------------------------------------------------------------
  1145. !BOC
  1146. character(len=*), parameter :: rname = mname//'/TM5_Prism_Done'
  1147. integer :: ireg
  1148. integer :: ivar
  1149. ! --- begin -----------------------------------
  1150. deallocate( part_id )
  1151. deallocate( var_shape )
  1152. deallocate( sp_var_shape )
  1153. ! clear descriptions:
  1154. do ivar = 1, ncplvar
  1155. deallocate( CplVar(ivar)%var_id )
  1156. deallocate( CplVar(ivar)%var_name )
  1157. if ( associated(CplVar(ivar)%cache_data) ) deallocate( CplVar(ivar)%cache_data )
  1158. end do
  1159. status = 0
  1160. END SUBROUTINE TM5_PRISM_DONE
  1161. !EOC
  1162. !--------------------------------------------------------------------------
  1163. ! TM5 !
  1164. !--------------------------------------------------------------------------
  1165. !BOP
  1166. !
  1167. ! !IROUTINE: InqCplVar
  1168. !
  1169. ! !DESCRIPTION: Inquire if there is a coupled variable with 'name'.
  1170. !\\
  1171. !\\
  1172. ! !INTERFACE:
  1173. !
  1174. SUBROUTINE InqCplVar( name, status, ivar, var_id, var_name, nlev )
  1175. !
  1176. ! !INPUT PARAMETERS:
  1177. !
  1178. character(len=*), intent(in) :: name
  1179. !
  1180. ! !OUTPUT PARAMETERS:
  1181. !
  1182. integer, intent(out) :: status
  1183. integer, intent(out), optional :: ivar
  1184. integer, intent(out), optional :: var_id(:)
  1185. character(len=*), intent(out), optional :: var_name(:)
  1186. integer, intent(out), optional :: nlev
  1187. !
  1188. ! !REVISION HISTORY:
  1189. ! 10 Sep 2013 - Ph. Le Sager -
  1190. !
  1191. ! !REMARKS:
  1192. !
  1193. !EOP
  1194. !------------------------------------------------------------------------
  1195. !BOC
  1196. character(len=*), parameter :: rname = mname//'/InqCplVar'
  1197. integer :: i, iv
  1198. ! --- begin -----------------------------------
  1199. ! loop over defined variables:
  1200. iv = -1
  1201. do i = 1, ncplvar
  1202. ! check name:
  1203. if ( trim(name) == trim(CplVar(i)%name) ) then
  1204. iv = i
  1205. exit
  1206. end if
  1207. end do
  1208. ! not found ?
  1209. if ( iv < 0 ) then
  1210. write (gol,'("name of cplvar not found : ",a)') trim(name) ; call goErr
  1211. write (gol,'(" available names : ")' ) ; call goErr
  1212. do i = 1, ncplvar
  1213. write (gol,'(" ",i4," ",a)') i, trim(CplVar(i)%name) ; call goErr
  1214. end do
  1215. end if
  1216. ! fill requested arguments:
  1217. if ( present(ivar ) ) ivar = iv
  1218. if ( present(var_id ) ) var_id = CplVar(iv)%var_id
  1219. if ( present(var_name) ) var_name = CplVar(iv)%var_name
  1220. if ( present(nlev ) ) nlev = CplVar(iv)%nlev
  1221. ! ok
  1222. status = 0
  1223. END SUBROUTINE InqCplVar
  1224. !EOC
  1225. ! **************************************************************************
  1226. ! ***
  1227. ! *** spectral field remapping routines
  1228. ! ***
  1229. ! **************************************************************************
  1230. !--------------------------------------------------------------------------
  1231. ! TM5 !
  1232. !--------------------------------------------------------------------------
  1233. !BOP
  1234. !
  1235. ! !IROUTINE: SHREMAP_INIT
  1236. !
  1237. ! !DESCRIPTION: Init TshRemap object
  1238. !\\
  1239. !\\
  1240. ! !INTERFACE:
  1241. !
  1242. SUBROUTINE SHREMAP_INIT( shR, status )
  1243. !
  1244. ! !USES:
  1245. !
  1246. use GO, only : NewDate
  1247. !
  1248. ! !OUTPUT PARAMETERS:
  1249. !
  1250. type(TshRemap), intent(out) :: shR
  1251. integer, intent(out) :: status
  1252. !
  1253. ! !REMARKS:
  1254. !
  1255. !EOP
  1256. !------------------------------------------------------------------------
  1257. !BOC
  1258. character(len=*), parameter :: rname = mname//'/shremap_Init'
  1259. ! --- begin ---------------------------------------
  1260. ! no time stored yet:
  1261. shR%t = NewDate(9999,9,9)
  1262. ! safety:
  1263. nullify( shR%remap )
  1264. ! nu truncation determined yet:
  1265. shR%shT = 0
  1266. status = 0
  1267. END SUBROUTINE SHREMAP_INIT
  1268. !EOC
  1269. !--------------------------------------------------------------------------
  1270. ! TM5 !
  1271. !--------------------------------------------------------------------------
  1272. !BOP
  1273. !
  1274. ! !IROUTINE: SHREMAP_DONE
  1275. !
  1276. ! !DESCRIPTION: deallocate var
  1277. !\\
  1278. !\\
  1279. ! !INTERFACE:
  1280. !
  1281. SUBROUTINE SHREMAP_DONE( shR, status )
  1282. !
  1283. ! !INPUT/OUTPUT PARAMETERS:
  1284. !
  1285. type(TshRemap), intent(inout) :: shR
  1286. !
  1287. ! !OUTPUT PARAMETERS:
  1288. !
  1289. integer, intent(out) :: status
  1290. !
  1291. ! !REMARKS:
  1292. !
  1293. !EOP
  1294. !------------------------------------------------------------------------
  1295. !BOC
  1296. character(len=*), parameter :: rname = mname//'/shremap_Done'
  1297. ! --- begin ---------------------------------------
  1298. if ( associated(shR%remap) ) deallocate( shR%remap )
  1299. status = 0
  1300. END SUBROUTINE SHREMAP_DONE
  1301. !EOC
  1302. !--------------------------------------------------------------------------
  1303. ! TM5 !
  1304. !--------------------------------------------------------------------------
  1305. !BOP
  1306. !
  1307. ! !IROUTINE: SHREMAP_SETUP
  1308. !
  1309. ! !DESCRIPTION:
  1310. !\\
  1311. !\\
  1312. ! !INTERFACE:
  1313. !
  1314. SUBROUTINE SHREMAP_SETUP( shR, spinf, spinf_nan, status )
  1315. !
  1316. ! !INPUT/OUTPUT PARAMETERS:
  1317. !
  1318. type(TshRemap), intent(inout) :: shR
  1319. !
  1320. ! !INPUT PARAMETERS:
  1321. !
  1322. real, intent(in) :: spinf(:,:,:) ! spectral info field
  1323. real, intent(in) :: spinf_nan ! not-a-number in spinf
  1324. !
  1325. ! !OUTPUT PARAMETERS:
  1326. !
  1327. integer, intent(out) :: status
  1328. !
  1329. ! !REMARKS:
  1330. !
  1331. !EOP
  1332. !------------------------------------------------------------------------
  1333. !BOC
  1334. character(len=*), parameter :: rname = mname//'/shremap_Setup'
  1335. integer :: nlev
  1336. integer :: sh_tripos(0:ifs_shT,0:ifs_shT)
  1337. integer :: vri, vm, vn, vp, vk
  1338. logical, allocatable :: sh_field(:,:,:)
  1339. integer :: i, j, k
  1340. real :: val
  1341. integer :: nzero, nerr
  1342. ! --- begin ---------------------------------------
  1343. ! number of levels:
  1344. nlev = size(spinf,2)
  1345. ! triangle position:
  1346. sh_tripos = 0
  1347. vp = 0
  1348. do vm = 0, ifs_shT
  1349. do vn = vm, ifs_shT
  1350. vp = vp + 1
  1351. sh_tripos(vm,vn) = vp
  1352. end do
  1353. end do
  1354. ! storage for mapping indices:
  1355. if ( .not. associated(shR%remap) ) then
  1356. allocate( shR%remap(ifs_shn_recv*2,nlev,3) )
  1357. end if
  1358. shR%remap = -1
  1359. ! flags for target values; not filled remains negative:
  1360. allocate( sh_field(2,ifs_shn,nlev) )
  1361. sh_field = .false.
  1362. ! loop over levels:
  1363. do k = 1, nlev
  1364. ! no zero's detected yet ...
  1365. nzero = 0
  1366. ! loop over spectral elements in layer:
  1367. do i = 1, ifs_shn_recv*2
  1368. !if (k==1) then
  1369. !write (gol,'(" k, i, spinf : ",2i6,f16.4)') k, i, spinf(i,k,1); call goPr
  1370. !endif
  1371. ! not a number ? then this is an extra element due to the partitioning
  1372. if ( spinf(i,k,1) == spinf_nan ) cycle
  1373. ! extract m, n, and level:
  1374. !
  1375. ! OLD : mmmnnnkk.0 real part
  1376. ! mmmnnnkk.5 imag part
  1377. !
  1378. !vri = 1 ! real part
  1379. !if ( spinf(i,k,1)-floor(spinf(i,k,1)) == 0.5 ) vri = 2 ! imaginary part
  1380. !vk = modulo( floor(spinf(i,k,1)), 100 ) ! level
  1381. !vn = modulo( floor((spinf(i,k,1)-vk)/100.0), 1000 ) ! n
  1382. !vm = floor(spinf(i,k,1)/100000.0) ! m
  1383. !
  1384. ! NEW : mmmnnn.kk real part
  1385. ! -mmmnnn.kk imag part
  1386. !
  1387. ! Note that real and imag for m=0,n=0 are both 000000.00 for nlev=1;
  1388. ! for nlev > 1, the values are both 000000.01
  1389. !
  1390. val = spinf(i,k,1)
  1391. if ( val > 0.0 ) then
  1392. ! positive value means real part:
  1393. vri = 1
  1394. else if ( val < 0.0 ) then
  1395. ! negative value means imag part:
  1396. vri = 2
  1397. else
  1398. ! zero values for both real and imag part of (0,0)
  1399. nzero = nzero + 1 ! counter for number of zero values found
  1400. !--
  1401. !! test number of zero values:
  1402. !select case ( nzero )
  1403. ! case ( 1 ) ; vri = 1 ! real part of (0,0)
  1404. ! case ( 2 ) ; vri = 2 ! imag part of (0,0)
  1405. ! case default
  1406. ! write (gol,'("found more than 2 zero values in spectral info, ")'); call goErr
  1407. ! write (gol,'("while only expected for real and imag part of (0,0)")'); call goErr
  1408. ! TRACEBACK; status=1; return
  1409. ! cycle ! next value from received array
  1410. !end select
  1411. !--
  1412. ! assume that the extra zero's are the imaginary part for m=0,
  1413. ! which is zero anyway and does not need to be received:
  1414. if ( (nzero == 1) .and. (nlev == 1) ) then
  1415. vri = 1 ! real part of (0,0) in spinf2d
  1416. else
  1417. cycle ! next value from received array
  1418. end if
  1419. end if
  1420. val = abs(val)
  1421. vk = nint( ( val - floor(val) )*100.0 ) ! level is fractional part
  1422. vn = modulo( floor(val), 1000 ) ! last 3 numbers is n
  1423. vm = floor( val/1000.0 ) ! first 3 numbers is m
  1424. ! trap surface level ...
  1425. if ( nlev == 1 ) vk = vk + 1
  1426. ! check ...
  1427. if ( (vri < 1) .or. (vri > 2) .or. &
  1428. (vm < 0) .or. (vm > ifs_shT) .or. (vn < vm) .or. (vn > ifs_shT) .or. &
  1429. ((nlev==1) .and. (vk/=1)) .or. &
  1430. ((nlev>1) .and. ((vk < 1) .or. (vk > nlev))) ) then
  1431. write (gol,'("strange values extracted from spectral info value:")') ; call goErr
  1432. write (gol,'(" spinf : ",f16.4)' ) spinf(i,k,1) ; call goErr
  1433. write (gol,'(" ri : ",i4," (1=real,2=imag)")' ) vri ; call goErr
  1434. write (gol,'(" m : ",i4," (0 .. ",i4,")")' ) vm, ifs_shT ; call goErr
  1435. write (gol,'(" n : ",i4," (m .. ",i4,")")' ) vn, ifs_shT ; call goErr
  1436. write (gol,'(" k : ",i4," (1 .. ",i4,")")' ) vk, nlev ; call goErr
  1437. write (gol,'(" nzero : ",i4)' ) nzero ; call goErr
  1438. TRACEBACK; status=1; return
  1439. end if
  1440. ! position in triangle:
  1441. vp = sh_tripos(vm,vn)
  1442. ! check ...
  1443. if ( (vp < 1) .or. (vp > ifs_shn) ) then
  1444. write (gol,'("strange triangle position:")' ) ; call goErr
  1445. write (gol,'(" ifs T : ",i4)' ) ifs_shT ; call goErr
  1446. write (gol,'(" m : ",i4)' ) vm ; call goErr
  1447. write (gol,'(" n : ",i4)' ) vm ; call goErr
  1448. write (gol,'(" p : ",i8)' ) vp ; call goErr
  1449. TRACEBACK; status=1; return
  1450. end if
  1451. ! store:
  1452. shR%remap(i,k,1) = vri
  1453. shR%remap(i,k,2) = vp
  1454. shR%remap(i,k,3) = vk
  1455. ! maximum truncation:
  1456. shR%shT = max( shR%shT, max( vm, vn ) )
  1457. ! flag ...
  1458. sh_field(shR%remap(i,k,1),shR%remap(i,k,2),shR%remap(i,k,3)) = .true.
  1459. end do ! received coeff i
  1460. end do ! level k
  1461. ! check on missing target values:
  1462. if ( .not. all(sh_field) ) then
  1463. ! error counter:
  1464. nerr = 0
  1465. ! loop over levels:
  1466. do k = 1, nlev
  1467. ! init triangle position counter:
  1468. vp = 0
  1469. ! loop over spectral triangle:
  1470. do vm = 0, ifs_shT
  1471. do vn = vm, ifs_shT
  1472. ! increase triangle position counter:
  1473. vp = vp + 1
  1474. ! negative values at unexpected places ?
  1475. if ( ( (vm==0) .and. (.not. sh_field(1,vp,k) ) ) .or. &
  1476. ( (vm> 0) .and. (.not. all(sh_field(:,vp,k))) ) ) then
  1477. ! increase error counter:
  1478. nerr = nerr + 1
  1479. ! intro message ?
  1480. if ( nerr == 1 ) then
  1481. write (gol,'("not all sh target values filled :")'); call goErr
  1482. write (gol,'(" ifs T : ",i4)') ifs_shT; call goErr
  1483. end if
  1484. ! show error:
  1485. write (gol,'(" (m, n) p, ; k ; real imag : (",2i5,") ",i8," ; ",i4," ; ",2l2)') vm, vn, vp, k, sh_field(:,vp,k); call goErr
  1486. end if ! negatives ?
  1487. end do ! n
  1488. end do ! m
  1489. end do ! level
  1490. ! leave ?
  1491. if (nerr>0) then
  1492. TRACEBACK; status=1; return
  1493. end if
  1494. end if ! some negatives ?
  1495. ! done
  1496. deallocate( sh_field )
  1497. status = 0
  1498. END SUBROUTINE SHREMAP_SETUP
  1499. !EOC
  1500. !--------------------------------------------------------------------------
  1501. ! TM5 !
  1502. !--------------------------------------------------------------------------
  1503. !BOP
  1504. !
  1505. ! !IROUTINE: SHREMAP_REMAP
  1506. !
  1507. ! !DESCRIPTION:
  1508. !\\
  1509. !\\
  1510. ! !INTERFACE:
  1511. !
  1512. SUBROUTINE SHREMAP_REMAP( shR, sh_recv, shi, sh_ri, status )
  1513. !
  1514. ! !USES:
  1515. !
  1516. use grid, only : TshGridInfo
  1517. !
  1518. ! !INPUT/OUTPUT PARAMETERS:
  1519. !
  1520. type(TshRemap), intent(inout) :: shR
  1521. !
  1522. ! !INPUT PARAMETERS:
  1523. !
  1524. real, intent(in) :: sh_recv(:,:,:)
  1525. type(TshGridInfo), intent(in) :: shi
  1526. !
  1527. ! !OUTPUT PARAMETERS:
  1528. !
  1529. real, intent(out) :: sh_ri(:,:,:)
  1530. integer, intent(out) :: status
  1531. !
  1532. ! !REMARKS:
  1533. !
  1534. !EOP
  1535. !------------------------------------------------------------------------
  1536. !BOC
  1537. character(len=*), parameter :: rname = mname//'/shremap_Remap'
  1538. integer :: nlev, i, k
  1539. ! --- begin ---------------------------------------
  1540. ! number of levels:
  1541. nlev = size(sh_recv,2)
  1542. ! check shape of input array ...
  1543. if ( any( shape(sh_recv) /= (/ifs_shn_recv*2,nlev,1/) ) ) then
  1544. write (gol,'("strange input shape :")' ) ; call goErr
  1545. write (gol,'(" sh_recv : ",3i6)' ) shape(sh_recv) ; call goErr
  1546. write (gol,'(" expected : ",3i6)' ) ifs_shn_recv*2,nlev,1 ; call goErr
  1547. TRACEBACK; status=1; return
  1548. end if
  1549. ! check shape of output array ...
  1550. if ( any( shape(sh_ri) /= (/2,ifs_shn,nlev/) ) ) then
  1551. write (gol,'("strange input shape :")' ) ; call goErr
  1552. write (gol,'(" sh : ",3i6)' ) shape(sh_ri) ; call goErr
  1553. write (gol,'(" expected : ",3i6)' ) 2,ifs_shn,nlev ; call goErr
  1554. TRACEBACK; status=1; return
  1555. end if
  1556. ! initial zero:
  1557. sh_ri = 0.0
  1558. ! loop over all elements of received array:
  1559. do k = 1, nlev
  1560. do i = 1, ifs_shn_recv*2
  1561. ! the triplet shR%remap(i,k,:) defines (/ 1=real/2=imag, traingle-position, level /)
  1562. ! all negative ?
  1563. ! o this is a dummy element due to partitioning
  1564. ! o this is the imaginary part for m=0, which should remain zero
  1565. if ( all( shR%remap(i,k,:) < 0 ) ) cycle
  1566. ! any negative ? should not be possible...
  1567. if ( any( shR%remap(i,k,:) < 0 ) ) then
  1568. write (gol,'("found strange mapping:")' ) ; call goErr
  1569. write (gol,'(" triangle point : ",i6)' ) i ; call goErr
  1570. write (gol,'(" level : ",i6)' ) k ; call goErr
  1571. write (gol,'(" mapping : ",3i6)' ) shR%remap(i,k,:) ; call goErr
  1572. end if
  1573. ! copy value from received array into spectral field:
  1574. sh_ri(shR%remap(i,k,1),shR%remap(i,k,2),shR%remap(i,k,3)) = sh_recv(i,k,1)
  1575. end do
  1576. end do
  1577. status = 0
  1578. END SUBROUTINE SHREMAP_REMAP
  1579. !EOC
  1580. !--------------------------------------------------------------------------
  1581. ! TM5 !
  1582. !--------------------------------------------------------------------------
  1583. !BOP
  1584. !
  1585. ! !IROUTINE: SetPrismTime_date
  1586. !
  1587. ! !DESCRIPTION: returns current time/date into prism format (seconds from
  1588. ! prism reference start.
  1589. !\\
  1590. !\\
  1591. ! !INTERFACE:
  1592. !
  1593. subroutine SetPrismTime_date( prism_t, t, status )
  1594. !
  1595. ! !USES:
  1596. !
  1597. use GO, only : TDate, NewDate, iTotal, operator(-)
  1598. !
  1599. ! !OUTPUT PARAMETERS:
  1600. !
  1601. integer, intent(out) :: prism_t ! seconds from start
  1602. !
  1603. ! !INPUT PARAMETERS:
  1604. !
  1605. type(TDate), intent(in) :: t
  1606. integer, intent(out) :: status
  1607. !
  1608. ! !REMARKS:
  1609. !
  1610. !EOP
  1611. !------------------------------------------------------------------------
  1612. !BOC
  1613. character(len=*), parameter :: rname = mname//'/SetPrismTime_date'
  1614. ! --- begin ----------------------------------------
  1615. ! seconds since start
  1616. prism_t = iTotal( t - NewDate(time6=PRISM_start_date), 'sec' )
  1617. status = 0
  1618. end subroutine SetPrismTime_date
  1619. !EOC
  1620. END MODULE TM5_PRISM