module climat_forcage_mod ! Forcage climatique avec climat GCM et fichier de forcage temporel index ! lecture de fichiers Snapshots climatiques à temps t ! lecture d'un fichier de forcage (temp carotte de glace) ! possibilite d'utiliser autant de snapshots que l'on veut (ntr) ! fichiers de forcage, Snapshot et valeurs de lapse rate definis dans fichier param 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 use netcdf use io_netcdf_grisli implicit none character(len=80) :: filin integer :: nft ! nft est le nombre de lignes a lire dans le fichier ! contenant le forcage climatique real,dimension(:),allocatable :: tdate ! time for climate forcing real,dimension(:),allocatable :: alphat real,dimension(:),allocatable :: alphap real,dimension(:),allocatable :: spert real :: mincoefbmelt ! butoirs mini real :: maxcoefbmelt ! butoirs maxi de coefbmshelf character(len=100) :: clim_ref_file ! climat de reference character(len=100) :: forcage_file1 ! fichier de forcage 1 (climat+topo) character(len=100) :: forcage_file2 ! fichier de forcage 2 (climat+topo) character(len=80) :: filforc ! nom du fichier forcage integer,parameter :: ntr=2 ! nb of snapshot files now explicitely specified character(len=200) :: filtr(ntr) ! fichiers de forcage real,dimension(ntr) :: ttr !< date des tranches (snapshots) real,dimension(nx,ny,ntr) :: delta, deltj, rapact real,dimension(nx,ny) :: delTatime, delTjtime, rapactime real,dimension(nx,ny) :: Zs !< surface topography above sea level integer :: typerun ! 1 for steady state using snap 1, ! 2 for steady state using snap 2 ... ! 0 for transient ! S0 topo de ref correspondant a la climato est lu dans lect_topo ! variables pour climato mensuelle : real,dimension(nx,ny,12) :: t2m0 !< initial monthly air temperature (climato) real,dimension(nx,ny,12) :: pr0 !< initial monthly precip (climato) ! climat des snapshots : real,dimension(nx,ny,12,ntr) :: t2m_ss !< t2m monthly for every GCM snapshot real,dimension(nx,ny,12,ntr) :: pr_ss !< precip monthly for every GCM snapshot real,dimension(nx,ny,ntr) :: orog_ss !< topo for every GCM snapshot real,dimension(nx,ny,12,ntr) :: t2m_sstopo0 !< t2m GCM snapshot on S0 real,dimension(nx,ny,12,ntr) :: pr_sstopo0 !< pr GCM snapshot on S0 real :: PYY !< constante pour calcul de la fraction solide des precips real :: psolid=2. !< temp limit between liquid and solid precip real,dimension(nx,ny) :: FT !< proportion of year with air temp below psolid real :: lapserate !< gradient vertical de temp annuel et juillet real :: rappact !< coef pour calcul du rapport accumulation logical :: annual !< True si forcage annuel, False si forcage mensuel contains !------------------------------------------------------------------------------------------ subroutine input_clim !routine qui permet d'initialiser les variables climatiques implicit none character(len=8) :: control !< label to check clim. forc. file (filin) is usable integer :: l !< in snapshot files:the first column is the mask, read but not used integer :: err !< recuperation erreur integer :: i,j,k,m real*8, dimension(:,:), pointer :: tab2d => null() !< tableau 2d pointer real*8, dimension(:,:,:), pointer :: tab3d => null() !< tableau 3d pointer deltatime(:,:)=0. deltjtime(:,:)=0. rapactime(:,:)=1. annual=.true. ! forcage avec fichier Tann et Tjuly ! lecture du climat de reference : climato mensuelle print*,'lecture du climat de reference: ',trim(DIRNAMEINP)//TRIM(clim_ref_file) call Read_Ncdf_var('precip',trim(DIRNAMEINP)//TRIM(clim_ref_file),tab3d) pr0(:,:,:)=tab3d(:,:,:) call Read_Ncdf_var('t2m',trim(DIRNAMEINP)//TRIM(clim_ref_file),tab3d) t2m0(:,:,:)=tab3d(:,:,:) filtr(1)=TRIM(DIRNAMEINP)//'Snapshots-GCM/'//TRIM(forcage_file1) filtr(2)=TRIM(DIRNAMEINP)//'Snapshots-GCM/'//TRIM(forcage_file2) ! fichiers donnant l'evolution temporelle ! ----------------------------------------- ! Pour rappel hemin40 : ! fichier forcage : signal-cycle-hemin.dat ! Pour anteis : ! fichier forcage : forcmodif4cycles.dat filin=trim(dirforcage)//trim(filforc) ! lecture des fichiers snapshots pour n'importe quel geoplace ! ------------------------------------------------------------- do k=1,ntr write(6,*) 'Read climate file :',trim(filtr(k)) call Read_Ncdf_var('pr',trim(filtr(k)),tab3d) pr_ss(:,:,:,k)=tab3d(:,:,:) call Read_Ncdf_var('tas',trim(filtr(k)),tab3d) t2m_ss(:,:,:,k)=tab3d(:,:,:) call Read_Ncdf_var('orog',trim(filtr(k)),tab2d) orog_ss(:,:,k)=tab2d(:,:) !~ read(num_forc,*) ttr(k) !~ read(num_forc,*) !~ read(num_forc,*) !~ do j=1,ny !~ do i=1,nx !~ read(num_forc,*) l,rapact(i,j,k),delta(i,j,k),deltj(i,j,k) !~ end do !~ end do !~ close(num_forc) end do ! calcul de t2m et pr sur la topo des obs : do m=1,12 do k=1,ntr t2m_sstopo0(:,:,m,k)=t2m_ss(:,:,m,k)-lapserate*(S0(:,:)-orog_ss(:,:,k)) pr_sstopo0(:,:,m,k)=pr_ss(:,:,m,k)*exp(rappact*(t2m_sstopo0(:,:,m,k)-t2m_ss(:,:,m,k))) enddo enddo ! lecture des fichiers d'evolution pour tout geoplace ! ---------------------------------------------------- open(num_forc,file=filin,status='old') ! print*,nft read(num_forc,*) control,nft print*,'control',control,nft ! determination of file size (line nb), allocation of perturbation array if (control.ne.'nb_lines') then write(6,*) filin,'indiquer le nb de ligne en debut de fichier:' write(6,*) 'le nb de lignes et le label de control nb_lines' stop endif if (.not.allocated(tdate)) THEN allocate(tdate(nft),stat=err) if (err/=0) then print *,"erreur à l'allocation du tableau tdate",err stop 4 end if end if if (.not.allocated(alphat)) THEN allocate(alphat(nft),stat=err) if (err/=0) then print *,"erreur à l'allocation du tableau tpert",err stop 4 end if end if if (.not.allocated(alphap)) THEN allocate(alphap(nft),stat=err) if (err/=0) then print *,"erreur à l'allocation du tableau tpert",err stop 4 end if end if if (.not.allocated(spert)) THEN allocate(spert(nft),stat=err) if (err/=0) then print *,"erreur à l'allocation du tableau spert",err stop 4 end if end if do i=1,nft read(num_forc,*) tdate(i),alphat(i),alphap(i),spert(i) end do close(num_forc) end subroutine input_clim !-------------------------------------------------------------------------------- subroutine init_forclim namelist/clim_forcage/clim_ref_file,ttr,forcage_file1,forcage_file2,typerun,rappact,mincoefbmelt,maxcoefbmelt,filforc,lapserate rewind(num_param) ! pour revenir au debut du fichier param_list.dat read(num_param,clim_forcage) ! formats pour les ecritures dans 42 428 format(A) write(num_rep_42,428)'!___________________________________________________________' write(num_rep_42,428) '&clim_forcage ! nom du bloc ' write(num_rep_42,*) write(num_rep_42,'(A,A)') 'clim_ref_file = ', clim_ref_file write(num_rep_42,'(A,A)') 'ttr = ', ttr(:) write(num_rep_42,'(A,A)') 'forcage_file1 = ', forcage_file1 write(num_rep_42,'(A,A)') 'forcage_file2 = ', forcage_file2 write(num_rep_42,'(A,A)') 'typerun = ', typerun write(num_rep_42,'(A,A)') 'rappact = ', rappact write(num_rep_42,*) 'mincoefbmelt = ', mincoefbmelt write(num_rep_42,*) 'maxcoefbmelt = ', maxcoefbmelt write(num_rep_42,'(A,A)') ' filforc = ', filforc write(num_rep_42,*) 'lapserate = ', lapserate write(num_rep_42,*)'/' write(num_rep_42,*) PYY=2.*PI/50. end subroutine init_forclim !--------------------------------------------------------------------- !forcage climatique au cours du temps subroutine forclim !$ USE OMP_LIB implicit none real :: coeft,coefp,coeftime !< coeftime is not used integer :: l !< dumm index for loops on snapsots files l=itr,ntr-1 integer :: itr !< nb of the current snapshot file (change with time) integer :: i,j,k,m integer :: ift !< indice correspondant au pas de temps du fichier de forcage real :: TEMP !< calcul du nbr de jour < psolid real,dimension(nx,ny) :: deltaZs ! diff temp par rapport a topo S0 real,dimension(nx,ny,12) :: prmois ! precip sur topo Zs real (kind=kind(0.d0)) :: timeloc if (typerun.eq.0) then timeloc = time else if ( (typerun .gt. 0) .and. (typerun .le. ubound(ttr,dim=1)) ) then timeloc = ttr(typerun) else write(*,*) "Unknown type of simulation for climat_forcage, abort" STOP end if if(timeloc.le.tdate(1)) then ! time avant le debut du fichier forcage sealevel=spert(1) coeft=alphat(1) coefp=alphap(1) ift=1 else if (timeloc.ge.tdate(nft)) then ! time apres la fin du fichier forcage sealevel=spert(nft) coeft=alphat(nft) coefp=alphap(nft) ift=nft else ! time en general ift = 1 do i=ift,nft-1 ! interpolation entre i et i+1 if((timeloc.ge.tdate(i)).and.(timeloc.lt.tdate(i+1))) then ! cas general sealevel=spert(i)+(spert(i+1)-spert(i))* & (timeloc-tdate(i))/(tdate(i+1)-tdate(i)) coeft=alphat(i)+(alphat(i+1)-alphat(i))* & (timeloc-tdate(i))/(tdate(i+1)-tdate(i)) coefp=alphap(i)+(alphap(i+1)-alphap(i))* & (timeloc-tdate(i))/(tdate(i+1)-tdate(i)) ift=i goto 100 endif end do endif 100 continue slv(:,:)=sealevel if(timeloc.le.ttr(1)) then ! time en dehors des limites du fichier forcage coeftime=0. !~ do j=1,ny !~ do i=1,nx !~ deltatime(i,j)=delta(i,j,1) !~ deltjtime(i,j)=deltj(i,j,1) !~ rapactime(i,j)=rapact(i,j,1) !~ end do !~ end do itr=1 else if (timeloc.ge.ttr(ntr)) then ! mis a 0 coeftime=1. itr=ntr !~ do j=1,ny !~ do i=1,nx !~ deltatime(i,j)=delta(i,j,ntr) !~ deltjtime(i,j)=deltj(i,j,ntr) !~ rapactime(i,j)=rapact(i,j,ntr) !~ end do !~ end do else ! interpolation entre l et l+1 : cas general itr=1 do l=itr,ntr-1 ! interpolation entre i et i+1 if((timeloc.ge.ttr(l)).and.(timeloc.lt.ttr(l+1))) then ! test tranche coeftime= (timeloc-ttr(l))/(ttr(l+1)-ttr(l)) ! do j=1,ny ! do i=1,nx ! delTatime(i,j)=delta(i,j,l)+ & ! (delta(i,j,l+1)-delta(i,j,l))*coeft ! delTjtime(i,j)=deltj(i,j,l)+ & ! (deltj(i,j,l+1)-deltj(i,j,l))*coeft ! rapactime(i,j)=rapact(i,j,l)+ & ! (rapact(i,j,l+1)-rapact(i,j,l))*coefp ! end do ! end do itr=l endif ! fin du test sur la tranche end do endif ! fin du test avant,apres,milieu ! test avec coefbmshelf fonction de coeft (1 a l'actuel et proche de 0 en glaciaire) coefbmshelf=coeft !coefbmshelf=1. ! modif pour test de sensibilite (tof 20 avril 01) coefbmshelf=max(coefbmshelf,mincoefbmelt) coefbmshelf=min(coefbmshelf,maxcoefbmelt) !if (geoplace(1:5).eq.'hemin') then !$OMP PARALLEL ! calcul de Tjuly et et Tann !$OMP WORKSHARE Zs(:,:)=max(slv(:,:),S(:,:)) deltaZs(:,:)= -lapserate*(Zs(:,:)-S0(:,:)) ! difference de temperature entre Zs et S0 !$OMP END WORKSHARE ! correction d'altitude do m=1,12 ! Temp et precip fonction de coeftime : if (itr.LT.ntr) then !$OMP WORKSHARE Tmois(:,:,m)= (t2m_sstopo0(:,:,m,itr) - t2m_sstopo0(:,:,m,itr+1))*(1-coeftime) + t2m0(:,:,m) + deltaZs(:,:) prmois(:,:,m)= pr0(:,:,m)*(coeftime+(1-coeftime)*(pr_sstopo0(:,:,m,itr)/pr_sstopo0(:,:,m,itr+1))) * exp(rappact*(deltaZs(:,:))) !$OMP END WORKSHARE else ! itr=ntr (time>tsnapmax) !$OMP WORKSHARE Tmois(:,:,m)= t2m0(:,:,m) + deltaZs(:,:) prmois(:,:,m)= pr0(:,:,m)*exp(rappact*(deltaZs(:,:))) !$OMP END WORKSHARE endif enddo !$OMP WORKSHARE Tann=sum(Tmois,dim=3)/12. Tjuly=Tmois(:,:,7) ! calcul de l'accumulation de neige : acc(:,:)=sum(prmois,dim=3,mask=Tmois < psolid) acc(:,:)=acc(:,:)*1000./ro ! conversion en glace !$OMP END WORKSHARE !$OMP END PARALLEL !debug !print*,'Tann, Tjuly',Tann(142,14),Tjuly(142,14) !print*,'acc',acc(142,14) end subroutine forclim end module climat_forcage_mod