source: codes/icosagcm/branches/SATURN_DYNAMICO/LMDZ.COMMON/libf/phystd/iniphysiq.F90 @ 222

Last change on this file since 222 was 222, checked in by ymipsl, 10 years ago

Creating temporary dynamico/lmdz/saturn branche

YM

File size: 3.2 KB
Line 
1subroutine iniphysiq(ngrid,nlayer, punjours, pdayref,ptimestep,          &
2                     plat,plon,parea,pcu,pcv,                            &
3                     prad,pg,pr,pcpp,iflag_phys)
4
5use dimphy, only : klev ! number of atmospheric levels
6use mod_grid_phy_lmdz, only : klon_glo ! number of atmospheric columns
7                                       ! (on full grid)
8use mod_phys_lmdz_para, only : klon_omp, & ! number of columns (on local omp grid)
9                               klon_omp_begin, & ! start index of local omp subgrid
10                               klon_omp_end, & ! end index of local omp subgrid
11                               klon_mpi_begin, & ! start indes of columns (on local mpi grid)
12                               klon_mpi
13use comgeomphy, only : airephy, & ! physics grid area (m2)
14                       cuphy, & ! cu coeff. (u_covariant = cu * u)
15                       cvphy, & ! cv coeff. (v_covariant = cv * v)
16                       rlond, & ! longitudes
17                       rlatd ! latitudes
18use infotrac, only : nqtot ! number of advected tracers
19
20implicit none
21
22real,intent(in) :: prad ! radius of the planet (m)
23real,intent(in) :: pg ! gravitational acceleration (m/s2)
24real,intent(in) :: pr ! ! reduced gas constant R/mu
25real,intent(in) :: pcpp ! specific heat Cp
26real,intent(in) :: punjours ! length (in s) of a standard day
27integer,intent(in) :: ngrid ! number of horizontal grid points in the physics (full grid)
28integer,intent(in) :: nlayer ! number of atmospheric layers
29real,intent(in) :: plat(klon_mpi) ! latitudes of the physics grid
30real,intent(in) :: plon(klon_mpi) ! longitudes of the physics grid
31real,intent(in) :: parea(klon_mpi) ! area (m2)
32real,intent(in) :: pcu(klon_mpi) ! cu coeff. (u_covariant = cu * u)
33real,intent(in) :: pcv(klon_mpi) ! cv coeff. (v_covariant = cv * v)
34integer,intent(in) :: pdayref ! reference day of for the simulation
35real,intent(in) :: ptimestep !physics time step (s)
36integer,intent(in) :: iflag_phys ! type of physics to be called
37
38integer :: ibegin,iend,offset
39character(len=20) :: modname='iniphysiq'
40character(len=80) :: abort_message
41
42IF (nlayer.NE.klev) THEN
43  write(*,*) 'STOP in ',trim(modname)
44  write(*,*) 'Problem with dimensions :'
45  write(*,*) 'nlayer     = ',nlayer
46  write(*,*) 'klev   = ',klev
47  abort_message = ''
48  CALL abort_gcm (modname,abort_message,1)
49ENDIF
50
51IF (ngrid.NE.klon_glo) THEN
52  write(*,*) 'STOP in ',trim(modname)
53  write(*,*) 'Problem with dimensions :'
54  write(*,*) 'ngrid     = ',ngrid
55  write(*,*) 'klon   = ',klon_glo
56  abort_message = ''
57  CALL abort_gcm (modname,abort_message,1)
58ENDIF
59
60!$OMP PARALLEL PRIVATE(ibegin,iend)
61!$OMP+         SHARED(parea,pcu,pcv,plon,plat)
62     
63offset=0
64airephy(1:klon_omp)=parea(offset+klon_omp_begin:offset+klon_omp_end)
65cuphy(1:klon_omp)=pcu(offset+klon_omp_begin:offset+klon_omp_end)
66cvphy(1:klon_omp)=pcv(offset+klon_omp_begin:offset+klon_omp_end)
67rlond(1:klon_omp)=plon(offset+klon_omp_begin:offset+klon_omp_end)
68rlatd(1:klon_omp)=plat(offset+klon_omp_begin:offset+klon_omp_end)
69
70! copy some fundamental parameters to physics
71! and do some initializations
72call inifis(klon_omp,nlayer,nqtot,pdayref,punjours,ptimestep, &
73            rlatd,rlond,airephy,prad,pg,pr,pcpp)
74
75!$OMP END PARALLEL
76
77
78end subroutine iniphysiq
Note: See TracBrowser for help on using the repository browser.