Changeset 2225 for LMDZ5/trunk
- Timestamp:
- Mar 11, 2015, 3:55:23 PM (10 years ago)
- Location:
- LMDZ5/trunk/libf
- Files:
-
- 6 edited
- 1 moved
Legend:
- Unmodified
- Added
- Removed
-
LMDZ5/trunk/libf/dyn3d/gcm.F
r2222 r2225 26 26 ! Only INCA needs these informations (from the Earth's physics) 27 27 USE indice_sol_mod 28 USE mod_phys_lmdz_para, ONLY : klon_mpi_para_nb 28 29 #endif 29 30 … … 33 34 ! dynamique -> physique pour l'initialisation 34 35 #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 38 38 #endif 39 39 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! … … 119 119 LOGICAL first 120 120 121 LOGICAL call_iniphys122 data call_iniphys/.true./121 ! LOGICAL call_iniphys 122 ! data call_iniphys/.true./ 123 123 124 124 c+jld variables test conservation energie … … 141 141 REAL :: heure 142 142 143 144 c-----------------------------------------------------------------------145 c variables pour l'initialisation de la physique :146 c ------------------------------------------------147 INTEGER ngridmx148 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, airefi153 154 143 c----------------------------------------------------------------------- 155 144 c Initialisations: … … 188 177 #ifdef CPP_PHYS 189 178 CALL Init_Phys_lmdz(iim,jjp1,llm,1,(/(jjm-1)*iim+2/)) 190 call InitComgeomphy 179 ! call InitComgeomphy ! now done in iniphysiq 191 180 #endif 192 181 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! … … 417 406 c ------------------------------- 418 407 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 442 409 ! Physics: 443 410 #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, 446 413 & iflag_phys) 447 414 #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)) 450 416 451 417 c numero de stockage pour les fichiers de redemarrage: -
LMDZ5/trunk/libf/dyn3dmem/gcm.F
r2222 r2225 1 1 ! 2 ! $Id $2 ! $Id: $ 3 3 ! 4 4 c … … 23 23 ! Only INCA needs these informations (from the Earth's physics) 24 24 USE indice_sol_mod 25 USE mod_phys_lmdz_omp_data, ONLY: klon_omp 25 26 #endif 26 27 27 28 #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 33 33 #endif 34 34 IMPLICIT NONE … … 104 104 105 105 LOGICAL lafin 106 c INTEGER ij,iq,l,i,j107 INTEGER i,j108 109 106 110 107 real time_step, t_wrt, t_ops 111 112 113 LOGICAL call_iniphys114 data call_iniphys/.true./115 108 116 109 c+jld variables test conservation energie … … 135 128 136 129 c----------------------------------------------------------------------- 137 c variables pour l'initialisation de la physique :138 c ------------------------------------------------139 INTEGER ngridmx140 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, airefi145 146 INTEGER :: ierr147 148 149 c-----------------------------------------------------------------------150 130 c Initialisations: 151 131 c ---------------- … … 164 144 c --------------------------------------- 165 145 c 166 ! Ehouarn: dump possibility of using defrun167 146 CALL conf_gcm( 99, .TRUE. ) 168 147 if (mod(iphysiq, iperiod) /= 0) call abort_gcm("conf_gcm", … … 181 160 #ifdef CPP_PHYS 182 161 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 184 169 CALL set_bands 185 #ifdef CPP_PHYS186 CALL Init_interface_dyn_phys187 #endif188 CALL barrier189 190 170 if (mpi_rank==0) call WriteBands 191 171 call Set_Distrib(distrib_caldyn) … … 195 175 c$OMP END PARALLEL 196 176 197 #ifdef CPP_PHYS198 c$OMP PARALLEL199 call InitComgeomphy 200 c$OMP END PARALLEL201 #endif177 !#ifdef CPP_PHYS 178 !c$OMP PARALLEL 179 ! call InitComgeomphy ! now done in iniphysiq 180 !c$OMP END PARALLEL 181 !#endif 202 182 203 183 c----------------------------------------------------------------------- … … 420 400 c Initialisation de la physique : 421 401 c ------------------------------- 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 445 403 ! Physics: 446 404 #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)) 453 410 454 411 -
LMDZ5/trunk/libf/dyn3dpar/gcm.F
r2222 r2225 24 24 ! Only INCA needs these informations (from the Earth's physics) 25 25 USE indice_sol_mod 26 USE mod_phys_lmdz_omp_data, ONLY: klon_omp 26 27 #endif 27 28 28 29 #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 34 34 #endif 35 35 IMPLICIT NONE … … 105 105 106 106 LOGICAL lafin 107 c INTEGER ij,iq,l,i,j108 INTEGER i,j109 110 107 111 108 real time_step, t_wrt, t_ops 112 109 113 114 LOGICAL call_iniphys115 data call_iniphys/.true./116 110 117 111 c+jld variables test conservation energie … … 136 130 137 131 c----------------------------------------------------------------------- 138 c variables pour l'initialisation de la physique :139 c ------------------------------------------------140 INTEGER ngridmx141 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, airefi146 147 INTEGER :: ierr148 149 150 c-----------------------------------------------------------------------151 132 c Initialisations: 152 133 c ---------------- … … 181 162 #ifdef CPP_PHYS 182 163 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 184 171 CALL set_bands 185 #ifdef CPP_PHYS186 CALL Init_interface_dyn_phys187 #endif188 CALL barrier189 190 172 if (mpi_rank==0) call WriteBands 191 173 call SetDistrib(jj_Nb_Caldyn) … … 195 177 c$OMP END PARALLEL 196 178 197 #ifdef CPP_PHYS198 c$OMP PARALLEL199 call InitComgeomphy 200 c$OMP END PARALLEL201 #endif179 !#ifdef CPP_PHYS 180 !c$OMP PARALLEL 181 ! call InitComgeomphy ! now done in iniphysiq 182 !c$OMP END PARALLEL 183 !#endif 202 184 203 185 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! … … 419 401 c Initialisation de la physique : 420 402 c ------------------------------- 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 444 404 ! Physics: 445 405 #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, 448 408 & iflag_phys) 449 409 #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)) 452 411 453 412 -
LMDZ5/trunk/libf/phydev/iniphysiq.F90
r1994 r2225 2 2 ! $Id: iniphysiq.F 1403 2010-07-01 09:02:53Z fairhead $ 3 3 ! 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 4 SUBROUTINE 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 12 20 USE phyaqua_mod, ONLY: iniaqua 13 21 IMPLICIT NONE 14 22 ! 15 23 !======================================================================= 16 !17 24 ! Initialisation of the physical constants and some positional and 18 25 ! 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 simulation25 !26 26 !======================================================================= 27 27 … … 34 34 REAL,INTENT(IN) :: pcpp ! specific heat Cp 35 35 REAL,INTENT(IN) :: punjours ! length (in s) of a standard day 36 INTEGER, INTENT(IN) :: ngrid ! number of horizontal grid points in the physics37 INTEGER, INTENT(IN) :: nlayer ! number of atmospheric layers38 REAL,INTENT(IN) :: plat(ngrid) ! latitudes of the physics grid39 REAL, INTENT(IN) :: plon(ngrid) ! longitudes of the physics grid40 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 simulation36 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 44 REAL,INTENT(IN) :: ptimestep !physics time step (s) 45 45 INTEGER,INTENT(IN) :: iflag_phys ! type of physics to be called 46 46 47 47 INTEGER :: ibegin,iend,offset 48 INTEGER :: i,j 48 49 CHARACTER (LEN=20) :: modname='iniphysiq' 49 50 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 51 61 IF (nlayer.NE.klev) THEN 52 62 WRITE(lunout,*) 'STOP in ',trim(modname) … … 58 68 ENDIF 59 69 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)) 68 75 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 72 119 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) 78 125 79 126 ! copy some fundamental parameters to physics … … 83 130 rcpp=pcpp 84 131 85 132 !$OMP END PARALLEL 86 133 87 134 ! Additional initializations for aquaplanets 88 135 !$OMP PARALLEL 89 136 IF (iflag_phys>=100) THEN 90 137 CALL iniaqua(klon_omp,rlatd,rlond,iflag_phys) 91 138 ENDIF 92 139 !$OMP END PARALLEL 93 140 94 141 END SUBROUTINE iniphysiq -
LMDZ5/trunk/libf/phylmd/iniphysiq.F90
r1993 r2225 3 3 4 4 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 5 SUBROUTINE 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 13 21 USE phyaqua_mod, ONLY: iniaqua 14 22 IMPLICIT NONE 15 23 16 24 ! ======================================================================= 17 18 25 ! Initialisation of the physical constants and some positional and 19 26 ! 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 simulation26 27 27 ! ======================================================================= 28 28 29 ! ym#include "dimensions.h"30 ! ym#include "dimphy.h"31 ! ym#include "comgeomphy.h"32 29 include "YOMCST.h" 33 30 include "iniprint.h" … … 38 35 REAL, INTENT (IN) :: pcpp ! specific heat Cp 39 36 REAL, INTENT (IN) :: punjours ! length (in s) of a standard day 40 INTEGER, INTENT (IN) :: ngrid ! number of horizontal grid points in the physics41 37 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) 47 45 INTEGER, INTENT (IN) :: pdayref ! reference day of for the simulation 48 46 REAL, INTENT (IN) :: ptimestep !physics time step (s) … … 50 48 51 49 INTEGER :: ibegin, iend, offset 50 INTEGER :: i,j 52 51 CHARACTER (LEN=20) :: modname = 'iniphysiq' 53 52 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(:) 54 62 55 63 IF (nlayer/=klev) THEN … … 62 70 END IF 63 71 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 75 129 76 130 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) 82 136 83 137 ! suphel => initialize some physical constants (orbital parameters, … … 86 140 CALL suphel 87 141 88 89 90 91 142 !$OMP END PARALLEL 143 144 ! check that physical constants set in 'suphel' are coherent 145 ! with values set in the dynamics: 92 146 IF (rday/=punjours) THEN 93 147 WRITE (lunout, *) 'iniphysiq: length of day discrepancy!!!' … … 142 196 143 197 ! Additional initializations for aquaplanets 144 198 !$OMP PARALLEL 145 199 IF (iflag_phys>=100) THEN 146 200 CALL iniaqua(klon_omp, rlatd, rlond, iflag_phys) 147 201 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 155 203 156 204 END SUBROUTINE iniphysiq -
LMDZ5/trunk/libf/phylmd/lmdz1d.F90
r2221 r2225 10 10 USE ioipsl, only: ju2ymds, ymds2ju, ioconf_calendar 11 11 use phys_state_var_mod 12 use comgeomphy13 12 use dimphy 14 13 use surface_data, only : type_ocean,ok_veget … … 471 470 call init_phys_lmdz(1,1,llm,1,(/1/)) 472 471 call suphel 473 call initcomgeomphy474 472 call infotrac_init 475 473 … … 606 604 rlon_rad(:)=rlon(:)*rpi/180. 607 605 608 call iniphysiq( ngrid,llm,rday,day_ini,timestep, &606 call iniphysiq(iim,jjm,llm,rday,day_ini,timestep, & 609 607 & rlat_rad,rlon_rad,airefi,zcufi,zcvfi,ra,rg,rd,rcpd,(/1/)) 610 608 print*,'apres iniphysiq' -
LMDZ5/trunk/libf/phymar/iniphysiq.F90
r2221 r2225 2 2 ! $Id: iniphysiq.F 1403 2010-07-01 09:02:53Z fairhead $ 3 3 ! 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======================================================================= 4 SUBROUTINE 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 !======================================================================= 30 27 31 28 32 cym#include "dimensions.h" 33 cym#include "dimphy.h" 34 cym#include "comgeomphy.h" 35 #include "iniprint.h" 29 include "iniprint.h" 36 30 37 38 39 40 41 42 INTEGER,INTENT(IN) :: ngrid ! number of horizontal grid points in the physics43 INTEGER,INTENT(IN) :: nlayer ! number of atmospheric layers44 REAL,INTENT(IN) :: plat(ngrid) ! latitudes of the physics grid45 REAL,INTENT(IN) :: plon(ngrid) ! longitudes of the physics grid46 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 simulation50 51 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 52 46 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(:) 56 60 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) 64 112 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) 65 122 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 77 126 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) 85 133 86 134 ! copy some fundamental parameters to physics 87 88 89 90 135 rradius=prad 136 rg=pg 137 rr=pr 138 rcpp=pcpp 91 139 92 140 !$OMP END PARALLEL … … 102 150 !$OMP END PARALLEL 103 151 104 ! RETURN105 !9999 CONTINUE106 ! abort_message ='Cette version demande les fichier rnatur.dat107 ! & et surf.def'108 ! CALL abort_gcm (modname,abort_message,1)109 110 152 END
Note: See TracChangeset
for help on using the changeset viewer.