Ignore:
Timestamp:
Aug 25, 2015, 5:14:59 PM (10 years ago)
Author:
Ehouarn Millour
Message:

More on physics/dynamics separation and cleanup:

  • Set things up so that all physics-related initializations are done via iniphysiq.
  • Created a "geometry_mod.F90" module in phy_common to store information on the loacl grid (i.e. replaces comgeomphy) and moreover give these variables more obvious names (e.g.: rlond => longitude, rlatd => latitude, airephy => cell_area).
  • removed obsolete comgeomphy.h and comgeomphy.F90

EM

Location:
LMDZ5/trunk/libf/dynlonlat_phylonlat
Files:
7 edited

Legend:

Unmodified
Added
Removed
  • LMDZ5/trunk/libf/dynlonlat_phylonlat/calfis_loc.F

    r2333 r2351  
    3232c   .........
    3333      USE dimphy
    34       USE mod_phys_lmdz_para, mpi_root_xx=>mpi_root
     34      USE mod_phys_lmdz_mpi_data, mpi_root_xx=>mpi_master
     35      USE mod_phys_lmdz_omp_data, ONLY: klon_omp, klon_omp_begin
     36      USE mod_const_mpi, ONLY: COMM_LMDZ
    3537      USE mod_interface_dyn_phys
    3638      USE IOPHY
  • LMDZ5/trunk/libf/dynlonlat_phylonlat/calfis_p.F

    r2239 r2351  
    2929! If using physics
    3030      USE dimphy
    31       USE mod_phys_lmdz_para, mpi_root_xx=>mpi_root
     31      USE mod_phys_lmdz_mpi_data, mpi_root_xx=>mpi_master
     32      USE mod_phys_lmdz_omp_data, ONLY: klon_omp, klon_omp_begin
     33      USE mod_const_mpi, ONLY: COMM_LMDZ
    3234      USE mod_interface_dyn_phys
    3335      USE IOPHY
  • LMDZ5/trunk/libf/dynlonlat_phylonlat/mod_interface_dyn_phys.F90

    r2239 r2351  
    77 
    88 
     9CONTAINS
     10 
    911#ifdef CPP_PARA
    1012! Interface with parallel physics,
    11 CONTAINS
    12  
    1313  SUBROUTINE Init_interface_dyn_phys
    1414    USE mod_phys_lmdz_mpi_data
     
    5454 
    5555  END SUBROUTINE Init_interface_dyn_phys
     56#else
     57  SUBROUTINE Init_interface_dyn_phys
     58  ! dummy routine for seq case
     59  END SUBROUTINE Init_interface_dyn_phys
    5660#endif
    5761! of #ifdef CPP_PARA
  • LMDZ5/trunk/libf/dynlonlat_phylonlat/phydev/iniphysiq_mod.F90

    r2347 r2351  
    66CONTAINS
    77
    8 SUBROUTINE iniphysiq(iim,jjm,nlayer,punjours, pdayref,ptimestep, &
     8SUBROUTINE iniphysiq(iim,jjm,nlayer, &
     9                     nbp, communicator, &
     10                     punjours, pdayref,ptimestep, &
    911                     rlatu,rlatv,rlonu,rlonv,aire,cu,cv, &
    1012                     prad,pg,pr,pcpp,iflag_phys)
    11   USE dimphy, ONLY: klev ! number of atmospheric levels
    12   USE mod_grid_phy_lmdz, ONLY: klon_glo ! number of atmospheric columns
    13                                         ! (on full grid)
     13  USE dimphy, ONLY: init_dimphy
     14  USE mod_grid_phy_lmdz, ONLY: klon_glo,  & ! number of atmospheric columns (on full grid)
     15                               regular_lonlat  ! regular longitude-latitude grid type
    1416  USE mod_phys_lmdz_para, ONLY: klon_omp, & ! number of columns (on local omp grid)
    1517                                klon_omp_begin, & ! start index of local omp subgrid
    1618                                klon_omp_end, & ! end index of local omp subgrid
    1719                                klon_mpi_begin ! start indes of columns (on local mpi grid)
    18   USE comgeomphy, ONLY: initcomgeomphy, &
    19                         airephy, & ! physics grid area (m2)
    20                         cuphy, & ! cu coeff. (u_covariant = cu * u)
    21                         cvphy, & ! cv coeff. (v_covariant = cv * v)
    22                         rlond, & ! longitudes
    23                         rlatd ! latitudes
     20  USE geometry_mod, ONLY : init_geometry
    2421  USE infotrac, ONLY: nqtot, type_trac
    2522  USE infotrac_phy, ONLY: init_infotrac_phy
     
    3027  USE inifis_mod, ONLY: inifis
    3128  USE phyaqua_mod, ONLY: iniaqua
     29  USE physics_distribution_mod, ONLY : init_physics_distribution
    3230  USE regular_lonlat_mod, ONLY : init_regular_lonlat, &
    3331                                 east, west, north, south, &
    3432                                 north_east, north_west, &
    3533                                 south_west, south_east
     34  USE mod_interface_dyn_phys, ONLY :  init_interface_dyn_phys
    3635  USE nrtype, ONLY: pi
    3736  IMPLICIT NONE
     
    5352  INTEGER, INTENT (IN) :: iim ! number of atmospheric coulumns along longitudes
    5453  INTEGER, INTENT (IN) :: jjm  ! number of atompsheric columns along latitudes
     54  INTEGER, INTENT(IN) :: nbp ! number of physics columns for this MPI process
     55  INTEGER, INTENT(IN) :: communicator ! MPI communicator
    5556  REAL, INTENT (IN) :: rlatu(jjm+1) ! latitudes of the physics grid
    5657  REAL, INTENT (IN) :: rlatv(jjm) ! latitude boundaries of the physics grid
     
    6566
    6667  INTEGER :: ibegin,iend,offset
    67   INTEGER :: i,j
     68  INTEGER :: i,j,k
    6869  CHARACTER (LEN=20) :: modname='iniphysiq'
    6970  CHARACTER (LEN=80) :: abort_message
     
    7576
    7677  ! global array, on full physics grid:
    77   REAL,ALLOCATABLE :: latfi(:)
    78   REAL,ALLOCATABLE :: lonfi(:)
    79   REAL,ALLOCATABLE :: cufi(:)
    80   REAL,ALLOCATABLE :: cvfi(:)
    81   REAL,ALLOCATABLE :: airefi(:)
    82 
    83   IF (nlayer.NE.klev) THEN
    84     WRITE(lunout,*) 'STOP in ',trim(modname)
    85     WRITE(lunout,*) 'Problem with dimensions :'
    86     WRITE(lunout,*) 'nlayer     = ',nlayer
    87     WRITE(lunout,*) 'klev   = ',klev
    88     abort_message = ''
    89     CALL abort_gcm (modname,abort_message,1)
    90   ENDIF
    91 
    92   !call init_phys_lmdz(iim,jjm+1,llm,1,(/(jjm-1)*iim+2/))
     78  REAL,ALLOCATABLE :: latfi_glo(:)
     79  REAL,ALLOCATABLE :: lonfi_glo(:)
     80  REAL,ALLOCATABLE :: cufi_glo(:)
     81  REAL,ALLOCATABLE :: cvfi_glo(:)
     82  REAL,ALLOCATABLE :: airefi_glo(:)
     83  REAL,ALLOCATABLE :: boundslonfi_glo(:,:)
     84  REAL,ALLOCATABLE :: boundslatfi_glo(:,:)
     85
     86  ! local arrays, on given MPI/OpenMP domain:
     87  REAL,ALLOCATABLE,SAVE :: latfi(:)
     88  REAL,ALLOCATABLE,SAVE :: lonfi(:)
     89  REAL,ALLOCATABLE,SAVE :: cufi(:)
     90  REAL,ALLOCATABLE,SAVE :: cvfi(:)
     91  REAL,ALLOCATABLE,SAVE :: airefi(:)
     92  REAL,ALLOCATABLE,SAVE :: boundslonfi(:,:)
     93  REAL,ALLOCATABLE,SAVE :: boundslatfi(:,:)
     94!$OMP THREADPRIVATE (latfi,lonfi,cufi,cvfi,airefi,boundslonfi,boundslatfi)
     95
     96  ! Initialize Physics distibution and parameters and interface with dynamics
     97  CALL init_physics_distribution(regular_lonlat,4, &
     98                                 nbp,iim,jjm+1,nlayer,communicator)
     99  CALL init_interface_dyn_phys
    93100 
    94101  ! init regular global longitude-latitude grid points and boundaries
     
    115122
    116123  ! Generate global arrays on full physics grid
    117   ALLOCATE(latfi(klon_glo),lonfi(klon_glo),cufi(klon_glo),cvfi(klon_glo))
    118   ALLOCATE(airefi(klon_glo))
     124  ALLOCATE(latfi_glo(klon_glo),lonfi_glo(klon_glo))
     125  ALLOCATE(cufi_glo(klon_glo),cvfi_glo(klon_glo))
     126  ALLOCATE(airefi_glo(klon_glo))
     127  ALLOCATE(boundslonfi_glo(klon_glo,4))
     128  ALLOCATE(boundslatfi_glo(klon_glo,4))
    119129
    120130    ! North pole
    121     latfi(1)=rlatu(1)
    122     lonfi(1)=0.
    123     cufi(1) = cu(1)
    124     cvfi(1) = cv(1)
     131    latfi_glo(1)=rlatu(1)
     132    lonfi_glo(1)=0.
     133    cufi_glo(1) = cu(1)
     134    cvfi_glo(1) = cv(1)
     135    boundslonfi_glo(1,north_east)=0
     136    boundslatfi_glo(1,north_east)=PI/2
     137    boundslonfi_glo(1,north_west)=2*PI
     138    boundslatfi_glo(1,north_west)=PI/2
     139    boundslonfi_glo(1,south_west)=2*PI
     140    boundslatfi_glo(1,south_west)=rlatv(1)
     141    boundslonfi_glo(1,south_east)=0
     142    boundslatfi_glo(1,south_east)=rlatv(1)
    125143    DO j=2,jjm
    126144      DO i=1,iim
    127         latfi((j-2)*iim+1+i)= rlatu(j)
    128         lonfi((j-2)*iim+1+i)= rlonv(i)
    129         cufi((j-2)*iim+1+i) = cu((j-1)*iim+1+i)
    130         cvfi((j-2)*iim+1+i) = cv((j-1)*iim+1+i)
     145        k=(j-2)*iim+1+i
     146        latfi_glo(k)= rlatu(j)
     147        lonfi_glo(k)= rlonv(i)
     148        cufi_glo(k) = cu((j-1)*iim+1+i)
     149        cvfi_glo(k) = cv((j-1)*iim+1+i)
     150        boundslonfi_glo(k,north_east)=rlonu(i)
     151        boundslatfi_glo(k,north_east)=rlatv(j-1)
     152        boundslonfi_glo(k,north_west)=rlonu(i+1)
     153        boundslatfi_glo(k,north_west)=rlatv(j-1)
     154        boundslonfi_glo(k,south_west)=rlonu(i+1)
     155        boundslatfi_glo(k,south_west)=rlatv(j)
     156        boundslonfi_glo(k,south_east)=rlonu(i)
     157        boundslatfi_glo(k,south_east)=rlatv(j)
    131158      ENDDO
    132159    ENDDO
    133160    ! South pole
    134     latfi(klon_glo)= rlatu(jjm+1)
    135     lonfi(klon_glo)= 0.
    136     cufi(klon_glo) = cu((iim+1)*jjm+1)
    137     cvfi(klon_glo) = cv((iim+1)*jjm-iim)
     161    latfi_glo(klon_glo)= rlatu(jjm+1)
     162    lonfi_glo(klon_glo)= 0.
     163    cufi_glo(klon_glo) = cu((iim+1)*jjm+1)
     164    cvfi_glo(klon_glo) = cv((iim+1)*jjm-iim)
     165    boundslonfi_glo(klon_glo,north_east)= 0
     166    boundslatfi_glo(klon_glo,north_east)= rlatv(jjm)
     167    boundslonfi_glo(klon_glo,north_west)= 2*PI
     168    boundslatfi_glo(klon_glo,north_west)= rlatv(jjm)
     169    boundslonfi_glo(klon_glo,south_west)= 2*PI
     170    boundslatfi_glo(klon_glo,south_west)= -PI/2
     171    boundslonfi_glo(klon_glo,south_east)= 0
     172    boundslatfi_glo(klon_glo,south_east)= -Pi/2
    138173
    139174    ! build airefi(), mesh area on physics grid
    140     CALL gr_dyn_fi(1,iim+1,jjm+1,klon_glo,aire,airefi)
     175    CALL gr_dyn_fi(1,iim+1,jjm+1,klon_glo,aire,airefi_glo)
    141176    ! Poles are single points on physics grid
    142     airefi(1)=sum(aire(1:iim,1))
    143     airefi(klon_glo)=sum(aire(1:iim,jjm+1))
     177    airefi_glo(1)=sum(aire(1:iim,1))
     178    airefi_glo(klon_glo)=sum(aire(1:iim,jjm+1))
    144179
    145180    ! Sanity check: do total planet area match between physics and dynamics?
    146181    total_area_dyn=sum(aire(1:iim,1:jjm+1))
    147     total_area_phy=sum(airefi(1:klon_glo))
     182    total_area_phy=sum(airefi_glo(1:klon_glo))
    148183    IF (total_area_dyn/=total_area_phy) THEN
    149184      WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
     
    158193
    159194!$OMP PARALLEL
     195  ! Now generate local lon/lat/cu/cv/area/bounds arrays
     196  ALLOCATE(latfi(klon_omp),lonfi(klon_omp),cufi(klon_omp),cvfi(klon_omp))
     197  ALLOCATE(airefi(klon_omp))
     198  ALLOCATE(boundslonfi(klon_omp,4))
     199  ALLOCATE(boundslatfi(klon_omp,4))
     200
     201
     202  offset = klon_mpi_begin - 1
     203  airefi(1:klon_omp) = airefi_glo(offset+klon_omp_begin:offset+klon_omp_end)
     204  cufi(1:klon_omp) = cufi_glo(offset+klon_omp_begin:offset+klon_omp_end)
     205  cvfi(1:klon_omp) = cvfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
     206  lonfi(1:klon_omp) = lonfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
     207  latfi(1:klon_omp) = latfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
     208  boundslonfi(1:klon_omp,:) = boundslonfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
     209  boundslatfi(1:klon_omp,:) = boundslatfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
     210
     211  ! copy over local grid longitudes and latitudes
     212  CALL init_geometry(klon_omp,lonfi,latfi,boundslonfi,boundslatfi, &
     213                     airefi,cufi,cvfi)
     214
    160215  ! Initialize physical constants in physics:
    161216  CALL inifis(prad,pg,pr,pcpp)
     
    164219  CALL init_infotrac_phy(nqtot,type_trac)
    165220
    166   ! Now generate local lon/lat/cu/cv/area arrays
    167   CALL initcomgeomphy
    168 
    169   offset = klon_mpi_begin - 1
    170   airephy(1:klon_omp) = airefi(offset+klon_omp_begin:offset+klon_omp_end)
    171   cuphy(1:klon_omp) = cufi(offset+klon_omp_begin:offset+klon_omp_end)
    172   cvphy(1:klon_omp) = cvfi(offset+klon_omp_begin:offset+klon_omp_end)
    173   rlond(1:klon_omp) = lonfi(offset+klon_omp_begin:offset+klon_omp_end)
    174   rlatd(1:klon_omp) = latfi(offset+klon_omp_begin:offset+klon_omp_end)
    175 
    176221  ! Additional initializations for aquaplanets
    177222  IF (iflag_phys>=100) THEN
    178     CALL iniaqua(klon_omp,rlatd,rlond,iflag_phys)
     223    CALL iniaqua(klon_omp,iflag_phys)
    179224  ENDIF
    180225!$OMP END PARALLEL
  • LMDZ5/trunk/libf/dynlonlat_phylonlat/phylmd/ce0l.F90

    r2349 r2351  
    2626  USE filtreg_mod,    ONLY: inifilr
    2727  USE iniphysiq_mod,  ONLY: iniphysiq
     28  USE mod_const_mpi,  ONLY: comm_lmdz
    2829#ifdef inca
    2930  USE indice_sol_mod, ONLY: nbsrf, is_oce, is_sic, is_ter, is_lic
     
    106107  CALL read_distrib()
    107108  CALL init_mod_hallo()
    108   CALL Init_Phys_lmdz(iim,jjp1,llm,mpi_size,distrib_phys)
    109   CALL init_interface_dyn_phys()
    110 #else
    111   CALL Init_Phys_lmdz(iim,jjp1,llm,1,(/(jjm-1)*iim+2/))
    112109#endif
    113110  WRITE(lunout,*)'---> klon=',klon
     
    126123
    127124  CALL inifilr()
    128   CALL iniphysiq(iim,jjm,llm,daysec,day_ini,dtphys/nsplit_phys, &
     125  CALL iniphysiq(iim,jjm,llm, &
     126                 distrib_phys(mpi_rank),comm_lmdz, &
     127                 daysec,day_ini,dtphys/nsplit_phys, &
    129128                 rlatu,rlatv,rlonu,rlonv,aire,cu,cv,rad,g,r,cpp,iflag_phys)
    130129  IF(pressure_exner) CALL test_disvert
  • LMDZ5/trunk/libf/dynlonlat_phylonlat/phylmd/iniphysiq_mod.F90

    r2347 r2351  
    66CONTAINS
    77
    8 SUBROUTINE iniphysiq(ii,jj,nlayer,punjours, pdayref,ptimestep, &
     8SUBROUTINE iniphysiq(ii,jj,nlayer, &
     9                     nbp, communicator, &
     10                     punjours, pdayref,ptimestep, &
    911                     rlatu,rlatv,rlonu,rlonv,aire,cu,cv,       &
    1012                     prad,pg,pr,pcpp,iflag_phys)
    11   USE dimphy, ONLY: klev ! number of atmospheric levels
    12   USE mod_grid_phy_lmdz, ONLY: klon_glo ! number of atmospheric columns
    13                                         ! (on full grid)
     13  USE dimphy, ONLY: init_dimphy
     14  USE mod_grid_phy_lmdz, ONLY: klon_glo,  & ! number of atmospheric columns (on full grid)
     15                               regular_lonlat  ! regular longitude-latitude grid type
    1416  USE mod_phys_lmdz_para, ONLY: klon_omp, & ! number of columns (on local omp grid)
    1517                                klon_omp_begin, & ! start index of local omp subgrid
    1618                                klon_omp_end, & ! end index of local omp subgrid
    1719                                klon_mpi_begin ! start indes of columns (on local mpi grid)
     20  USE geometry_mod, ONLY : init_geometry
    1821  USE vertical_layers_mod, ONLY : init_vertical_layers
    1922  USE infotrac, ONLY: nqtot,nqo,nbtr,tname,ttext,type_trac,&
     
    2629                      indnum_fn_num,index_trac,&
    2730                      niso,ntraceurs_zone,ntraciso
     31#ifdef REPROBUS
     32  USE CHEM_REP, ONLY : Init_chem_rep_phys
     33#endif
    2834  USE control_mod, ONLY: dayref,anneeref,day_step,nday,offline
    29   USE comgeomphy, ONLY: initcomgeomphy, &
    30                         airephy, & ! physics grid area (m2)
    31                         cuphy, & ! cu coeff. (u_covariant = cu * u)
    32                         cvphy, & ! cv coeff. (v_covariant = cv * v)
    33                         rlond, & ! longitudes
    34                         rlatd ! latitudes
    3535  USE inifis_mod, ONLY: inifis
    3636  USE time_phylmdz_mod, ONLY: init_time
     
    3838  USE phystokenc_mod, ONLY: init_phystokenc
    3939  USE phyaqua_mod, ONLY: iniaqua
     40  USE physics_distribution_mod, ONLY : init_physics_distribution
    4041  USE regular_lonlat_mod, ONLY : init_regular_lonlat, &
    4142                                 east, west, north, south, &
    4243                                 north_east, north_west, &
    4344                                 south_west, south_east
     45  USE mod_interface_dyn_phys, ONLY :  init_interface_dyn_phys
    4446  IMPLICIT NONE
    4547
     
    6466  INTEGER, INTENT (IN) :: ii ! number of atmospheric columns along longitudes
    6567  INTEGER, INTENT (IN) :: jj ! number of atompsheric columns along latitudes
     68  INTEGER, INTENT(IN) :: nbp ! number of physics columns for this MPI process
     69  INTEGER, INTENT(IN) :: communicator ! MPI communicator
    6670  REAL, INTENT (IN) :: rlatu(jj+1) ! latitudes of the physics grid
    6771  REAL, INTENT (IN) :: rlatv(jj) ! latitude boundaries of the physics grid
     
    7680
    7781  INTEGER :: ibegin, iend, offset
    78   INTEGER :: i,j
     82  INTEGER :: i,j,k
    7983  CHARACTER (LEN=20) :: modname = 'iniphysiq'
    8084  CHARACTER (LEN=80) :: abort_message
     
    8690
    8791  ! global array, on full physics grid:
    88   REAL,ALLOCATABLE :: latfi(:)
    89   REAL,ALLOCATABLE :: lonfi(:)
    90   REAL,ALLOCATABLE :: cufi(:)
    91   REAL,ALLOCATABLE :: cvfi(:)
    92   REAL,ALLOCATABLE :: airefi(:)
    93 
    94   IF (nlayer/=klev) THEN
    95     WRITE (lunout, *) 'nlayer     = ', nlayer
    96     WRITE (lunout, *) 'klev   = ', klev
    97     CALL abort_gcm(modname, 'Problem with dimensions', 1)
    98   END IF
    99 
    100   !call init_phys_lmdz(ii,jj+1,llm,1,(/(jj-1)*ii+2/))
     92  REAL,ALLOCATABLE :: latfi_glo(:)
     93  REAL,ALLOCATABLE :: lonfi_glo(:)
     94  REAL,ALLOCATABLE :: cufi_glo(:)
     95  REAL,ALLOCATABLE :: cvfi_glo(:)
     96  REAL,ALLOCATABLE :: airefi_glo(:)
     97  REAL,ALLOCATABLE :: boundslonfi_glo(:,:)
     98  REAL,ALLOCATABLE :: boundslatfi_glo(:,:)
     99
     100  ! local arrays, on given MPI/OpenMP domain:
     101  REAL,ALLOCATABLE,SAVE :: latfi(:)
     102  REAL,ALLOCATABLE,SAVE :: lonfi(:)
     103  REAL,ALLOCATABLE,SAVE :: cufi(:)
     104  REAL,ALLOCATABLE,SAVE :: cvfi(:)
     105  REAL,ALLOCATABLE,SAVE :: airefi(:)
     106  REAL,ALLOCATABLE,SAVE :: boundslonfi(:,:)
     107  REAL,ALLOCATABLE,SAVE :: boundslatfi(:,:)
     108!$OMP THREADPRIVATE (latfi,lonfi,cufi,cvfi,airefi,boundslonfi,boundslatfi)
     109
     110  ! Initialize Physics distibution and parameters and interface with dynamics
     111  CALL init_physics_distribution(regular_lonlat,4, &
     112                                 nbp,ii,jj+1,nlayer,communicator)
     113  CALL init_interface_dyn_phys
    101114 
    102115  ! init regular global longitude-latitude grid points and boundaries
     
    123136
    124137  ! Generate global arrays on full physics grid
    125   ALLOCATE(latfi(klon_glo),lonfi(klon_glo),cufi(klon_glo),cvfi(klon_glo))
    126   ALLOCATE(airefi(klon_glo))
     138  ALLOCATE(latfi_glo(klon_glo),lonfi_glo(klon_glo))
     139  ALLOCATE(cufi_glo(klon_glo),cvfi_glo(klon_glo))
     140  ALLOCATE(airefi_glo(klon_glo))
     141  ALLOCATE(boundslonfi_glo(klon_glo,4))
     142  ALLOCATE(boundslatfi_glo(klon_glo,4))
     143 
    127144
    128145  IF (klon_glo>1) THEN ! general case
    129146    ! North pole
    130     latfi(1)=rlatu(1)
    131     lonfi(1)=0.
    132     cufi(1) = cu(1)
    133     cvfi(1) = cv(1)
     147    latfi_glo(1)=rlatu(1)
     148    lonfi_glo(1)=0.
     149    cufi_glo(1) = cu(1)
     150    cvfi_glo(1) = cv(1)
     151    boundslonfi_glo(1,north_east)=0
     152    boundslatfi_glo(1,north_east)=PI/2
     153    boundslonfi_glo(1,north_west)=2*PI
     154    boundslatfi_glo(1,north_west)=PI/2
     155    boundslonfi_glo(1,south_west)=2*PI
     156    boundslatfi_glo(1,south_west)=rlatv(1)
     157    boundslonfi_glo(1,south_east)=0
     158    boundslatfi_glo(1,south_east)=rlatv(1)
    134159    DO j=2,jj
    135160      DO i=1,ii
    136         latfi((j-2)*ii+1+i)= rlatu(j)
    137         lonfi((j-2)*ii+1+i)= rlonv(i)
    138         cufi((j-2)*ii+1+i) = cu((j-1)*(ii+1)+i)
    139         cvfi((j-2)*ii+1+i) = cv((j-1)*(ii+1)+i)
     161        k=(j-2)*ii+1+i
     162        latfi_glo(k)= rlatu(j)
     163        lonfi_glo(k)= rlonv(i)
     164        cufi_glo(k) = cu((j-1)*ii+1+i)
     165        cvfi_glo(k) = cv((j-1)*ii+1+i)
     166        boundslonfi_glo(k,north_east)=rlonu(i)
     167        boundslatfi_glo(k,north_east)=rlatv(j-1)
     168        boundslonfi_glo(k,north_west)=rlonu(i+1)
     169        boundslatfi_glo(k,north_west)=rlatv(j-1)
     170        boundslonfi_glo(k,south_west)=rlonu(i+1)
     171        boundslatfi_glo(k,south_west)=rlatv(j)
     172        boundslonfi_glo(k,south_east)=rlonu(i)
     173        boundslatfi_glo(k,south_east)=rlatv(j)
    140174      ENDDO
    141175    ENDDO
    142176    ! South pole
    143     latfi(klon_glo)= rlatu(jj+1)
    144     lonfi(klon_glo)= 0.
    145     cufi(klon_glo) = cu((ii+1)*jj+1)
    146     cvfi(klon_glo) = cv((ii+1)*jj-ii)
    147 
    148     ! build airefi(), mesh area on physics grid
    149     CALL gr_dyn_fi(1,ii+1,jj+1,klon_glo,aire,airefi)
     177    latfi_glo(klon_glo)= rlatu(jj+1)
     178    lonfi_glo(klon_glo)= 0.
     179    cufi_glo(klon_glo) = cu((ii+1)*jj+1)
     180    cvfi_glo(klon_glo) = cv((ii+1)*jj-ii)
     181    boundslonfi_glo(klon_glo,north_east)= 0
     182    boundslatfi_glo(klon_glo,north_east)= rlatv(jj)
     183    boundslonfi_glo(klon_glo,north_west)= 2*PI
     184    boundslatfi_glo(klon_glo,north_west)= rlatv(jj)
     185    boundslonfi_glo(klon_glo,south_west)= 2*PI
     186    boundslatfi_glo(klon_glo,south_west)= -PI/2
     187    boundslonfi_glo(klon_glo,south_east)= 0
     188    boundslatfi_glo(klon_glo,south_east)= -Pi/2
     189
     190    ! build airefi_glo(), mesh area on physics grid
     191    CALL gr_dyn_fi(1,ii+1,jj+1,klon_glo,aire,airefi_glo)
    150192    ! Poles are single points on physics grid
    151     airefi(1)=sum(aire(1:ii,1))
    152     airefi(klon_glo)=sum(aire(1:ii,jj+1))
     193    airefi_glo(1)=sum(aire(1:ii,1))
     194    airefi_glo(klon_glo)=sum(aire(1:ii,jj+1))
    153195
    154196    ! Sanity check: do total planet area match between physics and dynamics?
    155197    total_area_dyn=sum(aire(1:ii,1:jj+1))
    156     total_area_phy=sum(airefi(1:klon_glo))
     198    total_area_phy=sum(airefi_glo(1:klon_glo))
    157199    IF (total_area_dyn/=total_area_phy) THEN
    158200      WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
     
    167209  ELSE ! klon_glo==1, running the 1D model
    168210    ! just copy over input values
    169     latfi(1)=rlatu(1)
    170     lonfi(1)=rlonv(1)
    171     cufi(1)=cu(1)
    172     cvfi(1)=cv(1)
    173     airefi(1)=aire(1,1)
     211    latfi_glo(1)=rlatu(1)
     212    lonfi_glo(1)=rlonv(1)
     213    cufi_glo(1)=cu(1)
     214    cvfi_glo(1)=cv(1)
     215    airefi_glo(1)=aire(1,1)
     216    boundslonfi_glo(1,north_east)=rlonu(1)
     217    boundslatfi_glo(1,north_east)=PI/2
     218    boundslonfi_glo(1,north_west)=rlonu(2)
     219    boundslatfi_glo(1,north_west)=PI/2
     220    boundslonfi_glo(1,south_west)=rlonu(2)
     221    boundslatfi_glo(1,south_west)=rlatv(1)
     222    boundslonfi_glo(1,south_east)=rlonu(1)
     223    boundslatfi_glo(1,south_east)=rlatv(1)
    174224  ENDIF ! of IF (klon_glo>1)
    175225
    176226!$OMP PARALLEL DEFAULT(SHARED) COPYIN(/temps/)
    177   ! Initialize physical constants in physics:
    178   CALL inifis(punjours,prad,pg,pr,pcpp)
    179   CALL init_time(annee_ref,day_ref,day_ini,start_time,nday,ptimestep)
    180 
    181   ! Copy over "offline" settings
    182   CALL init_phystokenc(offline,istphy)
     227  ! Now generate local lon/lat/cu/cv/area/bounds arrays
     228  ALLOCATE(latfi(klon_omp),lonfi(klon_omp),cufi(klon_omp),cvfi(klon_omp))
     229  ALLOCATE(airefi(klon_omp))
     230  ALLOCATE(boundslonfi(klon_omp,4))
     231  ALLOCATE(boundslatfi(klon_omp,4))
     232
     233
     234  offset = klon_mpi_begin - 1
     235  airefi(1:klon_omp) = airefi_glo(offset+klon_omp_begin:offset+klon_omp_end)
     236  cufi(1:klon_omp) = cufi_glo(offset+klon_omp_begin:offset+klon_omp_end)
     237  cvfi(1:klon_omp) = cvfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
     238  lonfi(1:klon_omp) = lonfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
     239  latfi(1:klon_omp) = latfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
     240  boundslonfi(1:klon_omp,:) = boundslonfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
     241  boundslatfi(1:klon_omp,:) = boundslatfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
     242
     243  ! copy over local grid longitudes and latitudes
     244  CALL init_geometry(klon_omp,lonfi,latfi,boundslonfi,boundslatfi, &
     245                     airefi,cufi,cvfi)
    183246
    184247  ! copy over preff , ap(), bp(), etc
    185248  CALL init_vertical_layers(nlayer,preff,scaleheight, &
    186249                            ap,bp,presnivs,pseudoalt)
     250
     251  ! Initialize physical constants in physics:
     252  CALL inifis(punjours,prad,pg,pr,pcpp)
     253
     254  CALL init_time(annee_ref,day_ref,day_ini,start_time,nday,ptimestep)
     255
     256  ! Initialize dimphy module
     257  CALL Init_dimphy(klon_omp,nlayer)
     258
     259  ! Copy over "offline" settings
     260  CALL init_phystokenc(offline,istphy)
    187261
    188262  ! Initialize tracer names, numbers, etc. for physics
     
    197271                         niso,ntraceurs_zone,ntraciso)
    198272
    199   ! Now generate local lon/lat/cu/cv/area arrays
    200   CALL initcomgeomphy
    201 
    202   offset = klon_mpi_begin - 1
    203   airephy(1:klon_omp) = airefi(offset+klon_omp_begin:offset+klon_omp_end)
    204   cuphy(1:klon_omp) = cufi(offset+klon_omp_begin:offset+klon_omp_end)
    205   cvphy(1:klon_omp) = cvfi(offset+klon_omp_begin:offset+klon_omp_end)
    206   rlond(1:klon_omp) = lonfi(offset+klon_omp_begin:offset+klon_omp_end)
    207   rlatd(1:klon_omp) = latfi(offset+klon_omp_begin:offset+klon_omp_end)
     273  ! Initializations for Reprobus
     274  IF (type_trac == 'repr') THEN
     275#ifdef REPROBUS
     276    CALL Init_chem_rep_phys(klon_omp,nlayer)
     277#endif
     278  ENDIF
    208279
    209280  ! Additional initializations for aquaplanets
    210281  IF (iflag_phys>=100) THEN
    211     CALL iniaqua(klon_omp, rlatd, rlond, iflag_phys)
     282    CALL iniaqua(klon_omp,iflag_phys)
    212283  END IF
    213284!$OMP END PARALLEL
  • LMDZ5/trunk/libf/dynlonlat_phylonlat/phymar/iniphysiq_mod.F90

    r2347 r2351  
    66CONTAINS
    77
    8 SUBROUTINE iniphysiq(iim,jjm,nlayer,punjours, pdayref,ptimestep, &
     8SUBROUTINE iniphysiq(iim,jjm,nlayer,
     9                     nbp, communicator, &
     10                     punjours, pdayref,ptimestep, &
    911                     rlatu,rlatv,rlonu,rlonv,aire,cu,cv,         &
    1012                     prad,pg,pr,pcpp,iflag_phys)
    11   USE dimphy, ONLY: klev ! number of atmospheric levels
    12   USE mod_grid_phy_lmdz, ONLY: klon_glo ! number of atmospheric columns
    13                                         ! (on full grid)
     13  USE dimphy, ONLY: init_dimphy
     14  USE mod_grid_phy_lmdz, ONLY: klon_glo,  & ! number of atmospheric columns (on full grid)
     15                               regular_lonlat  ! regular longitude-latitude grid type
    1416  USE mod_phys_lmdz_para, ONLY: klon_omp, & ! number of columns (on local omp grid)
    1517                                klon_omp_begin, & ! start index of local omp subgrid
    1618                                klon_omp_end, & ! end index of local omp subgrid
    1719                                klon_mpi_begin ! start indes of columns (on local mpi grid)
    18   USE comgeomphy, ONLY: initcomgeomphy, &
    19                         airephy, & ! physics grid area (m2)
    20                         cuphy, & ! cu coeff. (u_covariant = cu * u)
    21                         cvphy, & ! cv coeff. (v_covariant = cv * v)
    22                         rlond, & ! longitudes
    23                         rlatd ! latitudes
     20  USE geometry_mod, ONLY : init_geometry
    2421  USE comcstphy, ONLY: rradius, & ! planet radius (m)
    2522                       rr, & ! recuced gas constant: R/molar mass of atm
     
    2724                       rcpp  ! specific heat of the atmosphere
    2825!  USE phyaqua_mod, ONLY: iniaqua
     26  USE physics_distribution_mod, ONLY : init_physics_distribution
    2927  USE regular_lonlat_mod, ONLY : init_regular_lonlat, &
    3028                                 east, west, north, south, &
    3129                                 north_east, north_west, &
    3230                                 south_west, south_east
     31  USE mod_interface_dyn_phys, ONLY :  init_interface_dyn_phys
    3332  USE nrtype, ONLY: pi
    3433  IMPLICIT NONE
     
    5049  INTEGER, INTENT (IN) :: iim ! number of atmospheric coulumns along longitudes
    5150  INTEGER, INTENT (IN) :: jjm  ! number of atompsheric columns along latitudes
     51  INTEGER, INTENT(IN) :: nbp ! number of physics columns for this MPI process
     52  INTEGER, INTENT(IN) :: communicator ! MPI communicator
    5253  REAL, INTENT (IN) :: rlatu(jjm+1) ! latitudes of the physics grid
    5354  REAL, INTENT (IN) :: rlatv(jjm) ! latitude boundaries of the physics grid
     
    6263
    6364  INTEGER :: ibegin,iend,offset
    64   INTEGER :: i,j
     65  INTEGER :: i,j,k
    6566  CHARACTER (LEN=20) :: modname='iniphysiq'
    6667  CHARACTER (LEN=80) :: abort_message
     
    7273
    7374  ! global array, on full physics grid:
    74   REAL,ALLOCATABLE :: latfi(:)
    75   REAL,ALLOCATABLE :: lonfi(:)
    76   REAL,ALLOCATABLE :: cufi(:)
    77   REAL,ALLOCATABLE :: cvfi(:)
    78   REAL,ALLOCATABLE :: airefi(:)
    79  
    80   IF (nlayer.NE.klev) THEN
    81     WRITE(lunout,*) 'STOP in ',trim(modname)
    82     WRITE(lunout,*) 'Problem with dimensions :'
    83     WRITE(lunout,*) 'nlayer     = ',nlayer
    84     WRITE(lunout,*) 'klev   = ',klev
    85     abort_message = ''
    86     CALL abort_gcm (modname,abort_message,1)
    87   ENDIF
    88 
     75  REAL,ALLOCATABLE :: latfi_glo(:)
     76  REAL,ALLOCATABLE :: lonfi_glo(:)
     77  REAL,ALLOCATABLE :: cufi_glo(:)
     78  REAL,ALLOCATABLE :: cvfi_glo(:)
     79  REAL,ALLOCATABLE :: airefi_glo(:)
     80  REAL,ALLOCATABLE :: boundslonfi_glo(:,:)
     81  REAL,ALLOCATABLE :: boundslatfi_glo(:,:)
     82
     83  ! local arrays, on given MPI/OpenMP domain:
     84  REAL,ALLOCATABLE,SAVE :: latfi(:)
     85  REAL,ALLOCATABLE,SAVE :: lonfi(:)
     86  REAL,ALLOCATABLE,SAVE :: cufi(:)
     87  REAL,ALLOCATABLE,SAVE :: cvfi(:)
     88  REAL,ALLOCATABLE,SAVE :: airefi(:)
     89  REAL,ALLOCATABLE,SAVE :: boundslonfi(:,:)
     90  REAL,ALLOCATABLE,SAVE :: boundslatfi(:,:)
     91!$OMP THREADPRIVATE (latfi,lonfi,cufi,cvfi,airefi,boundslonfi,boundslatfi)
     92
     93  ! Initialize Physics distibution and parameters and interface with dynamics
     94  CALL init_physics_distribution(regular_lonlat,4, &
     95                                 nbp,ii,jj+1,nlayer,communicator)
     96  CALL init_interface_dyn_phys
     97 
    8998  ! init regular global longitude-latitude grid points and boundaries
    9099  ALLOCATE(boundslon_reg(iim,2))
     
    110119
    111120  ! Generate global arrays on full physics grid
    112   ALLOCATE(latfi(klon_glo),lonfi(klon_glo),cufi(klon_glo),cvfi(klon_glo))
    113   ALLOCATE(airefi(klon_glo))
     121  ALLOCATE(latfi_glo(klon_glo),lonfi_glo(klon_glo))
     122  ALLOCATE(cufi_glo(klon_glo),cvfi_glo(klon_glo))
     123  ALLOCATE(airefi_glo(klon_glo))
     124  ALLOCATE(boundslonfi_glo(klon_glo,4))
     125  ALLOCATE(boundslatfi_glo(klon_glo,4))
    114126
    115127  IF (klon_glo>1) THEN ! general case
    116128    ! North pole
    117     latfi(1)=rlatu(1)
    118     lonfi(1)=0.
    119     cufi(1) = cu(1)
    120     cvfi(1) = cv(1)
     129    latfi_glo(1)=rlatu(1)
     130    lonfi_glo(1)=0.
     131    cufi_glo(1) = cu(1)
     132    cvfi_glo(1) = cv(1)
     133    boundslonfi_glo(1,north_east)=0
     134    boundslatfi_glo(1,north_east)=PI/2
     135    boundslonfi_glo(1,north_west)=2*PI
     136    boundslatfi_glo(1,north_west)=PI/2
     137    boundslonfi_glo(1,south_west)=2*PI
     138    boundslatfi_glo(1,south_west)=rlatv(1)
     139    boundslonfi_glo(1,south_east)=0
     140    boundslatfi_glo(1,south_east)=rlatv(1)
    121141    DO j=2,jjm
    122142      DO i=1,iim
    123         latfi((j-2)*iim+1+i)= rlatu(j)
    124         lonfi((j-2)*iim+1+i)= rlonv(i)
    125         cufi((j-2)*iim+1+i) = cu((j-1)*iim+1+i)
    126         cvfi((j-2)*iim+1+i) = cv((j-1)*iim+1+i)
     143        k=(j-2)*iim+1+i
     144        latfi_glo((j-2)*iim+1+i)= rlatu(j)
     145        lonfi_glo((j-2)*iim+1+i)= rlonv(i)
     146        cufi_glo((j-2)*iim+1+i) = cu((j-1)*iim+1+i)
     147        cvfi_glo((j-2)*iim+1+i) = cv((j-1)*iim+1+i)
     148        boundslonfi_glo(k,north_east)=rlonu(i)
     149        boundslatfi_glo(k,north_east)=rlatv(j-1)
     150        boundslonfi_glo(k,north_west)=rlonu(i+1)
     151        boundslatfi_glo(k,north_west)=rlatv(j-1)
     152        boundslonfi_glo(k,south_west)=rlonu(i+1)
     153        boundslatfi_glo(k,south_west)=rlatv(j)
     154        boundslonfi_glo(k,south_east)=rlonu(i)
     155        boundslatfi_glo(k,south_east)=rlatv(j)
    127156      ENDDO
    128157    ENDDO
    129158    ! South pole
    130     latfi(klon_glo)= rlatu(jjm+1)
    131     lonfi(klon_glo)= 0.
    132     cufi(klon_glo) = cu((iim+1)*jjm+1)
    133     cvfi(klon_glo) = cv((iim+1)*jjm-iim)
     159    latfi_glo(klon_glo)= rlatu(jjm+1)
     160    lonfi_glo(klon_glo)= 0.
     161    cufi_glo(klon_glo) = cu((iim+1)*jjm+1)
     162    cvfi_glo(klon_glo) = cv((iim+1)*jjm-iim)
     163    boundslonfi_glo(klon_glo,north_east)= 0
     164    boundslatfi_glo(klon_glo,north_east)= rlatv(jjm)
     165    boundslonfi_glo(klon_glo,north_west)= 2*PI
     166    boundslatfi_glo(klon_glo,north_west)= rlatv(jjm)
     167    boundslonfi_glo(klon_glo,south_west)= 2*PI
     168    boundslatfi_glo(klon_glo,south_west)= -PI/2
     169    boundslonfi_glo(klon_glo,south_east)= 0
     170    boundslatfi_glo(klon_glo,south_east)= -Pi/2
    134171
    135172    ! build airefi(), mesh area on physics grid
    136     CALL gr_dyn_fi(1,iim+1,jjm+1,klon_glo,aire,airefi)
     173    CALL gr_dyn_fi(1,iim+1,jjm+1,klon_glo,aire,airefi_glo)
    137174    ! Poles are single points on physics grid
    138     airefi(1)=sum(aire(1:iim,1))
    139     airefi(klon_glo)=sum(aire(1:iim,jjm+1))
     175    airefi_glo(1)=sum(aire(1:iim,1))
     176    airefi_glo(klon_glo)=sum(aire(1:iim,jjm+1))
    140177
    141178    ! Sanity check: do total planet area match between physics and dynamics?
    142179    total_area_dyn=sum(aire(1:iim,1:jjm+1))
    143     total_area_phy=sum(airefi(1:klon_glo))
     180    total_area_phy=sum(airefi_glo(1:klon_glo))
    144181    IF (total_area_dyn/=total_area_phy) THEN
    145182      WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
     
    154191  ELSE ! klon_glo==1, running the 1D model
    155192    ! just copy over input values
    156     latfi(1)=rlatu(1)
    157     lonfi(1)=rlonv(1)
    158     cufi(1)=cu(1)
    159     cvfi(1)=cv(1)
    160     airefi(1)=aire(1,1)
     193    latfi_glo(1)=rlatu(1)
     194    lonfi_glo(1)=rlonv(1)
     195    cufi_glo(1)=cu(1)
     196    cvfi_glo(1)=cv(1)
     197    airefi_glo(1)=aire(1,1)
     198    boundslonfi_glo(1,north_east)=rlonu(1)
     199    boundslatfi_glo(1,north_east)=PI/2
     200    boundslonfi_glo(1,north_west)=rlonu(2)
     201    boundslatfi_glo(1,north_west)=PI/2
     202    boundslonfi_glo(1,south_west)=rlonu(2)
     203    boundslatfi_glo(1,south_west)=rlatv(1)
     204    boundslonfi_glo(1,south_east)=rlonu(1)
     205    boundslatfi_glo(1,south_east)=rlatv(1)
    161206  ENDIF ! of IF (klon_glo>1)
    162207
    163208!$OMP PARALLEL
    164   ! Now generate local lon/lat/cu/cv/area arrays
    165   CALL initcomgeomphy
    166      
     209  ! Now generate local lon/lat/cu/cv/area/bounds arrays
     210  ALLOCATE(latfi(klon_omp),lonfi(klon_omp),cufi(klon_omp),cvfi(klon_omp))
     211  ALLOCATE(airefi(klon_omp))
     212  ALLOCATE(boundslonfi(klon_omp,4))
     213  ALLOCATE(boundslatfi(klon_omp,4))
     214
     215
    167216  offset = klon_mpi_begin - 1
    168   airephy(1:klon_omp) = airefi(offset+klon_omp_begin:offset+klon_omp_end)
    169   cuphy(1:klon_omp) = cufi(offset+klon_omp_begin:offset+klon_omp_end)
    170   cvphy(1:klon_omp) = cvfi(offset+klon_omp_begin:offset+klon_omp_end)
    171   rlond(1:klon_omp) = lonfi(offset+klon_omp_begin:offset+klon_omp_end)
    172   rlatd(1:klon_omp) = latfi(offset+klon_omp_begin:offset+klon_omp_end)
     217  airefi(1:klon_omp) = airefi_glo(offset+klon_omp_begin:offset+klon_omp_end)
     218  cufi(1:klon_omp) = cufi_glo(offset+klon_omp_begin:offset+klon_omp_end)
     219  cvfi(1:klon_omp) = cvfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
     220  lonfi(1:klon_omp) = lonfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
     221  latfi(1:klon_omp) = latfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
     222  boundslonfi(1:klon_omp,:) = boundslonfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
     223  boundslatfi(1:klon_omp,:) = boundslatfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
     224
     225  ! copy over local grid longitudes and latitudes
     226  CALL init_geometry(klon_omp,lonfi,latfi,boundslonfi,boundslatfi, &
     227                     airefi,cufi,cvfi)
     228
    173229
    174230! copy some fundamental parameters to physics
Note: See TracChangeset for help on using the changeset viewer.