source: trunk/SOURCES/climat_forcage_mod.f90 @ 224

Last change on this file since 224 was 224, checked in by dumas, 5 years ago

format of the parameter file created by the code, possibility to use it in reading

File size: 12.8 KB
Line 
1module  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
40  integer :: typerun ! 1 for steady state using snap 1,
41                     ! 2 for steady state using snap 2 ...
42                     ! 0 for transient
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
65contains
66
67!------------------------------------------------------------------------------------------
68subroutine 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
127      t2m_sstopo0(:,:,m,k)=t2m_ss(:,:,m,k)-lapserate*(S0(:,:)-orog_ss(:,:,k))
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
186end subroutine input_clim
187
188
189!--------------------------------------------------------------------------------
190subroutine init_forclim
191
192
193  namelist/clim_forcage/clim_ref_file,ttr,forcage_file1,forcage_file2,typerun,rappact,mincoefbmelt,maxcoefbmelt,filforc,lapserate
194
195  rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
196  read(num_param,clim_forcage)
197
198  write(num_rep_42,'(A)')'!___________________________________________________________' 
199  write(num_rep_42,'(A)') '&clim_forcage                                ! nom du bloc '
200  write(num_rep_42,*)
201  write(num_rep_42,'(A,A,A)') 'clim_ref_file = "',trim(clim_ref_file),'"'
202  write(num_rep_42,*) 'ttr = ', ttr
203  write(num_rep_42,'(A,A,A)') 'forcage_file1 = "',trim(forcage_file1),'"'
204  write(num_rep_42,'(A,A,A)') 'forcage_file2 = "',trim(forcage_file2),'"'
205  write(num_rep_42,*) 'typerun = ', typerun
206  write(num_rep_42,*) 'rappact = ', rappact
207  write(num_rep_42,*) 'mincoefbmelt = ', mincoefbmelt
208  write(num_rep_42,*) 'maxcoefbmelt = ', maxcoefbmelt
209  write(num_rep_42,'(A,A,A)') ' filforc      = "',trim(filforc),'"'
210  write(num_rep_42,*) 'lapserate = ', lapserate
211  write(num_rep_42,*)'/'                           
212  write(num_rep_42,*)
213
214  PYY=2.*PI/50.
215
216end subroutine init_forclim
217
218
219!---------------------------------------------------------------------
220!forcage climatique au cours du temps
221
222subroutine forclim
223
224!$ USE OMP_LIB
225
226  implicit none
227  real :: coeft,coefp,coeftime       !< coeftime is not used
228  integer :: l                       !< dumm index for loops on snapsots files  l=itr,ntr-1
229  integer :: itr                     !< nb of the current snapshot file (change with time)
230  integer :: i,j,k,m
231  integer :: ift                     !< indice correspondant au pas de temps du fichier de forcage
232  real :: TEMP                       !< calcul du nbr de jour < psolid
233  real,dimension(nx,ny) :: deltaZs   ! diff temp par rapport a topo S0
234  real,dimension(nx,ny,12) :: prmois ! precip sur topo Zs
235 
236  real (kind=kind(0.d0)) ::  timeloc
237 
238  if (typerun.eq.0) then
239     timeloc = time
240  else if ( (typerun .gt. 0) .and. (typerun .le. ubound(ttr,dim=1)) ) then
241     timeloc = ttr(typerun)
242  else
243     write(*,*) "Unknown type of simulation for climat_forcage, abort"
244     STOP
245  end if
246
247  if(timeloc.le.tdate(1)) then   ! time avant le debut du fichier forcage
248     sealevel=spert(1)
249     coeft=alphat(1)
250     coefp=alphap(1)
251     ift=1
252
253  else if (timeloc.ge.tdate(nft)) then     ! time apres la fin du fichier forcage
254     sealevel=spert(nft)
255     coeft=alphat(nft)
256     coefp=alphap(nft)
257     ift=nft
258
259  else                         ! time en general
260
261     ift = 1     
262     do i=ift,nft-1  ! interpolation entre i et i+1
263
264        if((timeloc.ge.tdate(i)).and.(timeloc.lt.tdate(i+1))) then  ! cas general
265           sealevel=spert(i)+(spert(i+1)-spert(i))*    &
266                (timeloc-tdate(i))/(tdate(i+1)-tdate(i))
267
268           coeft=alphat(i)+(alphat(i+1)-alphat(i))*      &
269                (timeloc-tdate(i))/(tdate(i+1)-tdate(i))
270
271           coefp=alphap(i)+(alphap(i+1)-alphap(i))*      &
272                (timeloc-tdate(i))/(tdate(i+1)-tdate(i))
273               
274           ift=i
275           goto 100
276        endif
277     end do
278  endif
279100 continue
280
281  slv(:,:)=sealevel
282
283  if(timeloc.le.ttr(1)) then   ! time en dehors des limites du fichier forcage
284    coeftime=0.
285!~     do j=1,ny
286!~       do i=1,nx
287!~         deltatime(i,j)=delta(i,j,1)
288!~         deltjtime(i,j)=deltj(i,j,1)
289!~         rapactime(i,j)=rapact(i,j,1)
290!~       end do
291!~     end do
292    itr=1
293  else if (timeloc.ge.ttr(ntr)) then        ! mis a 0
294    coeftime=1.
295    itr=ntr
296!~     do j=1,ny
297!~       do i=1,nx
298!~         deltatime(i,j)=delta(i,j,ntr)
299!~         deltjtime(i,j)=deltj(i,j,ntr)
300!~         rapactime(i,j)=rapact(i,j,ntr)
301!~       end do
302!~     end do
303  else               !  interpolation entre l et l+1 : cas general
304     itr=1
305     do l=itr,ntr-1    ! interpolation entre i et i+1
306        if((timeloc.ge.ttr(l)).and.(timeloc.lt.ttr(l+1))) then  ! test tranche
307           coeftime= (timeloc-ttr(l))/(ttr(l+1)-ttr(l))
308!           do j=1,ny
309!              do i=1,nx
310!                 delTatime(i,j)=delta(i,j,l)+             &
311!                      (delta(i,j,l+1)-delta(i,j,l))*coeft
312!                 delTjtime(i,j)=deltj(i,j,l)+             &
313!                      (deltj(i,j,l+1)-deltj(i,j,l))*coeft
314!                 rapactime(i,j)=rapact(i,j,l)+            &
315!                      (rapact(i,j,l+1)-rapact(i,j,l))*coefp
316!              end do
317!           end do
318           itr=l
319        endif        ! fin du test sur la tranche
320     end do
321  endif   ! fin du test avant,apres,milieu
322
323
324! test avec coefbmshelf fonction de coeft (1 a l'actuel et proche de 0 en glaciaire)
325  coefbmshelf=coeft
326!coefbmshelf=1. ! modif pour test de sensibilite (tof 20 avril 01)
327  coefbmshelf=max(coefbmshelf,mincoefbmelt)
328  coefbmshelf=min(coefbmshelf,maxcoefbmelt)
329
330!if (geoplace(1:5).eq.'hemin') then
331
332  !$OMP PARALLEL
333!     calcul de Tjuly et et Tann
334  !$OMP WORKSHARE
335  Zs(:,:)=max(slv(:,:),S(:,:))
336  deltaZs(:,:)= -lapserate*(Zs(:,:)-S0(:,:)) ! difference de temperature entre Zs et S0
337  !$OMP END WORKSHARE
338! correction d'altitude
339do m=1,12 
340! Temp et precip fonction de coeftime :
341  if (itr.LT.ntr) then
342    !$OMP WORKSHARE
343    Tmois(:,:,m)=  (t2m_sstopo0(:,:,m,itr) - t2m_sstopo0(:,:,m,itr+1))*(1-coeftime) + t2m0(:,:,m) + deltaZs(:,:)
344    prmois(:,:,m)= pr0(:,:,m)*(coeftime+(1-coeftime)*(pr_sstopo0(:,:,m,itr)/pr_sstopo0(:,:,m,itr+1))) * exp(rappact*(deltaZs(:,:)))
345    !$OMP END WORKSHARE
346  else ! itr=ntr (time>tsnapmax)
347    !$OMP WORKSHARE
348    Tmois(:,:,m)=  t2m0(:,:,m) + deltaZs(:,:)
349    prmois(:,:,m)= pr0(:,:,m)*exp(rappact*(deltaZs(:,:)))
350    !$OMP END WORKSHARE
351  endif           
352enddo 
353
354
355
356  !$OMP WORKSHARE
357  Tann=sum(Tmois,dim=3)/12.
358  Tjuly=Tmois(:,:,7)
359 
360
361
362
363! calcul de l'accumulation de neige :
364  acc(:,:)=sum(prmois,dim=3,mask=Tmois < psolid)
365  acc(:,:)=acc(:,:)*1000./ro  ! conversion en glace
366  !$OMP END WORKSHARE
367  !$OMP END PARALLEL
368!debug
369!print*,'Tann, Tjuly',Tann(142,14),Tjuly(142,14) 
370!print*,'acc',acc(142,14)
371
372end subroutine forclim
373
374end module  climat_forcage_mod
Note: See TracBrowser for help on using the repository browser.