[215] | 1 | module climat_forcage_mod |
---|
| 2 | |
---|
| 3 | ! Forcage climatique avec climat GCM et fichier de forcage temporel index |
---|
| 4 | ! lecture de fichiers Snapshots climatiques à temps t |
---|
| 5 | ! lecture d'un fichier de forcage (temp carotte de glace) |
---|
| 6 | ! possibilite d'utiliser autant de snapshots que l'on veut (ntr) |
---|
| 7 | ! fichiers de forcage, Snapshot et valeurs de lapse rate definis dans fichier param |
---|
| 8 | use module3d_phy,only:nx,ny,sealevel,slv,S,S0,Tmois,Tann,Tjuly,acc,num_forc,num_param,num_rep_42,coefbmshelf,PI,ro,time,dirforcage,dirnameinp |
---|
| 9 | use netcdf |
---|
| 10 | use io_netcdf_grisli |
---|
| 11 | |
---|
| 12 | |
---|
| 13 | implicit none |
---|
| 14 | character(len=80) :: filin |
---|
| 15 | integer :: nft ! nft est le nombre de lignes a lire dans le fichier |
---|
| 16 | ! contenant le forcage climatique |
---|
| 17 | real,dimension(:),allocatable :: tdate ! time for climate forcing |
---|
| 18 | real,dimension(:),allocatable :: alphat |
---|
| 19 | real,dimension(:),allocatable :: alphap |
---|
| 20 | real,dimension(:),allocatable :: spert |
---|
| 21 | |
---|
| 22 | |
---|
| 23 | real :: mincoefbmelt ! butoirs mini |
---|
| 24 | real :: maxcoefbmelt ! butoirs maxi de coefbmshelf |
---|
| 25 | character(len=100) :: clim_ref_file ! climat de reference |
---|
| 26 | character(len=100) :: forcage_file1 ! fichier de forcage 1 (climat+topo) |
---|
| 27 | character(len=100) :: forcage_file2 ! fichier de forcage 2 (climat+topo) |
---|
| 28 | |
---|
| 29 | character(len=80) :: filforc ! nom du fichier forcage |
---|
| 30 | |
---|
| 31 | |
---|
| 32 | integer,parameter :: ntr=2 ! nb of snapshot files now explicitely specified |
---|
| 33 | character(len=200) :: filtr(ntr) ! fichiers de forcage |
---|
| 34 | |
---|
| 35 | real,dimension(ntr) :: ttr !< date des tranches (snapshots) |
---|
| 36 | real,dimension(nx,ny,ntr) :: delta, deltj, rapact |
---|
| 37 | real,dimension(nx,ny) :: delTatime, delTjtime, rapactime |
---|
| 38 | real,dimension(nx,ny) :: Zs !< surface topography above sea level |
---|
| 39 | |
---|
[220] | 40 | integer :: typerun ! 1 for steady state using snap 1, |
---|
| 41 | ! 2 for steady state using snap 2 ... |
---|
| 42 | ! 0 for transient |
---|
[215] | 43 | |
---|
| 44 | ! S0 topo de ref correspondant a la climato est lu dans lect_topo |
---|
| 45 | |
---|
| 46 | ! variables pour climato mensuelle : |
---|
| 47 | real,dimension(nx,ny,12) :: t2m0 !< initial monthly air temperature (climato) |
---|
| 48 | real,dimension(nx,ny,12) :: pr0 !< initial monthly precip (climato) |
---|
| 49 | |
---|
| 50 | ! climat des snapshots : |
---|
| 51 | real,dimension(nx,ny,12,ntr) :: t2m_ss !< t2m monthly for every GCM snapshot |
---|
| 52 | real,dimension(nx,ny,12,ntr) :: pr_ss !< precip monthly for every GCM snapshot |
---|
| 53 | real,dimension(nx,ny,ntr) :: orog_ss !< topo for every GCM snapshot |
---|
| 54 | |
---|
| 55 | real,dimension(nx,ny,12,ntr) :: t2m_sstopo0 !< t2m GCM snapshot on S0 |
---|
| 56 | real,dimension(nx,ny,12,ntr) :: pr_sstopo0 !< pr GCM snapshot on S0 |
---|
| 57 | |
---|
| 58 | real :: PYY !< constante pour calcul de la fraction solide des precips |
---|
| 59 | real :: psolid=2. !< temp limit between liquid and solid precip |
---|
| 60 | real,dimension(nx,ny) :: FT !< proportion of year with air temp below psolid |
---|
| 61 | real :: lapserate !< gradient vertical de temp annuel et juillet |
---|
| 62 | real :: rappact !< coef pour calcul du rapport accumulation |
---|
| 63 | logical :: annual !< True si forcage annuel, False si forcage mensuel |
---|
| 64 | |
---|
| 65 | contains |
---|
| 66 | |
---|
| 67 | !------------------------------------------------------------------------------------------ |
---|
| 68 | subroutine input_clim !routine qui permet d'initialiser les variables climatiques |
---|
| 69 | |
---|
| 70 | implicit none |
---|
| 71 | character(len=8) :: control !< label to check clim. forc. file (filin) is usable |
---|
| 72 | integer :: l !< in snapshot files:the first column is the mask, read but not used |
---|
| 73 | integer :: err !< recuperation erreur |
---|
| 74 | integer :: i,j,k,m |
---|
| 75 | real*8, dimension(:,:), pointer :: tab2d => null() !< tableau 2d pointer |
---|
| 76 | real*8, dimension(:,:,:), pointer :: tab3d => null() !< tableau 3d pointer |
---|
| 77 | |
---|
| 78 | deltatime(:,:)=0. |
---|
| 79 | deltjtime(:,:)=0. |
---|
| 80 | rapactime(:,:)=1. |
---|
| 81 | annual=.true. ! forcage avec fichier Tann et Tjuly |
---|
| 82 | |
---|
| 83 | ! lecture du climat de reference : climato mensuelle |
---|
| 84 | print*,'lecture du climat de reference: ',trim(DIRNAMEINP)//TRIM(clim_ref_file) |
---|
| 85 | call Read_Ncdf_var('precip',trim(DIRNAMEINP)//TRIM(clim_ref_file),tab3d) |
---|
| 86 | pr0(:,:,:)=tab3d(:,:,:) |
---|
| 87 | call Read_Ncdf_var('t2m',trim(DIRNAMEINP)//TRIM(clim_ref_file),tab3d) |
---|
| 88 | t2m0(:,:,:)=tab3d(:,:,:) |
---|
| 89 | |
---|
| 90 | filtr(1)=TRIM(DIRNAMEINP)//'Snapshots-GCM/'//TRIM(forcage_file1) |
---|
| 91 | filtr(2)=TRIM(DIRNAMEINP)//'Snapshots-GCM/'//TRIM(forcage_file2) |
---|
| 92 | |
---|
| 93 | ! fichiers donnant l'evolution temporelle |
---|
| 94 | ! ----------------------------------------- |
---|
| 95 | ! Pour rappel hemin40 : |
---|
| 96 | ! fichier forcage : signal-cycle-hemin.dat |
---|
| 97 | ! Pour anteis : |
---|
| 98 | ! fichier forcage : forcmodif4cycles.dat |
---|
| 99 | filin=trim(dirforcage)//trim(filforc) |
---|
| 100 | |
---|
| 101 | |
---|
| 102 | ! lecture des fichiers snapshots pour n'importe quel geoplace |
---|
| 103 | ! ------------------------------------------------------------- |
---|
| 104 | do k=1,ntr |
---|
| 105 | write(6,*) 'Read climate file :',trim(filtr(k)) |
---|
| 106 | call Read_Ncdf_var('pr',trim(filtr(k)),tab3d) |
---|
| 107 | pr_ss(:,:,:,k)=tab3d(:,:,:) |
---|
| 108 | call Read_Ncdf_var('tas',trim(filtr(k)),tab3d) |
---|
| 109 | t2m_ss(:,:,:,k)=tab3d(:,:,:) |
---|
| 110 | call Read_Ncdf_var('orog',trim(filtr(k)),tab2d) |
---|
| 111 | orog_ss(:,:,k)=tab2d(:,:) |
---|
| 112 | !~ read(num_forc,*) ttr(k) |
---|
| 113 | !~ read(num_forc,*) |
---|
| 114 | !~ read(num_forc,*) |
---|
| 115 | !~ do j=1,ny |
---|
| 116 | !~ do i=1,nx |
---|
| 117 | !~ read(num_forc,*) l,rapact(i,j,k),delta(i,j,k),deltj(i,j,k) |
---|
| 118 | !~ end do |
---|
| 119 | !~ end do |
---|
| 120 | !~ close(num_forc) |
---|
| 121 | end do |
---|
| 122 | |
---|
| 123 | |
---|
| 124 | ! calcul de t2m et pr sur la topo des obs : |
---|
| 125 | do m=1,12 |
---|
| 126 | do k=1,ntr |
---|
[219] | 127 | t2m_sstopo0(:,:,m,k)=t2m_ss(:,:,m,k)-lapserate*(S0(:,:)-orog_ss(:,:,k)) |
---|
[215] | 128 | pr_sstopo0(:,:,m,k)=pr_ss(:,:,m,k)*exp(rappact*(t2m_sstopo0(:,:,m,k)-t2m_ss(:,:,m,k))) |
---|
| 129 | enddo |
---|
| 130 | enddo |
---|
| 131 | |
---|
| 132 | |
---|
| 133 | ! lecture des fichiers d'evolution pour tout geoplace |
---|
| 134 | ! ---------------------------------------------------- |
---|
| 135 | open(num_forc,file=filin,status='old') |
---|
| 136 | ! print*,nft |
---|
| 137 | read(num_forc,*) control,nft |
---|
| 138 | print*,'control',control,nft |
---|
| 139 | ! determination of file size (line nb), allocation of perturbation array |
---|
| 140 | |
---|
| 141 | if (control.ne.'nb_lines') then |
---|
| 142 | write(6,*) filin,'indiquer le nb de ligne en debut de fichier:' |
---|
| 143 | write(6,*) 'le nb de lignes et le label de control nb_lines' |
---|
| 144 | stop |
---|
| 145 | endif |
---|
| 146 | |
---|
| 147 | if (.not.allocated(tdate)) THEN |
---|
| 148 | allocate(tdate(nft),stat=err) |
---|
| 149 | if (err/=0) then |
---|
| 150 | print *,"erreur à l'allocation du tableau tdate",err |
---|
| 151 | stop 4 |
---|
| 152 | end if |
---|
| 153 | end if |
---|
| 154 | |
---|
| 155 | if (.not.allocated(alphat)) THEN |
---|
| 156 | allocate(alphat(nft),stat=err) |
---|
| 157 | if (err/=0) then |
---|
| 158 | print *,"erreur à l'allocation du tableau tpert",err |
---|
| 159 | stop 4 |
---|
| 160 | end if |
---|
| 161 | end if |
---|
| 162 | |
---|
| 163 | if (.not.allocated(alphap)) THEN |
---|
| 164 | allocate(alphap(nft),stat=err) |
---|
| 165 | if (err/=0) then |
---|
| 166 | print *,"erreur à l'allocation du tableau tpert",err |
---|
| 167 | stop 4 |
---|
| 168 | end if |
---|
| 169 | end if |
---|
| 170 | |
---|
| 171 | if (.not.allocated(spert)) THEN |
---|
| 172 | allocate(spert(nft),stat=err) |
---|
| 173 | if (err/=0) then |
---|
| 174 | print *,"erreur à l'allocation du tableau spert",err |
---|
| 175 | stop 4 |
---|
| 176 | end if |
---|
| 177 | end if |
---|
| 178 | |
---|
| 179 | do i=1,nft |
---|
| 180 | read(num_forc,*) tdate(i),alphat(i),alphap(i),spert(i) |
---|
| 181 | end do |
---|
| 182 | |
---|
| 183 | close(num_forc) |
---|
| 184 | |
---|
| 185 | |
---|
| 186 | end subroutine input_clim |
---|
| 187 | |
---|
| 188 | |
---|
| 189 | !-------------------------------------------------------------------------------- |
---|
| 190 | subroutine init_forclim |
---|
| 191 | |
---|
| 192 | |
---|
[220] | 193 | namelist/clim_forcage/clim_ref_file,ttr,forcage_file1,forcage_file2,typerun,rappact,mincoefbmelt,maxcoefbmelt,filforc,lapserate |
---|
[215] | 194 | |
---|
| 195 | rewind(num_param) ! pour revenir au debut du fichier param_list.dat |
---|
| 196 | read(num_param,clim_forcage) |
---|
| 197 | |
---|
| 198 | |
---|
| 199 | ! formats pour les ecritures dans 42 |
---|
| 200 | 428 format(A) |
---|
| 201 | |
---|
| 202 | write(num_rep_42,428)'!___________________________________________________________' |
---|
| 203 | write(num_rep_42,428) '&clim_forcage ! nom du bloc ' |
---|
| 204 | write(num_rep_42,*) |
---|
| 205 | write(num_rep_42,'(A,A)') 'clim_ref_file = ', clim_ref_file |
---|
| 206 | write(num_rep_42,'(A,A)') 'ttr = ', ttr(:) |
---|
| 207 | write(num_rep_42,'(A,A)') 'forcage_file1 = ', forcage_file1 |
---|
| 208 | write(num_rep_42,'(A,A)') 'forcage_file2 = ', forcage_file2 |
---|
[220] | 209 | write(num_rep_42,'(A,A)') 'typerun = ', typerun |
---|
[215] | 210 | write(num_rep_42,'(A,A)') 'rappact = ', rappact |
---|
| 211 | write(num_rep_42,*) 'mincoefbmelt = ', mincoefbmelt |
---|
| 212 | write(num_rep_42,*) 'maxcoefbmelt = ', maxcoefbmelt |
---|
| 213 | write(num_rep_42,'(A,A)') ' filforc = ', filforc |
---|
| 214 | write(num_rep_42,*) 'lapserate = ', lapserate |
---|
| 215 | write(num_rep_42,*)'/' |
---|
| 216 | write(num_rep_42,*) |
---|
| 217 | |
---|
| 218 | PYY=2.*PI/50. |
---|
| 219 | |
---|
| 220 | end subroutine init_forclim |
---|
| 221 | |
---|
| 222 | |
---|
| 223 | !--------------------------------------------------------------------- |
---|
| 224 | !forcage climatique au cours du temps |
---|
| 225 | |
---|
| 226 | subroutine forclim |
---|
| 227 | |
---|
| 228 | !$ USE OMP_LIB |
---|
| 229 | |
---|
| 230 | implicit none |
---|
| 231 | real :: coeft,coefp,coeftime !< coeftime is not used |
---|
| 232 | integer :: l !< dumm index for loops on snapsots files l=itr,ntr-1 |
---|
| 233 | integer :: itr !< nb of the current snapshot file (change with time) |
---|
| 234 | integer :: i,j,k,m |
---|
| 235 | integer :: ift !< indice correspondant au pas de temps du fichier de forcage |
---|
| 236 | real :: TEMP !< calcul du nbr de jour < psolid |
---|
| 237 | real,dimension(nx,ny) :: deltaZs ! diff temp par rapport a topo S0 |
---|
| 238 | real,dimension(nx,ny,12) :: prmois ! precip sur topo Zs |
---|
| 239 | |
---|
[220] | 240 | real (kind=kind(0.d0)) :: timeloc |
---|
[215] | 241 | |
---|
[220] | 242 | if (typerun.eq.0) then |
---|
| 243 | timeloc = time |
---|
| 244 | else if ( (typerun .gt. 0) .and. (typerun .le. ubound(ttr,dim=1)) ) then |
---|
| 245 | timeloc = ttr(typerun) |
---|
| 246 | else |
---|
| 247 | write(*,*) "Unknown type of simulation for climat_forcage, abort" |
---|
| 248 | STOP |
---|
| 249 | end if |
---|
| 250 | |
---|
| 251 | if(timeloc.le.tdate(1)) then ! time avant le debut du fichier forcage |
---|
[215] | 252 | sealevel=spert(1) |
---|
| 253 | coeft=alphat(1) |
---|
| 254 | coefp=alphap(1) |
---|
| 255 | ift=1 |
---|
| 256 | |
---|
[220] | 257 | else if (timeloc.ge.tdate(nft)) then ! time apres la fin du fichier forcage |
---|
[215] | 258 | sealevel=spert(nft) |
---|
| 259 | coeft=alphat(nft) |
---|
| 260 | coefp=alphap(nft) |
---|
| 261 | ift=nft |
---|
| 262 | |
---|
| 263 | else ! time en general |
---|
| 264 | |
---|
| 265 | ift = 1 |
---|
| 266 | do i=ift,nft-1 ! interpolation entre i et i+1 |
---|
| 267 | |
---|
[220] | 268 | if((timeloc.ge.tdate(i)).and.(timeloc.lt.tdate(i+1))) then ! cas general |
---|
[215] | 269 | sealevel=spert(i)+(spert(i+1)-spert(i))* & |
---|
[220] | 270 | (timeloc-tdate(i))/(tdate(i+1)-tdate(i)) |
---|
[215] | 271 | |
---|
| 272 | coeft=alphat(i)+(alphat(i+1)-alphat(i))* & |
---|
[220] | 273 | (timeloc-tdate(i))/(tdate(i+1)-tdate(i)) |
---|
[215] | 274 | |
---|
| 275 | coefp=alphap(i)+(alphap(i+1)-alphap(i))* & |
---|
[220] | 276 | (timeloc-tdate(i))/(tdate(i+1)-tdate(i)) |
---|
[215] | 277 | |
---|
| 278 | ift=i |
---|
| 279 | goto 100 |
---|
| 280 | endif |
---|
| 281 | end do |
---|
| 282 | endif |
---|
| 283 | 100 continue |
---|
| 284 | |
---|
| 285 | slv(:,:)=sealevel |
---|
| 286 | |
---|
[220] | 287 | if(timeloc.le.ttr(1)) then ! time en dehors des limites du fichier forcage |
---|
[215] | 288 | coeftime=0. |
---|
| 289 | !~ do j=1,ny |
---|
| 290 | !~ do i=1,nx |
---|
| 291 | !~ deltatime(i,j)=delta(i,j,1) |
---|
| 292 | !~ deltjtime(i,j)=deltj(i,j,1) |
---|
| 293 | !~ rapactime(i,j)=rapact(i,j,1) |
---|
| 294 | !~ end do |
---|
| 295 | !~ end do |
---|
| 296 | itr=1 |
---|
[220] | 297 | else if (timeloc.ge.ttr(ntr)) then ! mis a 0 |
---|
[215] | 298 | coeftime=1. |
---|
| 299 | itr=ntr |
---|
| 300 | !~ do j=1,ny |
---|
| 301 | !~ do i=1,nx |
---|
| 302 | !~ deltatime(i,j)=delta(i,j,ntr) |
---|
| 303 | !~ deltjtime(i,j)=deltj(i,j,ntr) |
---|
| 304 | !~ rapactime(i,j)=rapact(i,j,ntr) |
---|
| 305 | !~ end do |
---|
| 306 | !~ end do |
---|
| 307 | else ! interpolation entre l et l+1 : cas general |
---|
| 308 | itr=1 |
---|
| 309 | do l=itr,ntr-1 ! interpolation entre i et i+1 |
---|
[220] | 310 | if((timeloc.ge.ttr(l)).and.(timeloc.lt.ttr(l+1))) then ! test tranche |
---|
| 311 | coeftime= (timeloc-ttr(l))/(ttr(l+1)-ttr(l)) |
---|
[215] | 312 | ! do j=1,ny |
---|
| 313 | ! do i=1,nx |
---|
| 314 | ! delTatime(i,j)=delta(i,j,l)+ & |
---|
| 315 | ! (delta(i,j,l+1)-delta(i,j,l))*coeft |
---|
| 316 | ! delTjtime(i,j)=deltj(i,j,l)+ & |
---|
| 317 | ! (deltj(i,j,l+1)-deltj(i,j,l))*coeft |
---|
| 318 | ! rapactime(i,j)=rapact(i,j,l)+ & |
---|
| 319 | ! (rapact(i,j,l+1)-rapact(i,j,l))*coefp |
---|
| 320 | ! end do |
---|
| 321 | ! end do |
---|
| 322 | itr=l |
---|
| 323 | endif ! fin du test sur la tranche |
---|
| 324 | end do |
---|
| 325 | endif ! fin du test avant,apres,milieu |
---|
| 326 | |
---|
| 327 | |
---|
| 328 | ! test avec coefbmshelf fonction de coeft (1 a l'actuel et proche de 0 en glaciaire) |
---|
| 329 | coefbmshelf=coeft |
---|
| 330 | !coefbmshelf=1. ! modif pour test de sensibilite (tof 20 avril 01) |
---|
| 331 | coefbmshelf=max(coefbmshelf,mincoefbmelt) |
---|
| 332 | coefbmshelf=min(coefbmshelf,maxcoefbmelt) |
---|
| 333 | |
---|
| 334 | !if (geoplace(1:5).eq.'hemin') then |
---|
| 335 | |
---|
| 336 | !$OMP PARALLEL |
---|
| 337 | ! calcul de Tjuly et et Tann |
---|
| 338 | !$OMP WORKSHARE |
---|
| 339 | Zs(:,:)=max(slv(:,:),S(:,:)) |
---|
| 340 | deltaZs(:,:)= -lapserate*(Zs(:,:)-S0(:,:)) ! difference de temperature entre Zs et S0 |
---|
| 341 | !$OMP END WORKSHARE |
---|
| 342 | ! correction d'altitude |
---|
| 343 | do m=1,12 |
---|
| 344 | ! Temp et precip fonction de coeftime : |
---|
| 345 | if (itr.LT.ntr) then |
---|
| 346 | !$OMP WORKSHARE |
---|
| 347 | Tmois(:,:,m)= (t2m_sstopo0(:,:,m,itr) - t2m_sstopo0(:,:,m,itr+1))*(1-coeftime) + t2m0(:,:,m) + deltaZs(:,:) |
---|
| 348 | prmois(:,:,m)= pr0(:,:,m)*(coeftime+(1-coeftime)*(pr_sstopo0(:,:,m,itr)/pr_sstopo0(:,:,m,itr+1))) * exp(rappact*(deltaZs(:,:))) |
---|
| 349 | !$OMP END WORKSHARE |
---|
| 350 | else ! itr=ntr (time>tsnapmax) |
---|
| 351 | !$OMP WORKSHARE |
---|
| 352 | Tmois(:,:,m)= t2m0(:,:,m) + deltaZs(:,:) |
---|
| 353 | prmois(:,:,m)= pr0(:,:,m)*exp(rappact*(deltaZs(:,:))) |
---|
| 354 | !$OMP END WORKSHARE |
---|
| 355 | endif |
---|
| 356 | enddo |
---|
| 357 | |
---|
| 358 | |
---|
| 359 | |
---|
| 360 | !$OMP WORKSHARE |
---|
| 361 | Tann=sum(Tmois,dim=3)/12. |
---|
| 362 | Tjuly=Tmois(:,:,7) |
---|
| 363 | |
---|
| 364 | |
---|
| 365 | |
---|
| 366 | |
---|
| 367 | ! calcul de l'accumulation de neige : |
---|
| 368 | acc(:,:)=sum(prmois,dim=3,mask=Tmois < psolid) |
---|
| 369 | acc(:,:)=acc(:,:)*1000./ro ! conversion en glace |
---|
| 370 | !$OMP END WORKSHARE |
---|
| 371 | !$OMP END PARALLEL |
---|
| 372 | !debug |
---|
| 373 | !print*,'Tann, Tjuly',Tann(142,14),Tjuly(142,14) |
---|
| 374 | !print*,'acc',acc(142,14) |
---|
| 375 | |
---|
| 376 | end subroutine forclim |
---|
| 377 | |
---|
| 378 | end module climat_forcage_mod |
---|