chgres_cube 1.14.0
Loading...
Searching...
No Matches
utils.F90
Go to the documentation of this file.
1module utilities
2
3contains
7
12 subroutine error_handler(string, rc)
13
14 use mpi_f08
15
16 implicit none
17
18 character(len=*), intent(in) :: string
19
20 integer, intent(in) :: rc
21
22 integer :: ierr
23
24 print*,"- FATAL ERROR: ", trim(string)
25 print*,"- IOSTAT IS: ", rc
26 call mpi_abort(mpi_comm_world, 999, ierr)
27
28 end subroutine error_handler
29
34 subroutine netcdf_err( err, string )
35
36 use mpi_f08
37 use netcdf
38
39 implicit none
40 integer, intent(in) :: err
41 character(len=*), intent(in) :: string
42 character(len=256) :: errmsg
43 integer :: iret
44
45 if( err.EQ.nf90_noerr )return
46 errmsg = nf90_strerror(err)
47 print*,''
48 print*,'FATAL ERROR: ', trim(string), ': ', trim(errmsg)
49 print*,'STOP.'
50 call mpi_abort(mpi_comm_world, 999, iret)
51
52 return
53 end subroutine netcdf_err
54
62function to_upper(strIn) result(strOut)
63
64 implicit none
65
66 character(len=*), intent(in) :: strin
67 character(len=len(strIn)) :: strout
68 integer :: i,j
69
70 do i = 1, len(strin)
71 j = iachar(strin(i:i))
72 if (j>= iachar("a") .and. j<=iachar("z") ) then
73 strout(i:i) = achar(iachar(strin(i:i))-32)
74 else
75 strout(i:i) = strin(i:i)
76 end if
77 end do
78
79end function to_upper
80
87subroutine to_lower(strIn)
88
89 implicit none
90
91 character(len=*), intent(inout) :: strIn
92 character(len=len(strIn)) :: strOut
93 integer :: i,j
94
95 do i = 1, len(strin)
96 j = iachar(strin(i:i))
97 if (j>= iachar("A") .and. j<=iachar("Z") ) then
98 strout(i:i) = achar(iachar(strin(i:i))+32)
99 else
100 strout(i:i) = strin(i:i)
101 end if
102 end do
103 strin(:) = strout(:)
104end subroutine to_lower
105
120subroutine handle_grib_error(vname,lev,method,value,varnum,read_from_input, iret,var,var8,var3d)
121
122 use, intrinsic :: ieee_arithmetic
123 use esmf
124
125 implicit none
126
127 real(esmf_kind_r4), intent(in) :: value
128 logical, intent(inout) :: read_from_input(:)
129 real(esmf_kind_r4), intent(inout), optional :: var(:,:)
130 real(esmf_kind_r8), intent(inout), optional :: var8(:,:)
131 real(esmf_kind_r8), intent(inout), optional :: var3d(:,:,:)
132
133 character(len=20), intent(in) :: vname, lev, method
134 character(len=200) :: err_msg
135
136 integer, intent(in) :: varnum
137 integer, intent(inout) :: iret
138
139 iret = 0
140 if (varnum == 9999) then
141 print*, "WARNING: ", trim(vname), " NOT FOUND AT LEVEL ", lev, " IN EXTERNAL FILE ", &
142 "AND NO ENTRY EXISTS IN VARMAP TABLE. VARIABLE WILL NOT BE USED."
143 iret = 1
144
145 return
146 endif
147
148 if (trim(method) == "skip" ) then
149 print*, "WARNING: SKIPPING ", trim(vname), " IN FILE"
150 read_from_input(varnum) = .false.
151 iret = 1
152 elseif (trim(method) == "set_to_fill") then
153 print*, "WARNING: ,", trim(vname), " NOT AVAILABLE AT LEVEL ", trim(lev), &
154 ". SETTING EQUAL TO FILL VALUE OF ", value
155 if(present(var)) var(:,:) = value
156 if(present(var8)) var8(:,:) = value
157 if(present(var3d)) var3d(:,:,:) = value
158 elseif (trim(method) == "set_to_NaN") then
159 print*, "WARNING: ,", trim(vname), " NOT AVAILABLE AT LEVEL ", trim(lev), &
160 ". SETTING EQUAL TO NaNs"
161 if(present(var)) var(:,:) = ieee_value(var,ieee_quiet_nan)
162 if(present(var8)) var8(:,:) = ieee_value(var8,ieee_quiet_nan)
163 if(present(var3d)) var3d(:,:,:) = ieee_value(var3d,ieee_quiet_nan)
164 elseif (trim(method) == "stop") then
165 err_msg="READING " // trim(vname) // " at level " //lev// ". TO MAKE THIS NON-" // &
166 "FATAL, CHANGE STOP TO SKIP FOR THIS VARIABLE IN YOUR VARMAP FILE."
167 call error_handler(err_msg, iret)
168 elseif (trim(method) == "intrp") then
169 print*, "WARNING: ,"//trim(vname)//" NOT AVAILABLE AT LEVEL "//trim(lev)// &
170 ". WILL INTERPOLATE INTERSPERSED MISSING LEVELS AND/OR FILL MISSING"//&
171 " LEVELS AT EDGES."
172 else
173 err_msg="ERROR USING MISSING_VAR_METHOD. PLEASE SET VALUES IN" // &
174 " VARMAP TABLE TO ONE OF: set_to_fill, set_to_NaN, intrp, skip, or stop."
175 call error_handler(err_msg, 1)
176 endif
177
178end subroutine handle_grib_error
179
186recursive subroutine quicksort(a, first, last)
187 implicit none
188 real*8 a(*), x, t
189 integer first, last
190 integer i, j
191
192 x = a( (first+last) / 2 )
193 i = first
194 j = last
195 do
196 do while (a(i) < x)
197 i=i+1
198 end do
199 do while (x < a(j))
200 j=j-1
201 end do
202 if (i >= j) exit
203 t = a(i); a(i) = a(j); a(j) = t
204 i=i+1
205 j=j-1
206 end do
207 if (first < i-1) call quicksort(a, first, i-1)
208 if (j+1 < last) call quicksort(a, j+1, last)
209end subroutine quicksort
210
228
229subroutine check_soilt(soilt, landmask, skint,ICET_DEFAULT,i_input,j_input,lsoil_input)
230 use esmf
231 implicit none
232 integer, intent(in) :: i_input, j_input, lsoil_input
233 real(esmf_kind_r8), intent(inout) :: soilt(i_input,j_input,lsoil_input)
234 real(esmf_kind_r8), intent(in) :: skint(i_input,j_input)
235 real, intent(in) :: ICET_DEFAULT
236 integer(esmf_kind_i4), intent(in) :: landmask(i_input,j_input)
237
238 integer :: i, j, k
239
240 do k=1,lsoil_input
241 do j = 1, j_input
242 do i = 1, i_input
243 if (landmask(i,j) == 0_esmf_kind_i4 ) then
244 soilt(i,j,k) = skint(i,j)
245 else if (landmask(i,j) == 1_esmf_kind_i4 .and. soilt(i,j,k) > 350.0_esmf_kind_r8) then
246 soilt(i,j,k) = skint(i,j)
247 else if (landmask(i,j) == 2_esmf_kind_i4 ) then
248 soilt(i,j,k) = icet_default
249 endif
250 enddo
251 enddo
252 enddo
253end subroutine check_soilt
254
263
264subroutine check_cnwat(cnwat,i_input,j_input)
265 use esmf
266 implicit none
267 integer, intent(in) :: i_input, j_input
268 real(esmf_kind_r8), intent(inout) :: cnwat(i_input,j_input)
269
270 real(esmf_kind_r8) :: max_cnwat = 0.5
271
272 integer :: i, j
273
274 do i = 1,i_input
275 do j = 1,j_input
276 if (cnwat(i,j) .gt. max_cnwat) cnwat(i,j) = 0.0_esmf_kind_r8
277 enddo
278 enddo
279end subroutine check_cnwat
280
299SUBROUTINE dint2p(PPIN,XXIN,NPIN,PPOUT,XXOUT,NPOUT &
300 ,LINLOG,XMSG,IER)
301 IMPLICIT NONE
302
303! NCL code for pressure level interpolation
304!
305! This code was designed for one simple task. It has since
306! been mangled and abused for assorted reasons. For example,
307! early gfortran compilers had some issues with automatic arrays.
308! Hence, the C-Wrapper was used to create 'work' arrays which
309! were then passed to this code. The original focused (non-NCL)
310! task was to handle PPIN & PPOUT that had the same 'monotonicity.'
311! Extra code was added to handle the more general case.
312! Blah-Blah: Punch line: it is embarrassingly convoluted!!!
313!
314! ! input types
315 INTEGER NPIN,NPOUT,LINLOG,IER
316 real*8 ppin(npin),xxin(npin),ppout(npout),xmsg
317 ! output
318 real*8 xxout(npout)
319 ! work
320 real*8 pin(npin),xin(npin),p(npin),x(npin)
321 real*8 pout(npout),xout(npout)
322
323! local
324 INTEGER NP,NL,NLMAX,NLSAVE,NP1,NO1,N1,N2,LOGLIN, &
325 NLSTRT
326 real*8 slope,pa,pb,pc
327
328 loglin = abs(linlog)
329
330! error check: enough points: pressures consistency?
331
332 ier = 0
333 IF (npout.GT.0) THEN
334 DO np = 1,npout
335 xxout(np) = xmsg
336 END DO
337 END IF
338! Jili Dong input levels have to be the same as output levels:
339! we only interpolate for levels with missing variables
340! IF (.not. all(PPIN .eq. PPOUT)) IER = IER+1
341
342 IF (npin.LT.2 .OR. npout.LT.1) ier = ier + 1
343
344 IF (ier.NE.0) THEN
345! PRINT *,'INT2P: error exit: ier=',IER
346 RETURN
347 END IF
348
349! should *input arrays* be reordered? want p(1) > p(2) > p(3) etc
350! so that it will match order for which code was originally designed
351! copy to 'work' arrays
352
353 np1 = 0
354 no1 = 0
355 IF (ppin(1).LT.ppin(2)) THEN
356 np1 = npin + 1
357 END IF
358 IF (ppout(1).LT.ppout(2)) THEN
359 no1 = npout + 1
360 END IF
361
362 DO np = 1,npin
363 pin(np) = ppin(abs(np1-np))
364 xin(np) = xxin(abs(np1-np))
365 END DO
366
367 DO np = 1,npout
368 pout(np) = ppout(abs(no1-np))
369 END DO
370
371! eliminate XIN levels with missing data.
372! . This can happen with observational data.
373
374 nl = 0
375 DO np = 1,npin
376 IF (xin(np).NE.xmsg .AND. pin(np).NE.xmsg) THEN
377 nl = nl + 1
378 p(nl) = pin(np)
379 x(nl) = xin(np)
380 END IF
381 END DO
382 nlmax = nl
383
384 ! all missing data
385 IF (nlmax.LT.2) THEN
386 ier = ier + 1000
387 print *,'INT2P: ier=',ier
388 RETURN
389 END IF
390
391! ===============> pressure in decreasing order <================
392! perform the interpolation [pin(1)>pin(2)>...>pin(npin)]
393! ( p ,x)
394! ------------------------- p(nl+1), x(nl+1) example (200,5)
395! .
396! ------------------------- pout(np), xout(np) (250,?)
397! .
398! ------------------------- p(nl) , x(nl) (300,10)
399
400
401! exact p-level matches
402 nlstrt = 1
403 nlsave = 1
404 DO np = 1,npout
405 xout(np) = xmsg
406 DO nl = nlstrt,nlmax
407 IF (pout(np).EQ.p(nl)) THEN
408 xout(np) = x(nl)
409 nlsave = nl + 1
410 GO TO 10
411 END IF
412 END DO
413 10 nlstrt = nlsave
414 END DO
415
416 IF (loglin.EQ.1) THEN
417 DO np = 1,npout
418 DO nl = 1,nlmax - 1
419 IF (pout(np).LT.p(nl) .AND. pout(np).GT.p(nl+1)) THEN
420 slope = (x(nl)-x(nl+1))/ (p(nl)-p(nl+1))
421 xout(np) = x(nl+1) + slope* (pout(np)-p(nl+1))
422 END IF
423 END DO
424 END DO
425 ELSE
426 DO np = 1,npout
427 DO nl = 1,nlmax - 1
428 IF (pout(np).LT.p(nl) .AND. pout(np).GT.p(nl+1)) THEN
429 pa = log(p(nl))
430 pb = log(pout(np))
431! special case: In case someome inadvertently enter p=0.
432 if (p(nl+1).gt.0.d0) then
433 pc = log(p(nl+1))
434 else
435 pc = log(1.e-4)
436 end if
437
438 slope = (x(nl)-x(nl+1))/ (pa-pc)
439 xout(np) = x(nl+1) + slope* (pb-pc)
440 END IF
441 END DO
442 END DO
443 END IF
444
445! extrapolate?
446! . use the 'last' valid slope for extrapolating
447
448 IF (linlog.LT.0) THEN
449 DO np = 1,npout
450 DO nl = 1,nlmax
451 IF (pout(np).GT.p(1)) THEN
452 IF (loglin.EQ.1) THEN
453 slope = (x(2)-x(1))/ (p(2)-p(1))
454 xout(np) = x(1) + slope* (pout(np)-p(1))
455 ELSE
456 pa = log(p(2))
457 pb = log(pout(np))
458 pc = log(p(1))
459 slope = (x(2)-x(1))/ (pa-pc)
460 xout(np) = x(1) + slope* (pb-pc)
461 END IF
462 ELSE IF (pout(np).LT.p(nlmax)) THEN
463 n1 = nlmax
464 n2 = nlmax - 1
465 IF (loglin.EQ.1) THEN
466 slope = (x(n1)-x(n2))/ (p(n1)-p(n2))
467 xout(np) = x(n1) + slope* (pout(np)-p(n1))
468 ELSE
469 pa = log(p(n1))
470 pb = log(pout(np))
471 pc = log(p(n2))
472 slope = (x(n1)-x(n2))/ (pa-pc)
473 !XOUT(NP) = X(N1) + SLOPE* (PB-PC) !bug fixed below
474 xout(np) = x(n1) + slope* (pb-pa)
475 END IF
476 END IF
477 END DO
478 END DO
479 END IF
480
481! place results in the return array;
482! . possibly .... reverse to original order
483
484 if (no1.GT.0) THEN
485 DO np = 1,npout
486 n1 = abs(no1-np)
487 ppout(np) = pout(n1)
488 xxout(np) = xout(n1)
489 END DO
490 ELSE
491 DO np = 1,npout
492 ppout(np) = pout(np)
493 xxout(np) = xout(np)
494 END DO
495 END IF
496
497
498 RETURN
499 END SUBROUTINE dint2p
500end module utilities