Changeset 2225 for LMDZ5


Ignore:
Timestamp:
Mar 11, 2015, 3:55:23 PM (9 years ago)
Author:
Ehouarn Millour
Message:

Some cleanup and tidying up in the dynamics/physics interface.
EM

Location:
LMDZ5/trunk/libf
Files:
6 edited
1 moved

Legend:

Unmodified
Added
Removed
  • LMDZ5/trunk/libf/dyn3d/gcm.F

    r2222 r2225  
    2626! Only INCA needs these informations (from the Earth's physics)
    2727      USE indice_sol_mod
     28      USE mod_phys_lmdz_para, ONLY : klon_mpi_para_nb
    2829#endif
    2930
     
    3334! dynamique -> physique pour l'initialisation
    3435#ifdef CPP_PHYS
    35       USE dimphy
    36       USE comgeomphy
    37       USE mod_phys_lmdz_para, ONLY : klon_mpi_para_nb
     36!      USE dimphy
     37!      USE comgeomphy
    3838#endif
    3939!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
     
    119119      LOGICAL first
    120120
    121       LOGICAL call_iniphys
    122       data call_iniphys/.true./
     121!      LOGICAL call_iniphys
     122!      data call_iniphys/.true./
    123123
    124124c+jld variables test conservation energie
     
    141141      REAL :: heure
    142142
    143 
    144 c-----------------------------------------------------------------------
    145 c    variables pour l'initialisation de la physique :
    146 c    ------------------------------------------------
    147       INTEGER ngridmx
    148       PARAMETER( ngridmx = 2+(jjm-1)*iim - 1/jjm   )
    149       REAL zcufi(ngridmx),zcvfi(ngridmx)
    150       REAL latfi(ngridmx),lonfi(ngridmx)
    151       REAL airefi(ngridmx)
    152       SAVE latfi, lonfi, airefi
    153 
    154143c-----------------------------------------------------------------------
    155144c   Initialisations:
     
    188177#ifdef CPP_PHYS
    189178      CALL Init_Phys_lmdz(iim,jjp1,llm,1,(/(jjm-1)*iim+2/))
    190       call InitComgeomphy
     179!      call InitComgeomphy ! now done in iniphysiq
    191180#endif
    192181!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
     
    417406c   -------------------------------
    418407
    419       IF (call_iniphys.and.(iflag_phys==1.or.iflag_phys>=100)) THEN
    420          latfi(1)=rlatu(1)
    421          lonfi(1)=0.
    422          zcufi(1) = cu(1)
    423          zcvfi(1) = cv(1)
    424          DO j=2,jjm
    425             DO i=1,iim
    426                latfi((j-2)*iim+1+i)= rlatu(j)
    427                lonfi((j-2)*iim+1+i)= rlonv(i)
    428                zcufi((j-2)*iim+1+i) = cu((j-1)*iip1+i)
    429                zcvfi((j-2)*iim+1+i) = cv((j-1)*iip1+i)
    430             ENDDO
    431          ENDDO
    432          latfi(ngridmx)= rlatu(jjp1)
    433          lonfi(ngridmx)= 0.
    434          zcufi(ngridmx) = cu(ip1jm+1)
    435          zcvfi(ngridmx) = cv(ip1jm-iim)
    436          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,aire,airefi)
    437          ! Poles are single points on physics grid
    438          airefi(1)=sum(aire(1:iim))
    439          airefi(ngridmx)=sum(aire(ip1jmp1-(iim+1)+1:ip1jmp1-1))
    440 !         WRITE(lunout,*)
    441 !     .       'GCM: WARNING!!! vitesse verticale nulle dans la physique'
     408      IF ((iflag_phys==1).or.(iflag_phys>=100)) THEN
    442409! Physics:
    443410#ifdef CPP_PHYS
    444          CALL iniphysiq(ngridmx,llm,daysec,day_ini,dtphys/nsplit_phys,
    445      &                latfi,lonfi,airefi,zcufi,zcvfi,rad,g,r,cpp,
     411         CALL iniphysiq(iim,jjm,llm,daysec,day_ini,dtphys/nsplit_phys,
     412     &                rlatu,rlonv,aire,cu,cv,rad,g,r,cpp,
    446413     &                iflag_phys)
    447414#endif
    448          call_iniphys=.false.
    449       ENDIF ! of IF (call_iniphys.and.(iflag_phys.eq.1))
     415      ENDIF ! of IF ((iflag_phys==1).or.(iflag_phys>=100))
    450416
    451417c  numero de stockage pour les fichiers de redemarrage:
  • LMDZ5/trunk/libf/dyn3dmem/gcm.F

    r2222 r2225  
    11!
    2 ! $Id$
     2! $Id: $
    33!
    44c
     
    2323! Only INCA needs these informations (from the Earth's physics)
    2424      USE indice_sol_mod
     25      USE mod_phys_lmdz_omp_data, ONLY: klon_omp
    2526#endif
    2627
    2728#ifdef CPP_PHYS
    28       USE mod_grid_phy_lmdz
    29       USE mod_phys_lmdz_para, ONLY : klon_mpi_para_nb
    30       USE mod_phys_lmdz_omp_data, ONLY: klon_omp
    31       USE dimphy
    32       USE comgeomphy
     29!      USE mod_grid_phy_lmdz
     30!      USE mod_phys_lmdz_para, ONLY : klon_mpi_para_nb
     31!      USE dimphy
     32!      USE comgeomphy
    3333#endif
    3434      IMPLICIT NONE
     
    104104
    105105      LOGICAL lafin
    106 c      INTEGER ij,iq,l,i,j
    107       INTEGER i,j
    108 
    109106
    110107      real time_step, t_wrt, t_ops
    111 
    112 
    113       LOGICAL call_iniphys
    114       data call_iniphys/.true./
    115108
    116109c+jld variables test conservation energie
     
    135128
    136129c-----------------------------------------------------------------------
    137 c    variables pour l'initialisation de la physique :
    138 c    ------------------------------------------------
    139       INTEGER ngridmx
    140       PARAMETER( ngridmx = 2+(jjm-1)*iim - 1/jjm   )
    141       REAL zcufi(ngridmx),zcvfi(ngridmx)
    142       REAL latfi(ngridmx),lonfi(ngridmx)
    143       REAL airefi(ngridmx)
    144       SAVE latfi, lonfi, airefi
    145      
    146       INTEGER :: ierr
    147 
    148 
    149 c-----------------------------------------------------------------------
    150130c   Initialisations:
    151131c   ----------------
     
    164144c  ---------------------------------------
    165145c
    166 ! Ehouarn: dump possibility of using defrun
    167146      CALL conf_gcm( 99, .TRUE. )
    168147      if (mod(iphysiq, iperiod) /= 0) call abort_gcm("conf_gcm",
     
    181160#ifdef CPP_PHYS
    182161        CALL Init_Phys_lmdz(iim,jjp1,llm,mpi_size,distrib_phys)
    183 #endif
     162!#endif
     163!      CALL set_bands
     164!#ifdef CPP_PHYS
     165      CALL Init_interface_dyn_phys
     166#endif
     167      CALL barrier
     168
    184169      CALL set_bands
    185 #ifdef CPP_PHYS
    186       CALL Init_interface_dyn_phys
    187 #endif
    188       CALL barrier
    189 
    190170      if (mpi_rank==0) call WriteBands
    191171      call Set_Distrib(distrib_caldyn)
     
    195175c$OMP END PARALLEL
    196176
    197 #ifdef CPP_PHYS
    198 c$OMP PARALLEL
    199       call InitComgeomphy
    200 c$OMP END PARALLEL
    201 #endif
     177!#ifdef CPP_PHYS
     178!c$OMP PARALLEL
     179!      call InitComgeomphy ! now done in iniphysiq
     180!c$OMP END PARALLEL
     181!#endif
    202182
    203183c-----------------------------------------------------------------------
     
    420400c   Initialisation de la physique :
    421401c   -------------------------------
    422       IF (call_iniphys.and.(iflag_phys==1.or.iflag_phys>=100)) THEN
    423          latfi(1)=rlatu(1)
    424          lonfi(1)=0.
    425          zcufi(1) = cu(1)
    426          zcvfi(1) = cv(1)
    427          DO j=2,jjm
    428             DO i=1,iim
    429                latfi((j-2)*iim+1+i)= rlatu(j)
    430                lonfi((j-2)*iim+1+i)= rlonv(i)
    431                zcufi((j-2)*iim+1+i) = cu((j-1)*iip1+i)
    432                zcvfi((j-2)*iim+1+i) = cv((j-1)*iip1+i)
    433             ENDDO
    434          ENDDO
    435          latfi(ngridmx)= rlatu(jjp1)
    436          lonfi(ngridmx)= 0.
    437          zcufi(ngridmx) = cu(ip1jm+1)
    438          zcvfi(ngridmx) = cv(ip1jm-iim)
    439          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,aire,airefi)
    440          ! Poles are single points on physics grid
    441          airefi(1)=sum(aire(1:iim))
    442          airefi(ngridmx)=sum(aire(ip1jmp1-(iim+1)+1:ip1jmp1-1))
    443 !         WRITE(lunout,*)
    444 !     .       'GCM: WARNING!!! vitesse verticale nulle dans la physique'
     402      IF ((iflag_phys==1).or.(iflag_phys>=100)) THEN
    445403! Physics:
    446404#ifdef CPP_PHYS
    447          CALL iniphysiq(ngridmx,llm,daysec,day_ini,dtphys/nsplit_phys,
    448      &                latfi,lonfi,airefi,zcufi,zcvfi,rad,g,r,cpp,
    449      &                iflag_phys)
    450 #endif
    451          call_iniphys=.false.
    452       ENDIF ! of IF (call_iniphys.and.(iflag_phys==1.or.iflag_phys>=100))
     405         CALL iniphysiq(iim,jjm,llm,daysec,day_ini,dtphys/nsplit_phys,
     406     &                rlatu,rlonv,aire,cu,cv,rad,g,r,cpp,
     407     &                iflag_phys)
     408#endif
     409      ENDIF ! of IF ((iflag_phys==1).or.(iflag_phys>=100))
    453410
    454411
  • LMDZ5/trunk/libf/dyn3dpar/gcm.F

    r2222 r2225  
    2424! Only INCA needs these informations (from the Earth's physics)
    2525      USE indice_sol_mod
     26      USE mod_phys_lmdz_omp_data, ONLY: klon_omp
    2627#endif
    2728
    2829#ifdef CPP_PHYS
    29       USE mod_grid_phy_lmdz
    30       USE mod_phys_lmdz_para, ONLY : klon_mpi_para_nb
    31       USE mod_phys_lmdz_omp_data, ONLY: klon_omp
    32       USE dimphy
    33       USE comgeomphy
     30!      USE mod_grid_phy_lmdz
     31!      USE mod_phys_lmdz_para, ONLY : klon_mpi_para_nb
     32!      USE dimphy
     33!      USE comgeomphy
    3434#endif
    3535      IMPLICIT NONE
     
    105105
    106106      LOGICAL lafin
    107 c      INTEGER ij,iq,l,i,j
    108       INTEGER i,j
    109 
    110107
    111108      real time_step, t_wrt, t_ops
    112109
    113 
    114       LOGICAL call_iniphys
    115       data call_iniphys/.true./
    116110
    117111c+jld variables test conservation energie
     
    136130
    137131c-----------------------------------------------------------------------
    138 c    variables pour l'initialisation de la physique :
    139 c    ------------------------------------------------
    140       INTEGER ngridmx
    141       PARAMETER( ngridmx = 2+(jjm-1)*iim - 1/jjm   )
    142       REAL zcufi(ngridmx),zcvfi(ngridmx)
    143       REAL latfi(ngridmx),lonfi(ngridmx)
    144       REAL airefi(ngridmx)
    145       SAVE latfi, lonfi, airefi
    146      
    147       INTEGER :: ierr
    148 
    149 
    150 c-----------------------------------------------------------------------
    151132c   Initialisations:
    152133c   ----------------
     
    181162#ifdef CPP_PHYS
    182163        CALL Init_Phys_lmdz(iim,jjp1,llm,mpi_size,distrib_phys)
    183 #endif
     164!#endif
     165!      CALL set_bands
     166!#ifdef CPP_PHYS
     167      CALL Init_interface_dyn_phys
     168#endif
     169      CALL barrier
     170
    184171      CALL set_bands
    185 #ifdef CPP_PHYS
    186       CALL Init_interface_dyn_phys
    187 #endif
    188       CALL barrier
    189 
    190172      if (mpi_rank==0) call WriteBands
    191173      call SetDistrib(jj_Nb_Caldyn)
     
    195177c$OMP END PARALLEL
    196178
    197 #ifdef CPP_PHYS
    198 c$OMP PARALLEL
    199       call InitComgeomphy
    200 c$OMP END PARALLEL
    201 #endif
     179!#ifdef CPP_PHYS
     180!c$OMP PARALLEL
     181!      call InitComgeomphy ! now done in iniphysiq
     182!c$OMP END PARALLEL
     183!#endif
    202184
    203185!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
     
    419401c   Initialisation de la physique :
    420402c   -------------------------------
    421       IF (call_iniphys.and.(iflag_phys==1.or.iflag_phys>=100)) THEN
    422          latfi(1)=rlatu(1)
    423          lonfi(1)=0.
    424          zcufi(1) = cu(1)
    425          zcvfi(1) = cv(1)
    426          DO j=2,jjm
    427             DO i=1,iim
    428                latfi((j-2)*iim+1+i)= rlatu(j)
    429                lonfi((j-2)*iim+1+i)= rlonv(i)
    430                zcufi((j-2)*iim+1+i) = cu((j-1)*iip1+i)
    431                zcvfi((j-2)*iim+1+i) = cv((j-1)*iip1+i)
    432             ENDDO
    433          ENDDO
    434          latfi(ngridmx)= rlatu(jjp1)
    435          lonfi(ngridmx)= 0.
    436          zcufi(ngridmx) = cu(ip1jm+1)
    437          zcvfi(ngridmx) = cv(ip1jm-iim)
    438          CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,aire,airefi)
    439          ! Poles are single points on physics grid
    440          airefi(1)=sum(aire(1:iim))
    441          airefi(ngridmx)=sum(aire(ip1jmp1-(iim+1)+1:ip1jmp1-1))
    442 !         WRITE(lunout,*)
    443 !     .       'GCM: WARNING!!! vitesse verticale nulle dans la physique'
     403      IF ((iflag_phys==1).or.(iflag_phys>=100)) THEN
    444404! Physics:
    445405#ifdef CPP_PHYS
    446          CALL iniphysiq(ngridmx,llm,daysec,day_ini,dtphys/nsplit_phys,
    447      &                latfi,lonfi,airefi,zcufi,zcvfi,rad,g,r,cpp,
     406         CALL iniphysiq(iim,jjm,llm,daysec,day_ini,dtphys/nsplit_phys,
     407     &                rlatu,rlonv,aire,cu,cv,rad,g,r,cpp,
    448408     &                iflag_phys)
    449409#endif
    450          call_iniphys=.false.
    451       ENDIF ! of IF (call_iniphys.and.(iflag_phys==1.or.iflag_phys>=100))
     410      ENDIF ! of IF ((iflag_phys==1).or.(iflag_phys>=100))
    452411
    453412
  • LMDZ5/trunk/libf/phydev/iniphysiq.F90

    r1994 r2225  
    22! $Id: iniphysiq.F 1403 2010-07-01 09:02:53Z fairhead $
    33!
    4 SUBROUTINE iniphysiq(ngrid, nlayer, punjours, pdayref, ptimestep, plat, plon, &
    5     parea, pcu, pcv, prad, pg, pr, pcpp, iflag_phys)
    6   USE dimphy, ONLY: klev
    7   USE mod_grid_phy_lmdz, ONLY: klon_glo
    8   USE mod_phys_lmdz_para, ONLY: klon_omp, klon_omp_begin, klon_omp_end, &
    9     klon_mpi_begin
    10   USE comgeomphy, ONLY: airephy, cuphy, cvphy, rlond, rlatd
    11   USE comcstphy, ONLY: rradius, rg, rr, rcpp
     4SUBROUTINE iniphysiq(iim,jjm,nlayer,punjours, pdayref,ptimestep,         &
     5                     rlatu,rlonv,aire,cu,cv,                             &
     6                     prad,pg,pr,pcpp,iflag_phys)
     7  USE dimphy, ONLY: klev ! number of atmospheric levels
     8  USE mod_grid_phy_lmdz, ONLY: klon_glo ! number of atmospheric columns
     9                                        ! (on full grid)
     10  USE mod_phys_lmdz_para, ONLY: klon_omp, & ! number of columns (on local omp grid)
     11                                klon_omp_begin, & ! start index of local omp subgrid
     12                                klon_omp_end, & ! end index of local omp subgrid
     13                                klon_mpi_begin ! start indes of columns (on local mpi grid)
     14  USE comgeomphy, ONLY: initcomgeomphy, &
     15                        airephy, & ! physics grid area (m2)
     16                        cuphy, & ! cu coeff. (u_covariant = cu * u)
     17                        cvphy, & ! cv coeff. (v_covariant = cv * v)
     18                        rlond, & ! longitudes
     19                        rlatd ! latitudes
    1220  USE phyaqua_mod, ONLY: iniaqua
    1321  IMPLICIT NONE
    1422  !
    1523  !=======================================================================
    16   !
    1724  !   Initialisation of the physical constants and some positional and
    1825  !   geometrical arrays for the physics
    19   !
    20   !
    21   !    ngrid                 Size of the horizontal grid.
    22   !                          All internal loops are performed on that grid.
    23   !    nlayer                Number of vertical layers.
    24   !    pdayref               Day of reference for the simulation
    25   !
    2626  !=======================================================================
    2727 
     
    3434  REAL,INTENT(IN) :: pcpp ! specific heat Cp
    3535  REAL,INTENT(IN) :: punjours ! length (in s) of a standard day
    36   INTEGER,INTENT(IN) :: ngrid ! number of horizontal grid points in the physics
    37   INTEGER,INTENT(IN) :: nlayer ! number of atmospheric layers
    38   REAL,INTENT(IN) :: plat(ngrid) ! latitudes of the physics grid
    39   REAL,INTENT(IN) :: plon(ngrid) ! longitudes of the physics grid
    40   REAL,INTENT(IN) :: parea(klon_glo) ! area (m2)
    41   REAL,INTENT(IN) :: pcu(klon_glo) ! cu coeff. (u_covariant = cu * u)
    42   REAL,INTENT(IN) :: pcv(klon_glo) ! cv coeff. (v_covariant = cv * v)
    43   INTEGER,INTENT(IN) :: pdayref ! reference day of for the simulation
     36  INTEGER, INTENT (IN) :: nlayer ! number of atmospheric layers
     37  INTEGER, INTENT (IN) :: iim ! number of atmospheric coulumns along longitudes
     38  INTEGER, INTENT (IN) :: jjm  ! number of atompsheric columns along latitudes
     39  REAL, INTENT (IN) :: rlatu(jjm+1) ! latitudes of the physics grid
     40  REAL, INTENT (IN) :: rlonv(iim+1) ! longitudes of the physics grid
     41  REAL, INTENT (IN) :: aire(iim+1,jjm+1) ! area of the dynamics grid (m2)
     42  REAL, INTENT (IN) :: cu((iim+1)*(jjm+1)) ! cu coeff. (u_covariant = cu * u)
     43  REAL, INTENT (IN) :: cv((iim+1)*jjm) ! cv coeff. (v_covariant = cv * v)
    4444  REAL,INTENT(IN) :: ptimestep !physics time step (s)
    4545  INTEGER,INTENT(IN) :: iflag_phys ! type of physics to be called
    4646
    4747  INTEGER :: ibegin,iend,offset
     48  INTEGER :: i,j
    4849  CHARACTER (LEN=20) :: modname='iniphysiq'
    4950  CHARACTER (LEN=80) :: abort_message
    50  
     51  REAL :: total_area_phy, total_area_dyn
     52
     53
     54  ! global array, on full physics grid:
     55  REAL,ALLOCATABLE :: latfi(:)
     56  REAL,ALLOCATABLE :: lonfi(:)
     57  REAL,ALLOCATABLE :: cufi(:)
     58  REAL,ALLOCATABLE :: cvfi(:)
     59  REAL,ALLOCATABLE :: airefi(:)
     60
    5161  IF (nlayer.NE.klev) THEN
    5262    WRITE(lunout,*) 'STOP in ',trim(modname)
     
    5868  ENDIF
    5969
    60   IF (ngrid.NE.klon_glo) THEN
    61     WRITE(lunout,*) 'STOP in ',trim(modname)
    62     WRITE(lunout,*) 'Problem with dimensions :'
    63     WRITE(lunout,*) 'ngrid     = ',ngrid
    64     WRITE(lunout,*) 'klon   = ',klon_glo
    65     abort_message = ''
    66     CALL abort_gcm (modname,abort_message,1)
    67   ENDIF
     70  !call init_phys_lmdz(iim,jjm+1,llm,1,(/(jjm-1)*iim+2/))
     71 
     72  ! Generate global arrays on full physics grid
     73  ALLOCATE(latfi(klon_glo),lonfi(klon_glo),cufi(klon_glo),cvfi(klon_glo))
     74  ALLOCATE(airefi(klon_glo))
    6875
    69   !$OMP PARALLEL PRIVATE(ibegin,iend) &
    70   !$OMP          SHARED(parea,pcu,pcv,plon,plat)
    71      
     76    ! North pole
     77    latfi(1)=rlatu(1)
     78    lonfi(1)=0.
     79    cufi(1) = cu(1)
     80    cvfi(1) = cv(1)
     81    DO j=2,jjm
     82      DO i=1,iim
     83        latfi((j-2)*iim+1+i)= rlatu(j)
     84        lonfi((j-2)*iim+1+i)= rlonv(i)
     85        cufi((j-2)*iim+1+i) = cu((j-1)*iim+1+i)
     86        cvfi((j-2)*iim+1+i) = cv((j-1)*iim+1+i)
     87      ENDDO
     88    ENDDO
     89    ! South pole
     90    latfi(klon_glo)= rlatu(jjm+1)
     91    lonfi(klon_glo)= 0.
     92    cufi(klon_glo) = cu((iim+1)*jjm+1)
     93    cvfi(klon_glo) = cv((iim+1)*jjm-iim)
     94
     95    ! build airefi(), mesh area on physics grid
     96    CALL gr_dyn_fi(1,iim+1,jjm+1,klon_glo,aire,airefi)
     97    ! Poles are single points on physics grid
     98    airefi(1)=sum(aire(1:iim,1))
     99    airefi(klon_glo)=sum(aire(1:iim,jjm+1))
     100
     101    ! Sanity check: do total planet area match between physics and dynamics?
     102    total_area_dyn=sum(aire(1:iim,1:jjm+1))
     103    total_area_phy=sum(airefi(1:klon_glo))
     104    IF (total_area_dyn/=total_area_phy) THEN
     105      WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
     106      WRITE (lunout, *) '     in the dynamics total_area_dyn=', total_area_dyn
     107      WRITE (lunout, *) '  but in the physics total_area_phy=', total_area_phy
     108      IF (abs(total_area_dyn-total_area_phy)>0.00001*total_area_dyn) THEN
     109        ! stop here if the relative difference is more than 0.001%
     110        abort_message = 'planet total surface discrepancy'
     111        CALL abort_gcm(modname, abort_message, 1)
     112      ENDIF
     113    ENDIF
     114
     115!$OMP PARALLEL
     116  ! Now generate local lon/lat/cu/cv/area arrays
     117  CALL initcomgeomphy
     118
    72119  offset = klon_mpi_begin - 1
    73   airephy(1:klon_omp) = parea(offset+klon_omp_begin:offset+klon_omp_end)
    74   cuphy(1:klon_omp) = pcu(offset+klon_omp_begin:offset+klon_omp_end)
    75   cvphy(1:klon_omp) = pcv(offset+klon_omp_begin:offset+klon_omp_end)
    76   rlond(1:klon_omp) = plon(offset+klon_omp_begin:offset+klon_omp_end)
    77   rlatd(1:klon_omp) = plat(offset+klon_omp_begin:offset+klon_omp_end)
     120  airephy(1:klon_omp) = airefi(offset+klon_omp_begin:offset+klon_omp_end)
     121  cuphy(1:klon_omp) = cufi(offset+klon_omp_begin:offset+klon_omp_end)
     122  cvphy(1:klon_omp) = cvfi(offset+klon_omp_begin:offset+klon_omp_end)
     123  rlond(1:klon_omp) = lonfi(offset+klon_omp_begin:offset+klon_omp_end)
     124  rlatd(1:klon_omp) = latfi(offset+klon_omp_begin:offset+klon_omp_end)
    78125
    79126  ! copy some fundamental parameters to physics
     
    83130  rcpp=pcpp
    84131
    85   !$OMP END PARALLEL
     132!$OMP END PARALLEL
    86133
    87134  ! Additional initializations for aquaplanets
    88   !$OMP PARALLEL
     135!$OMP PARALLEL
    89136  IF (iflag_phys>=100) THEN
    90137    CALL iniaqua(klon_omp,rlatd,rlond,iflag_phys)
    91138  ENDIF
    92   !$OMP END PARALLEL
     139!$OMP END PARALLEL
    93140
    94141END SUBROUTINE iniphysiq
  • LMDZ5/trunk/libf/phylmd/iniphysiq.F90

    r1993 r2225  
    33
    44
    5 
    6 SUBROUTINE iniphysiq(ngrid, nlayer, punjours, pdayref, ptimestep, plat, plon, &
    7     parea, pcu, pcv, prad, pg, pr, pcpp, iflag_phys)
    8   USE dimphy, ONLY: klev
    9   USE mod_grid_phy_lmdz, ONLY: klon_glo
    10   USE mod_phys_lmdz_para, ONLY: klon_omp, klon_omp_begin, klon_omp_end, &
    11     klon_mpi_begin
    12   USE comgeomphy, ONLY: airephy, cuphy, cvphy, rlond, rlatd
     5SUBROUTINE iniphysiq(iim,jjm,nlayer,punjours, pdayref,ptimestep,         &
     6                     rlatu,rlonv,aire,cu,cv,                             &
     7                     prad,pg,pr,pcpp,iflag_phys)
     8  USE dimphy, ONLY: klev ! number of atmospheric levels
     9  USE mod_grid_phy_lmdz, ONLY: klon_glo ! number of atmospheric columns
     10                                        ! (on full grid)
     11  USE mod_phys_lmdz_para, ONLY: klon_omp, & ! number of columns (on local omp grid)
     12                                klon_omp_begin, & ! start index of local omp subgrid
     13                                klon_omp_end, & ! end index of local omp subgrid
     14                                klon_mpi_begin ! start indes of columns (on local mpi grid)
     15  USE comgeomphy, ONLY: initcomgeomphy, &
     16                        airephy, & ! physics grid area (m2)
     17                        cuphy, & ! cu coeff. (u_covariant = cu * u)
     18                        cvphy, & ! cv coeff. (v_covariant = cv * v)
     19                        rlond, & ! longitudes
     20                        rlatd ! latitudes
    1321  USE phyaqua_mod, ONLY: iniaqua
    1422  IMPLICIT NONE
    1523
    1624  ! =======================================================================
    17 
    1825  ! Initialisation of the physical constants and some positional and
    1926  ! geometrical arrays for the physics
    20 
    21 
    22   ! ngrid                 Size of the horizontal grid.
    23   ! All internal loops are performed on that grid.
    24   ! nlayer                Number of vertical layers.
    25   ! pdayref               Day of reference for the simulation
    26 
    2727  ! =======================================================================
    2828
    29   ! ym#include "dimensions.h"
    30   ! ym#include "dimphy.h"
    31   ! ym#include "comgeomphy.h"
    3229  include "YOMCST.h"
    3330  include "iniprint.h"
     
    3835  REAL, INTENT (IN) :: pcpp ! specific heat Cp
    3936  REAL, INTENT (IN) :: punjours ! length (in s) of a standard day
    40   INTEGER, INTENT (IN) :: ngrid ! number of horizontal grid points in the physics
    4137  INTEGER, INTENT (IN) :: nlayer ! number of atmospheric layers
    42   REAL, INTENT (IN) :: plat(ngrid) ! latitudes of the physics grid
    43   REAL, INTENT (IN) :: plon(ngrid) ! longitudes of the physics grid
    44   REAL, INTENT (IN) :: parea(klon_glo) ! area (m2)
    45   REAL, INTENT (IN) :: pcu(klon_glo) ! cu coeff. (u_covariant = cu * u)
    46   REAL, INTENT (IN) :: pcv(klon_glo) ! cv coeff. (v_covariant = cv * v)
     38  INTEGER, INTENT (IN) :: iim ! number of atmospheric columns along longitudes
     39  INTEGER, INTENT (IN) :: jjm ! number of atompsheric columns along latitudes
     40  REAL, INTENT (IN) :: rlatu(jjm+1) ! latitudes of the physics grid
     41  REAL, INTENT (IN) :: rlonv(iim+1) ! longitudes of the physics grid
     42  REAL, INTENT (IN) :: aire(iim+1,jjm+1) ! area of the dynamics grid (m2)
     43  REAL, INTENT (IN) :: cu((iim+1)*(jjm+1)) ! cu coeff. (u_covariant = cu * u)
     44  REAL, INTENT (IN) :: cv((iim+1)*jjm) ! cv coeff. (v_covariant = cv * v)
    4745  INTEGER, INTENT (IN) :: pdayref ! reference day of for the simulation
    4846  REAL, INTENT (IN) :: ptimestep !physics time step (s)
     
    5048
    5149  INTEGER :: ibegin, iend, offset
     50  INTEGER :: i,j
    5251  CHARACTER (LEN=20) :: modname = 'iniphysiq'
    5352  CHARACTER (LEN=80) :: abort_message
     53  REAL :: total_area_phy, total_area_dyn
     54
     55
     56  ! global array, on full physics grid:
     57  REAL,ALLOCATABLE :: latfi(:)
     58  REAL,ALLOCATABLE :: lonfi(:)
     59  REAL,ALLOCATABLE :: cufi(:)
     60  REAL,ALLOCATABLE :: cvfi(:)
     61  REAL,ALLOCATABLE :: airefi(:)
    5462
    5563  IF (nlayer/=klev) THEN
     
    6270  END IF
    6371
    64   IF (ngrid/=klon_glo) THEN
    65     WRITE (lunout, *) 'STOP in ', trim(modname)
    66     WRITE (lunout, *) 'Problem with dimensions :'
    67     WRITE (lunout, *) 'ngrid     = ', ngrid
    68     WRITE (lunout, *) 'klon   = ', klon_glo
    69     abort_message = ''
    70     CALL abort_gcm(modname, abort_message, 1)
    71   END IF
    72 
    73   !$OMP PARALLEL PRIVATE(ibegin,iend) &
    74   !$OMP          SHARED(parea,pcu,pcv,plon,plat)
     72  !call init_phys_lmdz(iim,jjm+1,llm,1,(/(jjm-1)*iim+2/))
     73 
     74  ! Generate global arrays on full physics grid
     75  ALLOCATE(latfi(klon_glo),lonfi(klon_glo),cufi(klon_glo),cvfi(klon_glo))
     76  ALLOCATE(airefi(klon_glo))
     77
     78  IF (klon_glo>1) THEN ! general case
     79    ! North pole
     80    latfi(1)=rlatu(1)
     81    lonfi(1)=0.
     82    cufi(1) = cu(1)
     83    cvfi(1) = cv(1)
     84    DO j=2,jjm
     85      DO i=1,iim
     86        latfi((j-2)*iim+1+i)= rlatu(j)
     87        lonfi((j-2)*iim+1+i)= rlonv(i)
     88        cufi((j-2)*iim+1+i) = cu((j-1)*iim+1+i)
     89        cvfi((j-2)*iim+1+i) = cv((j-1)*iim+1+i)
     90      ENDDO
     91    ENDDO
     92    ! South pole
     93    latfi(klon_glo)= rlatu(jjm+1)
     94    lonfi(klon_glo)= 0.
     95    cufi(klon_glo) = cu((iim+1)*jjm+1)
     96    cvfi(klon_glo) = cv((iim+1)*jjm-iim)
     97
     98    ! build airefi(), mesh area on physics grid
     99    CALL gr_dyn_fi(1,iim+1,jjm+1,klon_glo,aire,airefi)
     100    ! Poles are single points on physics grid
     101    airefi(1)=sum(aire(1:iim,1))
     102    airefi(klon_glo)=sum(aire(1:iim,jjm+1))
     103
     104    ! Sanity check: do total planet area match between physics and dynamics?
     105    total_area_dyn=sum(aire(1:iim,1:jjm+1))
     106    total_area_phy=sum(airefi(1:klon_glo))
     107    IF (total_area_dyn/=total_area_phy) THEN
     108      WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
     109      WRITE (lunout, *) '     in the dynamics total_area_dyn=', total_area_dyn
     110      WRITE (lunout, *) '  but in the physics total_area_phy=', total_area_phy
     111      IF (abs(total_area_dyn-total_area_phy)>0.00001*total_area_dyn) THEN
     112        ! stop here if the relative difference is more than 0.001%
     113        abort_message = 'planet total surface discrepancy'
     114        CALL abort_gcm(modname, abort_message, 1)
     115      ENDIF
     116    ENDIF
     117  ELSE ! klon_glo==1, running the 1D model
     118    ! just copy over input values
     119    latfi(1)=rlatu(1)
     120    lonfi(1)=rlonv(1)
     121    cufi(1)=cu(1)
     122    cvfi(1)=cv(1)
     123    airefi(1)=aire(1,1)
     124  ENDIF ! of IF (klon_glo>1)
     125
     126!$OMP PARALLEL
     127  ! Now generate local lon/lat/cu/cv/area arrays
     128  CALL initcomgeomphy
    75129
    76130  offset = klon_mpi_begin - 1
    77   airephy(1:klon_omp) = parea(offset+klon_omp_begin:offset+klon_omp_end)
    78   cuphy(1:klon_omp) = pcu(offset+klon_omp_begin:offset+klon_omp_end)
    79   cvphy(1:klon_omp) = pcv(offset+klon_omp_begin:offset+klon_omp_end)
    80   rlond(1:klon_omp) = plon(offset+klon_omp_begin:offset+klon_omp_end)
    81   rlatd(1:klon_omp) = plat(offset+klon_omp_begin:offset+klon_omp_end)
     131  airephy(1:klon_omp) = airefi(offset+klon_omp_begin:offset+klon_omp_end)
     132  cuphy(1:klon_omp) = cufi(offset+klon_omp_begin:offset+klon_omp_end)
     133  cvphy(1:klon_omp) = cvfi(offset+klon_omp_begin:offset+klon_omp_end)
     134  rlond(1:klon_omp) = lonfi(offset+klon_omp_begin:offset+klon_omp_end)
     135  rlatd(1:klon_omp) = latfi(offset+klon_omp_begin:offset+klon_omp_end)
    82136
    83137    ! suphel => initialize some physical constants (orbital parameters,
     
    86140  CALL suphel
    87141
    88   !$OMP END PARALLEL
    89 
    90     ! check that physical constants set in 'suphel' are coherent
    91     ! with values set in the dynamics:
     142!$OMP END PARALLEL
     143
     144  ! check that physical constants set in 'suphel' are coherent
     145  ! with values set in the dynamics:
    92146  IF (rday/=punjours) THEN
    93147    WRITE (lunout, *) 'iniphysiq: length of day discrepancy!!!'
     
    142196
    143197  ! Additional initializations for aquaplanets
    144   !$OMP PARALLEL
     198!$OMP PARALLEL
    145199  IF (iflag_phys>=100) THEN
    146200    CALL iniaqua(klon_omp, rlatd, rlond, iflag_phys)
    147201  END IF
    148   !$OMP END PARALLEL
    149 
    150   ! RETURN
    151   ! 9999  CONTINUE
    152   ! abort_message ='Cette version demande les fichier rnatur.dat
    153   ! & et surf.def'
    154   ! CALL abort_gcm (modname,abort_message,1)
     202!$OMP END PARALLEL
    155203
    156204END SUBROUTINE iniphysiq
  • LMDZ5/trunk/libf/phylmd/lmdz1d.F90

    r2221 r2225  
    1010      USE ioipsl, only: ju2ymds, ymds2ju, ioconf_calendar
    1111      use phys_state_var_mod
    12       use comgeomphy
    1312      use dimphy
    1413      use surface_data, only : type_ocean,ok_veget
     
    471470      call init_phys_lmdz(1,1,llm,1,(/1/))
    472471      call suphel
    473       call initcomgeomphy
    474472      call infotrac_init
    475473
     
    606604      rlon_rad(:)=rlon(:)*rpi/180.
    607605
    608       call iniphysiq(ngrid,llm,rday,day_ini,timestep,                        &
     606      call iniphysiq(iim,jjm,llm,rday,day_ini,timestep,                        &
    609607     &     rlat_rad,rlon_rad,airefi,zcufi,zcvfi,ra,rg,rd,rcpd,(/1/))
    610608      print*,'apres iniphysiq'
  • LMDZ5/trunk/libf/phymar/iniphysiq.F90

    r2221 r2225  
    22! $Id: iniphysiq.F 1403 2010-07-01 09:02:53Z fairhead $
    33!
    4       SUBROUTINE iniphysiq(ngrid,nlayer,
    5      $           punjours,
    6      $           pdayref,ptimestep,
    7      $           plat,plon,parea,pcu,pcv,
    8      $           prad,pg,pr,pcpp,iflag_phys)
    9       USE dimphy, only : klev
    10       USE mod_grid_phy_lmdz, only : klon_glo
    11       USE mod_phys_lmdz_para, only : klon_omp,klon_omp_begin,
    12      &                               klon_omp_end,klon_mpi_begin
    13       USE comgeomphy, only : airephy,cuphy,cvphy,rlond,rlatd
    14       USE comcstphy, only : rradius,rg,rr,rcpp
    15 
    16       IMPLICIT NONE
    17 c
    18 c=======================================================================
    19 c
    20 c   Initialisation of the physical constants and some positional and
    21 c   geometrical arrays for the physics
    22 c
    23 c
    24 c    ngrid                 Size of the horizontal grid.
    25 c                          All internal loops are performed on that grid.
    26 c    nlayer                Number of vertical layers.
    27 c    pdayref               Day of reference for the simulation
    28 c
    29 c=======================================================================
     4SUBROUTINE iniphysiq(iim,jjm,nlayer,punjours, pdayref,ptimestep,         &
     5                     rlatu,rlonv,aire,cu,cv,                             &
     6                     prad,pg,pr,pcpp,iflag_phys)
     7  USE dimphy, ONLY: klev ! number of atmospheric levels
     8  USE mod_grid_phy_lmdz, ONLY: klon_glo ! number of atmospheric columns
     9                                        ! (on full grid)
     10  USE mod_phys_lmdz_para, ONLY: klon_omp, & ! number of columns (on local omp grid)
     11                                klon_omp_begin, & ! start index of local omp subgrid
     12                                klon_omp_end, & ! end index of local omp subgrid
     13                                klon_mpi_begin ! start indes of columns (on local mpi grid)
     14  USE comgeomphy, ONLY: initcomgeomphy, &
     15                        airephy, & ! physics grid area (m2)
     16                        cuphy, & ! cu coeff. (u_covariant = cu * u)
     17                        cvphy, & ! cv coeff. (v_covariant = cv * v)
     18                        rlond, & ! longitudes
     19                        rlatd ! latitudes
     20  USE phyaqua_mod, ONLY: iniaqua
     21  IMPLICIT NONE
     22  !
     23  !=======================================================================
     24  !   Initialisation of the physical constants and some positional and
     25  !   geometrical arrays for the physics
     26  !=======================================================================
    3027 
    3128 
    32 cym#include "dimensions.h"
    33 cym#include "dimphy.h"
    34 cym#include "comgeomphy.h"
    35 #include "iniprint.h"
     29  include "iniprint.h"
    3630
    37       REAL,INTENT(IN) :: prad ! radius of the planet (m)
    38       REAL,INTENT(IN) :: pg ! gravitational acceleration (m/s2)
    39       REAL,INTENT(IN) :: pr ! ! reduced gas constant R/mu
    40       REAL,INTENT(IN) :: pcpp ! specific heat Cp
    41       REAL,INTENT(IN) :: punjours ! length (in s) of a standard day
    42       INTEGER,INTENT(IN) :: ngrid ! number of horizontal grid points in the physics
    43       INTEGER,INTENT(IN) :: nlayer ! number of atmospheric layers
    44       REAL,INTENT(IN) :: plat(ngrid) ! latitudes of the physics grid
    45       REAL,INTENT(IN) :: plon(ngrid) ! longitudes of the physics grid
    46       REAL,INTENT(IN) :: parea(klon_glo) ! area (m2)
    47       REAL,INTENT(IN) :: pcu(klon_glo) ! cu coeff. (u_covariant = cu * u)
    48       REAL,INTENT(IN) :: pcv(klon_glo) ! cv coeff. (v_covariant = cv * v)
    49       INTEGER,INTENT(IN) :: pdayref ! reference day of for the simulation
    50       REAL,INTENT(IN) :: ptimestep !physics time step (s)
    51       INTEGER,INTENT(IN) :: iflag_phys ! type of physics to be called
     31  REAL,INTENT(IN) :: prad ! radius of the planet (m)
     32  REAL,INTENT(IN) :: pg ! gravitational acceleration (m/s2)
     33  REAL,INTENT(IN) :: pr ! ! reduced gas constant R/mu
     34  REAL,INTENT(IN) :: pcpp ! specific heat Cp
     35  REAL,INTENT(IN) :: punjours ! length (in s) of a standard day
     36  INTEGER, INTENT (IN) :: nlayer ! number of atmospheric layers
     37  INTEGER, INTENT (IN) :: iim ! number of atmospheric coulumns along longitudes
     38  INTEGER, INTENT (IN) :: jjm  ! number of atompsheric columns along latitudes
     39  REAL, INTENT (IN) :: rlatu(jjm+1) ! latitudes of the physics grid
     40  REAL, INTENT (IN) :: rlonv(iim+1) ! longitudes of the physics grid
     41  REAL, INTENT (IN) :: aire(iim+1,jjm+1) ! area of the dynamics grid (m2)
     42  REAL, INTENT (IN) :: cu((iim+1)*(jjm+1)) ! cu coeff. (u_covariant = cu * u)
     43  REAL, INTENT (IN) :: cv((iim+1)*jjm) ! cv coeff. (v_covariant = cv * v)
     44  REAL,INTENT(IN) :: ptimestep !physics time step (s)
     45  INTEGER,INTENT(IN) :: iflag_phys ! type of physics to be called
    5246
    53       INTEGER :: ibegin,iend,offset
    54       CHARACTER (LEN=20) :: modname='iniphysiq'
    55       CHARACTER (LEN=80) :: abort_message
     47  INTEGER :: ibegin,iend,offset
     48  INTEGER :: i,j
     49  CHARACTER (LEN=20) :: modname='iniphysiq'
     50  CHARACTER (LEN=80) :: abort_message
     51  REAL :: total_area_phy, total_area_dyn
     52
     53
     54  ! global array, on full physics grid:
     55  REAL,ALLOCATABLE :: latfi(:)
     56  REAL,ALLOCATABLE :: lonfi(:)
     57  REAL,ALLOCATABLE :: cufi(:)
     58  REAL,ALLOCATABLE :: cvfi(:)
     59  REAL,ALLOCATABLE :: airefi(:)
    5660 
    57       IF (nlayer.NE.klev) THEN
    58          write(lunout,*) 'STOP in ',trim(modname)
    59          write(lunout,*) 'Problem with dimensions :'
    60          write(lunout,*) 'nlayer     = ',nlayer
    61          write(lunout,*) 'klev   = ',klev
    62          abort_message = ''
    63          CALL abort_gcm (modname,abort_message,1)
     61  IF (nlayer.NE.klev) THEN
     62    WRITE(lunout,*) 'STOP in ',trim(modname)
     63    WRITE(lunout,*) 'Problem with dimensions :'
     64    WRITE(lunout,*) 'nlayer     = ',nlayer
     65    WRITE(lunout,*) 'klev   = ',klev
     66    abort_message = ''
     67    CALL abort_gcm (modname,abort_message,1)
     68  ENDIF
     69
     70
     71  ! Generate global arrays on full physics grid
     72  ALLOCATE(latfi(klon_glo),lonfi(klon_glo),cufi(klon_glo),cvfi(klon_glo))
     73  ALLOCATE(airefi(klon_glo))
     74
     75  IF (klon_glo>1) THEN ! general case
     76    ! North pole
     77    latfi(1)=rlatu(1)
     78    lonfi(1)=0.
     79    cufi(1) = cu(1)
     80    cvfi(1) = cv(1)
     81    DO j=2,jjm
     82      DO i=1,iim
     83        latfi((j-2)*iim+1+i)= rlatu(j)
     84        lonfi((j-2)*iim+1+i)= rlonv(i)
     85        cufi((j-2)*iim+1+i) = cu((j-1)*iim+1+i)
     86        cvfi((j-2)*iim+1+i) = cv((j-1)*iim+1+i)
     87      ENDDO
     88    ENDDO
     89    ! South pole
     90    latfi(klon_glo)= rlatu(jjm+1)
     91    lonfi(klon_glo)= 0.
     92    cufi(klon_glo) = cu((iim+1)*jjm+1)
     93    cvfi(klon_glo) = cv((iim+1)*jjm-iim)
     94
     95    ! build airefi(), mesh area on physics grid
     96    CALL gr_dyn_fi(1,iim+1,jjm+1,klon_glo,aire,airefi)
     97    ! Poles are single points on physics grid
     98    airefi(1)=sum(aire(1:iim,1))
     99    airefi(klon_glo)=sum(aire(1:iim,jjm+1))
     100
     101    ! Sanity check: do total planet area match between physics and dynamics?
     102    total_area_dyn=sum(aire(1:iim,1:jjm+1))
     103    total_area_phy=sum(airefi(1:klon_glo))
     104    IF (total_area_dyn/=total_area_phy) THEN
     105      WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
     106      WRITE (lunout, *) '     in the dynamics total_area_dyn=', total_area_dyn
     107      WRITE (lunout, *) '  but in the physics total_area_phy=', total_area_phy
     108      IF (abs(total_area_dyn-total_area_phy)>0.00001*total_area_dyn) THEN
     109        ! stop here if the relative difference is more than 0.001%
     110        abort_message = 'planet total surface discrepancy'
     111        CALL abort_gcm(modname, abort_message, 1)
    64112      ENDIF
     113    ENDIF
     114  ELSE ! klon_glo==1, running the 1D model
     115    ! just copy over input values
     116    latfi(1)=rlatu(1)
     117    lonfi(1)=rlonv(1)
     118    cufi(1)=cu(1)
     119    cvfi(1)=cv(1)
     120    airefi(1)=aire(1,1)
     121  ENDIF ! of IF (klon_glo>1)
    65122
    66       IF (ngrid.NE.klon_glo) THEN
    67          write(lunout,*) 'STOP in ',trim(modname)
    68          write(lunout,*) 'Problem with dimensions :'
    69          write(lunout,*) 'ngrid     = ',ngrid
    70          write(lunout,*) 'klon   = ',klon_glo
    71          abort_message = ''
    72          CALL abort_gcm (modname,abort_message,1)
    73       ENDIF
    74 
    75 !$OMP PARALLEL PRIVATE(ibegin,iend)
    76 !$OMP+         SHARED(parea,pcu,pcv,plon,plat)
     123!$OMP PARALLEL
     124  ! Now generate local lon/lat/cu/cv/area arrays
     125  CALL initcomgeomphy
    77126     
    78       offset=klon_mpi_begin-1
    79       airephy(1:klon_omp)=parea(offset+klon_omp_begin:
    80      &                          offset+klon_omp_end)
    81       cuphy(1:klon_omp)=pcu(offset+klon_omp_begin:offset+klon_omp_end)
    82       cvphy(1:klon_omp)=pcv(offset+klon_omp_begin:offset+klon_omp_end)
    83       rlond(1:klon_omp)=plon(offset+klon_omp_begin:offset+klon_omp_end)
    84       rlatd(1:klon_omp)=plat(offset+klon_omp_begin:offset+klon_omp_end)
     127  offset = klon_mpi_begin - 1
     128  airephy(1:klon_omp) = airefi(offset+klon_omp_begin:offset+klon_omp_end)
     129  cuphy(1:klon_omp) = cufi(offset+klon_omp_begin:offset+klon_omp_end)
     130  cvphy(1:klon_omp) = cvfi(offset+klon_omp_begin:offset+klon_omp_end)
     131  rlond(1:klon_omp) = lonfi(offset+klon_omp_begin:offset+klon_omp_end)
     132  rlatd(1:klon_omp) = latfi(offset+klon_omp_begin:offset+klon_omp_end)
    85133
    86134! copy some fundamental parameters to physics
    87       rradius=prad
    88       rg=pg
    89       rr=pr
    90       rcpp=pcpp
     135  rradius=prad
     136  rg=pg
     137  rr=pr
     138  rcpp=pcpp
    91139
    92140!$OMP END PARALLEL
     
    102150!$OMP END PARALLEL
    103151
    104 !      RETURN
    105 !9999  CONTINUE
    106 !      abort_message ='Cette version demande les fichier rnatur.dat
    107 !     & et surf.def'
    108 !      CALL abort_gcm (modname,abort_message,1)
    109 
    110152      END
Note: See TracChangeset for help on using the changeset viewer.