Changeset 93 for trunk/SOURCES/climat-perturb_mod-0.4.f90
- Timestamp:
- 11/28/16 16:03:50 (8 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/SOURCES/climat-perturb_mod-0.4.f90
r65 r93 14 14 15 15 use module3d_phy,only:nx,ny,S,S0,Tann,Tjuly,precip,acc,Ylat,num_forc,num_param,num_rep_42,dirforcage,dirnameinp,tafor,time,sealevel,coefbmshelf 16 use netcdf 17 use io_netcdf_grisli 16 18 17 19 implicit none … … 41 43 !! Routine qui permet d'initialiser les variations temporelles des variables climatiques 42 44 !> 43 subroutine input_clim !routine qui permet d'initialiser les variations temporelles des variables climatiques 45 46 subroutine input_clim 44 47 45 48 implicit none 49 character(len=100) :: precip_file ! fichier precipitations 50 character(len=100) :: temp_annual_file ! fichier temperature annuelle 51 real :: coef_dens ! pour corriger si donnees en eq. eau 52 logical :: temp_param ! si utilisation de temperature parametree 53 real*8, dimension(:,:), pointer :: data_2D => null() ! donnees lues dans le netcdf 54 55 46 56 character(len=8) :: control !label to check clim. forc. file (filin) is usable 47 57 character(len=80):: filin 48 58 integer :: err !< pour l'allocation des tableaux 49 integer :: i 59 integer :: i,j 60 61 namelist/climat_acc_T_gen/precip_file,coef_dens,temp_annual_file 62 63 428 format(A) 64 rewind(num_param) ! pour revenir au debut du fichier param_list.dat 65 read(num_param,climat_acc_T_gen) 66 67 write(num_rep_42,428)'!___________________________________________________________' 68 write(num_rep_42,428)'! module lect_clim_acc_T_ant_gen ' 69 write(num_rep_42,climat_acc_T_gen) 70 write(num_rep_42,428)'!___________________________________________________________' 71 72 73 ! precipitation 74 precip_file = trim(dirnameinp)//trim(precip_file) 75 76 call Read_ncdf_var('precip',trim(precip_file),data_2D) 77 precip(:,:)=data_2D(:,:) 78 !call lect_datfile(nx,ny,precip,1,precip_file) 79 80 precip(:,:)=precip(:,:)*coef_dens 81 acc(:,:)=precip(:,:) 82 83 if ((trim(temp_annual_file).eq.'no').or.(trim(temp_annual_file).eq.'NO')) then 84 temp_param=.true. 85 else 86 temp_param=.false. 87 end if 88 89 ! temperature en surface 90 91 test_param: if (.not.temp_param) then 92 temp_annual_file = trim(dirnameinp)//trim(temp_annual_file) 93 94 95 call Read_ncdf_var('Tann',trim(temp_annual_file),data_2D) 96 Tann(:,:)=data_2D(:,:) 97 ! call lect_input(3,'Tann',1,Tann,temp_annual_file,trim(dirnameinp)//trim(runname)//'.nc') 98 !call lect_datfile(nx,ny,Tann,1,temp_annual_file) ! temperature annuelle 99 100 else ! parametrisation de Fortuin pour la temperature annuelle. 101 102 do j=1,ny 103 do i=1,nx 104 105 7 if (s0(i,j).le.200.) then ! shelfs 106 tann(i,j)=49.642-0.943*abs(ylat(i,j)) 107 else if ((s0(i,j).gt.200.).and.(s0(i,j).lt.1500.)) then ! pente 108 tann(i,j)=36.689-0.005102*s0(i,j)-0.725*abs(ylat(i,j)) 109 else if (s0(i,j).ge.1500.) then ! plateau 110 tann(i,j)=7.405-0.014285*s0(i,j)-0.180*abs(ylat(i,j)) 111 endif 112 end do 113 end do 114 end if test_param 115 116 ta0(:,:)=tann(:,:) 117 118 119 ! pour la temperature d'ete, idem parametrisation huybrechts 120 do j=1,ny 121 do i=1,nx 122 123 tjuly(i,j)=tann(i,j)-17.65+0.00222*s0(i,j)& 124 +0.40802*abs(ylat(i,j)) 125 end do 126 end do 127 128 129 130 50 131 51 132 ! Lecture du forcage … … 138 219 139 220 140 !!!!!!!! ATTENTION AJOUTE POUR TEST MAIS A REMETTRE AU PROPRE PLUS TARD C. DUMAS !!!!!!!!!! 141 !!!!!!!! ancien input_climat_ref de lect_clim_act_anteis 142 ! accumulation de Philippe 143 filin='accumHUY40km.dat' 144 call lect_eis(nx,ny,precip,filin,DIRNAMEINP) 145 !====================================== La reponse est 42 =========== 146 write(num_rep_42,*) 'fichier accum : ', filin 147 148 ! cas particulier de Vostok 149 ivo=101 150 jvo=62 151 do j=jvo-1,jvo+1 152 do i=ivo-1,ivo+1 153 precip(i,j)=0.02 ! valeur plus faible a Vostok. 154 end do 155 end do 156 acc(:,:)=precip(:,:) 157 158 ! temperature en surface : 159 ! parametrisation de Fortuin pour la temperature annuelle. 160 do j=1,ny 161 do i=1,nx 162 163 if (s0(i,j).le.200.) then ! shelfs 164 tann(i,j)=49.642-0.943*abs(ylat(i,j)) 165 else if ((s0(i,j).gt.200.).and.(s0(i,j).lt.1500.)) then ! pente 166 tann(i,j)=36.689-0.005102*s0(i,j)-0.725*abs(ylat(i,j)) 167 else if (s0(i,j).ge.1500.) then ! plateau 168 tann(i,j)=7.405-0.014285*s0(i,j)-0.180*abs(ylat(i,j)) 169 endif 170 171 ta0(i,j)=tann(i,j) 172 ! pour la temperature d'ete, idem parametrisation huybrechts 173 tjuly(i,j)=tann(i,j)-17.65+0.00222*s0(i,j)& 174 +0.40802*abs(ylat(i,j)) 175 end do 176 end do 177 !!!!!!!! FIN MODIF TEMPORAIRE !!!!!!!!!! 178 221 !cdc Commente pour être compatible avec lecture fichiers Cat Schoofing 222 223 !~ !!!!!!!! ATTENTION AJOUTE POUR TEST MAIS A REMETTRE AU PROPRE PLUS TARD C. DUMAS !!!!!!!!!! 224 !~ !!!!!!!! ancien input_climat_ref de lect_clim_act_anteis 225 !~ ! accumulation de Philippe 226 !~ filin='accumHUY40km.dat' 227 !~ call lect_eis(nx,ny,precip,filin,DIRNAMEINP) 228 !~ !====================================== La reponse est 42 =========== 229 !~ write(num_rep_42,*) 'fichier accum : ', filin 230 231 !~ ! cas particulier de Vostok 232 !~ ivo=101 233 !~ jvo=62 234 !~ do j=jvo-1,jvo+1 235 !~ do i=ivo-1,ivo+1 236 !~ precip(i,j)=0.02 ! valeur plus faible a Vostok. 237 !~ end do 238 !~ end do 239 !~ acc(:,:)=precip(:,:) 240 241 !~ ! temperature en surface : 242 !~ ! parametrisation de Fortuin pour la temperature annuelle. 243 !~ do j=1,ny 244 !~ do i=1,nx 245 246 !~ if (s0(i,j).le.200.) then ! shelfs 247 !~ tann(i,j)=49.642-0.943*abs(ylat(i,j)) 248 !~ else if ((s0(i,j).gt.200.).and.(s0(i,j).lt.1500.)) then ! pente 249 !~ tann(i,j)=36.689-0.005102*s0(i,j)-0.725*abs(ylat(i,j)) 250 !~ else if (s0(i,j).ge.1500.) then ! plateau 251 !~ tann(i,j)=7.405-0.014285*s0(i,j)-0.180*abs(ylat(i,j)) 252 !~ endif 253 254 !~ ta0(i,j)=tann(i,j) 255 !~ ! pour la temperature d'ete, idem parametrisation huybrechts 256 !~ tjuly(i,j)=tann(i,j)-17.65+0.00222*s0(i,j)& 257 !~ +0.40802*abs(ylat(i,j)) 258 !~ end do 259 !~ end do 260 !~ !!!!!!!! FIN MODIF TEMPORAIRE !!!!!!!!!! 261 !cdc fin Commente pour être compatible avec lecture fichiers Cat Schoofing 179 262 180 263 return
Note: See TracChangeset
for help on using the changeset viewer.