source: trunk/LMDZ.GENERIC/libf/phystd/iniphysiq.F90 @ 1303

Last change on this file since 1303 was 1216, checked in by emillour, 11 years ago

Generic model:
Major cleanup, in order to ease the use of LMDZ.GENERIC with (parallel) dynamics
in LMDZ.COMMON: (NB: this will break LMDZ.UNIVERSAL, which should be thrashed
in the near future)

  • Updated makegcm_* scripts (and makdim) and added the "-full" (to enforce full recomputation of the model) option
  • In dyn3d: converted control.h to module control_mod.F90 and converted iniadvtrac.F to module infotrac.F90
  • Added module mod_const_mpi.F90 in dyn3d (not used in serial mode)
  • Rearanged input/outputs routines everywhere to handle serial/MPI cases. physdem.F => phyredem.F90 , phyetat0.F => phyetat0.F90 ; all read/write routines for startfi files are gathered in module iostart.F90
  • added parallelism related routines init_phys_lmdz.F90, comgeomphy.F90, dimphy.F90, iniphysiq.F90, mod_grid_phy_lmdz.F90, mod_phys_lmdz_mpi_data.F90, mod_phys_lmdz_mpi_transfert.F90, mod_phys_lmdz_omp_data.F90, mod_phys_lmdz_omp_transfert.F90, mod_phys_lmdz_para.F90, mod_phys_lmdz_transfert_para.F90 in phymars and mod_const_mpi.F90 in dyn3d (for compliance with parallelism)
  • added created generic routines 'planetwide_maxval' and 'planetwide_minval', in module "planetwide_mod", that enable obtaining the max and min of a field over the whole planet. This should be further imroved with computation of means (possibly area weighed), etc.

EM

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