Changeset 1201 for LMDZ4/branches/LMDZ4-dev/libf/phylmd
- Timestamp:
- Jul 7, 2009, 4:01:00 PM (15 years ago)
- Location:
- LMDZ4/branches/LMDZ4-dev/libf/phylmd
- Files:
-
- 2 added
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
LMDZ4/branches/LMDZ4-dev/libf/phylmd/hgardfou.F
r1145 r1201 1 1 ! 2 ! $Id$ 2 3 SUBROUTINE hgardfou (t,tsol,text) 3 4 use dimphy … … 12 13 REAL t(klon,klev), tsol(klon,nbsrf) 13 14 CHARACTER*(*) text 15 character (len=20) :: modname = 'hgardfou' 16 character (len=80) :: abort_message 14 17 C 15 18 INTEGER i, k, nsrf … … 124 127 c 125 128 IF (.NOT. ok) THEN 126 PRINT*, 'hgardfou s arrete ',text127 CALL abort 129 abort_message= 'hgardfou s arrete '//text 130 CALL abort_gcm (modname,abort_message,1) 128 131 ENDIF 129 132 -
LMDZ4/branches/LMDZ4-dev/libf/phylmd/physiq.F
r1196 r1201 4 4 5 5 SUBROUTINE physiq (nlon,nlev, 6 . debut,lafin, rjourvrai,gmtime,pdtphys,6 . debut,lafin,jD_cur, jH_cur,pdtphys, 7 7 . paprs,pplay,pphi,pphis,presnivs,clesphy0, 8 8 . u,v,t,qx, … … 60 60 c debut---input-L-variable logique indiquant le premier passage 61 61 c lafin---input-L-variable logique indiquant le dernier passage 62 c rjour---input-R-numero du jour de l'experience63 c gmtime--input-R-temps universel dans la journee (0 a 86400 s)62 c jD_cur--input-R-jour courant a l'appel de la physique (jour julien) 63 c jH_cur--input-R-heure courante a l'appel de la physique (jour julien) 64 64 c pdtphys-input-R-pas d'integration pour la physique (seconde) 65 65 c paprs---input-R-pression pour chaque inter-couche (en Pa) … … 186 186 INTEGER nlon 187 187 INTEGER nlev 188 REAL rjourvrai189 REAL gmtime 188 REAL :: jD_cur, jH_cur 189 190 190 REAL pdtphys 191 191 LOGICAL debut, lafin … … 703 703 c Conditions aux limites 704 704 c 705 INTEGER julien 705 ! 706 ! Gestion calendrier 707 ! 708 REAL :: jD_1jan, jH_1jan 709 INTEGER :: year_cur, mth_cur, day_cur, days_elapsed 710 REAL :: hour, day_since_equinox 711 ! Date de l'equinoxe de printemps 712 INTEGER, parameter :: mth_eq=3, day_eq=21 713 REAL :: jD_eq 714 715 LOGICAL, parameter :: new_orbit = .true. 716 706 717 c 707 718 INTEGER lmt_pas … … 1087 1098 c$OMP THREADPRIVATE(first) 1088 1099 1100 integer iunit 1101 1089 1102 logical, save:: read_climoz ! read ozone climatology 1090 1103 integer, save:: ncid_climoz ! NetCDF file containing ozone climatology … … 1114 1127 write(lunout,*) 'DEBUT DE PHYSIQ !!!!!!!!!!!!!!!!!!!!' 1115 1128 write(lunout,*) 1116 s 'nlon,klev,nqtot,debut,lafin, rjourvrai,gmtime,pdtphys'1129 s 'nlon,klev,nqtot,debut,lafin, jD_cur, jH_cur,pdtphys' 1117 1130 write(lunout,*) 1118 s nlon,klev,nqtot,debut,lafin, rjourvrai,gmtime,pdtphys1131 s nlon,klev,nqtot,debut,lafin, jD_cur, jH_cur,pdtphys 1119 1132 1120 1133 write(lunout,*) 'papers, play, phi, u, v, t, omega' … … 1184 1197 1185 1198 c====================================================================== 1186 xjour = rjourvrai 1199 ! Gestion calendrier 1200 ! 1201 call ju2ymds(jD_cur+jH_cur, year_cur, mth_cur, day_cur, hour) 1202 call ymds2ju(year_cur, 1, 1, 0., jD_1jan) 1203 jH_1jan = jD_1jan - int (jD_1jan) 1204 jD_1jan = int (jD_1jan) 1205 xjour = jD_cur - jD_1jan 1206 days_elapsed = jD_cur - jD_1jan 1207 1187 1208 c 1188 1209 c Si c'est le debut, il faut initialiser plusieurs choses … … 1472 1493 1473 1494 cXXXPB Positionner date0 pour initialisation de ORCHIDEE 1474 date0 = zjulian 1475 C date0 = day_ini 1495 date0 = jD_ref 1476 1496 WRITE(*,*) 'physiq date0 : ',date0 1477 1497 c … … 1491 1511 CALL VTe(VTphysiq) 1492 1512 CALL VTb(VTinca) 1493 iii = MOD(NINT(xjour),360) 1494 calday = FLOAT(iii) + gmtime 1495 WRITE(lunout,*) 'initial time ', xjour, calday 1513 ! iii = MOD(NINT(xjour),360) 1514 ! calday = FLOAT(iii) + jH_cur 1515 calday = FLOAT(days_elapsed) + jH_cur 1516 WRITE(lunout,*) 'initial time chemini', days_elapsed, calday 1496 1517 1497 1518 CALL chemini( … … 1534 1555 ! 1535 1556 itap = itap + 1 1536 julien = MOD(NINT(xjour),360)1537 if (julien .eq. 0) julien = 3601538 1539 1557 ! 1540 1558 ! Update fraction of the sub-surfaces (pctsrf) and … … 1542 1560 ! on the surface fraction. 1543 1561 ! 1544 CALL change_srf_frac(itap, dtime, julien,1562 CALL change_srf_frac(itap, dtime, days_elapsed+1, 1545 1563 * pctsrf, falb1, falb2, ftsol, u10m, v10m, pbl_tke) 1546 1547 1564 1548 1565 ! Tendances bidons pour les processus qui n'affectent pas certaines … … 1697 1714 if (read_climoz) then 1698 1715 C Ozone climatology from a NetCDF file 1699 call regr_pr_av(ncid_climoz, "tro3", julien, press_climoz, 1716 call regr_pr_av(ncid_climoz, "tro3", days_elapsed+1, 1717 & press_climoz, 1700 1718 $ paprs, wo) 1701 1719 ! Convert from mole fraction of ozone to column density of ozone in a … … 1707 1725 C "zmasse" changes a little.) 1708 1726 else 1709 CALL ozonecm(real(julien), rlat, paprs, wo) 1727 CALL ozonecm(real(days_elapsed+1), rlat, paprs, wo) 1728 1710 1729 end if 1711 1730 ENDIF … … 1749 1768 ! doit donc etre placé avant radlwsw et pbl_surface 1750 1769 1770 ! calcul selon la routine utilisee pour les planetes 1771 if (new_orbit) then 1772 call ymds2ju(year_cur, mth_eq, day_eq,0., jD_eq) 1773 day_since_equinox = (jD_cur + jH_cur) - jD_eq 1774 ! day_since_equinox = (jD_cur) - jD_eq 1775 call solarlong(day_since_equinox, zlongi, dist) 1776 else 1777 ! calcul selon la routine utilisee pour l'AR4 1751 1778 ! choix entre calcul de la longitude solaire vraie ou valeur fixee a 1752 1779 ! solarlong0 1753 1754 if (solarlong0<-999.) then1755 CALL orbite(FLOAT(julien),zlongi,dist)1756 else1757 zlongi=solarlong0 ! longitude solaire vraie1758 dist=1. ! distance au soleil / moyenne1780 if (solarlong0<-999.) then 1781 CALL orbite(FLOAT(days_elapsed+1),zlongi,dist) 1782 else 1783 zlongi=solarlong0 ! longitude solaire vraie 1784 dist=1. ! distance au soleil / moyenne 1785 endif 1759 1786 endif 1760 1761 if(prt_level.ge.1) print*,'Longitude solaire ',zlongi,solarlong01787 if(prt_level.ge.1) & 1788 & write(lunout,*)'Longitude solaire ',zlongi,solarlong0,dist 1762 1789 1763 1790 ! Avec ou sans cycle diurne 1764 1791 IF (cycle_diurne) THEN 1765 1792 zdtime=dtime*FLOAT(radpas) ! pas de temps du rayonnement (s) 1766 CALL zenang(zlongi, gmtime,zdtime,rlat,rlon,rmu0,fract)1793 CALL zenang(zlongi,jH_cur,zdtime,rlat,rlon,rmu0,fract) 1767 1794 ELSE 1768 1795 CALL angle(zlongi, rlat, fract, rmu0) … … 1797 1824 1798 1825 CALL pbl_surface( 1799 e dtime, date0, itap, julien,1826 e dtime, date0, itap, days_elapsed+1, 1800 1827 e debut, lafin, 1801 1828 e rlon, rlat, rugoro, rmu0, … … 2094 2121 2095 2122 if (itop_con(i).gt.klev-3) then 2096 print*,'La convection monte trop haut ' 2097 print*,'itop_con(,',i,',)=',itop_con(i) 2123 if(prt_level >= 9) then 2124 write(lunout,*)'La convection monte trop haut ' 2125 write(lunout,*)'itop_con(,',i,',)=',itop_con(i) 2126 endif 2098 2127 endif 2099 2128 ENDDO … … 2635 2664 IF (.NOT. aerosol_couple) 2636 2665 & CALL readaerosol_optic( 2637 & debut, new_aod, flag_aerosol, rjourvrai, pdtphys,2666 & debut, new_aod, flag_aerosol, jD_cur-jD_ref, pdtphys, 2638 2667 & pplay, paprs, t_seri, rhcl, 2639 2668 & mass_solu_aero, mass_solu_aero_pi, … … 2764 2793 CALL VTe(VTphysiq) 2765 2794 CALL VTb(VTinca) 2766 calday = FLOAT( julien) + gmtime2795 calday = FLOAT(days_elapsed + 1) + jH_cur 2767 2796 2768 2797 IF (config_inca == 'aero') THEN … … 2775 2804 2776 2805 CALL chemhook_begin (calday, 2777 $ julien,2778 $ gmtime,2806 $ days_elapsed, 2807 $ jH_cur, 2779 2808 $ pctsrf(1,1), 2780 2809 $ rlat, … … 3106 3135 IF (is_sequential) THEN 3107 3136 3108 CALL aaam_bud (27,klon,klev, rjourvrai,gmtime,3137 CALL aaam_bud (27,klon,klev,jD_cur-jD_ref,jH_cur, 3109 3138 C ra,rg,romega, 3110 3139 C rlat,rlon,pphis, … … 3133 3162 3134 3163 call phytrac ( 3135 I itap, julien, gmtime, debut,3164 I itap, days_elapsed+1, jH_cur, debut, 3136 3165 I lafin, dtime, u, v, t, 3137 3166 I paprs, pplay, pmfu, pmfd, … … 3351 3380 write(lunout,*) 'FIN DE PHYSIQ !!!!!!!!!!!!!!!!!!!!' 3352 3381 write(lunout,*) 3353 s 'nlon,klev,nqtot,debut,lafin, rjourvrai,gmtime,pdtphys pct tlos'3382 s 'nlon,klev,nqtot,debut,lafin,jD_cur, jH_cur, pdtphys pct tlos' 3354 3383 write(lunout,*) 3355 s nlon,klev,nqtot,debut,lafin, rjourvrai,gmtime,pdtphys,3384 s nlon,klev,nqtot,debut,lafin, jD_cur, jH_cur ,pdtphys, 3356 3385 s pctsrf(igout,is_ter), pctsrf(igout,is_lic),pctsrf(igout,is_oce), 3357 3386 s pctsrf(igout,is_sic) … … 3447 3476 ENDIF 3448 3477 3478 ! first=.false. 3449 3479 3450 3480 RETURN
Note: See TracChangeset
for help on using the changeset viewer.