Changeset 1516 for trunk/LMDZ.GENERIC/libf/phystd/soil.F
- Timestamp:
- Mar 17, 2016, 9:45:36 AM (9 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.GENERIC/libf/phystd/soil.F
r1315 r1516 4 4 & capcal,fluxgrd) 5 5 6 use comsoil_h, only: layer, mlayer, volcapa 6 use comsoil_h, only: layer, mlayer, volcapa, inertiedat 7 use comcstfi_mod, only: pi, daysec 8 use planete_mod, only: year_day 9 use comgeomfi_h, only: long, lati 7 10 8 11 implicit none … … 49 52 ! local variables: 50 53 integer ig,ik 51 54 real :: inertia_min,inertia_max 55 real :: diurnal_skin ! diurnal skin depth (m) 56 real :: annual_skin ! anuual skin depth (m) 52 57 53 58 ! 0. Initialisations and preprocessing step … … 118 123 !write(*,*)'soil: ig=',ig,' capcal(ig)=',capcal(ig) 119 124 enddo ! of do ig=1,ngrid 120 125 126 ! Additional checks: is the vertical discretization sufficient 127 ! to resolve diurnal and annual waves? 128 do ig=1,ngrid 129 ! extreme inertia for this column 130 inertia_min=minval(inertiedat(ig,:)) 131 inertia_max=maxval(inertiedat(ig,:)) 132 ! diurnal and annual skin depth 133 diurnal_skin=(inertia_min/volcapa)*sqrt(daysec/pi) 134 annual_skin=(inertia_max/volcapa)*sqrt(year_day*daysec/pi) 135 if (0.5*diurnal_skin<layer(1)) then 136 ! one should have the fist layer be at least half of diurnal skin depth 137 write(*,*) "soil Error: grid point ig=",ig 138 write(*,*) " longitude=",long(ig)*(180./pi) 139 write(*,*) " latitude=",lati(ig)*(180./pi) 140 write(*,*) " first soil layer depth ",layer(1) 141 write(*,*) " not small enough for a diurnal skin depth of ", 142 & diurnal_skin 143 write(*,*) " change soil layer distribution (comsoil_h.F90)" 144 stop 145 endif 146 if (2.*annual_skin>layer(nsoil)) then 147 ! one should have the full soil be at least twice the diurnal skin depth 148 write(*,*) "soil Error: grid point ig=",ig 149 write(*,*) " longitude=",long(ig)*(180./pi) 150 write(*,*) " latitude=",lati(ig)*(180./pi) 151 write(*,*) " total soil layer depth ",layer(nsoil) 152 write(*,*) " not large enough for an annual skin depth of ", 153 & annual_skin 154 write(*,*) " change soil layer distribution (comsoil_h.F90)" 155 stop 156 endif 157 enddo ! of do ig=1,ngrid 158 121 159 else ! of if (firstcall) 122 160
Note: See TracChangeset
for help on using the changeset viewer.