Changeset 2570
- Timestamp:
- Oct 22, 2021, 5:23:36 PM (3 years ago)
- Location:
- trunk/LMDZ.MARS
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/README
r2566 r2570 3485 3485 else 3486 3486 zcondicea, zfallice 3487 3488 == 22/10/2021 == EM 3489 Improve tests by "check_physics_fields": also look for 3490 negative values of tracers. -
trunk/LMDZ.MARS/libf/phymars/check_fields.F90
r2541 r2570 11 11 contains 12 12 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 13 subroutine check_physics_fields(message,temp,u,v,pplev,q) 14 use dimphy, only: klon, klev 15 use tracer_mod, only: nqmx 16 implicit none 17 character(len=*),intent(in):: message ! text message for outputs 17 18 real,intent(in) :: temp(klon,klev) 18 19 real,intent(in) :: u(klon,klev) ! zonal wind (m/s) 19 20 real,intent(in) :: v(klon,klev) ! meridional wind (m/s) 20 21 real,intent(in) :: pplev(klon,klev+1) ! pressure at level interfaces (Pa) 22 real,intent(in) :: q(klon,klev,nqmx) ! tracer mixing ratio (.../kg_air) 21 23 22 24 character(len=50) :: name="check_physics_fields" 23 logical ok_t,ok_w,ok_ps25 logical :: ok_t,ok_w,ok_ps,ok_q 24 26 25 27 ! 1. Initialisations 26 ok_t=.true. 27 ok_w=.true. 28 ok_ps=.true. 28 ok_q=.true. 29 29 30 30 ! 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 31 call check_temperature(message,temp,ok_t) 32 call check_winds(message,u,v,ok_w) 33 call check_ps(message,pplev(:,1),ok_ps) 34 if (nqmx>=1) then 35 call check_tracers(message,q,ok_q) 36 endif 37 38 if ((.not.ok_t).or.(.not.ok_w).or.(.not.ok_ps).or.(.not.ok_q)) then 36 39 ! Something is wrong, might as well stop here 37 40 call abort_physic(trim(name),trim(message)//" Invalid field values",1) … … 41 44 42 45 43 subroutine check_temperature(temp,ok,temp_min,temp_max) 44 use dimphy, only: klon, klev 45 implicit none 46 subroutine check_temperature(message,temp,ok,temp_min,temp_max) 47 use dimphy, only: klon, klev 48 implicit none 49 character(len=*),intent(in):: message ! text message for outputs 46 50 real,intent(in) :: temp(klon,klev) 47 51 logical,intent(out) :: ok ! returns .true. if everything OK, .false. otherwise … … 75 79 if (temp(i,k).ne.temp(i,k)) then 76 80 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 78 83 endif 79 84 ! check for temperatures too low 80 85 if (temp(i,k).lt.tmin) then 81 86 ok=.false. 82 write(*,*)trim( name)//" temp(i,k)=",temp(i,k)," for i=",i," k=",k,&83 " <",tmin87 write(*,*)trim(message)//" "//trim(name)//" temp(i,k)=",temp(i,k),& 88 " for i=",i," k=",k," <",tmin 84 89 endif 85 90 ! check for temperatures too high 86 91 if (temp(i,k).gt.tmax) then 87 92 ok=.false. 88 write(*,*)trim( name)//" temp(i,k)=",temp(i,k)," for i=",i," k=",k,&89 " >",tmax93 write(*,*)trim(message)//" "//trim(name)//" temp(i,k)=",temp(i,k),& 94 " for i=",i," k=",k," >",tmax 90 95 endif 91 96 enddo … … 94 99 end subroutine check_temperature 95 100 96 subroutine check_winds(u,v,ok,wind_max) 97 use dimphy, only: klon, klev 98 implicit none 101 subroutine check_winds(message,u,v,ok,wind_max) 102 use dimphy, only: klon, klev 103 implicit none 104 character(len=*),intent(in):: message ! text message for outputs 99 105 real,intent(in) :: u(klon,klev) ! zonal wind (m/s) 100 106 real,intent(in) :: v(klon,klev) ! meridional wind (m/s) … … 123 129 if (u(i,k).ne.u(i,k)) then 124 130 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 126 133 endif 127 134 if (v(i,k).ne.v(i,k)) then 128 135 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 130 138 endif 131 139 ! check for magnitudes too high 132 140 if (abs(u(i,k)).gt.wmax) then 133 141 ok=.false. 134 write(*,*)trim( name)//" u(i,k)=",u(i,k)," for i=",i," k=",k,&135 " >",wmax142 write(*,*)trim(message)//" "//trim(name)//" u(i,k)=",u(i,k),& 143 " for i=",i," k=",k," >",wmax 136 144 endif 137 145 if (abs(v(i,k)).gt.wmax) then 138 146 ok=.false. 139 write(*,*)trim( name)//" v(i,k)=",v(i,k)," for i=",i," k=",k,&140 " >",wmax147 write(*,*)trim(message)//" "//trim(name)//" v(i,k)=",v(i,k),& 148 " for i=",i," k=",k," >",wmax 141 149 endif 142 150 enddo … … 145 153 end subroutine check_winds 146 154 147 subroutine check_ps( ps,ok,ps_min,ps_max)155 subroutine check_ps(message,ps,ok,ps_min,ps_max) 148 156 use dimphy, only: klon 149 157 implicit none 158 character(len=*),intent(in):: message ! text message for outputs 150 159 real,intent(in) :: ps(klon) 151 160 logical,intent(out) :: ok ! returns .true. if everything OK, .false. otherwise … … 178 187 if (ps(i).ne.ps(i)) then 179 188 ok=.false. 180 write(*,*)trim( name)//" ps(i)=",ps(i)," for i=",i189 write(*,*)trim(message)//" "//trim(name)//" ps(i)=",ps(i)," for i=",i 181 190 endif 182 191 ! check for pressures too low 183 192 if (ps(i).lt.pmin) then 184 193 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 186 196 endif 187 197 ! check for pressures too high 188 198 if (ps(i).gt.pmax) then 189 199 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 191 202 endif 192 203 enddo … … 194 205 end subroutine check_ps 195 206 207 subroutine check_tracers(message,q,ok) 208 use dimphy, only: klon, klev 209 use tracer_mod, only: nqmx,noms 210 implicit none 211 character(len=*),intent(in):: message ! text message for outputs 212 real,intent(in) :: q(klon,klev,nqmx) ! tracer mixing ratio 213 logical,intent(inout) :: ok ! set to .false. if something is wrong 214 215 character(len=50) :: name="check_tracers" 216 integer :: i,k,iq 217 integer :: nb_bad_neg ! number of problematic points (negative values) 218 integer :: nb_bad_nan ! number of problematic points (NaNs) 219 220 ! 1. initializations 221 nb_bad_nan=0 222 223 ! 2. do the checking 224 do 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 244 enddo !of do iq=1,nqmx 245 246 end subroutine check_tracers 247 196 248 end module check_fields_mod -
trunk/LMDZ.MARS/libf/phymars/physiq_mod.F
r2566 r2570 743 743 if (check_physics_inputs) then 744 744 ! 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) 746 746 endif 747 747 … … 3978 3978 if (check_physics_outputs) then 3979 3979 ! 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) 3981 3981 endif 3982 3982
Note: See TracChangeset
for help on using the changeset viewer.