Changeset 1145 for trunk/LMDZ.GENERIC/libf
- Timestamp:
- Jan 4, 2014, 4:51:31 PM (11 years ago)
- Location:
- trunk/LMDZ.GENERIC/libf/phystd
- Files:
-
- 4 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.GENERIC/libf/phystd/callcorrk.F90
r1133 r1145 646 646 if(tlevrad(k).lt.tgasmin)then 647 647 print*,'Minimum temperature is outside the radiative' 648 print*,'transfer kmatrix bounds , exiting.'648 print*,'transfer kmatrix bounds' 649 649 print*,"k=",k," tlevrad(k)=",tlevrad(k) 650 650 print*,"tgasmin=",tgasmin 651 !print*,'WARNING, OVERRIDING FOR TEST' 652 call abort 653 !tlevrad(k)=tgasmin 651 if (strictboundcorrk) then 652 call abort 653 else 654 print*,'***********************************************' 655 print*,'we allow model to continue with tlevrad=tgasmin' 656 print*,' ... we assume we know what you are doing ... ' 657 print*,' ... but do not let this happen too often ... ' 658 print*,'***********************************************' 659 tlevrad(k)=tgasmin 660 endif 654 661 elseif(tlevrad(k).gt.tgasmax)then 655 662 print*,'Maximum temperature is outside the radiative' 656 663 print*,'transfer kmatrix bounds, exiting.' 657 print*,'level,grid,T',k,ig,tlevrad(k) 658 !print*,'WARNING, OVERRIDING FOR TEST' 659 call abort 660 !tlevrad(k)=tgasmax 664 print*,"k=",k," tlevrad(k)=",tlevrad(k) 665 print*,"tgasmax=",tgasmax 666 if (strictboundcorrk) then 667 call abort 668 else 669 print*,'***********************************************' 670 print*,'we allow model to continue with tlevrad=tgasmax' 671 print*,' ... we assume we know what you are doing ... ' 672 print*,' ... but do not let this happen too often ... ' 673 print*,'***********************************************' 674 tlevrad(k)=tgasmax 675 endif 661 676 endif 662 677 enddo … … 665 680 print*,'Minimum temperature is outside the radiative' 666 681 print*,'transfer kmatrix bounds, exiting.' 667 print*,"k=",k," t levrad(k)=",tlevrad(k)682 print*,"k=",k," tmid(k)=",tmid(k) 668 683 print*,"tgasmin=",tgasmin 669 !print*,'WARNING, OVERRIDING FOR TEST' 670 call abort 671 !tmid(k)=tgasmin 684 if (strictboundcorrk) then 685 call abort 686 else 687 print*,'***********************************************' 688 print*,'we allow model to continue with tmid=tgasmin' 689 print*,' ... we assume we know what you are doing ... ' 690 print*,' ... but do not let this happen too often ... ' 691 print*,'***********************************************' 692 tmid(k)=tgasmin 693 endif 672 694 elseif(tmid(k).gt.tgasmax)then 673 695 print*,'Maximum temperature is outside the radiative' 674 696 print*,'transfer kmatrix bounds, exiting.' 675 print*,'level,grid,T',k,ig,tlevrad(k) 676 !print*,'WARNING, OVERRIDING FOR TEST' 677 call abort 678 !tmid(k)=tgasmax 697 print*,"k=",k," tmid(k)=",tmid(k) 698 print*,"tgasmax=",tgasmax 699 if (strictboundcorrk) then 700 call abort 701 else 702 print*,'***********************************************' 703 print*,'we allow model to continue with tmid=tgasmin' 704 print*,' ... we assume we know what you are doing ... ' 705 print*,' ... but do not let this happen too often ... ' 706 print*,'***********************************************' 707 tmid(k)=tgasmax 708 endif 679 709 endif 680 710 enddo -
trunk/LMDZ.GENERIC/libf/phystd/callkeys.h
r1133 r1145 31 31 & , cloudlvl & 32 32 & , pceil & 33 & , strictboundcorrk & 33 34 & , rings_shadow 34 35 … … 38 39 & , callstats,calleofdump & 39 40 & , callgasvis,continuum,H2Ocont_simple,graybody & 41 & , strictboundcorrk & 40 42 & , rings_shadow 41 43 -
trunk/LMDZ.GENERIC/libf/phystd/gfluxi.F
r1016 r1145 104 104 ! -- same results for most thin atmospheres 105 105 ! -- and stabilizes integrations 106 NT = int(TLEV(2*L+1)*NTfac) - NTstar+1106 !NT = int(TLEV(2*L+1)*NTfac) - NTstar+1 107 107 !! For deep, opaque, thick first layers (e.g. Saturn) 108 108 !! what is below works much better, not unstable, ... 109 109 !! ... and actually fully accurate because 1st layer temp (JL) 110 !NT = int(TLEV(2*L)*NTfac) - NTstar+1110 NT = int(TLEV(2*L)*NTfac) - NTstar+1 111 111 !! (or this one yields same results 112 112 !NT = int( (TLEV(2*L)+TLEV(2*L+1))*0.5*NTfac ) - NTstar+1 -
trunk/LMDZ.GENERIC/libf/phystd/inifis.F
r1133 r1145 199 199 call getin("corrk",corrk) 200 200 write(*,*) " corrk = ",corrk 201 202 write(*,*) "prohibit calculations outside corrk T grid?" 203 strictboundcorrk=.true. ! default value 204 call getin("strictboundcorrk",strictboundcorrk) 205 write(*,*) "strictboundcorrk = ",strictboundcorrk 201 206 202 207 write(*,*) "call gaseous absorption in the visible bands?",
Note: See TracChangeset
for help on using the changeset viewer.