source: trunk/SOURCES/climat-perturb_mod-0.4.f90 @ 23

Last change on this file since 23 was 4, checked in by dumas, 10 years ago

initial import GRISLI trunk

File size: 6.2 KB
Line 
1!> \file climat-perturb_mod-0.4.f90
2!! Module pour les variations temporelles des variables climatiques
3!<
4
5!> \namespace  climat_perturb_mod
6!!Module pour les variations temporelles des variables climatiques
7!! \author ...
8!! \date ...
9!! @note Used modules:
10!! @note   - use module3D_phy
11!<
12module  climat_perturb_mod
13
14
15use module3d_phy
16implicit none
17
18integer nft             !< nombre de lignes a lire dans le fichier forcage climatique
19real,dimension(:),allocatable :: tdate          !< time for climate forcing
20real,dimension(:),allocatable :: tpert          !< temperature for climate forcing
21real,dimension(:),allocatable :: spert          !< sea surface perturbation
22
23real :: coefT                   !< pour modifier l'amplitude de la perturb. T
24real :: rappact                 !< pour le calcul du rapport 'accumulation
25integer :: retroac              !< 1-> full retroactions accum
26real :: rapbmshelf              !< pour calcul coefbmshelf
27real :: mincoefbmelt            !< butoirs mini
28real :: maxcoefbmelt            !< butoirs maxi de coefbmelt
29character(len=80) :: filforc    !< nom du fichier forcage
30
31! Pour l'instant tafor est global meme si inutile si on utilise
32! un forcage  (variation spatiales)
33
34contains
35
36!--------------------------------------------------------------------------------
37!> SUBROUTINE: input_clim 
38!! Routine qui permet d'initialiser les variations temporelles des variables climatiques
39!>
40  subroutine input_clim !routine qui permet d'initialiser les variations temporelles des variables climatiques
41
42    implicit none
43    character(len=8) :: control      !label to check clim. forc. file (filin) is usable
44    character(len=80):: filin
45
46
47    ! Lecture du forcage
48    !-----------------------
49    ! Le fichier de forcage est lu dans le fichier entree parametres
50
51    filin=trim(dirforcage)//trim(filforc)
52
53    open(num_forc,file=filin,status='old')
54
55    read(num_forc,*) control,nft
56
57    ! Determination of file size (line nb), allocation of perturbation array
58
59    if (control.ne.'nb_lines') then
60       write(6,*) filin,'indiquer le nb de ligne en debut de fichier:'
61       write(6,*) 'le nb de lignes et le label de control nb_lines'
62       stop
63    endif
64
65    ! Dimensionnement des tableaux tdate, ....
66
67    if (.not.allocated(tdate)) then
68       allocate(tdate(nft),stat=err)
69       if (err/=0) then
70          print *,"erreur a l'allocation du tableau Tdate",err
71          stop 4
72       end if
73    end if
74
75    if (.not.allocated(spert)) then
76       allocate(spert(nft),stat=err)
77       if (err/=0) then
78          print *,"erreur a l'allocation du tableau Spert",err
79          stop 4
80       end if
81    end if
82
83
84    if (.not.allocated(tpert)) then
85       allocate(tpert(nft),stat=err)
86       if (err/=0) then
87          print *,"erreur a l'allocation du tableau Tpert",err
88          stop 4
89       end if
90    end if
91
92    do i=1,nft
93       read(num_forc,*) tdate(i),spert(i),tpert(i)
94    end do
95    close(num_forc)
96
97    tpert(:)=tpert(:)*coefT
98
99
100  end subroutine input_clim
101
102!--------------------------------------------------------------------------------
103!> SUBROUTINE: init_forclim
104!! Routine qui permet d'initialiser les variables climatiques au cours du temps
105!>
106subroutine init_forclim
107
108
109namelist/clim_pert/coefT,rappact,retroac,rapbmshelf,mincoefbmelt,maxcoefbmelt,filforc
110
111rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
112read(num_param,clim_pert)
113
114! formats pour les ecritures dans 42
115428 format(A)
116
117write(num_rep_42,428)'!___________________________________________________________' 
118write(num_rep_42,428) '&clim_pert                      ! module climat_perturb_mod '
119write(num_rep_42,*)
120write(num_rep_42,*) 'coefT        = ', coefT 
121write(num_rep_42,*) 'rappact      = ', rappact 
122write(num_rep_42,*) 'rretroac     = ', retroac 
123write(num_rep_42,*) 'rapbmshelf   = ', rapbmshelf 
124write(num_rep_42,*) 'mincoefbmelt = ', mincoefbmelt
125write(num_rep_42,*) 'maxcoefbmelt = ', maxcoefbmelt
126write(num_rep_42,'(A,A)') ' filforc      = ', filforc
127write(num_rep_42,*)'/'                           
128write(num_rep_42,*)
129
130return
131end subroutine init_forclim
132
133!--------------------------------------------------------------------------------
134!> SUBROUTINE: forclim
135!! 
136!!  Routine qui permet le calcule climatiques au cours du temps
137!!  @note Au temps considere (time) attribue les scalaires
138!!  @note   - tafor : forcage en temperature
139!!  @note   - sealevel : forcage niveau des mers
140!!  @note   - coefbmelt : forcage fusion basale ice shelves
141!>
142subroutine forclim               !  au temps considere (time) attribue les scalaires
143  !  tafor : forcage en temperature
144  !  sealevel : forcage niveau des mers
145  !  coefbmelt : forcage fusion basale ice shelves
146
147  use module3d_phy
148  implicit none
149
150  !       time en dehors des limites du fichier forcage
151  if(time.lt.tdate(1)) then
152     tafor=tpert(1)
153     sealevel=spert(1)
154     ift=1
155
156  else if (time.ge.tdate(nft)) then
157     tafor=tpert(nft)
158     sealevel=spert(nft)
159     ift=nft
160
161  else
162     do i=1,nft-1 
163        if((time.ge.tdate(i)).and.(time.lt.tdate(i+1))) then ! entre i et i+1 : cas general
164           tafor=tpert(i)+(tpert(i+1)-tpert(i))*       &
165                (time-tdate(i))/(tdate(i+1)-tdate(i))
166           sealevel=spert(i)+(spert(i+1)-spert(i))*    &
167                (time-tdate(i))/(tdate(i+1)-tdate(i))
168           ift=i
169           goto 100
170        endif
171     end do
172  endif
173100 continue
174
175
176  !  coefmshelf est un coefficient qui fait varier bmgrz et bmshelf
177  !  en fonction de tafor
178
179  !   coefbmshelf=(1.+tafor/10.)           ! coefbmshelf=0 pour tafor=-10deg
180
181  if ((tafor.le.0).and.(tafor.gt.-5.)) then
182     coefbmshelf=(1.+tafor/rapbmshelf)     ! coefbmshelf=0 pour tafor=-7 standard precedent
183
184  else if (tafor.le.-5) then               ! lineaire en 0 a -10 et la valeur a -5
185     coefbmshelf=(1.-5./rapbmshelf)/5.*(tafor+10)
186
187  else
188     coefbmshelf=(1.+5.*tafor/rapbmshelf)  ! 5 fois plus efficace vers le chaud
189  endif
190
191  coefbmshelf=max(coefbmshelf,mincoefbmelt) 
192  coefbmshelf=min(coefbmshelf,maxcoefbmelt)
193
194  !      call massb_anteis_perturb() 
195  !  massb_perturb_Tparam remplace desormais massb_anteis_perturb
196  !  le code est le meme pour l'Antarctique (simplement passe en format libre)
197  !  dans le fichier massb-ant_perturb_Tparam.f90
198
199
200  call massb_perturb_Tparam
201
202
203
204
205
206end subroutine forclim
207
208end module  climat_perturb_mod
Note: See TracBrowser for help on using the repository browser.