Changeset 2910
- Timestamp:
- Mar 10, 2023, 6:05:41 PM (22 months ago)
- Location:
- trunk/LMDZ.MARS
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/README
r2909 r2910 3931 3931 It can create start with 1, 5 or 7 subslope abd needs mola64.nc to do so. 3932 3932 Correction in compute_mesh_avg for nslope=1: no average is done only a copy. 3933 3934 == 10/03/2023 == EM 3935 Make subslope_mola a module and turn (extremely) large static arrays into 3936 allocated ones and assume that the "mola64.nc" file is in "datadir" rather 3937 than a hard coded location. -
trunk/LMDZ.MARS/libf/dynphy_lonlat/phymars/newstart.F
r2909 r2910 55 55 use comslope_mod, ONLY: nslope,def_slope,def_slope_mean, 56 56 & subslope_dist,end_comslope_h,ini_comslope_h 57 USE subslope_mola_mod, ONLY: subslope_mola 58 57 59 implicit none 58 60 -
trunk/LMDZ.MARS/libf/dynphy_lonlat/phymars/subslope_mola.F90
r2909 r2910 1 MODULE subslope_mola_mod 2 3 IMPLICIT NONE 4 5 CONTAINS 6 1 7 SUBROUTINE subslope_mola(ngridmx,nslope,def_slope,subslope_dist) 2 8 … … 6 12 north_east, north_west, & 7 13 south_west, south_east 14 USE datafile_mod, ONLY: datadir 8 15 9 16 IMPLICIT NONE … … 15 22 integer :: jjm_mola, iim_mola 16 23 parameter(jjm_mola=180*resol, iim_mola=2*jjm_mola) 17 character*(*) dset ! Path to MCD datafiles18 parameter(dset ='/ccc/scratch/cont003/gen10391/vandemer/test_pem2_watercap/datadir/')19 24 integer :: ierr 20 25 ! returned status code (==0 if OK) 21 real ::theta_mola(iim_mola,jjm_mola)26 real,allocatable :: theta_mola(:,:) !theta_mola(iim_mola,jjm_mola) 22 27 !ED18 : slopes inclination (see getslope.F90) 23 real ::psi_mola(iim_mola,jjm_mola)28 real,allocatable :: psi_mola(:,:) !psi_mola(iim_mola,jjm_mola) 24 29 !ED18 : slopes orientation (idem) 25 30 REAL :: lon1, lon2, lat1, lat2 !bounds of the square … … 39 44 LOGICAL,PARAMETER :: rfz = .false. 40 45 INTEGER,PARAMETER :: igps = 1000 !first ig for RFZ location 41 46 47 ! allocate big arrays 48 allocate(theta_mola(iim_mola,jjm_mola),stat=ierr) 49 if (ierr/=0) then 50 write(*,*)"subslope_mola: allocation of theta_mola(:,:) failed!" 51 stop 52 endif 53 allocate(psi_mola(iim_mola,jjm_mola),stat=ierr) 54 if (ierr/=0) then 55 write(*,*)"subslope_mola: allocation of psi_mola(:,:) failed!" 56 stop 57 endif 58 59 42 60 !-------------Building of theta_mola and psi_mola 43 CALL mola(dset,ierr,theta_mola,psi_mola,resol,iim_mola,jjm_mola) 61 ! Assume that the mola file is to ben found in "datadir" 62 CALL mola(trim(datadir)//"/",& 63 ierr,theta_mola,psi_mola,resol,iim_mola,jjm_mola) 44 64 45 65 !-------------Building of distribution … … 83 103 max_crit,max_lon,max_lat 84 104 105 ! deallocate big arrays (this routine is only called once) 106 deallocate(theta_mola,psi_mola) 107 85 108 END SUBROUTINE subslope_mola 86 109 … … 159 182 trim(molafile) 160 183 endif 161 ierr=15 !set appropriate error code 162 return 184 stop 163 185 endif 164 186 … … 174 196 write(*,*)"Error in mola: <alt> not found" 175 197 endif 176 ierr=16 ! set appropriate error code 177 return 198 stop 178 199 endif 179 200 … … 364 385 INTEGER :: i_lon1, i_lon2 365 386 INTEGER :: j_lat1, j_lat2 366 REAL :: crit !function387 ! REAL :: crit !function 367 388 REAL :: val_crit !intermediate for evaluating the criterium 368 389 REAL :: max_crit !indicator for maximum slope criterium computed … … 531 552 532 553 533 END SUBROUTINE slopes_stat 534 554 CONTAINS 535 555 536 556 !*********************************************************************** … … 543 563 544 564 REAL crit 545 REAL theta, psi565 REAL,INTENT(IN) :: theta, psi 546 566 real, parameter :: pi = acos(-1.) 547 567 … … 553 573 554 574 575 END SUBROUTINE slopes_stat 576 577 578 579 END MODULE subslope_mola_mod
Note: See TracChangeset
for help on using the changeset viewer.