source: trunk/SOURCES/Ant45_CISM_files/lect-Ant_clim_CISM_45_dat.f90 @ 111

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

initial import GRISLI trunk

File size: 2.0 KB
Line 
1!> \file lect-Ant_clim_CISM_45_dat.f90
2!!Module pour la lecture du climat
3!<
4
5!> \namespace lect_clim_act_ant_45
6!! Module pour la lecture du climat
7!! \author ...
8!! \date ...
9!! @note Used module
10!! @note   - use module3D_phy
11!! @note   - use ablation_ann
12!<
13
14module lect_clim_act_ant_45
15
16  use module3D_phy
17  use ablation_ann                     ! le module pdd base sur Tann et Tjuly
18  use interface_input
19  character(len=100) :: precip_file      ! precipitations
20
21contains
22
23  subroutine input_climat_ref()
24
25    namelist/climat_CISM_45/precip_file
26
27428 format(A)
28    rewind(num_param)                     ! pour revenir au debut du fichier param_list.dat
29    read(num_param,climat_CISM_45)
30
31    write(num_rep_42,428)'!___________________________________________________________' 
32    write(num_rep_42,428)'!  module  lect_topo_ant_CISM      '
33    write(num_rep_42,climat_CISM_45)
34    write(num_rep_42,428)'!___________________________________________________________' 
35
36
37    precip_file  = trim(dirnameinp)//trim(precip_file)
38
39
40    call lect_input(3,'precip',1,precip,precip_file,trim(dirnameinp)//trim(runname)//'.nc')   
41    !call lect_datfile(nx,ny,precip,1,precip_file)                  ! precipitation
42
43    precip(:,:)=precip(:,:)/0.91
44    acc(:,:)=precip(:,:)
45
46
47    !    temperature en surface :
48    !    parametrisation de Fortuin pour la temperature annuelle.
49    do j=1,ny
50       do i=1,nx
51
527         if (s0(i,j).le.200.) then   ! shelfs
53             tann(i,j)=49.642-0.943*abs(ylat(i,j))
54          else if ((s0(i,j).gt.200.).and.(s0(i,j).lt.1500.)) then ! pente
55             tann(i,j)=36.689-0.005102*s0(i,j)-0.725*abs(ylat(i,j))
56          else if (s0(i,j).ge.1500.) then        ! plateau
57             tann(i,j)=7.405-0.014285*s0(i,j)-0.180*abs(ylat(i,j))
58          endif
59
60          ta0(i,j)=tann(i,j)
61
62          !           pour la temperature d'ete, idem parametrisation huybrechts
63          tjuly(i,j)=tann(i,j)-17.65+0.00222*s0(i,j)&
64               +0.40802*abs(ylat(i,j))
65       end do
66    end do
67
68  end subroutine input_climat_ref
69
70end module  lect_clim_act_ant_45
Note: See TracBrowser for help on using the repository browser.