Changeset 2909 for trunk/LMDZ.MARS/libf
- Timestamp:
- Mar 9, 2023, 9:38:24 AM (21 months ago)
- Location:
- trunk/LMDZ.MARS/libf
- Files:
-
- 1 added
- 8 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/libf/dynphy_lonlat/phymars/newstart.F
r2832 r2909 29 29 use surfdat_h, only: phisfi, z0, zmea, zstd, zsig, zgam, zthe, 30 30 & albedodat, z0_default, qsurf, tsurf, 31 & emis, hmons, summit, base, watercap 32 use comsoil_h, only: inertiedat, layer, mlayer, nsoilmx, tsoil 31 & emis, hmons, summit, base, watercap, 32 & ini_surfdat_h_slope_var,end_surfdat_h_slope_var 33 use comsoil_h, only: inertiedat, layer, mlayer, nsoilmx, tsoil, 34 & ini_comsoil_h_slope_var, end_comsoil_h_slope_var 33 35 use control_mod, only: day_step, iphysiq, anneeref, planet_type 34 36 use geometry_mod, only: longitude,latitude,cell_area … … 36 38 use phyredem, only: physdem0, physdem1 37 39 use iostart, only: open_startphy 38 use dimradmars_mod, only: albedo 40 use dimradmars_mod, only: albedo, 41 & ini_dimradmars_mod_slope_var,end_dimradmars_mod_slope_var 39 42 use dust_param_mod, only: tauscaling 40 43 use turb_mod, only: q2, wstar … … 50 53 USE exner_hyb_m, ONLY: exner_hyb 51 54 USE inichim_newstart_mod, ONLY: inichim_newstart 52 55 use comslope_mod, ONLY: nslope,def_slope,def_slope_mean, 56 & subslope_dist,end_comslope_h,ini_comslope_h 53 57 implicit none 54 58 … … 180 184 real,allocatable :: coefvmr(:) ! Correction coefficient when changing composition 181 185 real :: maxq 182 integer :: iloc(1), iqmax 186 integer :: iloc(1), iqmax, islope 183 187 ! sub-grid cloud fraction 184 188 real :: totcloudfrac(ngridmx) 189 !Variable to change the number of subslope 190 REAL, ALLOCATABLE :: default_def_slope(:) 191 REAL,ALLOCATABLE :: tsurf_old_slope(:,:) ! Surface temperature (K) 192 REAL,ALLOCATABLE :: emis_old_slope(:,:) ! Thermal IR surface emissivity 193 REAL,ALLOCATABLE :: qsurf_old_slope(:,:,:) ! tracer on surface (e.g. kg.m-2) 194 REAL,ALLOCATABLE :: watercap_old_slope(:,:) ! Surface water ice (kg.m-2) 195 REAL,ALLOCATABLE :: tsoil_old_slope(:,:,:) 196 REAL,ALLOCATABLE :: albedo_old_slope(:,:,:) ! Surface albedo in each solar band 197 integer :: iflat 198 integer :: nslope_old, nslope_new 199 200 185 201 186 202 c sortie visu pour les champs dynamiques … … 443 459 & day_ini,time,tsurf,tsoil,albedo,emis, 444 460 & q2,qsurf,tauscaling,totcloudfrac, 445 & wstar,watercap )461 & wstar,watercap,def_slope,def_slope_mean,subslope_dist) 446 462 447 463 ! copy albedo and soil thermal inertia … … 536 552 write(*,*) 'mons_ice : put underground ice layer according 537 553 $ to MONS derived data' 554 write(*,*) 'nslope : Change the number of subgrid scale 555 $ slope' 538 556 539 557 write(*,*) … … 633 651 write(*,*) ' new value of the subsurface temperature:',tsud 634 652 c nouvelle temperature sous la calotte permanente 635 do l=2,nsoilmx 636 tsoil(ngridmx,l) = tsud 637 end do 653 do islope=1,nslope 654 do l=2,nsoilmx 655 tsoil(ngridmx,l,islope) = tsud 656 end do 657 enddo 638 658 639 659 … … 689 709 patm = patm + ps(i,j)*aire(i,j) 690 710 airetot= airetot + aire(i,j) 691 pcap = pcap + aire(i,j)*qsurf(ig,igcm_co2)*g 711 DO islope=1,nslope 712 pcap = pcap + aire(i,j)*qsurf(ig,igcm_co2,islope)*g 713 & *subslope_dist(ig,islope)/cos(def_slope_mean(islope)*pi/180.) 714 ENDDO 692 715 ENDDO 693 716 ENDDO … … 776 799 DO iq =1, nqtot 777 800 DO ig=1,ngridmx 778 qsurf(ig,iq)=0. 801 DO islope=1,nslope 802 qsurf(ig,iq,islope)=0. 803 ENDDO 779 804 ENDDO 780 805 ENDDO … … 798 823 799 824 q(1:iip1,1:jjp1,1:llm,iq)=q(1:iip1,1:jjp1,1:llm,iq)*val 800 qsurf(1:ngridmx,iq )=qsurf(1:ngridmx,iq)*val825 qsurf(1:ngridmx,iq,:)=qsurf(1:ngridmx,iq,:)*val 801 826 802 827 c q=x : initialise tracer manually … … 828 853 read(*,*) val 829 854 DO ig=1,ngridmx 830 qsurf(ig,iq)=val 855 DO islope=1,nslope 856 qsurf(ig,iq,islope)=val 857 ENDDO 831 858 ENDDO 832 859 … … 864 891 if (ierr.eq.0) then 865 892 ! initialize tracer values 866 qsurf(:,iq)=profile(1) 893 do islope=1,nslope 894 qsurf(:,iq,islope)=profile(1) 895 enddo 867 896 do l=1,llm 868 897 q(:,:,l,iq)=profile(l+1) … … 1013 1042 1014 1043 do ig=1,ngridmx 1015 qsurf(ig,igcm_h2o_ice)=max(0.,qsurf(ig,igcm_h2o_ice)) 1044 do islope=1,nslope 1045 qsurf(ig,igcm_h2o_ice,islope)= 1046 & max(0.,qsurf(ig,igcm_h2o_ice,islope)) 1047 enddo 1016 1048 end do 1017 1049 … … 1022 1054 & =q(1:iip1,1:jjp1,1:llm,igcm_h2o_ice)* DoverH 1023 1055 1024 qsurf(1:ngridmx,igcm_hdo_ice )1025 & =qsurf(1:ngridmx,igcm_h2o_ice )*DoverH1056 qsurf(1:ngridmx,igcm_hdo_ice,:) 1057 & =qsurf(1:ngridmx,igcm_h2o_ice,:)*DoverH 1026 1058 1027 1059 … … 1249 1281 write(*,*)'also set negative values of surf ice to 0' 1250 1282 do ig=1,ngridmx 1251 qsurf(ig,igcm_h2o_ice)=min(val,qsurf(ig,igcm_h2o_ice)) 1252 qsurf(ig,igcm_h2o_ice)=max(0.,qsurf(ig,igcm_h2o_ice)) 1283 do islope=1,nslope 1284 qsurf(ig,igcm_h2o_ice,islope)= 1285 & min(val,qsurf(ig,igcm_h2o_ice,islope)) 1286 qsurf(ig,igcm_h2o_ice,islope)= 1287 & max(0.,qsurf(ig,igcm_h2o_ice,islope)) 1288 enddo 1253 1289 end do 1254 1290 … … 1261 1297 write(*,*) 'OK: remove surface ice for |lat|<45' 1262 1298 if (abs(rlatu(j)*180./pi).lt.45.) then 1263 qsurf(ig,igcm_h2o_ice)=0. 1299 do islope=1,nslope 1300 qsurf(ig,igcm_h2o_ice,islope)=0. 1301 enddo 1264 1302 end if 1265 1303 end do … … 1273 1311 if(ig.eq.1) j=1 1274 1312 if (rlatu(j)*180./pi.gt.80.) then 1275 qsurf(ig,igcm_h2o_ice)=1.e5 1276 write(*,*) 'ig=',ig,' H2O ice mass (kg/m2)= ', 1277 & qsurf(ig,igcm_h2o_ice) 1313 do islope=1,nslope 1314 qsurf(ig,igcm_h2o_ice,islope)=1.e5 1315 write(*,*) 'ig=',ig,', islope=', islope, 1316 & ' H2O ice mass (kg/m2)= ', 1317 & qsurf(ig,igcm_h2o_ice,islope) 1278 1318 write(*,*)' ==> Ice mesh South boundary (deg)= ', 1279 1319 & rlatv(j)*180./pi 1320 enddo 1280 1321 end if 1281 1322 enddo … … 1288 1329 if(ig.eq.1) j=1 1289 1330 if (rlatu(j)*180./pi.lt.-80.) then 1290 qsurf(ig,igcm_h2o_ice)=1.e5 1291 write(*,*) 'ig=',ig,' H2O ice mass (kg/m2)= ', 1292 & qsurf(ig,igcm_h2o_ice) 1331 do islope=1,nslope 1332 qsurf(ig,igcm_h2o_ice,islope)=1.e5 1333 write(*,*) 'ig=',ig,', islope=', islope, 1334 & ' H2O ice mass (kg/m2)= ', 1335 & qsurf(ig,igcm_h2o_ice,islope) 1293 1336 write(*,*)' ==> Ice mesh North boundary (deg)= ', 1294 1337 & rlatv(j-1)*180./pi 1338 enddo 1295 1339 end if 1296 1340 enddo … … 1307 1351 1308 1352 do ig=1, ngridmx 1309 tsurf(ig) = Tiso 1353 do islope=1,nslope 1354 tsurf(ig,islope) = Tiso 1355 enddo 1310 1356 end do 1311 1357 do l=2,nsoilmx 1312 1358 do ig=1, ngridmx 1313 tsoil(ig,l) = Tiso 1359 do islope=1,nslope 1360 tsoil(ig,l,islope) = Tiso 1361 enddo 1314 1362 end do 1315 1363 end do … … 1323 1371 else if (trim(modif) .eq. 'co2ice=0') then 1324 1372 do ig=1,ngridmx 1325 qsurf(ig,igcm_co2)=0 1326 emis(ig)=emis(ngridmx/2) 1373 do islope=1,nslope 1374 qsurf(ig,igcm_co2,islope)=0 1375 emis(ig,islope)=emis(ngridmx/2,islope) 1376 enddo 1327 1377 end do 1328 1378 … … 1670 1720 ! recast ith() into ithfi() 1671 1721 CALL gr_dyn_fi(nsoilmx,iip1,jjp1,ngridmx,ith,ithfi) 1722 1723 else if (trim(modif) .eq. 'nslope') then 1724 write(*,*) 'set a new number of subgrid scale slope' 1725 write(*,*) 'Current value=', nslope 1726 write(*,*) 'Enter value for nslope (ex: 1,5,7)?' 1727 ierr=1 1728 do while (ierr.ne.0) 1729 read(*,*,iostat=ierr) nslope_new 1730 enddo 1731 1732 write(*,*) 'This can take some time...' 1733 write(*,*) 'You can go grab a coffee and relax a bit' 1734 1735 if(nslope.eq.nslope_new) then 1736 write(*,*) 'The number of subslope you entered is the same' 1737 write(*,*) 'as the number written in startfi.nc. ' 1738 write(*,*) 'Nothing will be done' 1739 else 1740 1741 nslope_old=nslope 1742 nslope=nslope_new 1743 1744 call end_comslope_h 1745 call ini_comslope_h(ngridmx,nslope_new) 1746 1747 allocate(default_def_slope(nslope_new+1)) 1748 !Sub-grid scale subslopes 1749 if (nslope_new.eq.7) then 1750 default_def_slope(1) = -43. 1751 default_def_slope(2) = -19. 1752 default_def_slope(3) = -9. 1753 default_def_slope(4) = -3. 1754 default_def_slope(5) = 3. 1755 default_def_slope(6) = 9. 1756 default_def_slope(7) = 19. 1757 default_def_slope(8) = 43. 1758 elseif (nslope_new.eq.5) then 1759 default_def_slope(1) = -43. 1760 default_def_slope(2) = -9. 1761 default_def_slope(3) = -3. 1762 default_def_slope(4) = 3. 1763 default_def_slope(5) = 9. 1764 default_def_slope(6) = 43. 1765 elseif (nslope_new.eq.1) then 1766 default_def_slope(1) = -50. 1767 default_def_slope(2) = 50. 1768 endif 1769 1770 do islope=1,nslope_new+1 1771 def_slope(islope) = default_def_slope(islope) 1772 enddo 1773 1774 do islope=1,nslope_new 1775 def_slope_mean(islope) =(def_slope(islope)+def_slope(islope+1))/2. 1776 enddo 1777 1778 iflat = 1 1779 DO islope=2,nslope_new 1780 IF(abs(def_slope_mean(islope)).lt. 1781 & abs(def_slope_mean(iflat)))THEN 1782 iflat = islope 1783 ENDIF 1784 ENDDO 1785 1786 call subslope_mola(ngridmx,nslope_new,def_slope,subslope_dist) 1787 1788 ! Surfdat related stuff 1672 1789 1790 allocate(tsurf_old_slope(ngridmx,nslope_old)) 1791 allocate(qsurf_old_slope(ngridmx,nqtot,nslope_old)) 1792 allocate(emis_old_slope(ngridmx,nslope_old)) 1793 allocate(watercap_old_slope(ngridmx,nslope_old)) 1794 1795 tsurf_old_slope(:,:)=tsurf(:,:) 1796 qsurf_old_slope(:,:,:)=qsurf(:,:,:) 1797 emis_old_slope(:,:)=emis(:,:) 1798 watercap_old_slope(:,:)=watercap(:,:) 1799 1800 call end_surfdat_h_slope_var 1801 call ini_surfdat_h_slope_var(ngridmx,nqtot,nslope_new) 1802 1803 ! Comsoil related stuff (tsoil) 1804 1805 allocate(tsoil_old_slope(ngridmx,nsoilmx,nslope_old)) 1806 1807 tsoil_old_slope(:,:,:)=tsoil(:,:,:) 1808 1809 call end_comsoil_h_slope_var 1810 call ini_comsoil_h_slope_var(ngridmx,nslope_new) 1811 1812 ! Dimradmars related stuff (albedo) 1813 1814 allocate(albedo_old_slope(ngridmx,2,nslope_old)) 1815 albedo_old_slope(:,:,:)=albedo(:,:,:) 1816 1817 call end_dimradmars_mod_slope_var 1818 call ini_dimradmars_mod_slope_var(ngridmx,nslope_new) 1819 1820 if(nslope_old.eq.1 .and. nslope_new.gt.1) then 1821 do islope=1,nslope_new 1822 tsurf(:,islope)=tsurf_old_slope(:,1) 1823 qsurf(:,:,islope)=qsurf_old_slope(:,:,1) 1824 emis(:,islope)=emis_old_slope(:,1) 1825 watercap(:,islope)=watercap_old_slope(:,1) 1826 tsoil(:,:,islope)=tsoil_old_slope(:,:,1) 1827 albedo(:,:,islope)=albedo_old_slope(:,:,1) 1828 enddo 1829 elseif(nslope_new.eq.1) then 1830 tsurf(:,1)=tsurf_old_slope(:,iflat) 1831 qsurf(:,:,1)=qsurf_old_slope(:,:,iflat) 1832 emis(:,1)=emis_old_slope(:,iflat) 1833 watercap(:,1)=watercap_old_slope(:,iflat) 1834 tsoil(:,:,1)=tsoil_old_slope(:,:,iflat) 1835 albedo(:,:,1)=albedo_old_slope(:,:,iflat) 1836 elseif(nslope_old.eq.5 .and. nslope_new.eq.7) then 1837 do islope=1,nslope_new 1838 tsurf(:,islope)=tsurf_old_slope(:,iflat) 1839 qsurf(:,:,islope)=qsurf_old_slope(:,:,iflat) 1840 emis(:,islope)=emis_old_slope(:,iflat) 1841 watercap(:,islope)=watercap_old_slope(:,iflat) 1842 tsoil(:,:,islope)=tsoil_old_slope(:,:,iflat) 1843 albedo(:,:,islope)=albedo_old_slope(:,:,iflat) 1844 enddo 1845 elseif(nslope_old.eq.7 .and. nslope_new.eq.5) then 1846 do islope=1,nslope_new 1847 tsurf(:,islope)=tsurf_old_slope(:,iflat) 1848 qsurf(:,:,islope)=qsurf_old_slope(:,:,iflat) 1849 emis(:,islope)=emis_old_slope(:,iflat) 1850 watercap(:,islope)=watercap_old_slope(:,iflat) 1851 tsoil(:,:,islope)=tsoil_old_slope(:,:,iflat) 1852 albedo(:,:,islope)=albedo_old_slope(:,:,iflat) 1853 enddo 1854 else 1855 write(*,*)' Problem choice of nslope' 1856 write(*,*)' Value not taken into account' 1857 CALL ABORT 1858 endif 1859 1860 endif !nslope=nslope_new 1861 1673 1862 else 1674 1863 write(*,*) ' Unknown (misspelled?) option!!!' … … 1786 1975 & nsoilmx,ngridmx,llm, 1787 1976 & nqtot,dtphys,real(day_ini),0.0,cell_area, 1788 & albfi,ithfi) 1977 & albfi,ithfi,def_slope, 1978 & subslope_dist) 1789 1979 call physdem1("restartfi.nc",nsoilmx,ngridmx,llm,nqtot, 1790 1980 & dtphys,hour_ini, -
trunk/LMDZ.MARS/libf/phymars/comslope_mod.F90
r2900 r2909 23 23 24 24 CONTAINS 25 SUBROUTINE ini_comslope_h(ngrid )25 SUBROUTINE ini_comslope_h(ngrid,nslope_in) 26 26 27 27 IMPLICIT NONE 28 28 INTEGER, INTENT(IN) :: ngrid 29 INTEGER, INTENT(IN) :: nslope_in 29 30 30 allocate(def_slope(nslope +1))31 allocate(def_slope_mean(nslope ))32 allocate(sky_slope(nslope ))33 allocate(subslope_dist(ngrid,nslope ))31 allocate(def_slope(nslope_in+1)) 32 allocate(def_slope_mean(nslope_in)) 33 allocate(sky_slope(nslope_in)) 34 allocate(subslope_dist(ngrid,nslope_in)) 34 35 allocate(major_slope(ngrid)) 35 36 END SUBROUTINE ini_comslope_h … … 70 71 qsurf_meshavg(:,:) = 0. 71 72 73 if(nslope.eq.1) then 74 albedo_meshavg(:,:) = albedo_slope(:,:,1) 75 emis_meshavg(:) = emis_slope(:,1) 76 tsurf_meshavg(:) = tsurf_slope(:,1) 77 qsurf_meshavg(:,:) = qsurf_slope(:,:,1) 78 else 79 72 80 DO ig = 1,ngrid 73 81 DO islope = 1,nslope … … 84 92 ENDDO 85 93 94 endif 95 86 96 END SUBROUTINE compute_meshgridavg 87 97 -
trunk/LMDZ.MARS/libf/phymars/comsoil_h.F90
r2900 r2909 70 70 end subroutine end_comsoil_h 71 71 72 subroutine ini_comsoil_h_slope_var(ngrid,nslope) 73 74 implicit none 75 integer,intent(in) :: ngrid ! number of atmospheric columns 76 integer,intent(in) :: nslope ! number of sub grid slopes 77 78 allocate(tsoil(ngrid,nsoilmx,nslope)) ! soil temperatures 79 allocate(mthermdiff(ngrid,0:nsoilmx-1,nslope)) 80 allocate(thermdiff(ngrid,nsoilmx-1,nslope)) 81 allocate(coefd(ngrid,nsoilmx-1,nslope)) 82 allocate(alph(ngrid,nsoilmx-1,nslope)) 83 allocate(beta(ngrid,nsoilmx-1,nslope)) 84 85 end subroutine ini_comsoil_h_slope_var 86 87 88 subroutine end_comsoil_h_slope_var 89 90 implicit none 91 92 if (allocated(tsoil)) deallocate(tsoil) 93 if (allocated(mthermdiff)) deallocate(mthermdiff) 94 if (allocated(thermdiff)) deallocate(thermdiff) 95 if (allocated(coefd)) deallocate(coefd) 96 if (allocated(alph)) deallocate(alph) 97 if (allocated(beta)) deallocate(beta) 98 99 end subroutine end_comsoil_h_slope_var 100 72 101 end module comsoil_h -
trunk/LMDZ.MARS/libf/phymars/dimradmars_mod.F90
r2900 r2909 216 216 217 217 end subroutine end_dimradmars_mod 218 219 subroutine ini_dimradmars_mod_slope_var(ngrid,nslope) 220 221 implicit none 222 223 integer,intent(in) :: ngrid ! number of atmospheric columns 224 integer,intent(in) :: nslope ! number of subgrid scale slopes 225 226 allocate(albedo(ngrid,2,nslope)) 227 allocate(fluxrad_sky(ngrid,nslope)) 228 allocate(fluxrad(ngrid,nslope)) 229 230 end subroutine ini_dimradmars_mod_slope_var 231 232 subroutine end_dimradmars_mod_slope_var 233 234 implicit none 235 236 if (allocated(albedo)) deallocate(albedo) 237 if (allocated(fluxrad_sky)) deallocate(fluxrad_sky) 238 if (allocated(fluxrad)) deallocate(fluxrad) 239 240 end subroutine end_dimradmars_mod_slope_var 218 241 219 242 -
trunk/LMDZ.MARS/libf/phymars/phyetat0_mod.F90
r2896 r2909 134 134 endif ! of if (startphy_file) 135 135 136 if(nslope.ne.1) then137 call abort_physic(modname, &138 "phyetat0: For now, nslope should be 1 (set in comslope_mod)",1)139 endif140 141 allocate(default_def_slope(nslope+1))142 !Sub-grid scale subslopes143 if (nslope.eq.7) then144 default_def_slope(1) = -43.145 default_def_slope(2) = -19.146 default_def_slope(3) = -9.147 default_def_slope(4) = -3.148 default_def_slope(5) = 3.149 default_def_slope(6) = 9.150 default_def_slope(7) = 19.151 default_def_slope(8) = 43.152 elseif (nslope.eq.5) then153 default_def_slope(1) = -43.154 default_def_slope(2) = -9.155 default_def_slope(3) = -3.156 default_def_slope(4) = 3.157 default_def_slope(5) = 9.158 default_def_slope(6) = 43.159 elseif (nslope.eq.1) then160 default_def_slope(1) = 0.161 default_def_slope(2) = 0.162 endif163 164 136 if (startphy_file) then 165 137 call get_var("def_slope",def_slope,found) … … 167 139 startphy_slope=.false. 168 140 write(*,*)'slope_settings: Problem while reading <def_slope>' 169 write(*,*)'default def_slope will be used' 170 do islope=1,nslope+1 171 def_slope(islope) = default_def_slope(islope) 172 enddo 173 write(*,*)'computing corresponding distribution <subslope_dist>' 174 write(*,*)'For now, woth nslope=1, subslope_dist is straigforward' 175 write(*,*)'Later this operation will be done by newstart with a specific routine' 176 subslope_dist(:,:)=1. 177 !call subslope_mola(ngrid,nslope,def_slope,subslope_dist) 141 write(*,*)'Checking that nslope=1' 142 if(nslope.eq.1) then 143 write(*,*)'We take default def_slope and subslope dist' 144 allocate(default_def_slope(nslope+1)) 145 default_def_slope(1) = 0. 146 default_def_slope(2) = 0. 147 subslope_dist(:,:)=1. 148 else 149 write(*,*)'Variable def_slope is not present in the start and' 150 write(*,*)'you are trying to run with nslope!=1' 151 write(*,*)'This is not possible, you have to go through newstart' 152 endif 178 153 else 179 154 startphy_slope=.true. … … 181 156 if(.not.found) then 182 157 write(*,*)'slope_settings: Problem while reading <subslope_dist>' 183 write(*,*)'computing a new distribution' 184 write(*,*)'For now, woth nslope=1, subslope_dist is straigforward' 185 write(*,*)'Later this operation will be done by newstart with a specific routine' 186 subslope_dist(:,:)=1. 187 !call subslope_mola(ngrid,nslope,def_slope,subslope_dist) 158 write(*,*)'We have to abort.' 159 write(*,*)'Please check that nslope is coherent, you can modify this parameter with newstart' 160 call abort_physic(modname, & 161 "phyetat0: Failed loading <subslope_dist>",1) 188 162 endif 189 163 endif 190 else ! startphy_file 191 do islope=1,nslope+1 192 def_slope(islope) = default_def_slope(islope) 193 enddo 194 write(*,*)'computing corresponding distribution <subslope_dist>' 195 write(*,*)'For now, woth nslope=1, subslope_dist is straigforward' 196 write(*,*)'Later this operation will be done by newstart with a specific routine' 197 subslope_dist(:,:)=1. 198 !call subslope_mola(ngrid,nslope,def_slope,subslope_dist) 164 else ! (no startphy_file) 165 if(nslope.eq.1) then 166 allocate(default_def_slope(nslope+1)) 167 default_def_slope(1) = 0. 168 default_def_slope(2) = 0. 169 subslope_dist(:,:)=1. 170 else 171 write(*,*)'Without starfi, lets first run with nslope=1' 172 call abort_physic(modname, & 173 "phyetat0: No startfi and nslope!=1",1) 174 endif 199 175 endif 200 176 -
trunk/LMDZ.MARS/libf/phymars/phys_state_var_init_mod.F90
r2900 r2909 61 61 use dust_rad_adjust_mod, only: ini_dust_rad_adjust_mod, & 62 62 end_dust_rad_adjust_mod 63 use comslope_mod, ONLY: nslope 63 use comslope_mod, ONLY: nslope,end_comslope_h,ini_comslope_h 64 64 use netcdf 65 65 USE mod_phys_lmdz_para, ONLY: is_master,bcast … … 182 182 call ini_dust_rad_adjust_mod(ngrid) 183 183 184 ! allocate arrays in "comslope_mod" 185 call end_comslope_h 186 call ini_comslope_h(ngrid,nslope) 187 184 188 END SUBROUTINE phys_state_var_init 185 189 -
trunk/LMDZ.MARS/libf/phymars/physiq_mod.F
r2908 r2909 577 577 c ~~~~~~~~~~~~ 578 578 #ifndef MESOSCALE 579 580 call ini_comslope_h(ngrid)581 579 582 580 ! GCM. Read netcdf initial physical parameters. -
trunk/LMDZ.MARS/libf/phymars/surfdat_h.F90
r2900 r2909 116 116 end subroutine end_surfdat_h 117 117 118 subroutine ini_surfdat_h_slope_var(ngrid,nq,nslope) 119 120 implicit none 121 integer,intent(in) :: ngrid ! number of atmospheric columns 122 integer,intent(in) :: nq ! number of tracers 123 integer,intent(in) :: nslope ! number of sub-grid scale slope 124 allocate(qsurf(ngrid,nq,nslope)) 125 allocate(tsurf(ngrid,nslope)) 126 allocate(watercap(ngrid,nslope)) 127 allocate(emis(ngrid,nslope)) 128 allocate(capcal(ngrid,nslope)) 129 allocate(fluxgrd(ngrid,nslope)) 130 131 end subroutine ini_surfdat_h_slope_var 132 133 subroutine end_surfdat_h_slope_var 134 135 implicit none 136 137 if (allocated(qsurf)) deallocate(qsurf) 138 if (allocated(tsurf)) deallocate(tsurf) 139 if (allocated(watercap)) deallocate(watercap) 140 if (allocated(emis)) deallocate(emis) 141 if (allocated(capcal)) deallocate(capcal) 142 if (allocated(fluxgrd)) deallocate(fluxgrd) 143 144 end subroutine end_surfdat_h_slope_var 145 118 146 end module surfdat_h
Note: See TracChangeset
for help on using the changeset viewer.