Changeset 2570


Ignore:
Timestamp:
Oct 22, 2021, 5:23:36 PM (3 years ago)
Author:
emillour
Message:

Mars GCM:
Improve tests by "check_physics_fields": also look for
negative values of tracers.
EM

Location:
trunk/LMDZ.MARS
Files:
3 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/README

    r2566 r2570  
    34853485else
    34863486  zcondicea, zfallice
     3487
     3488== 22/10/2021 == EM
     3489Improve tests by "check_physics_fields": also look for
     3490negative values of tracers.
  • trunk/LMDZ.MARS/libf/phymars/check_fields.F90

    r2541 r2570  
    1111contains
    1212
    13 subroutine check_physics_fields(message,temp,u,v,pplev)
    14 use dimphy, only: klon, klev
    15 implicit none
    16 character(len=*),intent(in):: message
     13subroutine check_physics_fields(message,temp,u,v,pplev,q)
     14use dimphy, only: klon, klev
     15use tracer_mod, only: nqmx
     16implicit none
     17character(len=*),intent(in):: message ! text message for outputs
    1718real,intent(in) :: temp(klon,klev)
    1819real,intent(in) :: u(klon,klev) ! zonal wind (m/s)
    1920real,intent(in) :: v(klon,klev) ! meridional wind (m/s)
    2021real,intent(in) :: pplev(klon,klev+1) ! pressure at level interfaces (Pa)
     22real,intent(in) :: q(klon,klev,nqmx) ! tracer mixing ratio (.../kg_air)
    2123
    2224character(len=50) :: name="check_physics_fields"
    23 logical ok_t,ok_w,ok_ps
     25logical :: ok_t,ok_w,ok_ps,ok_q
    2426
    2527! 1. Initialisations
    26 ok_t=.true.
    27 ok_w=.true.
    28 ok_ps=.true.
     28ok_q=.true.
    2929
    3030! 2. Check temperature, winds and surface pressure
    31 call check_temperature(temp,ok_t)
    32 call check_winds(u,v,ok_w)
    33 call check_ps(pplev(:,1),ok_ps)
    34 
    35 if ((.not.ok_t).or.(.not.ok_w).or.(.not.ok_ps)) then
     31call check_temperature(message,temp,ok_t)
     32call check_winds(message,u,v,ok_w)
     33call check_ps(message,pplev(:,1),ok_ps)
     34if (nqmx>=1) then
     35  call check_tracers(message,q,ok_q)
     36endif
     37
     38if ((.not.ok_t).or.(.not.ok_w).or.(.not.ok_ps).or.(.not.ok_q)) then
    3639  ! Something is wrong, might as well stop here
    3740  call abort_physic(trim(name),trim(message)//" Invalid field values",1)
     
    4144
    4245
    43 subroutine check_temperature(temp,ok,temp_min,temp_max)
    44 use dimphy, only: klon, klev
    45 implicit none
     46subroutine check_temperature(message,temp,ok,temp_min,temp_max)
     47use dimphy, only: klon, klev
     48implicit none
     49character(len=*),intent(in):: message ! text message for outputs
    4650real,intent(in) :: temp(klon,klev)
    4751logical,intent(out) :: ok ! returns .true. if everything OK, .false. otherwise
     
    7579    if (temp(i,k).ne.temp(i,k)) then
    7680      ok=.false.
    77       write(*,*)trim(name)//" temp(i,k)=",temp(i,k)," for i=",i," k=",k
     81      write(*,*)trim(message)//" "//trim(name)//" temp(i,k)=",temp(i,k),&
     82                " for i=",i," k=",k
    7883    endif
    7984    ! check for temperatures too low
    8085    if (temp(i,k).lt.tmin) then
    8186      ok=.false.
    82       write(*,*)trim(name)//" temp(i,k)=",temp(i,k)," for i=",i," k=",k,&
    83       "<",tmin
     87      write(*,*)trim(message)//" "//trim(name)//" temp(i,k)=",temp(i,k),&
     88      " for i=",i," k=",k," <",tmin
    8489    endif
    8590    ! check for temperatures too high
    8691    if (temp(i,k).gt.tmax) then
    8792      ok=.false.
    88       write(*,*)trim(name)//" temp(i,k)=",temp(i,k)," for i=",i," k=",k,&
    89       ">",tmax
     93      write(*,*)trim(message)//" "//trim(name)//" temp(i,k)=",temp(i,k),&
     94      " for i=",i," k=",k," >",tmax
    9095    endif
    9196  enddo
     
    9499end subroutine check_temperature
    95100
    96 subroutine check_winds(u,v,ok,wind_max)
    97 use dimphy, only: klon, klev
    98 implicit none
     101subroutine check_winds(message,u,v,ok,wind_max)
     102use dimphy, only: klon, klev
     103implicit none
     104character(len=*),intent(in):: message ! text message for outputs
    99105real,intent(in) :: u(klon,klev) ! zonal wind (m/s)
    100106real,intent(in) :: v(klon,klev) ! meridional wind (m/s)
     
    123129    if (u(i,k).ne.u(i,k)) then
    124130      ok=.false.
    125       write(*,*)trim(name)//" u(i,k)=",u(i,k)," for i=",i," k=",k
     131      write(*,*)trim(message)//" "//trim(name)//" u(i,k)=",u(i,k),&
     132      " for i=",i," k=",k
    126133    endif
    127134    if (v(i,k).ne.v(i,k)) then
    128135      ok=.false.
    129       write(*,*)trim(name)//" v(i,k)=",v(i,k)," for i=",i," k=",k
     136      write(*,*)trim(message)//" "//trim(name)//" v(i,k)=",v(i,k),&
     137      " for i=",i," k=",k
    130138    endif
    131139    ! check for magnitudes too high
    132140    if (abs(u(i,k)).gt.wmax) then
    133141      ok=.false.
    134       write(*,*)trim(name)//" u(i,k)=",u(i,k)," for i=",i," k=",k,&
    135       ">",wmax
     142      write(*,*)trim(message)//" "//trim(name)//" u(i,k)=",u(i,k),&
     143      " for i=",i," k=",k," >",wmax
    136144    endif
    137145    if (abs(v(i,k)).gt.wmax) then
    138146      ok=.false.
    139       write(*,*)trim(name)//" v(i,k)=",v(i,k)," for i=",i," k=",k,&
    140       ">",wmax
     147      write(*,*)trim(message)//" "//trim(name)//" v(i,k)=",v(i,k),&
     148      " for i=",i," k=",k," >",wmax
    141149    endif
    142150  enddo
     
    145153end subroutine check_winds
    146154
    147 subroutine check_ps(ps,ok,ps_min,ps_max)
     155subroutine check_ps(message,ps,ok,ps_min,ps_max)
    148156use dimphy, only: klon
    149157implicit none
     158character(len=*),intent(in):: message ! text message for outputs
    150159real,intent(in) :: ps(klon)
    151160logical,intent(out) :: ok ! returns .true. if everything OK, .false. otherwise
     
    178187  if (ps(i).ne.ps(i)) then
    179188    ok=.false.
    180     write(*,*)trim(name)//" ps(i)=",ps(i)," for i=",i
     189    write(*,*)trim(message)//" "//trim(name)//" ps(i)=",ps(i)," for i=",i
    181190  endif
    182191  ! check for pressures too low
    183192  if (ps(i).lt.pmin) then
    184193    ok=.false.
    185     write(*,*)trim(name)//" ps(i)=",ps(i)," for i=",i,"<",pmin
     194    write(*,*)trim(message)//" "//trim(name)//" ps(i)=",ps(i)," for i=",i,&
     195    " <",pmin
    186196  endif
    187197  ! check for pressures too high
    188198  if (ps(i).gt.pmax) then
    189199    ok=.false.
    190     write(*,*)trim(name)//" ps(i)=",ps(i)," for i=",i,">",pmax
     200    write(*,*)trim(message)//" "//trim(name)//" ps(i)=",ps(i)," for i=",i,&
     201    " >",pmax
    191202  endif
    192203enddo
     
    194205end subroutine check_ps
    195206
     207subroutine check_tracers(message,q,ok)
     208use dimphy, only: klon, klev
     209use tracer_mod, only: nqmx,noms
     210implicit none
     211character(len=*),intent(in):: message ! text message for outputs
     212real,intent(in) :: q(klon,klev,nqmx) ! tracer mixing ratio
     213logical,intent(inout) :: ok ! set to .false. if something is wrong
     214
     215character(len=50) :: name="check_tracers"
     216integer :: i,k,iq
     217integer :: nb_bad_neg ! number of problematic points (negative values)
     218integer :: nb_bad_nan ! number of problematic points (NaNs)
     219
     220! 1. initializations
     221nb_bad_nan=0
     222
     223! 2. do the checking
     224do iq=1,nqmx
     225 nb_bad_neg=0 ! initialize
     226 do i=1,klon
     227  do k=1,klev
     228    ! look for NaN
     229    if (q(i,k,iq).ne.q(i,k,iq)) then
     230      ok=.false.
     231      nb_bad_nan=nb_bad_nan+1
     232    endif
     233    ! look for negative values
     234    if (q(i,k,iq).lt.0) then
     235!      ok=.false.
     236      nb_bad_neg=nb_bad_neg+1
     237    endif
     238  enddo
     239 enddo
     240 if (nb_bad_neg>0) then
     241   write(*,*)trim(message)//" "//trim(name)//" tracer "//trim(noms(iq))//&
     242   " contains ",nb_bad_neg," negative values!"
     243 endif
     244enddo !of do iq=1,nqmx
     245
     246end subroutine check_tracers
     247
    196248end module check_fields_mod
  • trunk/LMDZ.MARS/libf/phymars/physiq_mod.F

    r2566 r2570  
    743743      if (check_physics_inputs) then
    744744        ! Check the validity of input fields coming from the dynamics
    745         call check_physics_fields("begin physiq:",pt,pu,pv,pplev)
     745        call check_physics_fields("begin physiq:",pt,pu,pv,pplev,pq)
    746746      endif
    747747
     
    39783978      if (check_physics_outputs) then
    39793979        ! Check the validity of updated fields at the end of the physics step
    3980         call check_physics_fields("end of physiq:",zt,zu,zv,zplev)
     3980        call check_physics_fields("end of physiq:",zt,zu,zv,zplev,zq)
    39813981      endif
    39823982
Note: See TracChangeset for help on using the changeset viewer.