My Project
Functions/Subroutines | Variables
mod_boundschk Module Reference

Functions/Subroutines

subroutine setup_boundschk
 
subroutine boundschk
 

Variables

type(infovar), dimension(:), allocatable vlist
 

Function/Subroutine Documentation

◆ boundschk()

subroutine mod_boundschk::boundschk ( )

Definition at line 146 of file mod_boundschk.f90.

146  use all_vars
147  use mod_utils, only : pstop
148  use mod_wd
149  use mod_ncdio, only : archive
150  use mod_par, only : egid,ngid
151  implicit none
152  type(infovar), allocatable :: vlist_global(:)
153  integer :: icnt,pt,lay,i,j,k,nviolations,printproc,ierr
154  integer :: iiint
155  real(sp) :: val
156  integer, allocatable :: violations(:),tmp(:)
157  character(len=20) :: vname
158 
159 
160  iiint=iint
161  if(.not.boundschk_on) return
162  if(mod(iiint,chk_interval)/= 0) return
163 ! if(.not.boundschk_on .or. (mod(iint,chk_interval)/= 0)) return
164  allocate(tmp(nprocs))
165  tmp = 0
166 
167  !------------------------------------------------------
168  !Check if variables are in user-defined bounds
169  !------------------------------------------------------
170  icnt = 0
171 
172  ! vert-averaged x-velocity
173  val = maxval(abs(ua(1:n)))
174  if(val > veloc_mag_max)then
175  icnt = icnt + 1
176  pt = 1
177  vname = 'vert-averaged u'
178  ualoop: do i=1,n
179  if(abs(ua(i)) > veloc_mag_max)then
180  pt = i
181  val = abs(ua(i))
182  exit ualoop
183  endif
184  end do ualoop
185  call set_infovar(vlist(icnt),vname,2,val,veloc_mag_max,xc(pt)+vxmin,yc(pt)+vymin,myid,pt,1)
186  endif
187 
188  ! vert-averaged y-velocity
189  val = maxval(abs(va(1:n)))
190  if(val > veloc_mag_max)then
191  icnt = icnt + 1
192  pt = 1
193  vname = 'vert-averaged v'
194  valoop: do i=1,n
195  if(abs(va(i)) > veloc_mag_max)then
196  pt = i
197  val = abs(va(i))
198  exit valoop
199  endif
200  end do valoop
201  call set_infovar(vlist(icnt),vname,2,val,veloc_mag_max,xc(pt)+vxmin,yc(pt)+vymin,myid,pt,1)
202  endif
203 
204  ! x-velocity
205  val = maxval(abs(u(1:n,1:kbm1)))
206  if(val > veloc_mag_max)then
207  icnt = icnt + 1
208  pt = 1
209  vname = 'u'
210  uloop: do k=1,kbm1
211  do i=1,n
212  if(abs(u(i,k)) > veloc_mag_max)then
213  pt = i
214  lay = k
215  val = abs(u(i,k))
216  exit uloop
217  endif
218  end do
219  end do uloop
220  call set_infovar(vlist(icnt),vname,2,val,veloc_mag_max,xc(pt)+vxmin,yc(pt)+vymin,myid,pt,lay)
221  endif
222 
223  !hydrographic vars
224  if(.not. barotropic)then
225 
226  ! salinity - max
227  val = maxval(s1(1:m,1:kbm1))
228  if(val > salt_max)then
229  icnt = icnt + 1
230  pt = 1
231  vname = 'salinity'
232  smaxloop: do k=1,kbm1
233  do i=1,m
234  if(s1(i,k) > salt_max)then
235  pt = i
236  lay = k
237  val = s1(i,k)
238  exit smaxloop
239  endif
240  end do
241  end do smaxloop
242  call set_infovar(vlist(icnt),vname,1,val,salt_max,vx(pt)+vxmin,vy(pt)+vymin,myid,pt,lay)
243  endif
244 
245  ! salinity - min
246  val = minval(s1(1:m,1:kbm1))
247  if(val < salt_min)then
248  icnt = icnt + 1
249  pt = 1
250  vname = 'salinity'
251  sminloop: do k=1,kbm1
252  do i=1,m
253  if(s1(i,k) < salt_min)then
254  pt = i
255  lay = k
256  val = s1(i,k)
257  exit sminloop
258  endif
259  end do
260  end do sminloop
261  call set_infovar(vlist(icnt),vname,1,val,salt_min,vx(pt)+vxmin,vy(pt)+vymin,myid,pt,lay)
262  endif
263 
264  ! temp - max
265  val = maxval(t1(1:m,1:kbm1))
266  if(val > temp_max)then
267  icnt = icnt + 1
268  pt = 1
269  vname = 'temperature'
270  tmaxloop: do k=1,kbm1
271  do i=1,m
272  if(t1(i,k) > temp_max)then
273  pt = i
274  lay = k
275  val = t1(i,k)
276  exit tmaxloop
277  endif
278  end do
279  end do tmaxloop
280  call set_infovar(vlist(icnt),vname,1,val,temp_max,vx(pt)+vxmin,vy(pt)+vymin,myid,pt,lay)
281  endif
282 
283  ! temperature - min
284  val = minval(t1(1:m,1:kbm1))
285  if(val < temp_min)then
286  icnt = icnt + 1
287  pt = 1
288  vname = 'temperature'
289  tminloop: do k=1,kbm1
290  do i=1,m
291  if(t1(i,k) < temp_min)then
292  pt = i
293  lay = k
294  val = t1(i,k)
295  exit tminloop
296  endif
297  end do
298  end do tminloop
299  call set_infovar(vlist(icnt),vname,1,val,temp_min,vx(pt)+vxmin,vy(pt)+vymin,myid,pt,lay)
300  endif
301 
302 
303  endif !.not. barotropic
304 
305 
306 
307 
308 
309  !-----------------------------------------------------------
310  !Collect number of violations from each proc into vector tmp
311  !-----------------------------------------------------------
312  nviolations = icnt
313  tmp(1) = icnt
314 
315 
316  !check total violations and return if none found
317  if(sum(tmp) == 0)then
318  deallocate(tmp)
319  return
320  endif
321 
322  !look at violations, just pick lowest procid with violation to print violation
323  printproc = 1
324  do i=1,nprocs
325  if(tmp(i) > 0)then
326  printproc = i
327  exit
328  endif
329  end do
330 
331  !write violations to screen
332  if(printproc==myid)then
333  write(ipt,*)'WARNING: Variable(s) have exceeded user-defined thresholds'
334  do i=1,nviolations
335  call print_infovar(vlist(i),ipt)
336  end do
337  write(ipt,*)'ARCHIVING FRAME AND HALTING'
338  endif
339 
340  !dump frame to netcdf file
341  force_archive = .true.
342  call archive
343 
344  !shutdown
345  call pstop
346 
347  deallocate(tmp)
348 
real(sp), dimension(:), allocatable, target va
Definition: mod_main.f90:1104
subroutine archive
Definition: mod_ncdio.f90:153
logical barotropic
Definition: mod_main.f90:381
logical force_archive
Definition: mod_main.f90:805
integer myid
Definition: mod_main.f90:67
integer chk_interval
Definition: mod_main.f90:807
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target t1
Definition: mod_main.f90:1307
real(sp) temp_min
Definition: mod_main.f90:811
integer, target nprocs
Definition: mod_main.f90:72
integer m
Definition: mod_main.f90:56
real(sp) vymin
Definition: mod_main.f90:989
real(sp) salt_min
= bounds checking
Definition: mod_main.f90:813
real(sp), dimension(:,:), allocatable, target u
Definition: mod_main.f90:1268
real(sp), dimension(:,:), allocatable, target s1
Definition: mod_main.f90:1308
integer(itime) iint
Definition: mod_main.f90:850
real(sp), dimension(:), allocatable, target vx
Definition: mod_main.f90:1001
subroutine pstop
Definition: mod_utils.f90:273
integer n
Definition: mod_main.f90:55
real(sp), dimension(:), allocatable, target vy
Definition: mod_main.f90:1002
real(sp) temp_max
Definition: mod_main.f90:810
real(sp) veloc_mag_max
Definition: mod_main.f90:808
real(sp), dimension(:), allocatable, target ua
Definition: mod_main.f90:1103
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
logical boundschk_on
Definition: mod_main.f90:806
integer, dimension(:), pointer ngid
Definition: mod_par.f90:61
integer ipt
Definition: mod_main.f90:922
integer kbm1
Definition: mod_main.f90:65
integer, dimension(:), pointer egid
Definition: mod_par.f90:60
real(sp) salt_max
Definition: mod_main.f90:812
real(sp) vxmin
Definition: mod_main.f90:989
Here is the call graph for this function:
Here is the caller graph for this function:

◆ setup_boundschk()

subroutine mod_boundschk::setup_boundschk ( )

Definition at line 118 of file mod_boundschk.f90.

118  use all_vars
119  integer :: nvars
120 
121  !report the bounds
122  if(boundschk_on)then
123  write(ipt,*)'bounds checking : on'
124  write(ipt,*)'checking interval: ',chk_interval
125  write(ipt,*)'veloc_mag_max : ',veloc_mag_max
126  write(ipt,*)'zeta_mag_max : ',zeta_mag_max
127  write(ipt,*)'temp_max : ',temp_max
128  write(ipt,*)'temp_min : ',temp_min
129  write(ipt,*)'salt_max : ',salt_max
130  write(ipt,*)'salt_min : ',salt_min
131 
132  !allocate vlist
133  nvars = 7
134  allocate(vlist(nvars))
135 
136  else
137  allocate(vlist(1))
138  endif
139 
real(sp) zeta_mag_max
Definition: mod_main.f90:809
integer chk_interval
Definition: mod_main.f90:807
real(sp) temp_min
Definition: mod_main.f90:811
real(sp) salt_min
= bounds checking
Definition: mod_main.f90:813
real(sp) temp_max
Definition: mod_main.f90:810
real(sp) veloc_mag_max
Definition: mod_main.f90:808
logical boundschk_on
Definition: mod_main.f90:806
integer ipt
Definition: mod_main.f90:922
real(sp) salt_max
Definition: mod_main.f90:812
Here is the caller graph for this function:

Variable Documentation

◆ vlist

type(infovar), dimension(:), allocatable mod_boundschk::vlist

Definition at line 112 of file mod_boundschk.f90.

112  type(infovar), allocatable :: vlist(:)