source: trunk/SOURCES/GrIce2sea_files/climat_GrIce2sea_years_perturb_mod.f90 @ 66

Last change on this file since 66 was 66, checked in by aquiquet, 8 years ago

Dragging neff parameters in param_list && distinction between actual ice streams and SSA velocity used as Ub in SIA (slow SSA zones) in netcdf outputs

File size: 22.6 KB
Line 
1!> \file climat_GrIce2sea_mod.f90
2! forcage avec BM
3!! Module ou les variations temporelles des variables climatiques
4!! sont directement imposees
5!<
6
7!> \namespace  climat_Grice2sea_mod
8!! Module ou les variations temporelles des variables climatiques
9!! sont directement imposees
10!! \author Cat
11!! \date 31 oct
12!! @note Used modules:
13!! @note   - use module3D_phy
14!<
15module climat_Grice2sea_years_perturb_mod
16
17
18use module3d_phy,only: nx,ny,S,S0,H,Bsoc,acc,abl,BM,Tann,Tjuly,Ts,time,dt,num_param,num_rep_42,num_forc,dirforcage,dirnameinp,tafor,time,sealevel
19!use lect_climref_Ice2sea
20use netcdf
21use io_netcdf_grisli
22implicit none
23
24
25real :: T_lapse_rate              !< pour la temperature
26
27
28! declarations specifiques year by year
29
30real,dimension(:),allocatable        :: time_snap         !> date des snapshots  indice : nb de nb_snap)
31real,dimension(:,:,:),allocatable    :: smb_snap          !> bilan de masse des snapshots indices : nx,ny,nb_snap
32real                                 :: ecart_snap        !> ecart temporel entre les snapshots
33real                                 :: time_depart_snaps !> temps du debut premier snapshot
34integer                              :: nb_snap           !> nombre de snapshots
35
36! declaration pour le bilan de masse
37real,dimension(nx,ny)                :: bm_0              !> bilan de masse initial equ. S_0ref de Tamsin
38real,dimension(nx,ny)                :: bm_time           !> bilan de masse interpole entre 2 snapshots (~ S_t^RCM)
39
40! calcul du gradient
41real,dimension(nx,ny)                :: bm_ref_t          !> bilan de masse de reference au temps t
42real,dimension(nx,ny,10)             :: bm_ref_10         !> tableau des bm_ref pour moyenne 10 ans
43real                                 :: time_prec         !> temps du precedent snapshot
44integer                              :: icum              !> pour le test stockage dans bm_ref_10
45integer                              :: i_moy             !> nombre de snapshots stockes
46
47real,dimension(nx,ny)                :: grad_bm           !> gradient du bilan de masse
48real,dimension(nx,ny)                :: S_ref             !> surface de reference
49
50real                                 :: grad_N_smb_neg    !> SMB < 0 North
51real                                 :: grad_N_smb_pos    !> SMB >= 0 North
52real                                 :: grad_S_smb_neg    !> SMB < 0 South
53real                                 :: grad_S_smb_pos    !> SMB >= 0 South
54real                                 :: coef_smb_unit     ! pour corriger l'unite
55
56real,dimension(nx,ny)                :: TA0               !> Temp annuelle sur S0
57
58
59! aurel, pour climat type perturb:
60integer nft             !> nombre de lignes a lire dans le fichier forcage climatique
61real,dimension(:),allocatable :: tdate          !> time for climate forcing
62real,dimension(:),allocatable :: tpert          !> temperature for climate forcing
63real,dimension(:),allocatable :: spert          !> sea surface perturbation
64real :: coefT                    !> pour modifier l'amplitude de la perturb. T
65character(len=120) :: filforc    !> nom du fichier forcage
66integer :: pertsmb               !> boolean: do we modify the smb?
67real    :: rapsmb                !> if we modify the smb this is the equivalent of rappact
68
69
70
71! pour les lectures ncdf
72real*8, dimension(:,:,:),  pointer   :: tab3D             !< tableau 3d real pointer
73real*8, dimension(:,:),   pointer    :: tab               !< tableau 2d real pointer
74Character(len=10), dimension(3)      :: dimtestname       !> pour la sortie test netcdf
75integer                              :: ncidloc           !> pour la sortie test netcdf
76integer                              :: status            !> pour la sortie test netcdf
77
78integer                              :: massb_time        !< pour selectionner le type de calcul de smb
79! On a deux routines : celle avec un seul fichier (donnees observees) : massb_Ice2sea_fixe
80! Celle avec le bilan de masse sur plusieurs snapshots (annuels par ex.) entre lesquels on interpole.
81
82contains
83
84!--------------------------------------------------------------------------------
85!> SUBROUTINE: input_clim 
86!! Routine qui permet d'initialiser les variations temporelles des variables climatiques
87!>
88!--------------------------------------------------------------------------------
89
90subroutine input_clim 
91
92  character(len=100) :: smb_file           ! fichier smb
93  character(len=100) :: temp_annual_file   ! temperature annuelles
94  character(len=100) :: file_smb_snap      !> nom du fichier dans lequel sont les snapshots
95
96  integer            :: err                ! recuperation d'erreur
97  integer            :: i,j,k
98
99  !aurel for Tafor:
100  character(len=8) :: control      !label to check clim. forc. file (filin) is usable
101  character(len=80):: filin
102
103  namelist/clim_smb_T_gen/smb_file,coef_smb_unit,temp_annual_file
104  namelist/clim_snap/nb_snap,time_depart_snaps,ecart_snap,file_smb_snap,massb_time
105
106428 format(A)
107  rewind(num_param)                     ! pour revenir au debut du fichier param_list.dat
108  read(num_param,clim_smb_T_gen)
109
110  write(num_rep_42,428)'!________________________________________________________________' 
111  write(num_rep_42,428)'!  module climat_Grice2sea_years_mod lecture climat ref          '
112  write(num_rep_42,clim_smb_T_gen)
113  write(num_rep_42,428)'! smb_file          = fichier SMB (kg/m2/an)                     '
114  write(num_rep_42,428)'! coef_smb_unit     = coef passage m glace/an  (1/910 ou 1/918)  '
115  write(num_rep_42,428)'! temp_annual_file  = Temp moy annuelle  (°C)                    '
116  write(num_rep_42,428)'!________________________________________________________________'
117
118 
119! smb : surface mass balance
120  smb_file  = trim(dirnameinp)//trim(smb_file)
121
122  call Read_Ncdf_var('smb',smb_file,tab)
123
124!  call lect_input(3,'smb',1,bm,smb_file,trim(dirnameinp)//trim(runname)//'.nc')
125
126  bm(:,:)  = tab(:,:) * coef_smb_unit
127
128  !where ((H(:,:).lt.1.).and.(Bsoc(:,:).gt.0.))
129  !   bm(:,:) = bm(:,:) - 5.                     ! pour faire un masque a l'exterieur du Groenland actuel
130  !end where
131  where ((H(:,:).lt.1.).and.(Bsoc(:,:).gt.0.).and.(bm(:,:).gt.0.))
132     bm(:,:) = -0.01
133  end where
134!cdc test debug Hemin15 et Greeneem15
135!  where (bm(:,:).lt.-1000) bm(:,:)=0.
136  where (bm(:,:).eq.0) bm(:,:)=-5. !afq
137 
138  acc(:,:) = 0.
139  abl(:,:) = 0.
140
141  where (bm(:,:).gt.0.)
142     acc(:,:) = bm(:,:)   ! accumulation quand positif
143  elsewhere
144     abl(:,:) = - bm(:,:) ! ablation quand negatif
145  end where
146
147
148! surface temperature  Tann
149
150  temp_annual_file = trim(dirnameinp)//trim(temp_annual_file)
151
152  call Read_Ncdf_var('Tann',temp_annual_file,tab)
153  Tann(:,:)  = tab(:,:)
154!  call lect_input(3,'Tann',1,Tann,temp_annual_file,trim(dirnameinp)//trim(runname)//'.nc')
155
156!cdc test debug Hemin15
157!  where (Tann(:,:).lt.-100.) Tann(:,:)=10.
158  where (Tann(:,:).lt.-100.) Tann(:,:)=10. !afq
159
160  ta0(:,:)   = Tann(:,:)
161  Tjuly(:,:) = Tann(:,:)
162
163  rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
164  read(num_param,clim_snap)
165
166! formats pour les ecritures dans 42
167  rewind(num_param)                     ! pour revenir au debut du fichier param_list.dat
168  read(num_param,clim_snap)
169
170  write(num_rep_42,428)'!_______________________________________________________________________' 
171  write(num_rep_42,428)'!  module climat_Grice2sea_years_mod                                    '
172  write(num_rep_42,clim_snap)
173  write(num_rep_42,428)'! nb_snap           = nombre de snapshots                               '
174  write(num_rep_42,428)'! time_depart_snaps = debut du forçage                                  '
175  write(num_rep_42,428)'! ecart_snap        = ecart entre les snapshots                         '
176  write(num_rep_42,428)'! file_smb_snap     = fichier serie temp anomalie SMB de GCM            '
177  write(num_rep_42,428)'! massb_time        = 0:fixe, 1:interp snapshots, 2:snapsh+interp vert  '
178  write(num_rep_42,428)'!_______________________________________________________________________' 
179
180  if (massb_time == 1) then ! lecture des snapshots
181! allocation dynamique de time_snap
182     if (allocated(time_snap)) then
183        deallocate(time_snap,stat=err)
184        if (err/=0) then
185           print *,"Erreur à la desallocation de time_snap",err
186           stop
187        end if
188     end if
189
190     allocate(time_snap(nb_snap),stat=err)
191     if (err/=0) then
192        print *,"erreur a l'allocation du tableau time_snap ",err
193        print *,"nb_snap = ",nb_snap
194        stop
195     end if
196
197! remplissage de time_snap
198!write(6,*) 'time_snap'
199     do i=1,nb_snap
200        time_snap(i) = time_depart_snaps + ecart_snap * (i-1)
201!   write(6,*) i,time_snap(i)
202     end do
203
204! allocation dynamique de smb_snap
205     if (allocated(smb_snap)) then
206        deallocate(smb_snap,stat=err)
207        if (err/=0) then
208           print *,"Erreur à la desallocation de smb_snap",err
209           stop
210        end if
211     end if
212
213     allocate(smb_snap(nx,ny,nb_snap),stat=err)
214     if (err/=0) then
215        print *,"erreur a l'allocation du tableau smb_snap ",err
216        print *,"nx,ny,nb_snap = ",nx,',',ny,',',nb_snap
217        stop
218     end if
219
220! lecture de smb_snap
221     file_smb_snap  = trim(dirnameinp)//trim(file_smb_snap)
222     call Read_Ncdf_var('z',file_smb_snap,tab3D)
223     smb_snap (:,:,:) = Tab3D(:,:,:) * coef_smb_unit
224
225! ce sont des anomalies : ajoute les valeurs de reference
226     do k = 1,nb_snap
227        do j = 1,ny
228           do i = 1,nx
229              smb_snap (i,j,k) = smb_snap(i,j,k) + bm(i,j)
230           end do
231        end do
232     end do
233
234! copie la valeur de reference dans bm_0
235!     bm_0(:,:) = bm(:,:)
236  endif
237! copie la valeur de reference dans bm_0
238  bm_0(:,:) = bm(:,:) !afq marion dufresne: we do it even if we don't read any snapshot
239 
240  filin=trim(dirforcage)//trim(filforc)
241
242  open(num_forc,file=filin,status='old')
243 
244  read(num_forc,*) control,nft
245 
246    ! Determination of file size (line nb), allocation of perturbation array
247 
248  if (control.ne.'nb_lines') then
249     write(6,*) filin,'indiquer le nb de ligne en debut de fichier:'
250     write(6,*) 'le nb de lignes et le label de control nb_lines'
251     stop
252  endif
253   
254    ! Dimensionnement des tableaux tdate, ....
255
256  if (.not.allocated(tdate)) then
257     allocate(tdate(nft),stat=err)
258     if (err/=0) then
259        print *,"erreur a l'allocation du tableau Tdate",err
260        stop 4
261     end if
262  end if
263
264  if (.not.allocated(spert)) then
265     allocate(spert(nft),stat=err)
266     if (err/=0) then
267        print *,"erreur a l'allocation du tableau Spert",err
268        stop 4
269     end if
270  end if
271 
272  if (.not.allocated(tpert)) then
273     allocate(tpert(nft),stat=err)
274     if (err/=0) then
275        print *,"erreur a l'allocation du tableau Tpert",err
276        stop 4
277     end if
278  end if
279 
280  do i=1,nft
281     read(num_forc,*) tdate(i),spert(i),tpert(i)
282  end do
283  close(num_forc)
284 
285  tpert(:)=tpert(:)*coefT
286 
287
288 
289! ecriture de verification
290!file_smb_snap = 'test_sortie_smb.nc'
291!file_smb_snap  = trim(dirnameout)//trim(file_smb_snap)
292
293!ncidloc = 501
294
295!call  lect_netcdf_type
296!write(6,*) 'ncdf_type'
297
298!if (ncdf_type.eq.32) then
299!   status  = nf90_create(TRIM(file_smb_snap),NF90_WRITE,ncidloc)    ! ouverture du fichier
300!else if (ncdf_type.eq.64) then
301!   status  = nf90_create(trim(file_smb_snap),and(nf90_write,nf90_64bit_offset),ncidloc) ! r2d2
302!end if
303
304!   status  = nf90_close(ncidloc)                                 ! fermeture
305
306! call write_ncdf_dim('x',trim(file_smb_snap),nx)               ! dimensions des tableaux
307! call write_ncdf_dim('y',trim(file_smb_snap),ny)
308! call write_ncdf_dim('time',trim(file_smb_snap),nb_snap)
309
310!Tab3D(:,:,:) = smb_snap (:,:,:)
311!dimtestname(1) = 'x'
312!dimtestname(2) = 'y'
313!dimtestname(3) = 'time'
314
315!call write_ncdf_var(trim('smb'),dimtestname,trim(file_smb_snap),tab3D,'double')
316   
317end subroutine input_clim
318
319!--------------------------------------------------------------------------------
320!> SUBROUTINE: init_forclim
321!! Routine qui permet d'initialiser les variables climatiques au cours du temps
322!>
323subroutine init_forclim
324
325  implicit none
326  namelist/lapse_rates/T_lapse_rate
327  namelist/clim_pert_massb/coefT,filforc,pertsmb,rapsmb
328 
329  rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
330  read(num_param,lapse_rates)
331
332! formats pour les ecritures dans 42
333428 format(A)
334
335  rewind(num_param)                     ! pour revenir au debut du fichier param_list.dat
336  read(num_param,lapse_rates)
337
338  write(num_rep_42,428)'!________________________________________________________________' 
339  write(num_rep_42,428)'!  module climat_Grice2sea_years_mod                             '
340  write(num_rep_42,lapse_rates)
341  write(num_rep_42,428)'!T_lapse_rate       = lapse rate temp annuelle                   '
342  write(num_rep_42,428)'!________________________________________________________________' 
343
344 
345  rewind(num_param) 
346  read(num_param,clim_pert_massb)
347 
348  write(num_rep_42,428)'!___________________________________________________________'
349  write(num_rep_42,428) '&clim_pert                      ! module climat_perturb_mod '
350  write(num_rep_42,*)
351  write(num_rep_42,*) 'coefT        = ', coefT
352  write(num_rep_42,'(A,A)') ' filforc      = ', filforc
353  write(num_rep_42,*) 'pertsmb        = ', pertsmb
354  write(num_rep_42,*) 'rapsmb         = ', rapsmb
355  write(num_rep_42,*)'/'
356  write(num_rep_42,*)
357
358
359! appelle la routine de lecture des smb annuels
360  call input_clim
361  if (massb_time == 1) then ! lecture gradients smb
362     call init_grad_smb
363  endif
364
365  return
366end subroutine init_forclim
367
368!--------------------------------------------------------------------------------
369!> SUBROUTINE: forclim
370!! 
371!!  Routine qui permet le calcul climatique au cours du temps
372!!  @note Au temps considere (time) attribue les scalaires
373!!  @note   - tafor : forcage en temperature
374!!  @note   - sealevel : forcage niveau des mers
375!!  @note   - coefbmelt : forcage fusion basale ice shelves
376!>
377subroutine forclim               !  au temps considere (time)
378
379  implicit none
380
381  integer i,ift
382 
383  select case (massb_time)
384     case(0)
385        ! smb fixe et Tann avec lapse rate + Tafor
386       
387        if(time.lt.tdate(1)) then
388           tafor=tpert(1)
389           sealevel=spert(1)
390           ift=1
391
392        else if (time.ge.tdate(nft)) then
393           tafor=tpert(nft)
394           sealevel=spert(nft)
395           ift=nft
396           
397        else
398           do i=1,nft-1
399              if((time.ge.tdate(i)).and.(time.lt.tdate(i+1))) then ! entre i et i+1 : cas general
400                 tafor=tpert(i)+(tpert(i+1)-tpert(i))*       &
401                      (time-tdate(i))/(tdate(i+1)-tdate(i))
402                 sealevel=spert(i)+(spert(i+1)-spert(i))*    &
403                      (time-tdate(i))/(tdate(i+1)-tdate(i))
404                 ift=i
405                 goto 100
406              endif
407           end do
408        endif
409100     continue
410
411       
412        Tann (:,:) = Ta0 (:,:) + T_lapse_rate * (S(:,:)-S0(:,:)) +Tafor
413        Ts(:,:)    = Tann(:,:)
414
415        ! aurel marion dufresne: we might want to decrease the SMB during glacials..?
416        if (pertsmb.eq.1) then
417           bm(:,:) = bm_0(:,:) * exp( rapsmb *(Tann(:,:)-Ta0(:,:)))
418           if (Tafor.lt.0.) then
419              where(bm(:,:).lt.0.) bm(:,:)=min(bm(:,:)-Tafor*0.05,1.) !10 degrees less give 0.5 meter more ?
420           end if
421        end if
422       
423     case(1)
424        call massb_Ice2sea_RCM
425     case default
426        print *, 'Numero de massb invalide dans climat_Grice2sea_years_mod'
427        stop
428     end select
429end subroutine forclim
430
431
432subroutine massb_Ice2sea_RCM                ! calcule le mass balance
433
434  implicit none
435  integer             :: k_snap               ! pour calculer les indices de temps
436  real                :: time_bis             ! pour repliquer les dernieres annees
437  integer             :: k
438
439! calcul smb a partir fichier snapshots massb_Ice2sea_RCM
440! Calcule le mass balance d'apres un fichier de snapshots
441! avec la temperature parametree
442!     surface temperature et accumulation
443  Tann (:,:) = Ta0 (:,:) + T_lapse_rate * (S(:,:)-S0(:,:)) 
444  Ts(:,:)    = Tann(:,:)
445! calcule bm_time par interpolation entre deux snapshots
446! avant prend la valeur de reference
447! apres prend la derniere valeur
448! en general les snapshots vont de 2000 a 2200
449  if(time.lt.time_snap(1)) then              ! time avant le forcage
450     bm_time(:,:) = bm_0(:,:)
451     k_snap       = 1
452     S_ref(:,:)   = S(:,:)                   ! du coup sera la surface de reference avant le forcage
453     icum         = 0
454     i_moy        = 0
455     bm_ref_t(:,:)=  bm_0(:,:)               ! bilan de masse de reference
456     time_prec    = time
457  else if (time.ge.time_snap(nb_snap)) then  ! time apres le forcage
458     k_snap       = nb_snap
459     bm_time(:,:) =  smb_snap (:,:,k_snap)
460     if (abs(time-time_prec-1.).lt.dt) then   !
461        time_prec = time_prec + 1
462        icum      = 1
463     else
464        icum      = 0
465     endif
466  else                                       ! cas general
467     do k = 1 , nb_snap-1
468        if((time.ge.time_snap(k)).and.(time.lt.time_snap(k+1))) then ! entre k et k+1
469           bm_time(:,:) = smb_snap(:,:,k)+(smb_snap(:,:,k+1)-smb_snap(:,:,k)) *   &
470                (time-time_snap(k))/(time_snap(k+1)-time_snap(k))
471! exactement sur le snapshot et avec un ecart 1 an par rapport au precedent stockage
472!           write(6,*) 'time,tests',k,time,time-time_snap(k),time-time_prec-1.
473           if ((abs(time-time_snap(k)).le.dt).and.(abs(time-time_prec-1.).lt.dt)) then   
474              k_snap    = k
475              icum      = 1
476              time_prec = time_snap(k)     ! time_prec est le temps du precedent
477                                                 ! snapshot garde
478           else
479              icum    = 0
480           endif
481           exit
482        endif
483     end do
484  endif
485  call  grad_smb              !-----------------------------> A faire
486
487  if (massb_time == 1) then ! pas d'interpolation verticale
488     bm(:,:) = bm_time(:,:)
489  else if (massb_time == 2) then ! interpolation verticale
490     ! ajuste bm en fonction du temps et du gradient
491     bm(:,:) = bm_time(:,:) + grad_bm(:,:) *(S(:,:) - S_ref(:,:))
492     write(6,897) time, time_prec, icum, i_moy
493897  format('test temps smb   ',2(f0.3,1x),2(i0,1x))
494  endif
495
496! garde les 10 dernieres annees et calcule la moyenne
497  if (icum.eq.1) then                           ! stockage dans le tableau bm_ref_10
498     do k = 9,1,-1
499        bm_ref_10(:,:,k+1) = bm_ref_10(:,:,k)   ! on decale tous les elements
500     end do
501     bm_ref_10(:,:,k+1) = bm(:,:)               ! le plus recent est en position 1
502     i_moy              = i_moy +1              ! compte combien il y en a pour la moyenne
503     i_moy              = min(i_moy,10)
504     bm_ref_t(:,:)   = 0.
505     do k = 1,i_moy
506        bm_ref_t(:,:) = bm_ref_t(:,:) + bm_ref_10(:,:,k)
507     end do
508     bm_ref_t(:,:)    = bm_ref_t(:,:)/i_moy
509     write(6,898) time, time_prec, icum, i_moy
510898  format('cumul pour gradient  ',2(f0.3,1x),2(i0,1x))
511  end if
512end subroutine massb_Ice2sea_RCM
513 
514!------------------------------------------------------------------------------
515!> initialise le calcul du gradient vertical de smb
516subroutine init_grad_smb 
517
518  use module3D_phy
519  implicit none
520
521  character(len=120)       ::  file_grad_smb ! nom du fichier gradients de Tamsin
522  character(len=40)        ::  fin_ligne     ! fin de la ligne
523  real                     ::  grad
524
525  namelist/grad_smb/file_grad_smb
526
527! Dans lequel sont :
528! grad_N_smb_pos,grad_N_smb_neg,grad_S_smb_pos,grad_S_smb_neg,lim_lat
529
530  rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
531  read(num_param,grad_smb)
532
533
534! formats pour les ecritures dans 42
535428 format(A)
536
537  rewind(num_param)                     ! pour revenir au debut du fichier param_list.dat
538  read(num_param,grad_smb)
539
540  write(num_rep_42,428)'!________________________________________________________________' 
541  write(num_rep_42,428)'!  gradient smb climat_Grice2sea_years_mod                       '
542  write(num_rep_42,grad_smb)
543  write(num_rep_42,428)'!grad_smb       = fichier gradient SMB                           '
544  write(num_rep_42,428)'!________________________________________________________________' 
545
546
547  file_grad_smb=trim(dirnameinp)//'SMB-H-Feedback/gradients_18_07_2012/'//trim(file_grad_smb)
548  open(622,file=file_grad_smb)
549  do i=1,4
550     read(622,'(f9.4,A)') grad,fin_ligne
551     write(6,*) grad,fin_ligne
552
553     if (index(fin_ligne, "North").ne.0) then       ! North est dans la ligne
554        if (index(fin_ligne, "<").ne.0)  then       ! smb negatif
555           grad_N_smb_neg = grad
556        else if (index(fin_ligne, ">=").ne.0) then  ! smb positif
557           grad_N_smb_pos = grad
558        else
559           write(6,*) 'probleme lecture North fichier smb', file_grad_smb
560           STOP
561        end if
562
563     else if (index(fin_ligne, "South").ne.0) then   !South est dans la ligne
564        if (index(fin_ligne, "<").ne.0)  then        ! smb negatif
565           grad_S_smb_neg = grad
566        else if (index(fin_ligne, ">=").ne.0) then   ! smb positif
567           grad_S_smb_pos = grad
568        else
569           write(6,*) 'probleme lecture South fichier smb  ', file_grad_smb
570           STOP
571        end if
572
573     else
574        write(6,*) 'probleme lecture ni North ni South fichier smb  ', file_grad_smb
575        write(6,*) 'fin_ligne',fin_ligne,'   index North', index(fin_ligne, "North"),'   index South', index(fin_ligne, "South")
576        STOP
577     end if
578  end do
579
580  write(6,*) 'coefficients lus '
581  write(6,*) 'grad_N_smb_pos', grad_N_smb_pos
582  write(6,*) 'grad_N_smb_neg', grad_N_smb_neg
583  write(6,*) 'grad_S_smb_pos', grad_S_smb_pos
584  write(6,*) 'grad_S_smb_neg', grad_S_smb_neg
585
586  grad_N_smb_pos = coef_smb_unit * grad_N_smb_pos
587  grad_N_smb_neg = coef_smb_unit * grad_N_smb_neg
588  grad_S_smb_pos = coef_smb_unit * grad_S_smb_pos
589  grad_S_smb_neg = coef_smb_unit * grad_S_smb_neg
590
591  write(6,*) 'coefficients *coef_smb_unit '
592  write(6,*) 'grad_N_smb_pos', grad_N_smb_pos
593  write(6,*) 'grad_N_smb_neg', grad_N_smb_neg
594  write(6,*) 'grad_S_smb_pos', grad_S_smb_pos
595  write(6,*) 'grad_S_smb_neg', grad_S_smb_neg
596
597
598end subroutine init_grad_smb
599
600
601!------------------------------------------------------------------------------
602!> Calcule le gradient vertical de smb
603
604subroutine grad_smb 
605
606  use module3D_phy
607  implicit none
608
609  do j = 1,ny
610     do i = 1,nx
611        if (Ylat(i,j).gt.77.) then                     ! region nord
612           if (bm_ref_t(i,j).lt. 0.) then              ! smb negatif
613              grad_bm (i,j) = grad_N_smb_neg
614
615           else if (bm_ref_t(i,j).ge. 0.) then         ! smb positif
616              grad_bm (i,j) = grad_N_smb_pos 
617
618           end if
619
620        else if (Ylat(i,j).le.77.) then                ! region sud
621           if (bm_ref_t(i,j).lt. 0.) then              ! smb negatif
622              grad_bm (i,j) = grad_S_smb_neg
623
624           else if (bm_ref_t(i,j).ge. 0.) then         ! smb positif
625              grad_bm (i,j) = grad_S_smb_pos 
626
627           end if
628        end if
629     end do
630  end do
631
632end subroutine grad_smb
633
634end module  climat_Grice2sea_years_perturb_mod
Note: See TracBrowser for help on using the repository browser.