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

Last change on this file since 180 was 180, checked in by aquiquet, 6 years ago

Cleaning-up: unused variables removed

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