FV3 Bundle
topography.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 
21 
22 ! <CONTACT EMAIL="Bruce.Wyman@noaa.gov">
23 ! Bruce Wyman
24 ! </CONTACT>
25 
26 ! <HISTORY SRC="http://www.gfdl.noaa.gov/fms-cgi-bin/cvsweb.cgi/FMS/"/>
27 
28 ! <OVERVIEW>
29 ! Routines for creating land surface topography fields and land-water masks
30 ! for latitude-longitude grids.
31 ! </OVERVIEW>
32 
33 ! <DESCRIPTION>
34 ! This module generates realistic mountains and land-water masks
35 ! on a specified latitude-longitude grid by interpolating from the
36 ! 1/6 degree Navy mean topography and percent water data sets.
37 ! The fields that can be generated are mean and standard deviation
38 ! of topography within the specified grid boxes; and land-ocean (or
39 ! water) mask and land-ocean (or water) fractional area.
40 !
41 ! The interpolation scheme conserves the area-weighted average
42 ! of the input data by using module horiz_interp.
43 !
44 ! The interfaces get_gaussian_topog and gaussian_topog_init are documented in <LINK SRC="gaussian_topog.html">gaussian_topog_mod</LINK>.
45 ! </DESCRIPTION>
46 
48 use horiz_interp_mod, only: horiz_interp_type, horiz_interp_new, &
50 
51 use fms_mod, only: file_exist, check_nml_error, &
52  open_namelist_file, close_file, stdlog, &
53  mpp_pe, mpp_root_pe, write_version_number, &
54  open_ieee32_file, error_mesg, fatal, note, &
55  mpp_error
56 use fms_io_mod, only: read_data
57 use constants_mod, only: pi
58 use mpp_mod, only: input_nml_file
59 
60 implicit none
61 private
62 
63 public :: topography_init, &
68 
69 interface get_topog_mean
70  module procedure get_topog_mean_1d, get_topog_mean_2d
71 end interface
72 interface get_topog_stdev
73  module procedure get_topog_stdev_1d, get_topog_stdev_2d
74 end interface
75 interface get_ocean_frac
76  module procedure get_ocean_frac_1d, get_ocean_frac_2d
77 end interface
78 interface get_ocean_mask
79  module procedure get_ocean_mask_1d, get_ocean_mask_2d
80 end interface
81 interface get_water_frac
82  module procedure get_water_frac_1d, get_water_frac_2d
83 end interface
84 interface get_water_mask
85  module procedure get_water_mask_1d, get_water_mask_2d
86 end interface
87 
88 !-----------------------------------------------------------------------
89 ! <NAMELIST NAME="topography_nml">
90 ! <DATA NAME="topog_file" TYPE="character" DEFAULT="DATA/navy_topography.data">
91 ! Name of topography file.
92 ! </DATA>
93 ! <DATA NAME="water_file" TYPE="character" DEFAULT="DATA/navy_pctwater.data">
94 ! Name of percent water file.
95 ! </DATA>
96 
97  character(len=128) :: topog_file = 'DATA/navy_topography.data', &
98  water_file = 'DATA/navy_pctwater.data'
99  namelist /topography_nml/ topog_file, water_file
100 ! </NAMELIST>
101 
102 !-----------------------------------------------------------------------
103 ! --- resolution of the topography data set ---
104 ! <DATASET NAME="">
105 ! This module uses the 1/6 degree U.S. Navy mean topography
106 ! and percent water data sets.
107 !
108 ! These data sets have been re-formatted to separate 32-bit IEEE files.
109 ! The names of these files is specified by the <LINK SRC="#NAMELIST">namelist</LINK> input.
110 !
111 !The format for both files is as follows:
112 ! <PRE>
113 ! record = 1 nlon, nlat
114 ! record = 2 blon, blat
115 ! record = 3 data
116 ! </PRE>
117 !where:
118 ! <PRE>
119 ! nlon, nlat = The number of longitude and latitude points
120 ! in the horizontal grid. For the 1/6 degree
121 ! data sets this is 2160 x 1080. [integer]
122 ! blon, blat = The longitude and latitude grid box boundaries in degrees.
123 ! [real :: blon(nlon+1), blat(nlat+1)]
124 !
125 ! data = The topography or percent water data.
126 ! [real :: data(nlon,nlat)]
127 ! </PRE>
128 ! </DATASET>
129  integer :: unit
130  integer :: ipts, jpts
131  integer, parameter :: compute_stdev = 123 ! use this flag to
132  ! compute st dev
133 
134 !-----------------------------------------------------------------------
135 
136 ! Include variable "version" to be written to log file.
137 #include<file_version.h>
138 
139  logical :: module_is_initialized = .false.
140 
141 !-----------------------------------------------------------------------
142 
143  contains
144 
145 !#######################################################################
146 
147  subroutine topography_init ()
149  if ( module_is_initialized ) return
150 
151  call write_version_number("TOPOGRAPHY_MOD", version)
152  call read_namelist
153  module_is_initialized = .true.
154 
155  end subroutine topography_init
156 
157 !#######################################################################
158 
159 ! <FUNCTION NAME="get_topog_mean">
160 
161 ! <OVERVIEW>
162 ! Returns a "realistic" mean surface height field.
163 ! </OVERVIEW>
164 ! <DESCRIPTION>
165 ! Returns realistic mountains on a latitude-longtude grid.
166 ! The returned field is the mean topography for the given grid boxes.
167 ! Computed using a conserving area-weighted interpolation.
168 ! The current input data set is the 1/6 degree Navy mean topography.
169 ! </DESCRIPTION>
170 ! <TEMPLATE>
171 ! flag = <B>get_topog_mean</B> ( blon, blat, zmean )
172 ! </TEMPLATE>
173 
174 ! <IN NAME="blon" TYPE="real" DIM="(:)">
175 ! The longitude (in radians) at grid box boundaries.
176 ! </IN>
177 ! <IN NAME="blat" TYPE="real" DIM="(:)">
178 ! The latitude (in radians) at grid box boundaries.
179 ! </IN>
180 ! <OUT NAME="zmean" UNIT=" meter" TYPE="real" DIM="(:,:)">
181 ! The mean surface height (meters).
182 ! The size of this field must be size(blon)-1 by size(blat)-1.
183 ! </OUT>
184 ! <OUT NAME="get_topog_mean" TYPE="logical">
185 ! A logical value of TRUE is returned if the surface height field
186 ! was successfully created. A value of FALSE may be returned if the
187 ! input topography data set was not readable.
188 ! </OUT>
189 
190 ! <ERROR MSG="shape(zmean) is not equal to (/size(blon)-1,size(blat)-1/))" STATUS="FATAL">
191 ! Check the input grid size and output field size.
192 ! </ERROR>
193 
194  function get_topog_mean_1d (blon, blat, zmean)
196  real, intent(in), dimension(:) :: blon, blat
197  real, intent(out), dimension(:,:) :: zmean
198  logical :: get_topog_mean_1d
199 
200 !-----------------------------------------------------------------------
201  if (.not. module_is_initialized) call topography_init()
202 
203  if ( any(shape(zmean(:,:)) /= (/size(blon(:))-1,size(blat(:))-1/)) ) &
204  call error_mesg('get_topog_mean_1d','shape(zmean) is not&
205  & equal to (/size(blon)-1,size(blat)-1/))', fatal)
206 
207  if ( open_topog_file(topog_file) ) then
208  call interp_topog_1d ( blon, blat, zmean )
209  get_topog_mean_1d = .true.
210  else
211  get_topog_mean_1d = .false.
212  endif
213 
214 !-----------------------------------------------------------------------
215 
216  end function get_topog_mean_1d
217 
218 !############################################################
219 
220  function get_topog_mean_2d (blon, blat, zmean)
222  real, intent(in), dimension(:,:) :: blon, blat
223  real, intent(out), dimension(:,:) :: zmean
224  logical :: get_topog_mean_2d
225 
226 !-----------------------------------------------------------------------
227  if (.not. module_is_initialized) call topography_init()
228 
229  if ( any(shape(zmean(:,:)) /= (/size(blon,1)-1,size(blon,2)-1/)) .or. &
230  any(shape(zmean(:,:)) /= (/size(blat,1)-1,size(blat,2)-1/)) ) &
231  call error_mesg('get_topog_mean_2d','shape(zmean) is not&
232  & equal to (/size(blon,1)-1,size(blon,2)-1/))', fatal)
233 
234  if ( open_topog_file(topog_file) ) then
235  call interp_topog_2d ( blon, blat, zmean )
236  get_topog_mean_2d = .true.
237  else
238  get_topog_mean_2d = .false.
239  endif
240 
241 !-----------------------------------------------------------------------
242 
243  end function get_topog_mean_2d
244 
245 ! </FUNCTION>
246 !#######################################################################
247 
248 ! <FUNCTION NAME="get_topog_stdev">
249 
250 ! <OVERVIEW>
251 ! Returns a standard deviation of higher resolution topography with
252 ! the given model grid boxes.
253 ! </OVERVIEW>
254 ! <DESCRIPTION>
255 ! Returns the standard deviation of the "finer" input topography data set,
256 ! currently the Navy 1/6 degree mean topography data, within the
257 ! boundaries of the given input grid.
258 ! </DESCRIPTION>
259 ! <TEMPLATE>
260 ! flag = <B>get_topog_stdev</B> ( blon, blat, stdev )
261 ! </TEMPLATE>
262 ! <IN NAME="blon" TYPE="real" DIM="(:)">
263 ! The longitude (in radians) at grid box boundaries.
264 ! </IN>
265 ! <IN NAME="blat" TYPE="real" DIM="(:)">
266 ! The latitude (in radians) at grid box boundaries.
267 ! </IN>
268 ! <OUT NAME="stdev" UNITS="meter" TYPE="real" DIM="(:,:)">
269 ! The standard deviation of surface height (in meters) within
270 ! given input model grid boxes.
271 ! The size of this field must be size(blon)-1 by size(blat)-1.
272 ! </OUT>
273 ! <OUT NAME="get_topog_stdev" TYPE="logical">
274 ! A logical value of TRUE is returned if the output field was
275 ! successfully created. A value of FALSE may be returned if the
276 ! input topography data set was not readable.
277 ! </OUT>
278 
279  function get_topog_stdev_1d (blon, blat, stdev)
281  real, intent(in), dimension(:) :: blon, blat
282  real, intent(out), dimension(:,:) :: stdev
283  logical :: get_topog_stdev_1d
284 
285 !-----------------------------------------------------------------------
286  if (.not. module_is_initialized) call topography_init()
287 
288  if ( any(shape(stdev(:,:)) /= (/size(blon(:))-1,size(blat(:))-1/)) ) &
289  call error_mesg('get_topog_stdev','shape(stdev) is not&
290  & equal to (/size(blon)-1,size(blat)-1/))', fatal)
291 
292  if ( open_topog_file(topog_file) ) then
293  call interp_topog_1d ( blon, blat, stdev, flag=compute_stdev )
294  get_topog_stdev_1d = .true.
295  else
296  get_topog_stdev_1d = .false.
297  endif
298 
299 !-----------------------------------------------------------------------
300 
301  end function get_topog_stdev_1d
302 
303 !#######################################################################
304 
305  function get_topog_stdev_2d (blon, blat, stdev)
307  real, intent(in), dimension(:,:) :: blon, blat
308  real, intent(out), dimension(:,:) :: stdev
309  logical :: get_topog_stdev_2d
310 
311 !-----------------------------------------------------------------------
312  if (.not. module_is_initialized) call topography_init()
313 
314  if ( any(shape(stdev(:,:)) /= (/size(blon,1)-1,size(blon,2)-1/)) .or. &
315  any(shape(stdev(:,:)) /= (/size(blat,1)-1,size(blat,2)-1/)) ) &
316  call error_mesg('get_topog_stdev_2d','shape(stdev) is not&
317  & equal to (/size(blon,1)-1,size(blon,2)-1/))', fatal)
318 
319  if ( open_topog_file(topog_file) ) then
320  call interp_topog_2d ( blon, blat, stdev, flag=compute_stdev )
321  get_topog_stdev_2d = .true.
322  else
323  get_topog_stdev_2d = .false.
324  endif
325 
326 !-----------------------------------------------------------------------
327 
328  end function get_topog_stdev_2d
329 
330 ! </FUNCTION>
331 !#######################################################################
332 
333 ! <FUNCTION NAME="get_ocean_frac">
334 
335 ! <OVERVIEW>
336 ! Returns fractional area covered by ocean in a grid box.
337 ! </OVERVIEW>
338 ! <DESCRIPTION>
339 ! Returns fractional area covered by ocean in the given model grid boxes.
340 ! </DESCRIPTION>
341 ! <TEMPLATE>
342 ! flag = <B>get_ocean_frac</B> ( blon, blat, ocean_frac )
343 ! </TEMPLATE>
344 
345 ! <IN NAME="blon" UNITS="radians" TYPE="real" DIM="(:)">
346 ! The longitude (in radians) at grid box boundaries.
347 ! </IN>
348 ! <IN NAME="blat" UNITS="radians" TYPE="real" DIM="(:)">
349 ! The latitude (in radians) at grid box boundaries.
350 ! </IN>
351 ! <OUT NAME="ocean_frac" TYPE="real" DIM="(:,:)">
352 ! The fractional amount (0 to 1) of ocean in a grid box.
353 ! The size of this field must be size(blon)-1 by size(blat)-1.
354 ! </OUT>
355 ! <OUT NAME="get_ocean_frac" TYPE="logical">
356 ! A logical value of TRUE is returned if the output field
357 ! was successfully created. A value of FALSE may be returned
358 ! if the Navy 1/6 degree percent water data set was not readable.
359 ! </OUT>
360 
361  function get_ocean_frac_1d (blon, blat, ocean_frac)
363  real, intent(in), dimension(:) :: blon, blat
364  real, intent(out), dimension(:,:) :: ocean_frac
365  logical :: get_ocean_frac_1d
366 
367 !-----------------------------------------------------------------------
368  if (.not. module_is_initialized) call topography_init()
369 
370  if ( any(shape(ocean_frac(:,:)) /= (/size(blon(:))-1,size(blat(:))-1/)) ) &
371  call error_mesg('get_ocean_frac','shape(ocean_frac) is not&
372  & equal to (/size(blon)-1,size(blat)-1/))', fatal)
373 
374  if ( open_topog_file(water_file) ) then
375  call interp_water_1d ( blon, blat, ocean_frac, do_ocean=.true. )
376  get_ocean_frac_1d = .true.
377  else
378  get_ocean_frac_1d = .false.
379  endif
380 
381 !-----------------------------------------------------------------------
382 
383  end function get_ocean_frac_1d
384 
385 !#######################################################################
386 
387  function get_ocean_frac_2d (blon, blat, ocean_frac)
389  real, intent(in), dimension(:,:) :: blon, blat
390  real, intent(out), dimension(:,:) :: ocean_frac
391  logical :: get_ocean_frac_2d
392 
393 !-----------------------------------------------------------------------
394  if (.not. module_is_initialized) call topography_init()
395 
396  if ( any(shape(ocean_frac(:,:)) /= (/size(blon,1)-1,size(blon,2)-1/)) .or. &
397  any(shape(ocean_frac(:,:)) /= (/size(blat,1)-1,size(blat,2)-1/)) ) &
398  call error_mesg('get_ocean_frac_2d','shape(ocean_frac) is not&
399  & equal to (/size(blon,1)-1,size(blon,2)-1/))', fatal)
400 
401  if ( open_topog_file(water_file) ) then
402  call interp_water_2d ( blon, blat, ocean_frac, do_ocean=.true. )
403  get_ocean_frac_2d = .true.
404  else
405  get_ocean_frac_2d = .false.
406  endif
407 
408 !-----------------------------------------------------------------------
409 
410  end function get_ocean_frac_2d
411 
412 ! </FUNCTION>
413 !#######################################################################
414 ! <FUNCTION NAME="get_ocean_mask">
415 
416 ! <OVERVIEW>
417 ! Returns a land-ocean mask in a grid box.
418 ! </OVERVIEW>
419 ! <DESCRIPTION>
420 ! Returns a land-ocean mask in the given model grid boxes.
421 ! </DESCRIPTION>
422 ! <TEMPLATE>
423 ! flag = <B>get_ocean_mask</B> ( blon, blat, ocean_mask )
424 ! </TEMPLATE>
425 
426 ! <IN NAME="blon" UNITS="radians" TYPE="real" DIM="(:)">
427 ! The longitude (in radians) at grid box boundaries.
428 ! </IN>
429 ! <IN NAME="blat" UNITS="radians" TYPE="real" DIM="(:)">
430 ! The latitude (in radians) at grid box boundaries.
431 ! </IN>
432 ! <OUT NAME="ocean_frac" TYPE="real" DIM="(:,:)">
433 ! The fractional amount (0 to 1) of ocean in a grid box.
434 ! The size of this field must be size(blon)-1 by size(blat)-1.
435 ! </OUT>
436 ! <OUT NAME="get_ocean_mask" TYPE="logical">
437 ! A logical value of TRUE is returned if the output field
438 ! was successfully created. A value of FALSE may be returned
439 ! if the Navy 1/6 degree percent water data set was not readable.
440 ! </OUT>
441 
442  function get_ocean_mask_1d (blon, blat, ocean_mask)
444  real , intent(in), dimension(:) :: blon, blat
445  logical, intent(out), dimension(:,:) :: ocean_mask
446  logical :: get_ocean_mask_1d
447 
448  real, dimension(size(ocean_mask,1),size(ocean_mask,2)) :: ocean_frac
449 !-----------------------------------------------------------------------
450  if (.not. module_is_initialized) call topography_init()
451 
452  if ( get_ocean_frac(blon, blat, ocean_frac) ) then
453  where (ocean_frac > 0.50)
454  ocean_mask = .true.
455  elsewhere
456  ocean_mask = .false.
457  end where
458  get_ocean_mask_1d = .true.
459  else
460  get_ocean_mask_1d = .false.
461  endif
462 
463 !-----------------------------------------------------------------------
464 
465  end function get_ocean_mask_1d
466 
467 !#######################################################################
468 
469  function get_ocean_mask_2d (blon, blat, ocean_mask)
471  real , intent(in), dimension(:,:) :: blon, blat
472  logical, intent(out), dimension(:,:) :: ocean_mask
473  logical :: get_ocean_mask_2d
474 
475  real, dimension(size(ocean_mask,1),size(ocean_mask,2)) :: ocean_frac
476 !-----------------------------------------------------------------------
477  if (.not. module_is_initialized) call topography_init()
478 
479  if ( get_ocean_frac(blon, blat, ocean_frac) ) then
480  where (ocean_frac > 0.50)
481  ocean_mask = .true.
482  elsewhere
483  ocean_mask = .false.
484  end where
485  get_ocean_mask_2d = .true.
486  else
487  get_ocean_mask_2d = .false.
488  endif
489 
490 !-----------------------------------------------------------------------
491 
492  end function get_ocean_mask_2d
493 
494 ! </FUNCTION>
495 !#######################################################################
496 ! <FUNCTION NAME="get_water_frac">
497 
498 ! <OVERVIEW>
499 ! Returns fractional area covered by water.
500 ! </OVERVIEW>
501 ! <DESCRIPTION>
502 ! Returns the percent of water in a grid box.
503 ! </DESCRIPTION>
504 ! <TEMPLATE>
505 ! flag = <B>get_water_frac</B> ( blon, blat, water_frac )
506 ! </TEMPLATE>
507 
508 ! <IN NAME="blon" UNITS="radians" TYPE="real" DIM="(:)">
509 ! The longitude (in radians) at grid box boundaries.
510 ! </IN>
511 ! <IN NAME="blat" UNITS="radians" TYPE="real" DIM="(:)">
512 ! The latitude (in radians) at grid box boundaries.
513 ! </IN>
514 ! <OUT NAME="water_frac" TYPE="real" DIM="(:,:)">
515 ! The fractional amount (0 to 1) of water in a grid box.
516 ! The size of this field must be size(blon)-1 by size(blat)-1.
517 ! </OUT>
518 ! <OUT NAME="get_water_frac" TYPE="logical">
519 ! A logical value of TRUE is returned if the output field
520 ! was successfully created. A value of FALSE may be returned
521 ! if the Navy 1/6 degree percent water data set was not readable.
522 ! </OUT>
523 ! <ERROR MSG="shape(water_frac) is not equal to (/size(blon)-1,size(blat)-1/))" STATUS="FATAL">
524 ! Check the input grid size and output field size.
525 ! </ERROR>
526 
527  function get_water_frac_1d (blon, blat, water_frac)
529  real, intent(in), dimension(:) :: blon, blat
530  real, intent(out), dimension(:,:) :: water_frac
531  logical :: get_water_frac_1d
532 
533 !-----------------------------------------------------------------------
534  if (.not. module_is_initialized) call topography_init()
535 
536  if ( any(shape(water_frac(:,:)) /= (/size(blon(:))-1,size(blat(:))-1/)) ) &
537  call error_mesg('get_water_frac_1d','shape(water_frac) is not&
538  & equal to (/size(blon)-1,size(blat)-1/))', fatal)
539 
540  if ( open_topog_file(water_file) ) then
541  call interp_water_1d ( blon, blat, water_frac )
542  get_water_frac_1d = .true.
543  else
544  get_water_frac_1d = .false.
545  endif
546 
547 !-----------------------------------------------------------------------
548 
549  end function get_water_frac_1d
550 
551 !#######################################################################
552 
553  function get_water_frac_2d (blon, blat, water_frac)
555  real, intent(in), dimension(:,:) :: blon, blat
556  real, intent(out), dimension(:,:) :: water_frac
557  logical :: get_water_frac_2d
558 
559 !-----------------------------------------------------------------------
560  if (.not. module_is_initialized) call topography_init()
561 
562  if ( any(shape(water_frac(:,:)) /= (/size(blon,1)-1,size(blon,2)-1/)) .or. &
563  any(shape(water_frac(:,:)) /= (/size(blat,1)-1,size(blat,2)-1/)) ) &
564  call error_mesg('get_water_frac_2d','shape(water_frac) is not&
565  & equal to (/size(blon,1)-1,size(blon,2)-1/))', fatal)
566 
567  if ( open_topog_file(water_file) ) then
568  call interp_water_2d ( blon, blat, water_frac )
569  get_water_frac_2d = .true.
570  else
571  get_water_frac_2d = .false.
572  endif
573 
574 !-----------------------------------------------------------------------
575 
576  end function get_water_frac_2d
577 
578 ! </FUNCTION>
579 !#######################################################################
580 ! <FUNCTION NAME="get_water_mask">
581 
582 ! <OVERVIEW>
583 ! Returns a land-water mask in a grid box.
584 ! </OVERVIEW>
585 ! <DESCRIPTION>
586 ! Returns a land-water mask in the given model grid boxes.
587 ! </DESCRIPTION>
588 ! <TEMPLATE>
589 ! flag = <B>get_water_mask</B> ( blon, blat, water_mask )
590 ! </TEMPLATE>
591 
592 ! <IN NAME="blon" UNITS="radians" TYPE="real" DIM="(:)">
593 ! The longitude (in radians) at grid box boundaries.
594 ! </IN>
595 ! <IN NAME="blat" UNITS="radians" TYPE="real" DIM="(:)">
596 ! The latitude (in radians) at grid box boundaries.
597 ! </IN>
598 ! <OUT NAME="water_mask" TYPE="real" DIM="(:,:)">
599 ! A binary mask for water (true) or land (false).
600 ! The size of this field must be size(blon)-1 by size(blat)-1.
601 ! </OUT>
602 ! <OUT NAME="get_water_mask" TYPE="logical">
603 ! A logical value of TRUE is returned if the output field
604 ! was successfully created. A value of FALSE may be returned
605 ! if the Navy 1/6 degree percent water data set was not readable.
606 ! </OUT>
607 
608  function get_water_mask_1d (blon, blat, water_mask)
610  real , intent(in), dimension(:) :: blon, blat
611  logical, intent(out), dimension(:,:) :: water_mask
612  logical :: get_water_mask_1d
613 
614  real, dimension(size(water_mask,1),size(water_mask,2)) :: water_frac
615 !-----------------------------------------------------------------------
616  if (.not. module_is_initialized) call topography_init()
617 
618  if ( get_water_frac(blon, blat, water_frac) ) then
619  where (water_frac > 0.50)
620  water_mask = .true.
621  elsewhere
622  water_mask = .false.
623  end where
624  get_water_mask_1d = .true.
625  else
626  get_water_mask_1d = .false.
627  endif
628 
629 !-----------------------------------------------------------------------
630 
631  end function get_water_mask_1d
632 
633 !#######################################################################
634 
635  function get_water_mask_2d (blon, blat, water_mask)
637  real , intent(in), dimension(:,:) :: blon, blat
638  logical, intent(out), dimension(:,:) :: water_mask
639  logical :: get_water_mask_2d
640 
641  real, dimension(size(water_mask,1),size(water_mask,2)) :: water_frac
642 !-----------------------------------------------------------------------
643  if (.not. module_is_initialized) call topography_init()
644 
645  if ( get_water_frac(blon, blat, water_frac) ) then
646  where (water_frac > 0.50)
647  water_mask = .true.
648  elsewhere
649  water_mask = .false.
650  end where
651  get_water_mask_2d = .true.
652  else
653  get_water_mask_2d = .false.
654  endif
655 
656 !-----------------------------------------------------------------------
657 
658  end function get_water_mask_2d
659 
660 ! </FUNCTION>
661 
662 !#######################################################################
663 !################## private interfaces below here ##################
664 !#######################################################################
665 
666  function open_topog_file ( filename )
667  character(len=*), intent(in) :: filename
668  logical :: open_topog_file
669  real :: r_ipts, r_jpts
670  integer :: namelen
671 
672  namelen = len(trim(filename))
673  if ( file_exist(filename) .AND. filename(namelen-2:namelen) == '.nc') then
674  if (mpp_pe() == mpp_root_pe()) call mpp_error ('topography_mod', &
675  'Reading NetCDF formatted input data file: '//filename, note)
676  call read_data(filename, 'ipts', r_ipts, no_domain=.true.)
677  call read_data(filename, 'jpts', r_jpts, no_domain=.true.)
678  ipts = nint(r_ipts)
679  jpts = nint(r_jpts)
680  open_topog_file = .true.
681  else
682  if ( file_exist(filename) ) then
683  if (mpp_pe() == mpp_root_pe()) call mpp_error ('topography_mod', &
684  'Reading native formatted input data file: '//filename, note)
685  unit = open_ieee32_file(trim(filename), 'read')
686  read (unit) ipts, jpts
687  open_topog_file = .true.
688  else
689  open_topog_file = .false.
690  endif
691  endif
692 
693  end function open_topog_file
694 
695 !#######################################################################
696 
697  subroutine interp_topog_1d ( blon, blat, zout, flag )
698  real , intent(in) :: blon(:), blat(:)
699  real , intent(out) :: zout(:,:)
700  integer, intent(in), optional :: flag
701 
702  real :: xdat(ipts+1), ydat(jpts+1)
703  real :: zdat(ipts,jpts)
704  real :: zout2(size(zout,1),size(zout,2))
705 
706  call input_data ( topog_file, xdat, ydat, zdat )
707 
708  call horiz_interp ( zdat, xdat, ydat, blon, blat, zout )
709 
710 ! compute standard deviation if necessary
711  if (present(flag)) then
712  if (flag == compute_stdev) then
713  zdat = zdat*zdat
714  call horiz_interp ( zdat, xdat, ydat, blon, blat, zout2 )
715  zout = zout2 - zout*zout
716  where (zout > 0.0)
717  zout = sqrt( zout )
718  elsewhere
719  zout = 0.0
720  endwhere
721  endif
722  endif
723 
724  end subroutine interp_topog_1d
725 
726 !#######################################################################
727 
728  subroutine interp_topog_2d ( blon, blat, zout, flag )
729  real , intent(in) :: blon(:,:), blat(:,:)
730  real , intent(out) :: zout(:,:)
731  integer, intent(in), optional :: flag
732 
733  real :: xdat(ipts+1), ydat(jpts+1)
734  real :: zdat(ipts,jpts)
735  real :: zout2(size(zout,1),size(zout,2))
736  integer :: js, je
737  type(horiz_interp_type) :: interp
738 
739  call input_data ( topog_file, xdat, ydat, zdat )
740  call find_indices ( minval(blat), maxval(blat), ydat, js, je )
741 
742  call horiz_interp_new ( interp, xdat, ydat(js:je+1), blon, blat )
743  call horiz_interp ( interp, zdat(:,js:je), zout )
744 
745 ! compute standard deviation if necessary
746  if (present(flag)) then
747  if (flag == compute_stdev) then
748  zdat = zdat*zdat
749  call horiz_interp ( interp, zdat(:,js:je), zout2 )
750  zout = zout2 - zout*zout
751  where (zout > 0.0)
752  zout = sqrt( zout )
753  elsewhere
754  zout = 0.0
755  endwhere
756  endif
757  endif
758 
759  call horiz_interp_del ( interp )
760 
761  end subroutine interp_topog_2d
762 
763 !#######################################################################
764 
765  subroutine find_indices ( ybeg, yend, ydat, js, je )
766  real, intent(in) :: ybeg, yend, ydat(:)
767  integer, intent(out) :: js, je
768  integer :: j
769 
770  js = 1
771  do j = 1, size(ydat(:))-1
772  if (ybeg >= ydat(j) .and. ybeg <= ydat(j+1)) then
773  js = j
774  exit
775  endif
776  enddo
777 
778  je = size(ydat(:))-1
779  do j = js, size(ydat(:))-1
780  if (yend >= ydat(j) .and. yend <= ydat(j+1)) then
781  je = j
782  exit
783  endif
784  enddo
785 
786  !print '(a,i2,2(a,f10.5),2(a,i4))', "PE=",mpp_pe()," phs=",ybeg," phn=",yend," js=",js," je=",je
787 
788  end subroutine find_indices
789 
790 !#######################################################################
791 
792  subroutine input_data ( ifile, xdat, ydat, zdat )
793  character(len=*), intent(in) :: ifile
794  real, intent(out) :: xdat(ipts+1), ydat(jpts+1), zdat(ipts,jpts)
795  integer :: nc
796 
797  nc = len_trim(ifile)
798 
799 ! note: ipts,jpts,unit are global
800 
801  if ( file_exist(trim(ifile)) .AND. ifile(nc-2:nc) == '.nc') then
802  call read_data(trim(ifile), 'xdat', xdat, no_domain=.true.)
803  call read_data(trim(ifile), 'ydat', ydat, no_domain=.true.)
804  call read_data(trim(ifile), 'zdat', zdat, no_domain=.true.)
805  else
806  read (unit) xdat, ydat ! read lon/lat edges in radians
807  read (unit) zdat ! read land surface height in meters
808  call close_file (unit)
809  endif
810 
811  end subroutine input_data
812 
813 !#######################################################################
814 
815  subroutine interp_water_1d ( blon, blat, zout, do_ocean )
816  real , intent(in) :: blon(:), blat(:)
817  real , intent(out) :: zout(:,:)
818  logical, intent(in), optional :: do_ocean
819 
820  real :: xdat(ipts+1), ydat(jpts+1), zdat(ipts,jpts)
821 
822  call input_data ( water_file, xdat, ydat, zdat )
823 
824 ! only use designated ocean points
825  if (present(do_ocean)) then
826  if (do_ocean) call determine_ocean_points (zdat)
827  endif
828 
829 ! interpolate onto output grid
830  call horiz_interp ( zdat, xdat, ydat, blon, blat, zout )
831 
832  end subroutine interp_water_1d
833 
834 !#######################################################################
835 
836  subroutine interp_water_2d ( blon, blat, zout, do_ocean )
837  real , intent(in) :: blon(:,:), blat(:,:)
838  real , intent(out) :: zout(:,:)
839  logical, intent(in), optional :: do_ocean
840 
841  real :: xdat(ipts+1), ydat(jpts+1), zdat(ipts,jpts)
842 
843  call input_data ( water_file, xdat, ydat, zdat )
844 
845 ! only use designated ocean points
846  if (present(do_ocean)) then
847  if (do_ocean) call determine_ocean_points (zdat)
848  endif
849 
850 ! interpolate onto output grid
851  call horiz_interp ( zdat, xdat, ydat, blon, blat, zout )
852 
853  end subroutine interp_water_2d
854 
855 !#######################################################################
856 
857  subroutine determine_ocean_points ( pctwater )
858  real, intent(inout) :: pctwater(:,:)
859  logical :: ocean(size(pctwater,1),size(pctwater,2))
860  integer :: i, j, m, n, im, ip, jm, jp, new
861 
862  real :: ocean_pct_crit = .500
863 
864  ! resolution of the grid
865  m = size(pctwater,1)
866  n = size(pctwater,2)
867 
868  ! the 1/6 degree navy percent water data set
869  ! designates ocean grid boxes as 100 percent water
870  ! all other grid boxes have <= 99 percent water
871 
872  ! set a mask for ocean grid boxes
873  ocean = (pctwater > .999)
874  new = count(ocean)
875 
876  ! set land grid boxes that have sufficient amount of water
877  ! to ocean grid boxes when they are adjacent to ocean points
878  ! iterate until there are no new ocean points
879  do
880  if (new == 0) exit
881  new = 0
882 
883  do j = 1, n
884  do i = 1, m
885  if (.not.ocean(i,j) .and. pctwater(i,j) > ocean_pct_crit) then
886  im = i-1; ip = i+1; jm = j-1; jp = j+1
887  if (im == 0) im = m
888  if (ip == m+1) ip = 1
889  if (jm == 0) jm = 1
890  if (jp == n+1) jp = n
891  ! check the 8 grid boxes that surround this grid box
892  if (ocean(im,j ) .or. ocean(ip,j ) .or. ocean(i ,jm) .or. ocean(i ,jp) .or. &
893  ocean(im,jm) .or. ocean(ip,jm) .or. ocean(ip,jp) .or. ocean(im,jp)) then
894  ocean(i,j) = .true.
895  new = new + 1
896  endif
897  endif
898  enddo
899  enddo
900  !print *, 'new=',new
901 
902  enddo
903 
904  ! final step is to elimate water percentage if land
905  where (.not.ocean) pctwater = 0.
906 
907  end subroutine determine_ocean_points
908 
909 !#######################################################################
910 ! reads the namelist file, write namelist to log file,
911 ! and initializes constants
912 
913 subroutine read_namelist
915  integer :: unit, ierr, io
916 
917 ! read namelist
918 
919 #ifdef INTERNAL_FILE_NML
920  read (input_nml_file, topography_nml, iostat=io)
921  ierr = check_nml_error(io,'topography_nml')
922 #else
923  if ( file_exist('input.nml')) then
924  unit = open_namelist_file( )
925  ierr=1; do while (ierr /= 0)
926  read (unit, nml=topography_nml, iostat=io, end=10)
927  ierr = check_nml_error(io,'topography_nml')
928  enddo
929  10 call close_file (unit)
930  endif
931 #endif
932 
933 ! write version and namelist to log file
934 
935  if (mpp_pe() == mpp_root_pe()) then
936  unit = stdlog()
937  write (unit, nml=topography_nml)
938  endif
939 
940 end subroutine read_namelist
941 
942 !#######################################################################
943 
944 end module topography_mod
945 
946 ! <INFO>
947 
948 ! <TESTPROGRAM NAME="">
949 !
950 ! To run this program you will need the topography and percent water
951 ! data sets and use the following namelist (in file input.nml).
952 !
953 ! &amp;gaussian_topog_nml
954 ! height = 5000., 3000., 3000., 3000.,
955 ! olon = 90., 255., 285., 0.,
956 ! olat = 45., 45., -15., -90.,
957 ! wlon = 15., 10., 5., 180.,
958 ! wlat = 15., 25., 25., 20., /
959 !
960 ! program test
961 !
962 ! ! test program for topography and gaussian_topog modules
963 ! <PRE>
964 ! use topography_mod
965 ! implicit none
966 !
967 ! integer, parameter :: nlon=24, nlat=18
968 ! real :: x(nlon), y(nlat), xb(nlon+1), yb(nlat+1), z(nlon,nlat)
969 ! real :: hpi, rtd
970 ! integer :: i,j
971 ! logical :: a
972 !
973 ! ! gaussian mountain parameters
974 ! real, parameter :: ht=4000.
975 ! real, parameter :: x0=90., y0=45. ! origin in degrees
976 ! real, parameter :: xw=15., yw=15. ! half-width in degees
977 ! real, parameter :: xr=30., yr= 0. ! ridge-width in degrees
978 !
979 ! ! create lat/lon grid in radians
980 ! hpi = acos(0.0)
981 ! rtd = 90./hpi ! rad to deg
982 ! do i=1,nlon
983 ! xb(i) = 4.*hpi*real(i-1)/real(nlon)
984 ! enddo
985 ! xb(nlon+1) = xb(1)+4.*hpi
986 ! yb(1) = -hpi
987 ! do j=2,nlat
988 ! yb(j) = yb(j-1) + 2.*hpi/real(nlat)
989 ! enddo
990 ! yb(nlat+1) = hpi
991 ! ! mid-point of grid boxes
992 ! x(1:nlon) = 0.5*(xb(1:nlon)+xb(2:nlon+1))
993 ! y(1:nlat) = 0.5*(yb(1:nlat)+yb(2:nlat+1))
994 ! ! test topography_mod routines
995 ! a = get_topog_mean(xb,yb,z)
996 ! call printz ('get_topog_mean')
997 !
998 ! a = get_water_frac(xb,yb,z)
999 ! z = z*100. ! in percent
1000 ! call printz ('get_water_frac')
1001 !
1002 ! a = get_ocean_frac(xb,yb,z)
1003 ! z = z*100. ! in percent
1004 ! call printz ('get_ocean_frac')
1005 !
1006 ! ! test gaussian_topog_mod routines
1007 ! a = .true.
1008 ! z = get_gaussian_topog(x,y,ht,x0,y0,xw,yw,xr,yr)
1009 ! call printz ('get_gaussian_topog')
1010 !
1011 ! call gaussian_topog_init (x,y,z)
1012 ! call printz ('gaussian_topog_init')
1013 !
1014 ! contains
1015 !
1016 ! ! simple printout of topog/water array
1017 ! subroutine printz (lab)
1018 ! character(len=*), intent(in) :: lab
1019 ! if (a) then
1020 ! print '(/a)', trim(lab)
1021 ! else
1022 ! print '(/a)', 'no data available: '//trim(lab)
1023 ! return
1024 ! endif
1025 ! ! print full grid
1026 ! print '(3x,25i5)', (nint(x(i)*rtd),i=1,nlon)
1027 ! do j=nlat,1,-1
1028 ! print '(i3,25i5)', nint(y(j)*rtd), (nint(z(i,j)),i=1,nlon)
1029 ! enddo
1030 ! end subroutine printz
1031 !
1032 ! end program test
1033 ! </PRE>
1034 ! </TESTPROGRAM>
1035 
1036 ! <BUG>
1037 ! Water mask produces some possible erroneous water points along
1038 ! the coast of Antarctic (at about 90W).
1039 ! </BUG>
1040 
1041 ! <FUTURE>Use of netcdf data sets. </FUTURE>
1042 ! <FUTURE>Incorporate other topography and ocean data sets. </FUTURE>
1043 !
1044 ! </INFO>
Definition: fms.F90:20
subroutine interp_water_1d(blon, blat, zout, do_ocean)
Definition: topography.F90:816
subroutine, public topography_init()
Definition: topography.F90:148
logical function get_water_mask_2d(blon, blat, water_mask)
Definition: topography.F90:636
logical function get_topog_stdev_1d(blon, blat, stdev)
Definition: topography.F90:280
logical function get_ocean_frac_1d(blon, blat, ocean_frac)
Definition: topography.F90:362
subroutine, public horiz_interp_del(Interp)
logical module_is_initialized
Definition: topography.F90:139
logical function get_water_frac_2d(blon, blat, water_frac)
Definition: topography.F90:554
logical function get_ocean_mask_1d(blon, blat, ocean_mask)
Definition: topography.F90:443
subroutine interp_topog_2d(blon, blat, zout, flag)
Definition: topography.F90:729
Definition: mpp.F90:39
integer function, public check_nml_error(IOSTAT, NML_NAME)
Definition: fms.F90:658
logical function get_water_frac_1d(blon, blat, water_frac)
Definition: topography.F90:528
logical function get_topog_mean_2d(blon, blat, zmean)
Definition: topography.F90:221
logical function get_ocean_mask_2d(blon, blat, ocean_mask)
Definition: topography.F90:470
character(len=input_str_length), dimension(:), allocatable, target, public input_nml_file
Definition: mpp.F90:1378
subroutine find_indices(ybeg, yend, ydat, js, je)
Definition: topography.F90:766
subroutine, public gaussian_topog_init(lon, lat, zsurf)
logical function get_water_mask_1d(blon, blat, water_mask)
Definition: topography.F90:609
logical function get_topog_stdev_2d(blon, blat, stdev)
Definition: topography.F90:306
logical function open_topog_file(filename)
Definition: topography.F90:667
logical function get_topog_mean_1d(blon, blat, zmean)
Definition: topography.F90:195
real function, dimension(size(lon, 1), size(lat, 1)), public get_gaussian_topog(lon, lat, height, olond, olatd, wlond, wlatd, rlond, rlatd)
subroutine read_namelist
Definition: topography.F90:914
character(len=128) topog_file
Definition: topography.F90:97
subroutine interp_topog_1d(blon, blat, zout, flag)
Definition: topography.F90:698
integer, parameter compute_stdev
Definition: topography.F90:131
subroutine input_data(ifile, xdat, ydat, zdat)
Definition: topography.F90:793
character(len=128) water_file
Definition: topography.F90:97
subroutine interp_water_2d(blon, blat, zout, do_ocean)
Definition: topography.F90:837
subroutine, public error_mesg(routine, message, level)
Definition: fms.F90:529
logical function get_ocean_frac_2d(blon, blat, ocean_frac)
Definition: topography.F90:388
subroutine determine_ocean_points(pctwater)
Definition: topography.F90:858
real(fp), parameter, public pi