Changeset 07c3e71 in flexpart.git
- Timestamp:
- Jul 24, 2019, 12:16:51 PM (5 years ago)
- Branches:
- master, 10.4.1_pesei, GFS_025, bugfixes+enhancements, dev, release-10, release-10.4.1, scaling-bug
- Children:
- 25b4532
- Parents:
- 95a45d3
- Location:
- src
- Files:
-
- 6 edited
Legend:
- Unmodified
- Added
- Removed
-
src/FLEXPART.f90
r5fc7b68 r07c3e71 62 62 implicit none 63 63 64 integer :: i,j,ix,jy,inest 64 integer :: i,j,ix,jy,inest, iopt 65 65 integer :: idummy = -320 66 66 character(len=256) :: inline_options !pathfile, flexversion, arg2 … … 80 80 ! FLEXPART version string 81 81 flexversion_major = '10' ! Major version number, also used for species file names 82 flexversion='Version '//trim(flexversion_major)//'.4 (2019-07- 16)'82 flexversion='Version '//trim(flexversion_major)//'.4 (2019-07-23)' 83 83 verbosity=0 84 84 … … 109 109 print*,'FLEXPART is free software released under the GNU General Public License.' 110 110 111 112 ! Ingest inline options 113 !******************************************************* 111 114 if (inline_options(1:1).eq.'-') then 112 if (trim(inline_options).eq.'-v'.or.trim(inline_options).eq.'-v1') then 113 print*, 'Verbose mode 1: display detailed information during run' 114 verbosity=1 115 endif 116 if (trim(inline_options).eq.'-v2') then 117 print*, 'Verbose mode 2: display more detailed information during run' 118 verbosity=2 119 endif 115 print*,'inline_options:',inline_options 116 !verbose mode 117 iopt=index(inline_options,'v') 118 if (iopt.gt.0) then 119 verbosity=1 120 !print*, iopt, inline_options(iopt+1:iopt+1) 121 if (trim(inline_options(iopt+1:iopt+1)).eq.'2') then 122 !print*, 'verbosity=2' 123 print*, 'Verbose mode 2: display more detailed information during run' 124 verbosity=2 125 endif 126 endif 127 128 129 !debug mode 130 iopt=index(inline_options,'d') 131 if (iopt.gt.0) then 132 debug_mode=.true. 133 endif 134 endif 135 136 137 138 ! stop 139 ! if (trim(inline_options).eq.'-v'.or.trim(inline_options).eq.'-v1') then 140 ! print*, 'Verbose mode 1: display detailed information during run' 141 ! verbosity=1 142 ! endif 143 ! if (trim(inline_options).eq.'-v2') then 144 ! print*, 'Verbose mode 2: display more detailed information during run' 145 ! verbosity=2 146 ! endif 147 120 148 if (trim(inline_options).eq.'-i') then 121 149 print*, 'Info mode: provide detailed run specific information and stop' … … 135 163 print*, 'nzmax=',nzmax 136 164 print*,'nxshift=',nxshift 165 endif 166 167 if (verbosity.gt.0) then 137 168 write(*,*) 'call readpaths' 138 169 endif … … 193 224 ! Detect metdata format 194 225 !********************** 226 if (verbosity.gt.0) then 227 write(*,*) 'call detectformat' 228 endif 195 229 196 230 metdata_format = detectformat() -
src/com_mod.f90
r2eefa58 r07c3e71 762 762 ! Verbosity, testing flags, namelist I/O 763 763 !******************** 764 logical :: debug_mode=.false. 764 765 integer :: verbosity=0 765 766 integer :: info_flag=0 -
src/getvdep.f90
r5f9d14a r07c3e71 170 170 if ((ra+rb(i)+rc(i)).gt.0.) then 171 171 vd=1./(ra+rb(i)+rc(i)) 172 ! XXXXXXXXXXXXXXXXXXXXXXXXXX TEST173 ! vd=1./rc(i)174 ! XXXXXXXXXXXXXXXXXXXXXXXXXX TEST175 172 else 176 173 vd=9.999 … … 188 185 call partdep(nspec,density,fract,schmi,vset,raquer,ust,nyl,vdepo) 189 186 187 !if (debug_mode) then 188 ! print*,'getvdep:188: vdepo=', vdepo 189 !stop 190 !endif 190 191 191 192 ! 7. If no detailed parameterization available, take constant deposition -
src/part0.f90
re200b7a r07c3e71 103 103 x01=alog(d01/dquer)/xdummy 104 104 x02=alog(d02/dquer)/xdummy 105 105 !print*,'part0:: d02=' , d02 , 'd01=', d01 106 106 107 107 ! Area under Gauss-function is calculated and gives mass fraction of interval … … 109 109 110 110 fract(i)=0.5*(erf(x01)-erf(x02)) 111 111 !print*,'part0:: fract(',i,')', fract(i) 112 !print*,'part0:: fract', fract(i), x01, x02, erf(x01), erf(x02) 112 113 113 114 ! Geometric mean diameter of interval in [m] … … 115 116 116 117 dmean=1.E-6*exp(0.5*alog(d01*d02)) 117 118 !print*,'part0:: dmean=', dmean 118 119 119 120 ! Calculation of time independent parameters of each interval … … 132 133 vsh(i)=ga*density*dmean*dmean*cun/(18.*myl) 133 134 135 !print*,'part0:: vsh(',i,')', vsh(i) 136 134 137 end do 135 138 139 !stop 'part0' 140 136 141 end subroutine part0 -
src/partdep.f90
r5f9d14a r07c3e71 72 72 73 73 use par_mod 74 use com_mod, only: debug_mode 74 75 75 76 implicit none … … 110 111 111 112 vdep(ic)=vdep(ic)+vdepj*fract(ic,j) 113 114 !print*, 'partdep:113: ic', ic, 'vdep', vdep 115 !stop 112 116 end do 113 117 endif 118 119 114 120 end do 115 121 122 !if (debug_mode) then 123 ! print*, 'partdep:122:' 124 ! write(*,*) (vdep(ic), ic=1,nc) 125 !stop 126 !endif 127 116 128 end subroutine partdep -
src/readspecies.f90
r46f93d5 r07c3e71 296 296 if (density(pos_spec) .gt. 0) then 297 297 write(*,'(a)') ' Dry deposition is turned : ON' 298 else 298 if (reldiff(pos_spec).gt.0) then 299 stop 'density>0 (SPECIES is a particle) implies reldiff <=0 ' 300 endif 301 else 302 if (reldiff(pos_spec).le.0) then 303 stop 'density<=0 (SPECIES is a gas) implies reldiff >0 ' 304 endif 299 305 write(*,'(a)') ' Dry deposition is (density<0) : OFF' 300 306 end if … … 344 350 end if 345 351 346 if (dsigma(i).eq.0.) dsigma(i)=1.0001 ! avoid floating exception 352 ! if (dsigma(i).eq.0.) dsigma(i)=1.0001 ! avoid floating exception 353 if (dsigma(i).le.1.) then !dsigma(i)=1.0001 ! avoid floating exception 354 !write(*,*) '#### FLEXPART MODEL ERROR! ####' 355 write(*,*) '#### FLEXPART MODEL WARNING ####' 356 write(*,*) '#### in SPECIES_',aspecnumb, ' ####' 357 write(*,*) '#### from v10.4 dsigma has to be larger than 1 ####' 358 write(*,*) '#### to adapt older SPECIES files, ####' 359 write(*,*) '#### if dsigma was < 1 ####' 360 write(*,*) '#### use the reciprocal of the old dsigma ####' 361 if (.not.debug_mode) then 362 stop 363 else 364 write(*,*) 'debug mode: continue' 365 endif 366 endif 347 367 348 368 if ((reldiff(i).gt.0.).and.(density(i).gt.0.)) then
Note: See TracChangeset
for help on using the changeset viewer.