source: trunk/SOURCES/lake_rsl_mod.f90 @ 334

Last change on this file since 334 was 329, checked in by dumas, 3 years ago

hi-resolution lake update in lake_rsl_mod.f90

File size: 33.1 KB
Line 
1!*************************************************************************************************************************************
2! This script creates runoff directions from a topographic dataset based on the algorithm of Wang and Liu (2006) described in:
3! An efficient method for identifying and filling surface depressions in digital elevation model for hydrologic analysis and modelling
4! International Journal of Geographical Information Science, Vol. 20, No. 2, Feb. 2006, 193-213
5
6!* It avoids endroheic basins and provide consistent runoff directions for use in the IPSL model.
7!  It compares favourably with directions created with the 'Fill Sinks' algorithm implemented in
8!  GIS softwares (like ArcGIS or QGIS) which is also based on the Wang and Liu algorithm.
9!  There are some differences in runoff directions related to the way I have interpreted some
10!  specific cases that occur. Basins generated by these directions are thus not identical but
11!  major basins are relatively similar.
12!  The corrected topography is however rigourously similar to that obtained from GIS softwares,
13!  at the exception of the Antarctic continent where I noted some inconsistencies in these
14!  softwares due to lack of east-west periodicity and of understanding the increasing longitudinal resolution.
15
16!* INPUTS:
17!  - a topographic dataset
18
19!* OUTPUTS:
20!  - rtm: runoff directions
21!  - topo: initial topographie
22!  - orog_lake : corrected topographic dataset to avoid endorheic basin
23!  - cmsk: coastal mask
24
25!*************************************************************************************************************************************
26
27!* C Dumas - 26.11.2020
28
29module lake_rsl
30
31use module3D_phy,only:nx,ny,S,BSOC,H,sealevel,sealevel_2d,num_param,num_rep_42,time,dtmin,ice,mk,debug_3D,ro,row,dirnameinp,xcc,ycc,Bsoc0,S0,H0
32use io_netcdf_grisli
33
34!$ USE OMP_LIB
35
36implicit none
37
38!~ !* General & discrete parameters
39!~ real :: undef=-9999. !undef=1.e20
40
41
42!~ !* Parameters and arrays used in the Initialization step
43real, dimension(nx,ny) :: orog, orog_lake, lake_level
44integer, dimension(nx,ny) :: omsk, cmsk
45integer, dimension(nx,ny) :: cont!, conttrack, contsave, runoff_msk, mask_search, exut
46integer, dimension(nx,ny) :: ocean
47integer, dimension(nx,ny) :: pot_lake
48integer :: hi_res_runoff ! select 1 : interpolation on hi resolution topo or 0 : no interpolation
49
50integer :: countcont ! nbr de continents
51!~ integer :: nbocean, countocean
52
53!~ !* Parameters and arrays of the runoff calculation
54integer, dimension(nx,ny) :: rtm
55
56!~ !* Debug mode & nc files
57character(len=100) :: file_topo_hi_res, coord_topo_hi_res
58!~ character(len=15) :: invar
59!~ integer :: ncid, xdimid, ydimid, varid
60!~ integer :: latid, lonid, rtmid, orogid, cmskid, contid, oroglakeid, lake_levelid
61
62!* Routing code parameter
63integer, dimension(8,2) :: inc
64
65!~ !cdc variable ajoutees
66
67!~ integer, dimension(8) :: loccont ! cont des 8 points autour du point ii jj
68
69!~ integer :: sumcont ! somme des cont des 8 points autour du point ii jj
70
71
72!
73!~ integer :: maxloop ! pour debug evite boucle infinies
74!~ integer,dimension(2) :: ijcoord ! coordonnees ij du point le plus bas
75!~ integer,dimension(2) :: ijcoord_cmsk ! coordonnees ij du point cotier le plus bas
76!~ real, dimension(8) :: locorog ! orog des 8 points voisins
77integer,dimension(nx,ny) :: lake ! 1 si lac, 0 sinon
78!~ integer,dimension(8) :: locexut ! 1 si exutoire => permet de gerer les lacs
79
80!~ !integer :: contprint=0 ! pour debug et arret du print
81!~ integer :: contexut ! compte les points sous le niveau du lac
82
83!~ logical :: stop_lake
84!~ integer, dimension(2) :: ij_lake
85!~ integer,dimension(8) :: locomsk ! omsk des 8 voisins
86!~ integer :: kloc
87
88integer,dimension(80,2) :: inclake ! 80 voisins (+4 -4) en i,j
89integer, dimension(nx,ny) :: mask_extlake ! extension du lac sous la glace
90real, dimension(nx,ny) :: orog_extlake ! niveau du lac sur les zones glace a proximité (rayon de 4 pts de grille)
91!~ integer, dimension(nx,ny) :: clake ! point lac avec voisin point englace
92real, dimension(nx,ny) :: sealevel_2D_loc
93
94
95!~ ! topo haute resolution (10km) pour calcul des lacs
96!~ ! nbr de points topo hi res
97integer, parameter :: nxhi = 964
98integer, parameter :: nyhi = 964
99integer, parameter :: ninterpo = 4 ! interpolation avec les 4 points autours
100
101real, dimension(nxhi,nyhi) :: Bsoc0_hi, S0_hi, H0_hi ! topo haute resolution actuelle (ref)
102!~ ! anomalie de topo sur la haute resolution
103!~ !real, dimension(nxhi,nyhi) :: anom_topo_hires
104!~ ! anomalie de topo
105real,dimension(nx,ny) :: anom_Bsoc, anom_S, anom_H ! anomalie de topo vs topo ref (low res)
106
107!~ ! matrice passage 40km vers 10km
108
109real, dimension(nxhi,nyhi) :: xcc_hi,ycc_hi,xlong_hi,ylat_hi
110real, dimension(nxhi,nyhi,4) :: poids  ! poids pour passage de low vers hi_res
111integer,dimension(nxhi,nyhi,4) :: low_2_hi_i, low_2_hi_j ! correspondance grille low vers grille hi_res
112
113real, dimension(nxhi,nyhi) :: anom_Bsoc_hi, anom_S_hi, anom_H_hi ! anomalie de topo resol std vs resol_hi
114real, dimension(nxhi,nyhi) :: Bsoc_hi, S_hi, H_hi ! topo haute resolution a l'instant t
115real, dimension(nxhi,nyhi) :: orog_hi ! orographie hi_res
116real, dimension(nxhi,nyhi) :: S_interp_hi, H_interp_hi ! interpolation de S et H low res vers hi_res sans anomalie (pour faire topo englacee non bruitee)
117integer, dimension(nxhi,nyhi) :: omsk_hi, cmsk_hi
118integer, dimension(nxhi,nyhi) :: cont_hi ! masque des points continent hi_res
119integer, dimension(nxhi,nyhi) :: ocean_hi ! masque des points ocean hi_res
120integer, dimension(nxhi,nyhi) :: pot_lake_hi
121integer, dimension(nxhi,nyhi) :: lake_hi ! 1 si lac, 0 sinon
122integer, dimension(nxhi,nyhi) :: rtm_hi ! runoff direction (1 to 8)
123real, dimension(nxhi,nyhi) :: lake_level_hi ! niveau du lac sur grille hi_res
124real, dimension(nxhi,nyhi) :: orog_lake_hi  ! orographie a la surface des lacs
125
126integer :: countcont_hi ! nbr de continents sur grille hi_res
127
128contains
129
130!> SUBROUTINE: input_rsl
131!! Routine permettant d'introduire des variations locales du niveau marin
132!! par la presence de lacs identifies via un algo de routage sur la topo de GRISLI
133!>
134subroutine input_rsl
135
136  implicit none
137
138  namelist/lake_rsl/hi_res_runoff, file_topo_hi_res, coord_topo_hi_res
139
140  integer :: ii, i, j, k
141  integer :: ios
142  real*8, dimension(:,:), pointer   :: tab2d => null()   !< tableau 2d pointer
143  real, dimension(ninterpo) :: dist_i, dist_j, dist_ij
144 
145  rewind(num_param)
146  read(num_param,lake_rsl)
147
148  write(num_rep_42,'(A)')'!___________________________________________________________'
149  write(num_rep_42,'(A)') '&lake_rsl                                    ! nom du bloc '
150  write(num_rep_42,*)
151  write(num_rep_42,*) 'hi_res_runoff = ', hi_res_runoff
152  write(num_rep_42,'(A,A,A)') 'file_topo_hi_res = ', trim(file_topo_hi_res)
153  write(num_rep_42,'(A,A,A)') 'coord_topo_hi_res = ', trim(coord_topo_hi_res)
154  write(num_rep_42,*)'/'
155  write(num_rep_42,*)
156 
157!~   print*,'Debut input_rsl'
158 
159!* Continent number (outfile)
160!ncfilecont = "continent_number_test-0.05.nc"
161!ncfilecmsk = "coastal_mask_test-0.05.nc"
162!* Runoff & corrected topo
163!ncfilernff = "runoff_20Ma_Fred.nc"
164!ncfiletopoout = "Topo.20Ma_Fred.OK_runoff.nc"
165!ncfilernff = "runoff_NorthAm_15ka_test-0.05.nc"
166!ncfiletopoout = "Topo.NorthAm_15ka_test.OK_runoff-0.05.nc"
167
168!***************************
169!*** Initialization step ***
170!***************************
171
172!~ write(*,*) "Start Initialization step"
173
174!* The routing code (1er indice : k=1,8, 2eme indice i=1, j=2)
175inc(1,1) = 0   ! indice i N
176inc(1,2) = 1   ! indice j N
177inc(2,1) = 1   ! indice i NE
178inc(2,2) = 1   ! indice j NE
179inc(3,1) = 1   ! indice i E
180inc(3,2) = 0   ! indice j E
181inc(4,1) = 1   ! indice i SE
182inc(4,2) = -1  ! indice j SE
183inc(5,1) = 0   ! indice i S
184inc(5,2) = -1  ! indice j S
185inc(6,1) = -1  ! indice i SW
186inc(6,2) = -1  ! indice j SW
187inc(7,1) = -1  ! indice i W
188inc(7,2) = 0   ! indice j W
189inc(8,1) = -1  ! indice i NW
190inc(8,2) = 1   ! indice j NW
191
192!~ print*
193!~ print*,'rtm directions :'
194!~ print*,'1 :  N'
195!~ print*,'2 : NE'
196!~ print*,'3 :  E'
197!~ print*,'4 : SE'
198!~ print*,'5 :  S'
199!~ print*,'6 : SW'
200!~ print*,'7 :  W'
201!~ print*,'8 : NW'
202!~ print*
203
204
205!do ii=1,80
206ii = 1
207do j = -4,4
208  do i = -4,4
209    if (.not.((i.eq.0).and.(j.eq.0))) then
210      inclake(ii,1)=i
211      inclake(ii,2)=j
212      ii=ii+1
213    endif
214  enddo
215enddo
216
217
218if (hi_res_runoff.eq.1) then
219! lecture topo hi resolution :
220  call Read_Ncdf_var('Bsoc',trim(dirnameinp)//trim(file_topo_hi_res),tab2d) ! socle
221  Bsoc0_hi(:,:)=tab2d(:,:)
222  call Read_Ncdf_var('S',trim(dirnameinp)//trim(file_topo_hi_res),tab2d)    ! surface
223  S0_hi(:,:)=tab2d(:,:)
224  call Read_Ncdf_var('H',trim(dirnameinp)//trim(file_topo_hi_res),tab2d)    ! epaisseur
225  H0_hi(:,:)=tab2d(:,:)
226 
227! lecture des coordonnées geographiques
228! les coordonnees sont calculees en °dec avec GMT,
229! les longitudes sont comprises entre -180 et +180 (negative a l'Ouest de
230! Greenwich et positive a l'Est)
231  open(unit=2004,file=trim(dirnameinp)//trim(coord_topo_hi_res),iostat=ios)
232  do k=1,nxhi*nyhi
233     read(2004,*) i,j,xcc_hi(i,j),ycc_hi(i,j),xlong_hi(i,j),ylat_hi(i,j)
234  enddo
235  close(2004)
236             
237!  xmin=xcc(1,1)/1000.
238!  ymin=ycc(1,1)/1000.
239!  xmax=xcc(nx,ny)/1000.
240!  ymax=ycc(nx,ny)/1000.
241
242
243! calcul de la matrice de passage de 40km => hi_res (10km)
244  do j=1,nyhi
245    do i=1,nxhi
246      low_2_hi_i(i,j,1) = max(nint((i-1)/4.),1)
247      low_2_hi_i(i,j,2) = min(low_2_hi_i(i,j,1) + 1,nx)
248      low_2_hi_i(i,j,3) = max(nint((i-1)/4.),1)
249      low_2_hi_i(i,j,4) = min(low_2_hi_i(i,j,1) + 1,nx)
250      low_2_hi_j(i,j,1) = max(nint((j-1)/4.),1)
251      low_2_hi_j(i,j,2) = min(low_2_hi_j(i,j,1) + 1,ny)
252      low_2_hi_j(i,j,3) = min(low_2_hi_j(i,j,1) + 1,ny)
253      low_2_hi_j(i,j,4) = max(nint((j-1)/4.),1)
254! calcul de la distance entre le pt hi_res et les 4 points low_res autours
255      do k=1,ninterpo
256        dist_i(k)=abs(xcc(low_2_hi_i(i,j,k),low_2_hi_j(i,j,k))-xcc_hi(i,j))
257        dist_j(k)=abs(ycc(low_2_hi_i(i,j,k),low_2_hi_j(i,j,k))-ycc_hi(i,j))
258        dist_ij(k)=sqrt(dist_i(k)**2+dist_j(k)**2)
259      enddo
260      do k=1,ninterpo
261        poids(i,j,k)=(sum(dist_ij(:))-dist_ij(k))/(3*sum(dist_ij(:)))
262      enddo
263    enddo
264  enddo
265
266! calcul de l'anomalie de topo vs la topo de reference :
267  anom_Bsoc(:,:) = Bsoc(:,:) - Bsoc0(:,:)
268  anom_S(:,:) = S(:,:) - S0(:,:)
269  anom_H(:,:) = H(:,:) - H0(:,:)
270
271! interpolation 40km => hi_res
272  anom_Bsoc_hi(:,:)=0.
273  anom_S_hi(:,:)=0.
274  anom_H_hi(:,:)=0.
275
276  do k=1,ninterpo
277    do j=1,nyhi
278      do i=1,nxhi
279        anom_Bsoc_hi(i,j) = anom_Bsoc_hi(i,j) + poids(i,j,k)*anom_Bsoc(low_2_hi_i(i,j,k),low_2_hi_j(i,j,k))
280        anom_S_hi(i,j) = anom_S_hi(i,j) + poids(i,j,k)*anom_S(low_2_hi_i(i,j,k),low_2_hi_j(i,j,k))
281        anom_H_hi(i,j) = anom_H_hi(i,j) + poids(i,j,k)*anom_H(low_2_hi_i(i,j,k),low_2_hi_j(i,j,k))
282      enddo
283    enddo
284  enddo
285endif
286
287end subroutine input_rsl
288
289
290!------------------------------------------------------------------------
291! subroutine rsl : caclcul du niveau marin local
292!------------------------------------------------------------------------
293subroutine rsl
294
295integer :: i, j, k, l
296integer :: ii, jj
297! pour les lectures ncdf
298real*8, dimension(:,:),   pointer    :: tab => null()       !< tableau 2d real pointer
299
300real, dimension(nx,ny) :: orog_lake_tmp ! pour calculer orog_lake a partir de orog_lake_hi
301! pour debug : ecriture resultats hi_res
302character(len=100) :: ncfilecont, ncfilernff, ncfiletopoout,ncfilecmsk    ! nom fichier netcdf de sortie
303character(len=1),dimension(2) :: dimname ! nom des dimensions du fichier netcdf de sortie
304integer :: ncid, xdimid, ydimid, S_hivarid, H_hivarid, Bsoc_hivarid, rtm_hivarid, cont_hivarid, ocean_hivarid, cmsk_hivarid, orog_lake_hivarid, lake_hivarid
305integer :: status ! erreur operations netcdf
306
307ncfilecont = "cont-topo_hi_res.nc"
308dimname(1)='x'
309dimname(2)='y'
310
311!~ print*,'Debut rsl time = ', time
312! orography sur grille std :
313
314
315where ((BSOC(:,:)+H(:,:)*ro/row -sealevel).LT.0.)
316  orog(:,:) = BSOC(:,:) ! le point flotte et socle sous le nivx de la mer
317elsewhere
318  orog(:,:) = S(:,:)    ! point terre
319endwhere
320
321orog_lake(:,:)=orog(:,:) ! orographie a la surface des lacs
322
323where (orog(:,:).GT.sealevel) ! les depressions continentales seront ajoutees a ce masque plus loin
324  omsk(:,:) = 1
325elsewhere
326  omsk(:,:) = 0
327endwhere
328
329!~ print*,'1/ rsl'
330
331if (hi_res_runoff.eq.1) then
332! calcul de l'anomalie de topo vs la topo de reference :
333  anom_Bsoc(:,:) = Bsoc(:,:) - Bsoc0(:,:)
334  anom_S(:,:) = S(:,:) - S0(:,:)
335  anom_H(:,:) = H(:,:) - H0(:,:)
336 
337  anom_Bsoc_hi(:,:) = 0.
338  anom_S_hi(:,:) = 0.
339  S_interp_hi(:,:) = 0.
340  H_interp_hi(:,:) = 0.
341  anom_H_hi(:,:) = 0.
342 
343  ! calcul de la topo hi_res :
344  do k=1,ninterpo
345    do j=1,nyhi
346      do i=1,nxhi
347        anom_Bsoc_hi(i,j) = anom_Bsoc_hi(i,j) + poids(i,j,k)*anom_Bsoc(low_2_hi_i(i,j,k),low_2_hi_j(i,j,k))
348        anom_S_hi(i,j) = anom_S_hi(i,j) + poids(i,j,k)*anom_S(low_2_hi_i(i,j,k),low_2_hi_j(i,j,k))
349        S_interp_hi(i,j) = S_interp_hi(i,j) + poids(i,j,k)*S(low_2_hi_i(i,j,k),low_2_hi_j(i,j,k))
350        H_interp_hi(i,j) = H_interp_hi(i,j) + poids(i,j,k)*H(low_2_hi_i(i,j,k),low_2_hi_j(i,j,k))
351!        anom_H_hi(i,j) = anom_H_hi(i,j) + poids(i,j,k)*anom_H(low_2_hi_i(i,j,k),low_2_hi_j(i,j,k))
352      enddo
353    enddo
354  enddo
355 
356  Bsoc_hi(:,:) = Bsoc0_hi(:,:) + anom_Bsoc_hi(:,:)
357 
358  where (H_interp_hi(:,:).GT.1.) ! zones englacees
359    where (Bsoc_hi(:,:).GT.S_interp_hi(:,:)) ! socle qui sort de la glace
360      S_hi(:,:) = Bsoc_hi(:,:)
361      H_hi(:,:) = 0.
362    elsewhere ! socle sous la glace
363      S_hi(:,:) = S_interp_hi(:,:)
364      H_hi(:,:) = S_interp_hi(:,:) - Bsoc_hi(:,:)
365    endwhere
366  elsewhere ! zones non englacees
367    S_hi(:,:) = max(Bsoc_hi(:,:),sealevel)
368    H_hi(:,:) = 0.
369  endwhere
370 
371! orography sur grille hi_res :
372  where ((Bsoc_hi(:,:)+H_hi(:,:)*ro/row -sealevel).LT.0.)
373    orog_hi(:,:) = Bsoc_hi(:,:) ! le point flotte et socle sous le nivx de la mer
374  elsewhere
375    orog_hi(:,:) = S_hi(:,:)    ! point terre
376  endwhere
377 
378  orog_lake_hi(:,:)=orog_hi(:,:) ! orographie a la surface des lacs
379
380  where (orog_hi(:,:).GT.sealevel) ! les depressions continentales seront ajoutees a ce masque plus loin
381    omsk_hi(:,:) = 1
382  elsewhere
383    omsk_hi(:,:) = 0
384  endwhere
385
386endif
387
388!~ print*,'2/ rsl'
389! mise à jour du niveau des lacs et du routage :
390if (mod(abs(TIME),500.).lt.dtmin) then
391  call mask_orog_ocean(nx,ny,orog,omsk,cmsk,cont,ocean,countcont,pot_lake)
392
393!~   print*,'3/ rsl'
394  if (hi_res_runoff.eq.1) then
395    call mask_orog_ocean(nxhi,nyhi,orog_hi,omsk_hi,cmsk_hi,cont_hi,ocean_hi,countcont_hi,pot_lake_hi)
396    call runoff(nxhi,nyhi,countcont_hi,cont_hi,cmsk_hi,omsk_hi,pot_lake_hi,lake_hi,rtm_hi,lake_level_hi,orog_lake_hi)
397! interpolation des variables sur grille std :
398    orog_lake(:,:) = 0.
399    orog_lake_tmp(:,:) = 0.
400    lake(:,:) = 0
401    do j=1,ny
402      do i=1,nx
403        orog_lake_tmp(i,j) = max(sum(orog_lake_hi(i*4-3:i*4,j*4-3:j*4)) / 16., orog(i,j))
404!        lake(i,j) = nint(sum(lake_hi(i*4-3:i*4,j*4-3:j*4)) / 16.)  ! lac si au moins 50% des points sont lacs
405        lake(i,j) = floor(sum(lake_hi(i*4-3:i*4,j*4-3:j*4)) / 16.)  ! lac si 100% des points hi_res sont lacs
406      enddo
407    enddo
408    where (lake(:,:).EQ.1) orog_lake(:,:)=orog_lake_tmp(:,:)
409 
410  else
411    call runoff(nx,ny,countcont,cont,cmsk,omsk,pot_lake,lake,rtm,lake_level,orog_lake)
412  endif
413
414
415!~   print*,'3/ rsl'
416  call extlake(nx,ny,lake,ice,mk,orog_lake,inclake,mask_extlake,orog_extlake)
417
418!~   print*,'4/ rsl'
419
420  sealevel_2D_loc(:,:) = sealevel_2D(:,:) ! pour debug : A SUPPRIMER QUAND ON REPASSERA SUR SEALEVEL_2D
421  where (lake(:,:).eq.1)
422    sealevel_2D_loc(:,:)=orog_lake(:,:)
423  elsewhere
424    sealevel_2D_loc(:,:)=sealevel
425  endwhere
426
427  where ((mask_extlake.LT.9999).AND.(ice.EQ.1).AND.(mk.EQ.0).AND.(lake.EQ.0).AND.(BSOC.LT.orog_extlake))
428    sealevel_2D_loc(:,:) = orog_extlake(:,:)
429  endwhere
430
431  if (hi_res_runoff.eq.1) then
432
433!* Save mode
434    status = nf90_create(trim(ncfilecont),NF90_CLOBBER, ncid)
435!  status = nf90_close(ncid)
436    status = nf90_def_dim(ncid, "x", nxhi, xdimid)
437    status = nf90_def_dim(ncid, "y", nyhi, ydimid)
438    status = nf90_def_var(ncid, "S_hi", nf90_float, (/ xdimid, ydimid /), S_hivarid)
439    status = nf90_def_var(ncid, "H_hi", nf90_float, (/ xdimid, ydimid /), H_hivarid)
440    status = nf90_def_var(ncid, "Bsoc_hi", nf90_float, (/ xdimid, ydimid /), Bsoc_hivarid)
441    status = nf90_def_var(ncid, "rtm_hi", nf90_float, (/ xdimid, ydimid /), rtm_hivarid)
442    status = nf90_def_var(ncid, "cont_hi", nf90_float, (/ xdimid, ydimid /), cont_hivarid)
443    status = nf90_def_var(ncid, "ocean_hi", nf90_float, (/ xdimid, ydimid /), ocean_hivarid)
444    status = nf90_def_var(ncid, "cmsk_hi", nf90_float, (/ xdimid, ydimid /), cmsk_hivarid)
445    status = nf90_def_var(ncid, "orog_lake_hi", nf90_float, (/ xdimid, ydimid /), orog_lake_hivarid)
446    status = nf90_def_var(ncid, "lake_hi", nf90_float, (/ xdimid, ydimid /), lake_hivarid)
447 
448    status = nf90_enddef(ncid)
449 
450    status = nf90_put_var(ncid, S_hivarid, S_hi)
451    status = nf90_put_var(ncid, H_hivarid, H_hi)
452    status = nf90_put_var(ncid, Bsoc_hivarid, Bsoc_hi)
453    status = nf90_put_var(ncid, rtm_hivarid, rtm_hi)
454    status = nf90_put_var(ncid, cont_hivarid, cont_hi)
455    status = nf90_put_var(ncid, ocean_hivarid, ocean_hi)
456    status = nf90_put_var(ncid, cmsk_hivarid, cmsk_hi)
457    status = nf90_put_var(ncid, orog_lake_hivarid, orog_lake_hi)
458    status = nf90_put_var(ncid, lake_hivarid, lake_hi)
459 
460    status = nf90_close(ncid)
461  endif
462
463
464endif
465
466sealevel_2D(:,:) = sealevel_2D_loc(:,:)
467
468debug_3d(:,:,126)=sealevel_2D_loc(:,:)
469debug_3d(:,:,127)=cont(:,:)
470debug_3d(:,:,128)=lake(:,:)
471
472!~ print*,'Fin rsl'
473
474end subroutine rsl
475
476
477
478!----------------------------------------------------------------------------
479! Subroutine de recherche des points terre, cotiers et ocean
480!----------------------------------------------------------------------------
481subroutine mask_orog_ocean(nx_loc,ny_loc,orog_loc,omsk_loc,cmsk_loc,cont_loc,ocean_loc,countcont_loc,pot_lake_loc)
482
483implicit none
484
485integer,intent(in) :: nx_loc,ny_loc
486real,dimension(nx_loc,ny_loc), intent(in) :: orog_loc
487integer,dimension(nx_loc,ny_loc), intent(inout) :: omsk_loc
488integer,dimension(nx_loc,ny_loc), intent(out) :: cmsk_loc, cont_loc
489integer, dimension(nx_loc,ny_loc), intent(out) :: ocean_loc
490integer, intent(out) :: countcont_loc
491integer, dimension(nx_loc,ny_loc), intent(out) :: pot_lake_loc
492
493
494integer :: undefint=0
495integer, dimension(nx_loc,ny_loc) :: conttrack
496integer :: sumcont ! somme des cont des 8 points autour du point ii jj
497integer :: maxloop ! pour debug evite boucle infinies
498integer, dimension(8) :: ipinctab, jpinctab ! indice des points a scanner autour du point ii jj
499integer, dimension(8) :: loccont ! cont des 8 points autour du point ii jj
500integer, dimension(nx_loc,ny_loc) :: contsave
501
502integer :: nbocean, countocean
503integer :: countpot_lake ! nbr de dépressions (potential lakes)
504logical :: stop_lake
505integer, dimension(2) :: ij_lake
506
507integer,dimension(8) :: locomsk ! omsk des 8 voisins
508integer :: maxcont ! valeur max de cont parmi les 8 voisins
509integer :: sumocean ! somme des ocean des 8 points autour du point ii jj
510integer, dimension(8) :: lococean ! ocean des 8 points autour du point ii jj
511integer :: maxocean ! valeur max de ocean des 8 points autour du point ii jj
512integer,dimension(nx_loc,ny_loc) :: oceansave, true_ocean
513integer, parameter :: oceanmax=3000 ! nbr max de zones ocean ou lac
514logical, dimension(oceanmax) :: k_ocean ! permet d'identifier les ocean ouverts des lacs
515integer :: nbcont
516integer :: kloc
517integer :: i,j,ii,jj,k,l
518
519!* Initialize the coastal mask (cmsk_loc)
520cmsk_loc = 0
521!~ print*,'1/ mask_orog_ocean'
522
523!* Discriminate continents
524!~ write(*,*) "Number the different continents"
525cont_loc = undefint
526nbcont = 1
527conttrack = 0
528
529do j = 1, ny_loc
530  do i = 1, nx_loc
531!do j = 110,111
532!  do i = 182,183
533
534    !* Locate current grid point
535    ii = i; jj = j
536    sumcont=0
537    maxloop=0
538
539    !* Only on continental points
540    if ( omsk_loc(ii,jj) .eq. 1) then
541      do k=1,8
542        ipinctab(k) = max(1,min(nx_loc,ii+inc(k,1)))
543! cas cyclicité selon i :
544!        ipinctab(k) = ii+inc(k,1)
545!        if (ipinctab(k).eq.(nx + 1)) ipinctab(k) = 1
546!        if (ipinctab(k).eq.(0)) ipinctab(k) = nx
547        jpinctab(k) = max(1,min(ny_loc,jj+inc(k,2)))
548        loccont(k)=cont_loc(ipinctab(k),jpinctab(k))
549        sumcont = sumcont + cont_loc(ipinctab(k),jpinctab(k))
550      enddo
551 
552      if (sumcont.GT.0) then
553        cont_loc(ii,jj)=minval(loccont,mask=loccont>0)
554        do while (maxval(loccont,mask=loccont>0).GT.cont_loc(ii,jj).and.maxloop.lt.8) ! test si dans les 8 voisins certains ont un loccont > minval
555          maxcont=maxval(loccont,mask=loccont>0)
556          where (cont_loc(:,:).EQ.maxcont) cont_loc(:,:) = cont_loc(ii,jj)
557          where (loccont(:).EQ.maxcont) loccont(:) = cont_loc(ii,jj)
558          maxloop = maxloop + 1
559        enddo
560      else
561        cont_loc(ii,jj) = nbcont
562        nbcont = nbcont + 1
563      endif
564    conttrack(ii,jj) = 1
565    endif
566  enddo
567enddo
568!~ print*,'2/ mask_orog_ocean'   
569
570!* Reordering to get an incremental list
571countcont_loc = 1
572contsave(:,:)=cont_loc(:,:) ! copie cont for reordonation
573
574if (nbcont.gt.0) then ! on a trouve des continents !!!!
575  do k=1,nbcont
576    if (count(contsave(:,:).eq.k).gt.0) then ! there is continents associated with number k
577      where (cont_loc(:,:).eq.k) cont_loc(:,:) = countcont_loc
578      countcont_loc = countcont_loc + 1
579    endif
580  enddo
581else
582  print*,'CONTINENT NOT FOUND'
583endif
584!~ print*,'3/ mask_orog_ocean'
585
586! recherche des lacs et oceans sous le niveau de la mer
587!~ write(*,*) "Number the different lakes & oceans"
588ocean_loc = undefint
589nbocean = 1
590
591do j = 1, ny_loc
592  do i = 1, nx_loc
593    !* Locate current grid point
594    ii = i; jj = j
595    sumocean=0
596    maxloop=0
597
598    !* Only on ocean points
599    if ( omsk_loc(ii,jj) .eq. 0) then  ! omask=0 => orog < sealevel
600      do k=1,8
601        ipinctab(k) = max(1,min(nx_loc,ii+inc(k,1)))
602! cas cyclicité selon i :
603!        ipinctab(k) = ii+inc(k,1)
604!        if (ipinctab(k).eq.(nx + 1)) ipinctab(k) = 1
605!        if (ipinctab(k).eq.(0)) ipinctab(k) = nx
606        jpinctab(k) = max(1,min(ny_loc,jj+inc(k,2)))
607        lococean(k)=ocean_loc(ipinctab(k),jpinctab(k))
608        sumocean = sumocean + ocean_loc(ipinctab(k),jpinctab(k))
609      enddo
610
611      if (sumocean.GT.0) then
612        ocean_loc(ii,jj)=minval(lococean,mask=lococean>0)
613        do while (maxval(lococean,mask=lococean>0).GT.ocean_loc(ii,jj).and.maxloop.lt.8) ! test si dans les 8 voisins certains ont un lococean > minval
614          maxocean=maxval(lococean,mask=lococean>0)
615          where (ocean_loc(:,:).EQ.maxocean) ocean_loc(:,:) = ocean_loc(ii,jj)
616          where (lococean(:).EQ.maxocean) lococean(:) = ocean_loc(ii,jj)
617          maxloop = maxloop + 1
618        enddo
619      else
620        ocean_loc(ii,jj) = nbocean
621        nbocean = nbocean + 1
622      endif
623    endif
624  enddo
625enddo
626
627!~ print*,'4/ mask_orog_ocean' 
628
629!* Reordering to get an incremental list
630countocean = 1
631countpot_lake = 1 ! nombre de dépression intérieur (lacs potentiels
632oceansave(:,:)=ocean_loc(:,:) ! copie ocean for reordonation
633true_ocean(:,:)=0 ! point vraiment ocean
634pot_lake_loc(:,:)=0
635
636k_ocean(:)=.false.
637 
638
639if (nbocean.gt.0) then ! on a trouve des oceans !!!!
640! si l'ocean a un point qui touche un des bords de grille c'est bien un ocean, sinon c'est un lac
641  do i=1,nx_loc
642    if (oceansave(i,1).GT.0) k_ocean(oceansave(i,1))=.true.
643    if (oceansave(i,ny_loc).GT.0) k_ocean(oceansave(i,ny_loc))=.true.
644  enddo
645  do j=1,ny_loc
646    if (oceansave(1,j).GT.0) k_ocean(oceansave(1,j))=.true.
647    if (oceansave(nx_loc,j).GT.0) k_ocean(oceansave(nx_loc,j))=.true.
648  enddo
649 
650  do k=1,nbocean
651    if (count(oceansave(:,:).eq.k).gt.0) then ! there is oceans points associated with number k
652      if (k_ocean(k)) then  ! point is true ocean
653        where (oceansave(:,:).eq.k)
654          ocean_loc(:,:) = countocean
655        endwhere
656        countocean = countocean + 1
657      else ! point is a depresion under sea level
658        where (oceansave(:,:).eq.k)
659          pot_lake_loc(:,:) = countpot_lake
660          ocean_loc(:,:) = undefint
661        endwhere
662        countpot_lake = countpot_lake + 1
663      endif
664    endif
665  enddo
666else
667  print*,'OCEAN NOT FOUND'
668endif
669!~ print*,'5/ mask_orog_ocean countpot_lake', countpot_lake
670! les points pot_lake sont continent => on leur attribue le no du continent adjacent.
671
672do l=1,countpot_lake
673  stop_lake = .true.
674  ij_lake=minloc(pot_lake_loc(:,:),mask=(pot_lake_loc.EQ.l))
675  do while (stop_lake)
676    do k=1,8
677      ipinctab(k) = max(1,min(nx_loc,ij_lake(1)+inc(k,1)))
678! si cyclicite selon i :
679!      ipinctab(k) = ijcoord(1)+inc(k,1)
680!      if (ipinctab(k).eq.(nx + 1)) ipinctab(k) = 1
681!      if (ipinctab(k).eq.(0)) ipinctab(k) = nx
682      jpinctab(k) = max(1,min(ny_loc,ij_lake(2)+inc(k,2)))
683      locomsk(k)=omsk_loc(ipinctab(k),jpinctab(k))
684      loccont(k)=cont_loc(ipinctab(k),jpinctab(k))
685    enddo
686    if (maxval(loccont).GT.0) then ! voisin point terre
687      ! le lac appartient a ce continent
688      kloc=maxloc(loccont,dim=1)
689      where (pot_lake_loc(:,:).EQ.l)
690        cont_loc(:,:)=cont_loc(ipinctab(kloc),jpinctab(kloc))
691        omsk_loc(:,:)=1
692      endwhere
693      stop_lake=.false.
694    else
695      ij_lake(1) = ij_lake(1)+1
696    endif
697  enddo
698enddo
699!~ print*,'6/ mask_orog_ocean'
700
701! recherche des points cotiers
702where ((omsk_loc(:,:).EQ.1).AND.((eoshift(omsk_loc,shift=-1,dim=1).EQ.0) .OR. &
703  (eoshift(omsk_loc,shift=1,dim=1).EQ.0) .OR. (eoshift(omsk_loc,shift=-1,dim=2).EQ.0) .OR. (eoshift(omsk_loc,shift=1,dim=2).EQ.0)))
704  cmsk_loc(:,:) = 1
705endwhere
706!~ print*,'7/ mask_orog_ocean'
707
708end subroutine mask_orog_ocean
709
710
711
712!---------------------------------------------------------
713! subroutine qui calcul l'écoulement (runoff) et les lacs
714!---------------------------------------------------------
715subroutine runoff(nx_loc,ny_loc,countcont_loc,cont_loc,cmsk_loc,omsk_loc,pot_lake_loc,lake_loc,rtm_loc,lake_level_loc,orog_lake_loc)
716
717implicit none
718
719integer, intent(in) :: nx_loc, ny_loc
720integer, intent(in) :: countcont_loc ! nbr de continents
721
722integer,dimension(nx_loc,ny_loc), intent(in) :: cont_loc
723integer,dimension(nx_loc,ny_loc), intent(in) :: cmsk_loc
724integer,dimension(nx_loc,ny_loc), intent(in) :: omsk_loc
725integer,dimension(nx_loc,ny_loc), intent(in) :: pot_lake_loc
726
727integer,dimension(nx_loc,ny_loc), intent(out) :: lake_loc ! 1 si lac, 0 sinon
728integer,dimension(nx_loc,ny_loc), intent(out) :: rtm_loc ! runoff direction (1 to 8)
729real,dimension(nx_loc,ny_loc), intent(out) :: lake_level_loc ! lake level
730real,dimension(nx_loc,ny_loc), intent(inout) :: orog_lake_loc  ! orographie a la surface des lacs
731
732real,dimension(nx_loc,ny_loc) :: orog_extlake ! niveau du lac sur les zones glace a proximité (rayon de 4 pts de grille)
733
734integer,dimension(nx_loc,ny_loc) :: exut ! > 0 si exutoire, 0 sinon
735
736integer,dimension(nx_loc,ny_loc) :: mask_search
737integer,dimension(nx_loc,ny_loc) :: runoff_msk
738integer :: nbcont ! no du continent étudié
739integer,dimension(2) :: ijcoord ! coordonnees ij du point le plus bas
740
741integer,dimension(2) :: ijcoord_cmsk ! coordonnees ij du point cotier le plus bas
742integer :: contexut ! compte les points sous le niveau du lac
743real,dimension(8) :: locorog ! orog des 8 points voisins
744integer,dimension(8) :: locexut ! 1 si exutoire => permet de gerer les lacs
745integer,dimension(8) :: ipinctab,jpinctab ! indice des points a scanner autour du point ii jj
746integer :: countpot_lake ! nbr de dépressions (potential lakes)
747
748
749integer :: k,l
750
751!**************************
752!*** Runoff calculation ***
753!**************************
754
755write(*,*) "Start runoff calculation"
756
757
758!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
759! version inspiree article
760!$OMP PARALLEL
761!$OMP WORKSHARE
762mask_search(:,:)=0 ! masque des points a prendre en compte a l'iteration n : 0 pas traité, 1 en cours, 2 traité
763runoff_msk(:,:)=0  ! initialisation du mask points non traités => 0
764lake_loc(:,:)=0
765exut(:,:)=0 ! point exutoire d'un lac => 1
766!output(:,:)=0 ! initialisation des points de sortie vers la mer (0 pas sortie, 1 pt de sortie)
767!* We work on each continent separately for convenience
768!do nbcont = 1,1
769!$OMP END WORKSHARE
770!$OMP END PARALLEL
771
772do nbcont = 1,countcont_loc
773  write(*,*) "Working on continent nb :", nbcont
774! on demarre des points cotiers en partant du plus bas puis en remontant en altitude
775! mise a jour du masque des points a etudier :
776
777! demarrage par les points cotier du continent scanné :
778  where ((cont_loc.eq.nbcont).and.(cmsk_loc.eq.1)) mask_search(:,:)=1 ! masque de scan : points cote pour demarrer
779 
780! debut de la boucle de recherche du point a scanner
781  do while (minval(runoff_msk(:,:),mask=mask_search(:,:).EQ.1).EQ.0) ! test si il reste des points a scanner
782    ijcoord(:)=minloc(orog_lake_loc,mask=(mask_search(:,:).EQ.1))   ! position du point le plus bas etudié
783   
784    if (count(mask=((mask_search.eq.1).and.(cmsk_loc.eq.1))).GT.0) then! il reste des points cotiers non traités
785    ! position du point cotier le plus bas :
786      ijcoord_cmsk(:)=minloc(orog_lake_loc,mask=((mask_search(:,:).EQ.1).and.cmsk_loc(:,:).eq.1))
787      ! si altitude egale entre 2 points on traite en priorité le point cotier
788      if (orog_lake_loc(ijcoord(1),ijcoord(2)).EQ.orog_lake_loc(ijcoord_cmsk(1),ijcoord_cmsk(2))) then
789        ijcoord(:) = ijcoord_cmsk(:)
790      endif
791    endif
792    contexut = 0
793    do k=1,8
794      ipinctab(k) = max(1,min(nx_loc,ijcoord(1)+inc(k,1)))
795! si cyclicite selon i :
796!      ipinctab(k) = ijcoord(1)+inc(k,1)
797!      if (ipinctab(k).eq.(nx + 1)) ipinctab(k) = 1
798!      if (ipinctab(k).eq.(0)) ipinctab(k) = nx
799      jpinctab(k) = max(1,min(ny_loc,ijcoord(2)+inc(k,2)))
800     
801      ! on ajoute les points adjacents aux points scannés pour le prochain tour sauf si ce sont des points ocean ou dejà traités:
802      if ((omsk_loc(ipinctab(k),jpinctab(k)).EQ.1).AND.(mask_search(ipinctab(k),jpinctab(k)).EQ.0)) then
803        mask_search(ipinctab(k),jpinctab(k))=1
804! cas des depressions :
805! si l'un des points est plus bas que le point de base et n'est pas cotier => on remonte son orog => orog_lake_loc
806        if ((orog_lake_loc(ipinctab(k),jpinctab(k))).LE.(orog_lake_loc(ijcoord(1),ijcoord(2))) &
807          .AND.(cmsk_loc(ipinctab(k),jpinctab(k)).NE.1)) then
808          orog_lake_loc(ipinctab(k),jpinctab(k)) = orog_lake_loc(ijcoord(1),ijcoord(2))
809          lake_loc(ipinctab(k),jpinctab(k)) = 1
810          contexut= contexut + 1
811        endif
812      endif
813      locorog(k) = orog_lake_loc(ipinctab(k),jpinctab(k))
814      locexut(k) = exut(ipinctab(k),jpinctab(k))
815    enddo
816    if ((contexut.ge.1).and.(exut(ijcoord(1),ijcoord(2)).EQ.0)) exut(ijcoord(1),ijcoord(2)) = 1
817    rtm_loc(ijcoord(1),ijcoord(2))=minloc(locorog,dim=1) ! sens ecoulement donne par minloc de la topo des 8 voisins
818    ! cas d'un point lac
819    if (lake_loc(ijcoord(1),ijcoord(2)).eq.1) then
820      exut(ijcoord(1),ijcoord(2)) = maxval(locexut) + 1
821    ! ecoulement vers l'exutoire du lac
822      if (maxval(locexut).GT.0) then
823        rtm_loc(ijcoord(1),ijcoord(2))=minloc(locexut,dim=1,mask=(locexut.ge.1)) ! ecoulement vers l'exutoire
824      endif
825    endif
826    ! on supprime le point attribué de mask_search :
827    mask_search(ijcoord(1),ijcoord(2))=2
828    ! le point a ete traité :
829    runoff_msk(ijcoord(1),ijcoord(2))=1
830  enddo
831enddo
832
833! pour chaque lac calcul du niveau max :
834lake_level_loc(:,:)=-9999. ! initialisation
835do l=1,countpot_lake
836  where (pot_lake_loc(:,:).EQ.l)
837    lake_level_loc(:,:)=orog_lake_loc(:,:)
838  endwhere
839enddo
840
841end subroutine runoff
842
843!--------------------------------------------------------------------
844!  subroutine extlake : extension des lacs sous la glace
845!--------------------------------------------------------------------
846
847subroutine extlake(nx_loc,ny_loc,lake_loc,ice_loc,mk_loc,orog_lake_loc,inclake_loc,mask_extlake_loc,orog_extlake_loc)
848
849  implicit none
850
851  integer, intent(in) :: nx_loc, ny_loc
852  integer,dimension(nx_loc,ny_loc), intent(in) :: lake_loc ! 1 si lac, 0 sinon
853  integer,dimension(nx_loc,ny_loc), intent(in) :: ice_loc ! 1 si glace, 0 sinon
854  integer,dimension(nx_loc,ny_loc), intent(in) :: mk_loc ! 0 si terre, 1 sinon
855  real,dimension(nx_loc,ny_loc), intent(in) :: orog_lake_loc ! orographie a la surface des lacs
856  integer,dimension(80,2),intent(in) :: inclake_loc ! 80 voisins (+4 -4) en i,j
857  integer, dimension(nx_loc,ny_loc),intent(out) :: mask_extlake_loc ! extension du lac sous la glace
858  real, dimension(nx_loc,ny_loc),intent(out) :: orog_extlake_loc ! niveau du lac sur les zones glace a proximité (rayon de 4 pts de grille)
859
860  integer, dimension(8) :: ipinctab,jpinctab ! indice des points a scanner autour du point ii jj
861  integer, dimension(nx_loc,ny_loc) :: clake ! point lac avec voisin point englace
862  integer :: i,j,k
863
864! initialisation
865  clake(:,:)=0
866  orog_extlake_loc(:,:)=0.
867  mask_extlake_loc(:,:)=9999
868
869! recherche des points lacs en bordure de calotte => voisin pas lac englace et terre
870  where ((lake_loc(:,:).EQ.1).AND.(((eoshift(lake_loc,shift=-1,dim=1).EQ.0).AND.(eoshift(ice_loc,shift=-1,dim=1).EQ.1).AND.(eoshift(mk_loc,shift=-1,dim=1).EQ.0)) .OR. &
871    ((eoshift(lake_loc,shift=1,dim=1).EQ.0).AND.(eoshift(ice_loc,shift=1,dim=1).EQ.1).AND.(eoshift(mk_loc,shift=1,dim=1).EQ.0)) .OR. &
872    ((eoshift(lake_loc,shift=-1,dim=2).EQ.0).AND.(eoshift(ice_loc,shift=-1,dim=2).EQ.1).AND.(eoshift(mk_loc,shift=-1,dim=2).EQ.0)) .OR. &
873    ((eoshift(lake_loc,shift=1,dim=2).EQ.0).AND.(eoshift(ice_loc,shift=1,dim=2).EQ.1).AND.(eoshift(mk_loc,shift=1,dim=2).EQ.0))))
874    clake(:,:) = 1
875  endwhere
876
877
878
879! point terre et englace : mk=0 et H>0
880  do j=1,ny
881    do i=1,nx
882      if (clake(i,j).EQ.1) then
883        do k=1,80
884          ipinctab(k) = max(1,min(nx,i + inclake_loc(k,1)))
885          jpinctab(k) = max(1,min(ny,j + inclake_loc(k,2)))
886          mask_extlake_loc(ipinctab(k),jpinctab(k)) = 1
887          orog_extlake_loc(ipinctab(k),jpinctab(k))=orog_lake_loc(i,j)
888        enddo
889      endif
890    enddo
891  enddo
892
893end subroutine extlake
894
895end module lake_rsl
896
897
898
899
900
901
902
903
904
Note: See TracBrowser for help on using the repository browser.