Changeset 2131 for trunk/LMDZ.GENERIC/libf
- Timestamp:
- Apr 29, 2019, 2:52:42 PM (6 years ago)
- Location:
- trunk/LMDZ.GENERIC/libf/phystd
- Files:
-
- 4 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.GENERIC/libf/phystd/callkeys_mod.F90
r2062 r2131 20 20 logical,save :: specOLR 21 21 logical,save :: kastprof 22 !$OMP THREADPRIVATE(enertest,nonideal,meanOLR,kastprof) 22 logical,save :: diagdtau 23 !$OMP THREADPRIVATE(enertest,nonideal,meanOLR,kastprof,diagdtau) 23 24 logical,save :: newtonian 24 25 logical,save :: check_cpp_match -
trunk/LMDZ.GENERIC/libf/phystd/inifis_mod.F90
r2078 r2131 327 327 write(*,*)" specOLR = ",specOLR 328 328 329 write(*,*)"Output diagnostic optical thickness?" 330 diagdtau=.false. 331 call getin_p("diagdtau",diagdtau) 332 write(*,*)" diagdtau = ",diagdtau 333 ! sanity check 334 if (diagdtau.and.(ngrid.ne.1)) then 335 print*,"Diagnostic optical thickness can be output in 1D only !" 336 print*,"Start again with diagdtau=.false." 337 stop 338 endif 339 329 340 write(*,*)"Operate in kastprof mode?" 330 341 kastprof=.false. -
trunk/LMDZ.GENERIC/libf/phystd/optci.F90
r2032 r2131 11 11 use radinc_h, only: L_LEVELS, L_NLAYRAD, L_NSPECTI, L_NGAUSS, & 12 12 L_NLEVRAD, L_REFVAR, naerkind 13 use radcommon_h, only: gasi,tlimit,wrefVAR,Cmk,tgasref,pfgasref,wnoi,scalep,indi,glat_ig 13 use radcommon_h, only: gasi,tlimit,wrefVAR,Cmk,tgasref,pfgasref,wnoi,scalep,indi,glat_ig,gweight 14 14 use gases_h, only: gfrac, ngasmx, igas_N2, igas_He, igas_H2O, igas_H2 15 15 use comcstfi_mod, only: g, r, mugaz 16 use callkeys_mod, only: kastprof,continuum,graybody,H2Ocont_simple 16 use callkeys_mod, only: kastprof,continuum,graybody,H2Ocont_simple,diagdtau 17 17 implicit none 18 18 … … 84 84 !real*8 rho !! see test below 85 85 86 ! variables for diags of optical depth 87 real *8 INT_DTAU(L_NLAYRAD,L_NSPECTI) 88 character*2 str2 89 86 90 integer igas, jgas 87 91 … … 449 453 ! print*,DTAUI 450 454 ! call abort 455 456 457 ! Optical thickness for 1D diagnostics (added by JVO) 458 if (diagdtau) then ! diagtau can be true only if 1D 459 do l=1,L_NLAYRAD 460 do nw=1,L_NSPECTI 461 INT_DTAU(L,NW) = 0.0d+0 462 DO NG=1,L_NGAUSS 463 INT_DTAU(L,NW)= INT_DTAU(L,NW) + dtaui(L,nw,ng)*gweight(NG) 464 enddo 465 enddo 466 enddo 467 do nw=1,L_NSPECTI 468 write(str2,'(i2.2)') nw 469 call writediagfi(1,'dtaui'//str2,'Layer optical thickness in IR band '//str2,'',1,int_dtau(L_NLAYRAD:1:-1,nw)) 470 enddo 471 endif 472 451 473 452 474 -
trunk/LMDZ.GENERIC/libf/phystd/optcv.F90
r2032 r2131 10 10 11 11 use radinc_h, only: L_NLAYRAD, L_NLEVRAD, L_LEVELS, L_NSPECTV, L_NGAUSS, L_REFVAR, NAERKIND 12 use radcommon_h, only: gasv, tlimit, wrefVAR, Cmk, tgasref, pfgasref,wnov,scalep,indv,glat_ig 12 use radcommon_h, only: gasv, tlimit, wrefVAR, Cmk, tgasref, pfgasref,wnov,scalep,indv,glat_ig,gweight 13 13 use gases_h, only: gfrac, ngasmx, igas_H2, igas_H2O, igas_He, igas_N2 14 14 use comcstfi_mod, only: g, r, mugaz 15 use callkeys_mod, only: kastprof,continuum,graybody,H2Ocont_simple,callgasvis 15 use callkeys_mod, only: kastprof,continuum,graybody,H2Ocont_simple,callgasvis,diagdtau 16 16 17 17 implicit none … … 91 91 ! variables for k in units m^-1 92 92 real*8 dz(L_LEVELS) 93 94 ! variables for diags of optical depth 95 real*8 INT_DTAU(L_NLAYRAD,L_NSPECTV) 96 character*2 str2 93 97 94 98 integer igas, jgas … … 375 379 376 380 381 ! Optical thickness for 1D diagnostics (added by JVO) 382 if (diagdtau) then ! diagtau can be true only if 1D 383 do l=1,L_NLAYRAD 384 do nw=1,L_NSPECTV 385 INT_DTAU(L,NW) = 0.0d+0 386 DO NG=1,L_NGAUSS 387 INT_DTAU(L,NW)= INT_DTAU(L,NW) + dtauv(L,nw,ng)*gweight(NG) 388 enddo 389 enddo 390 enddo 391 do nw=1,L_NSPECTV 392 write(str2,'(i2.2)') nw 393 call writediagfi(1,'dtauv'//str2,'Layer optical thickness in VI band '//str2,'',1,int_dtau(L_NLAYRAD:1:-1,nw)) 394 enddo 395 endif 396 397 377 398 end subroutine optcv 378 399
Note: See TracChangeset
for help on using the changeset viewer.