Changeset 3698 for trunk/LMDZ.PLUTO/libf
- Timestamp:
- Mar 26, 2025, 11:05:58 AM (10 months ago)
- Location:
- trunk/LMDZ.PLUTO/libf/phypluto
- Files:
-
- 4 edited
-
callcorrk_pluto_mod.F90 (modified) (3 diffs)
-
physiq_mod.F90 (modified) (3 diffs)
-
surfprop.F90 (modified) (2 diffs)
-
writediagsoil.F90 (modified) (2 diffs)
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.PLUTO/libf/phypluto/callcorrk_pluto_mod.F90
r3685 r3698 28 28 strobel,vmrch4_proffix,specOLR,vmrch4fix,& 29 29 haze_radproffix,callmufi 30 use callkeys_mod, only: strictboundcorrk 30 31 use optcv_pluto_mod, only: optcv_pluto 31 32 use optci_pluto_mod, only: optci_pluto … … 210 211 real,save :: levdat(Nfine),vmrdat(Nfine) 211 212 real :: vmrch4(ngrid,nlayer) ! vmr ch4 from vmrch4_proffix 213 character(len=100) :: message 214 character(len=15),parameter :: subname="callcorrk_pluto" 212 215 213 216 !======================================================================= … … 634 637 if(tlevrad(k).lt.tgasmin)then 635 638 print*,'Minimum temperature is outside the radiative' 636 print*,'transfer kmatrix bounds , exiting.'639 print*,'transfer kmatrix bounds' 637 640 print*,'t(',k,')=',tlevrad(k),' < ',tgasmin 638 ! call abort 641 if (strictboundcorrk) then 642 message="Minimum temperature outside of kmatrix bounds" 643 call abort_physic(subname,message,1) 644 else 645 print*,'***********************************************' 646 print*,'we allow model to continue with tlevrad<tgasmin' 647 print*,' ... we assume we know what you are doing ... ' 648 print*,' ... but do not let this happen too often ... ' 649 print*,'***********************************************' 650 endif 639 651 elseif(tlevrad(k).gt.tgasmax)then 640 652 print*,'Maximum temperature is outside the radiative' 641 print*,'transfer kmatrix bounds , exiting.'653 print*,'transfer kmatrix bounds' 642 654 print*,'t(',k,')=',tlevrad(k),' > ',tgasmax 643 ! call abort 655 if (strictboundcorrk) then 656 message="Maximum temperature outside of kmatrix bounds" 657 call abort_physic(subname,message,1) 658 else 659 print*,'***********************************************' 660 print*,'we allow model to continue with tlevrad>tgasmax' 661 print*,' ... we assume we know what you are doing ... ' 662 print*,' ... but do not let this happen too often ... ' 663 print*,'***********************************************' 664 endif 644 665 endif 645 666 enddo -
trunk/LMDZ.PLUTO/libf/phypluto/physiq_mod.F90
r3695 r3698 70 70 use callcorrk_mod, only: callcorrk 71 71 use callcorrk_pluto_mod, only: callcorrk_pluto 72 use surfprop_mod, only: surfprop 72 73 use vdifc_mod, only: vdifc 73 74 use vdifc_pluto_mod, only: vdifc_pluto … … 1958 1959 ! Info about Ls, declin... 1959 1960 IF (fast) THEN 1960 if (is_master) write(*,*) ,'Ls=',zls*180./pi,' dec=',declin*180./pi1961 if (is_master) write(*,*) ,'zday=',zday,' ps=',globave1961 if (is_master) write(*,*)'Ls=',zls*180./pi,' dec=',declin*180./pi 1962 if (is_master) write(*,*)'zday=',zday,' ps=',globave 1962 1963 IF (lastcall) then 1963 if (is_master) write(*,*) ,'lastcall'1964 if (is_master) write(*,*)'lastcall' 1964 1965 ENDIF 1965 1966 ELSE 1966 if (is_master) write(*,*) ,'Ls=',zls*180./pi,'decli=',declin*180./pi,'zday=',zday1967 if (is_master) write(*,*)'Ls=',zls*180./pi,'decli=',declin*180./pi,'zday=',zday 1967 1968 ENDIF 1968 1969 … … 2139 2140 if (fast) then 2140 2141 call write_output("fluxrad","fluxrad","W m-2",fluxrad) 2141 call write_output("fluxgrd","fluxgrd","W m-2",fluxgrd)2142 2142 ! call write_output("dplanck","dplanck","W.s m-2 K-1",dplanck) 2143 endif 2144 2145 if (callsoil) then 2143 2146 ! "soil" variables 2144 call write_output("capcal"," capcal","W.s m-2 K-1",capcal)2147 call write_output("capcal","Surface Heat Capacity","W.s m-2 K-1",capcal) 2145 2148 call write_output("tsoil","tsoil","K",tsoil) 2146 2149 call write_output("therm_inertia","therm_inertia","S.I.",therm_inertia) -
trunk/LMDZ.PLUTO/libf/phypluto/surfprop.F90
r3539 r3698 1 module surfprop_mod 2 3 implicit none 4 5 contains 6 1 7 subroutine surfprop(ngrid,nq,pfra,qsurf,tsurface, & 2 8 capcal,adjust,dist,fluold,ptimestep,zls, & … … 507 513 end subroutine surfprop 508 514 515 end module surfprop_mod -
trunk/LMDZ.PLUTO/libf/phypluto/writediagsoil.F90
r3506 r3698 23 23 use mod_phys_lmdz_para, only : is_mpi_root, is_master, gather 24 24 use mod_grid_phy_lmdz, only : klon_glo, Grid1Dto2D_glo, & 25 nbp_lon, nbp_lat 25 nbp_lon, nbp_lat, & 26 grid_type, unstructured 26 27 27 28 implicit none … … 78 79 real px2(ngrid) 79 80 #endif 81 82 ! 0. This routine shoul only be used in lon-lat case 83 if (grid_type==unstructured) then 84 return 85 endif 80 86 81 87 ! 1. Initialization step
Note: See TracChangeset
for help on using the changeset viewer.
