source: branches/iLoveclim/SOURCES/GrIce2sea_files/climat_GrIce2sea_years_perturb_mod.f90

Last change on this file was 244, checked in by aquiquet, 5 years ago

Grisli-iloveclim branch merged to trunk at revision 243

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