Changeset 563 for trunk/LMDZ.MARS/libf/dyn3d
- Timestamp:
- Mar 6, 2012, 11:36:52 AM (13 years ago)
- Location:
- trunk/LMDZ.MARS/libf/dyn3d
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/libf/dyn3d/lect_start_archive.F
r311 r563 770 770 read(*,*) tmpval 771 771 qsurfold(1:imold+1,1:jmold+1,iq)=tmpval 772 E NDIF773 #ifdef NC_DOUBLE 774 ierr = NF_GET_VARA_DOUBLE(nid, nvarid,start,count,772 ELSE ! tracer exists in file, load it 773 #ifdef NC_DOUBLE 774 ierr = NF_GET_VARA_DOUBLE(nid, nvarid,start,count, 775 775 & qsurfold(1,1,iq)) 776 776 #else 777 ierr = NF_GET_VARA_REAL(nid, nvarid,start,count,777 ierr = NF_GET_VARA_REAL(nid, nvarid,start,count, 778 778 & qsurfold(1,1,iq)) 779 779 #endif 780 IF (ierr .NE. NF_NOERR) THEN781 PRINT*, "lect_start_archive: ",780 IF (ierr .NE. NF_NOERR) THEN 781 PRINT*, "lect_start_archive: ", 782 782 & " Failed loading <",trim(txt),">" 783 print*, "which (constant) value should it be initialized to?" 784 read(*,*) tmpval 785 qsurfold(1:imold+1,1:jmold+1,iq)=tmpval 783 stop 784 ENDIF 786 785 ENDIF 787 786 … … 952 951 read(*,*) tmpval 953 952 qold(1:imold+1,1:jmold+1,1:lmold,iq)=tmpval 954 E NDIF955 #ifdef NC_DOUBLE 956 ierr=NF_GET_VARA_DOUBLE(nid,nvarid,start,count,qold(1,1,1,iq))957 #else 958 ierr=NF_GET_VARA_REAL(nid,nvarid,start,count,qold(1,1,1,iq))959 #endif 960 IF (ierr .NE. NF_NOERR) THEN961 PRINT*, "lect_start_archive: ",953 ELSE ! tracer exists in file, load it 954 #ifdef NC_DOUBLE 955 ierr=NF_GET_VARA_DOUBLE(nid,nvarid,start,count,qold(1,1,1,iq)) 956 #else 957 ierr=NF_GET_VARA_REAL(nid,nvarid,start,count,qold(1,1,1,iq)) 958 #endif 959 IF (ierr .NE. NF_NOERR) THEN 960 PRINT*, "lect_start_archive: ", 962 961 & " Failed loading <",trim(txt),">" 963 print*, "which (constant) value should it be initialized to?" 964 read(*,*) tmpval 965 qold(1:imold+1,1:jmold+1,1:lmold,iq)=tmpval 962 stop 963 ENDIF 966 964 ENDIF 967 965 … … 1044 1042 write(*,*)'Old grid: mass of CO2 ice:',co2icetotalold 1045 1043 write(*,*)'New grid: mass of CO2 ice:',co2icetotal 1046 write(*,*)'Ratio new ice / old ice =',co2icetotal/co2icetotalold 1044 if (co2icetotalold.ne.0.) then 1045 write(*,*)'Ratio new ice / old ice =',co2icetotal/co2icetotalold 1046 endif 1047 1047 write(*,*) 1048 1048 -
trunk/LMDZ.MARS/libf/dyn3d/newstart.F
r480 r563 155 155 156 156 INTEGER :: nq,numvanle 157 character(len= 20) :: txt ! to store some text157 character(len=50) :: txt ! to store some text 158 158 integer :: count 159 real :: profile(llm+1) ! to store an atmospheric profile + surface value 159 160 160 161 ! MONS data: … … 482 483 write(*,*) 'q=0 : ALL tracer =zero' 483 484 write(*,*) 'q=x : give a specific uniform value to one tracer' 485 write(*,*) 'q=profile : specify a profile for a tracer' 484 486 write(*,*) 'ini_q : tracers initialisation for chemistry, water an 485 487 $d ice ' … … 779 781 qsurf(ig,iq)=val 780 782 ENDDO 783 784 c q=profile : initialize tracer with a given profile 785 c -------------------------------------------------- 786 else if (trim(modif) .eq. 'q=profile') then 787 write(*,*) 'Tracer profile will be sought in ASCII file' 788 write(*,*) "'profile_tracer' where 'tracer' is tracer name" 789 write(*,*) "(one value per line in file; starting with" 790 write(*,*) "surface value, the 1st atmospheric layer" 791 write(*,*) "followed by 2nd, etc. up to top of atmosphere)" 792 write(*,*) 'Which tracer do you want to set?' 793 do iq=1,nqmx 794 write(*,*)iq,' : ',trim(tnom(iq)) 795 enddo 796 write(*,*) '(choose between 1 and ',nqmx,')' 797 read(*,*) iq 798 if ((iq.lt.1).or.(iq.gt.nqmx)) then 799 ! wrong value for iq, go back to menu 800 write(*,*) "wrong input value:",iq 801 cycle 802 endif 803 ! look for input file 'profile_tracer' 804 txt="profile_"//trim(tnom(iq)) 805 open(41,file=trim(txt),status='old',form='formatted', 806 & iostat=ierr) 807 if (ierr.eq.0) then 808 ! OK, found file 'profile_...', load the profile 809 do l=1,llm+1 810 read(41,*,iostat=ierr) profile(l) 811 if (ierr.ne.0) then ! something went wrong 812 exit ! quit loop 813 endif 814 enddo 815 if (ierr.eq.0) then 816 ! initialize tracer values 817 qsurf(:,iq)=profile(1) 818 do l=1,llm 819 q(:,:,l,iq)=profile(l+1) 820 enddo 821 write(*,*)'OK, tracer ',trim(tnom(iq)),' initialized ', 822 & 'using values from file ',trim(txt) 823 else 824 write(*,*)'problem reading file ',trim(txt),' !' 825 write(*,*)'No modifications to tracer ',trim(tnom(iq)) 826 endif 827 else 828 write(*,*)'Could not find file ',trim(txt),' !' 829 write(*,*)'No modifications to tracer ',trim(tnom(iq)) 830 endif 831 781 832 782 833 c ini_q : Initialize tracers for chemistry
Note: See TracChangeset
for help on using the changeset viewer.