Changeset 2395 for LMDZ5/trunk/libf/phydev/phyetat0.F90
- Timestamp:
- Nov 18, 2015, 1:25:35 PM (9 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
LMDZ5/trunk/libf/phydev/phyetat0.F90
r1907 r2395 2 2 ! $Id $ 3 3 ! 4 subroutinephyetat0(fichnom)4 SUBROUTINE phyetat0(fichnom) 5 5 ! Load initial state for the physics 6 6 ! and do some resulting initializations 7 7 8 use iostart, only : open_startphy,get_field,close_startphy 9 use iophy, only : init_iophy_new 10 use phys_state_var_mod, only : rlat,rlon 8 USE dimphy, only: klon 9 USE iostart, ONLY : open_startphy,get_field,close_startphy 10 USE iophy, ONLY : init_iophy_new 11 USE geometry_mod, ONLY : longitude_deg, latitude_deg 11 12 12 implicit none 13 IMPLICIT NONE 13 14 14 character(len=*),intent(in) :: fichnom ! input file name15 CHARACTER(len=*),INTENT(in) :: fichnom ! input file name 15 16 16 ! open physics initial state file: 17 call open_startphy(fichnom) 17 REAL :: lon_startphy(klon), lat_startphy(klon) 18 INTEGER :: i 18 19 19 ! read latitudes 20 call get_field("latitude",rlat)20 ! open physics initial state file: 21 CALL open_startphy(fichnom) 21 22 22 ! read longitudes 23 call get_field("longitude",rlon) 23 ! read latitudes and make a sanity check (because already known from dyn) 24 CALL get_field("latitude",lat_startphy) 25 DO i=1,klon 26 IF (ABS(lat_startphy(i)-latitude_deg(i))>=1) THEN 27 WRITE(*,*) "phyetat0: Error! Latitude discrepancy wrt startphy file:",& 28 " i=",i," lat_startphy(i)=",lat_startphy(i),& 29 " latitude_deg(i)=",latitude_deg(i) 30 ! This is presumably serious enough to abort run 31 CALL abort_physic("phyetat0","discrepancy in latitudes!",1) 32 ENDIF 33 IF (ABS(lat_startphy(i)-latitude_deg(i))>=0.0001) THEN 34 WRITE(*,*) "phyetat0: Warning! Latitude discrepancy wrt startphy file:",& 35 " i=",i," lat_startphy(i)=",lat_startphy(i),& 36 " latitude_deg(i)=",latitude_deg(i) 37 ENDIF 38 ENDDO 24 39 25 ! read in other variables here ... 40 ! read longitudes and make a sanity check (because already known from dyn) 41 CALL get_field("longitude",lon_startphy) 42 DO i=1,klon 43 IF (ABS(lon_startphy(i)-longitude_deg(i))>=1) THEN 44 WRITE(*,*) "phyetat0: Error! Longitude discrepancy wrt startphy file:",& 45 " i=",i," lon_startphy(i)=",lon_startphy(i),& 46 " longitude_deg(i)=",longitude_deg(i) 47 ! This is presumably serious enough to abort run 48 CALL abort_physic("phyetat0","discrepancy in longitudes!",1) 49 ENDIF 50 IF (ABS(lon_startphy(i)-longitude_deg(i))>=0.0001) THEN 51 WRITE(*,*) "phyetat0: Warning! Longitude discrepancy wrt startphy file:",& 52 " i=",i," lon_startphy(i)=",lon_startphy(i),& 53 " longitude_deg(i)=",longitude_deg(i) 54 ENDIF 55 ENDDO 26 56 27 ! close file 28 call close_startphy 57 ! read in other variables here ... 29 58 30 ! do some more initializations 31 call init_iophy_new(rlat,rlon) 59 ! close file 60 CALL close_startphy 32 61 33 end subroutine phyetat0 62 ! do some more initializations 63 CALL init_iophy_new(latitude_deg,longitude_deg) 64 65 END SUBROUTINE phyetat0
Note: See TracChangeset
for help on using the changeset viewer.