FE-Project
Loading...
Searching...
No Matches
mod_mktopo.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
10!-------------------------------------------------------------------------------
11#include "scalelib.h"
13 !-----------------------------------------------------------------------------
14 !
15 !++ used modules
16 !
17 use scale_precision
18 use scale_io
19 use scale_prof
20 use scale_prc
21
22 use scale_const, only: &
23 pi => const_pi, &
24 rplanet => const_radius
25
26 use scale_element_base, only: &
36
37 use mod_atmos_component, only: &
39 use mod_atmos_mesh, only: atmosmesh
40
41 !-----------------------------------------------------------------------------
42 implicit none
43 private
44 !-----------------------------------------------------------------------------
45 !
46 !++ Public procedure
47 !
48 public :: mktopo_setup
49 public :: mktopo
50 public :: mktopo_write
51
52 !-----------------------------------------------------------------------------
53 !
54 !++ Public parameters & variables
55 !
56
57 integer, public :: mktopo_type = -1
58 integer, parameter, public :: i_ignore = 0
59 integer, parameter, public :: i_input_file = 1
60 integer, parameter, public :: i_flat = 2
61 integer, parameter, public :: i_bellshape = 3
62 integer, parameter, public :: i_scaer = 4
63 integer, parameter, public :: i_bellshape_global = 5
64 integer, parameter, public :: i_scaher_global = 6
65 integer, parameter, public :: i_barocwave_global_jw2006 = 7
66
67 !-----------------------------------------------------------------------------
68 !
69 !++ Private procedure
70 !
71
72 private :: mktopo_flat
73 private :: mktopo_bellshape
74 private :: mktopo_schaer_type_mountain
75 private :: mktopo_bellshape_global
76 private :: mktopo_barocwave_global_jw2006
77
78 !-----------------------------------------------------------------------------
79 !
80 !++ Private variables
81 !
82 character(len=H_LONG) :: out_basename = ''
83 character(len=H_MID) :: out_varname = 'topo'
84 character(len=H_MID) :: out_title = 'SCALE-DG TOPOGRAPHY'
85 character(len=H_SHORT) :: out_dtype = 'DEFAULT'
86 !-----------------------------------------------------------------------------
87contains
88 !-----------------------------------------------------------------------------
90!OCL SERIAL
91 subroutine mktopo_setup
92 implicit none
93
94 character(len=H_SHORT) :: toponame = 'NONE'
95
96 namelist / param_mktopo / &
97 toponame, &
98 out_basename, out_varname, &
99 out_title, out_dtype
100
101
102 integer :: ierr
103 !---------------------------------------------------------------------------
104
105 log_newline
106 log_info("MKTOPO_setup",*) 'Setup'
107
108 !--- read namelist
109 rewind(io_fid_conf)
110 read(io_fid_conf,nml=param_mktopo,iostat=ierr)
111 if( ierr < 0 ) then !--- missing
112 log_info("MKTOPO_setup",*) 'Not found namelist. Default used.'
113 elseif( ierr > 0 ) then !--- fatal error
114 log_error("MKTOPO_setup",*) 'Not appropriate names in namelist PARAM_MKTOPO. Check!'
115 call prc_abort
116 endif
117 log_nml(param_mktopo)
118
119 select case(trim(toponame))
120 case('NONE')
122 case('INPUT_FILE')
124 case('FLAT')
126 case('BELLSHAPE')
128 case('SCHAER')
130 case('BELLSHAPE_GLOBAL')
132 case('SCHAER_GLOBAL')
134 case('BAROCWAVE_GLOBAL_JW2006')
136 case default
137 log_error("MKTOPO_setup",*) 'Not appropriate toponame. Check!', toponame
138 call prc_abort
139 end select
140
141 return
142 end subroutine mktopo_setup
143
144 !-----------------------------------------------------------------------------
146!OCL SERIAL
147 subroutine mktopo( output, model_mesh, topography )
149 implicit none
150
151 logical, intent(out) :: output
152 class(atmosmesh), target, intent(in) :: model_mesh
153 class(meshtopography), intent(inout) :: topography
154
155 integer :: n
156 integer :: ke
157 class(localmesh3d), pointer :: lcmesh3d
158 class(meshbase3d), pointer :: mesh
159 class(localmesh2d), pointer :: lcmesh2d
160 class(meshbase2d), pointer :: mesh2d
161 !---------------------------------------------------------------------------
162
163 mesh => model_mesh%ptr_mesh
164
165 if ( mktopo_type == i_ignore ) then
166 log_newline
167 log_progress(*) 'skip making topography data'
168 output = .false.
169 else
170 log_newline
171 log_progress(*) 'start making topography data'
172
173 call prof_rapstart('_MkTOPO_main', 3)
174
175 call mesh%GetMesh2D( mesh2d )
176
177 select case( mktopo_type )
178 case ( i_input_file )
179 call mktopo_input_file( mesh2d, topography%topo )
180 case ( i_flat )
181 call mktopo_flat( mesh2d, topography%topo )
182 case ( i_bellshape )
183 call mktopo_bellshape( mesh2d, topography%topo )
184 case ( i_scaer )
185 call mktopo_schaer_type_mountain( mesh2d, topography%topo )
186 case ( i_bellshape_global )
187 call mktopo_bellshape_global( mesh2d, topography%topo )
188 case ( i_scaher_global )
189 call mktopo_schaer_type_global( mesh2d, topography%topo )
191 call mktopo_barocwave_global_jw2006( mesh2d, topography%topo )
192 end select
193
194 call prof_rapend ('_MkTOPO_main', 3)
195 log_progress(*) 'end making topography data'
196
197 output = .true.
198 end if
199
200 return
201 end subroutine mktopo
202
203 !-----------------------------------------------------------------------------
205!OCL SERIAL
206 subroutine mktopo_write( model_mesh, topography )
207
208 use scale_time, only: &
209 nowdate => time_nowdate, &
210 nowsubsec => time_nowsubsec, &
211 nowdaysec => time_nowdaysec
212 use scale_prc, only: prc_myrank
213
217 use scale_mesh_base2d, only: &
218 mftype2d_xy => meshbase2d_dimtypeid_xy
219
220 implicit none
221 class(atmosmesh), target, intent(in) :: model_mesh
222 class(meshtopography), intent(inout) :: topography
223
224 type(file_base_meshfield) :: file
225 class(meshbase2d), pointer :: mesh
226
227 logical :: file_existed
228 integer, parameter :: vid_topo = 1
229 !--------------------------------------------------
230
231 nullify( mesh )
232 call model_mesh%ptr_mesh%GetMesh2D( mesh )
233 select type( mesh )
234 class is (meshcubedspheredom2d)
235 call file%Init( 1, meshcubedsphere2d=mesh )
236 class is (meshrectdom2d)
237 call file%Init( 1, mesh2d=mesh )
238 end select
239
240 call file%Create( out_basename, out_title, out_dtype, &
241 file_existed, &
242 myrank=prc_myrank )
243
244 if ( .not. file_existed ) then
245 call file%Put_GlobalAttribute_time( nowdate, nowsubsec )
246 end if
247
248 call file%Def_Var( topography%topo, "TOPOGRAPHY", vid_topo, mftype2d_xy, 'XY' )
249 call file%End_def()
250 call file%Write_var2D( vid_topo, topography%topo, nowdaysec, nowdaysec )
251 call file%Close()
252 call file%Final()
253
254 return
255 end subroutine mktopo_write
256
257 !-- private---------------------------------------------------------------------
258
259 !-----------------------------------------------------------------------------
262!OCL SERIAL
263 subroutine mktopo_input_file( mesh, topo )
267 use scale_mesh_base2d, only: &
268 mftype2d_xy => meshbase2d_dimtypeid_xy
269
270 implicit none
271
272 class(meshbase2d), intent(in), target :: mesh
273 class(meshfield2d), intent(inout) :: topo
274
275 character(len=H_LONG) :: in_basename
276 character(len=H_MID) :: varname
277
278 namelist / param_mktopo_input_file / &
279 in_basename, &
280 varname
281
282 integer :: ierr
283 integer :: n
284 integer :: ke2d
285 type(localmesh2d), pointer :: lmesh2d
286
287 type(file_base_meshfield) :: file
288
289 logical :: file_existed
290 !--------------------------------
291
292 log_info("MKTOPO_input_file",*) 'Setup'
293
294 !--- read namelist
295 rewind(io_fid_conf)
296 read(io_fid_conf,nml=param_mktopo_input_file,iostat=ierr)
297 if( ierr < 0 ) then !--- missing
298 log_info("MKTOPO_input_file",*) 'Not found namelist. Default used.'
299 elseif( ierr > 0 ) then !--- fatal error
300 log_error("MKTOPO_input_file",*) 'Not appropriate names in namelist PARAM_MKTOPO_INPUT_FILE. Check!'
301 call prc_abort
302 endif
303 log_nml(param_mktopo_input_file)
304
305 !-- Read topography data
306
307 select type( mesh )
308 class is (meshcubedspheredom2d)
309 call file%Init( 1, meshcubedsphere2d=mesh )
310 class is (meshrectdom2d)
311 call file%Init( 1, mesh2d=mesh )
312 end select
313
314 call file%Open( in_basename, myrank=prc_myrank)
315 call file%Read_Var( mftype2d_xy, varname, topo )
316
317 call file%Close()
318 call file%Final()
319
320 return
321 end subroutine mktopo_input_file
322
323
324 !-----------------------------------------------------------------------------
327!OCL SERIAL
328 subroutine mktopo_flat( mesh, topo )
329 implicit none
330
331 class(meshbase2d), intent(in), target :: mesh
332 class(meshfield2d), intent(inout) :: topo
333
334 ! flat mountain parameter
335 real(rp) :: flat_height = 100.0_rp
336
337 namelist / param_mktopo_flat / &
338 flat_height
339
340 integer :: ierr
341 integer :: n
342 integer :: ke2d
343 type(localmesh2d), pointer :: lmesh2d
344 !--------------------------------
345
346 log_info("MKTOPO_flat",*) 'Setup'
347
348 !--- read namelist
349 rewind(io_fid_conf)
350 read(io_fid_conf,nml=param_mktopo_flat,iostat=ierr)
351 if( ierr < 0 ) then !--- missing
352 log_info("MKTOPO_flat",*) 'Not found namelist. Default used.'
353 elseif( ierr > 0 ) then !--- fatal error
354 log_error("MKTOPO_flat",*) 'Not appropriate names in namelist PARAM_MKTOPO_FLAT. Check!'
355 call prc_abort
356 endif
357 log_nml(param_mktopo_flat)
358
359 do n=1, mesh%LOCAL_MESH_NUM
360 lmesh2d => mesh%lcmesh_list(n)
361 !$omp parallel do
362 do ke2d=lmesh2d%NeS, lmesh2d%NeE
363 topo%local(n)%val(:,ke2d) = flat_height
364 end do
365 end do
366
367 return
368 end subroutine mktopo_flat
369
370
371 !-----------------------------------------------------------------------------
374!OCL SERIAL
375 subroutine mktopo_bellshape( mesh, topo )
376 implicit none
377
378 class(meshbase2d), intent(in), target :: mesh
379 class(meshfield2d), intent(inout) :: topo
380
381 ! bell-shaped mountain parameter
382 real(rp) :: bell_cx = 0.0e0_rp
383 real(rp) :: bell_cy = 0.0e3_rp
384 real(rp) :: bell_r = 2.e3_rp
385 real(rp) :: bell_height = 100.0_rp
386 logical :: bell_quasi_2d = .false.
387 namelist / param_mktopo_bellshape / &
388 bell_cx, &
389 bell_cy, &
390 bell_r, &
391 bell_height, &
392 bell_quasi_2d
393
394 integer :: ierr
395 integer :: n
396 integer :: ke2d
397 type(localmesh2d), pointer :: lmesh2d
398 class(elementbase2d), pointer :: elem
399
400 real(rp), allocatable :: r2(:)
401 real(rp) :: fac_y_quasi2d
402 !--------------------------------
403
404 log_info("MKTOPO_bellshape",*) 'Setup'
405
406 !--- read namelist
407 rewind(io_fid_conf)
408 read(io_fid_conf,nml=param_mktopo_bellshape,iostat=ierr)
409 if( ierr < 0 ) then !--- missing
410 log_info("MKTOPO_bellshape",*) 'Not found namelist. Default used.'
411 elseif( ierr > 0 ) then !--- fatal error
412 log_error("MKTOPO_bellshape",*) 'Not appropriate names in namelist PARAM_MKTOPO_BELLSHAPE. Check!'
413 call prc_abort
414 endif
415 log_nml(param_mktopo_bellshape)
416
417 if ( bell_quasi_2d ) then
418 fac_y_quasi2d = 0.0_rp
419 else
420 fac_y_quasi2d = 1.0_rp
421 end if
422
423 do n=1, mesh%LOCAL_MESH_NUM
424 lmesh2d => mesh%lcmesh_list(n)
425 elem => lmesh2d%refElem2D
426
427 allocate( r2(elem%Np) )
428 !$omp parallel do private(r2)
429 do ke2d=lmesh2d%NeS, lmesh2d%NeE
430 r2(:) = ( ( lmesh2d%pos_en(:,ke2d,1) - bell_cx )**2 &
431 + fac_y_quasi2d * ( lmesh2d%pos_en(:,ke2d,2) - bell_cy )**2 ) &
432 / bell_r**2
433 topo%local(n)%val(:,ke2d) = bell_height / ( 1.0_rp + r2(:) )
434 end do
435 deallocate( r2 )
436 end do
437
438 return
439 end subroutine mktopo_bellshape
440
441 !-----------------------------------------------------------------------------
447!OCL SERIAL
448 subroutine mktopo_schaer_type_mountain( mesh, topo )
449 implicit none
450
451 class(meshbase2d), intent(in), target :: mesh
452 class(meshfield2d), intent(inout) :: topo
453
454 ! Schaer-type mountain parameter
455 real(rp) :: schaer_cx = 25.e3_rp
456 real(rp) :: schaer_rx = 5.e3_rp
457 real(rp) :: schaer_lambda = 4.e3_rp
458 real(rp) :: schaer_height = 250.0_rp
459 namelist / param_mktopo_schaer / &
460 schaer_cx, &
461 schaer_rx, &
462 schaer_lambda, &
463 schaer_height
464
465 integer :: ierr
466 integer :: n
467 integer :: ke2d
468 type(localmesh2d), pointer :: lmesh2d
469 class(elementbase2d), pointer :: elem
470
471 real(rp), allocatable :: dist(:)
472 real(rp) :: fac_y_quasi2d
473 !--------------------------------
474
475 log_info("MKTOPO_Schaer",*) 'Setup'
476
477 !--- read namelist
478 rewind(io_fid_conf)
479 read(io_fid_conf,nml=param_mktopo_schaer,iostat=ierr)
480 if( ierr < 0 ) then !--- missing
481 log_info("MKTOPO_Schaer",*) 'Not found namelist. Default used.'
482 elseif( ierr > 0 ) then !--- fatal error
483 log_error("MKTOPO_Schaer",*) 'Not appropriate names in namelist PARAM_MKTOPO_SCAHER. Check!'
484 call prc_abort
485 endif
486 log_nml(param_mktopo_schaer)
487
488 do n=1, mesh%LOCAL_MESH_NUM
489 lmesh2d => mesh%lcmesh_list(n)
490 elem => lmesh2d%refElem2D
491
492 allocate( dist(elem%Np) )
493 !$omp parallel do private(dist)
494 do ke2d=lmesh2d%NeS, lmesh2d%NeE
495 dist(:) = exp( - ( lmesh2d%pos_en(:,ke2d,1) - schaer_cx )**2 / schaer_rx**2 )
496 topo%local(n)%val(:,ke2d) = schaer_height * dist(:) * ( cos( pi * ( lmesh2d%pos_en(:,ke2d,1) - schaer_cx ) / schaer_lambda ) )**2
497 end do
498 deallocate( dist )
499 end do
500
501 return
502 end subroutine mktopo_schaer_type_mountain
503
504 !-----------------------------------------------------------------------------
507 subroutine mktopo_bellshape_global( mesh, topo )
508 implicit none
509
510 class(meshbase2d), intent(in), target :: mesh
511 class(meshfield2d), intent(inout) :: topo
512
513 ! bell-shaped mountain parameter
514 real(rp) :: bell_clon = 0.0e0_rp
515 real(rp) :: bell_clat = 0.0e0_rp
516 real(rp) :: bell_r = 2.e3_rp
517 real(rp) :: bell_height = 100.0_rp
518 namelist / param_mktopo_bellshape_global / &
519 bell_clon, &
520 bell_clat, &
521 bell_r, &
522 bell_height
523
524 integer :: ierr
525 integer :: n
526 integer :: ke2d
527 type(localmesh2d), pointer :: lmesh2d
528 class(elementbase2d), pointer :: elem
529
530 real(rp), allocatable :: r(:)
531 !--------------------------------
532
533 log_info("MKTOPO_bellshape_global",*) 'Setup'
534
535 !--- read namelist
536 rewind(io_fid_conf)
537 read(io_fid_conf,nml=param_mktopo_bellshape_global,iostat=ierr)
538 if( ierr < 0 ) then !--- missing
539 log_info("MKTOPO_bellshape_global",*) 'Not found namelist. Default used.'
540 elseif( ierr > 0 ) then !--- fatal error
541 log_error("MKTOPO_bellshape_global",*) 'Not appropriate names in namelist PARAM_MKTOPO_BELLSHAPE_GLOBAL. Check!'
542 call prc_abort
543 endif
544 log_nml(param_mktopo_bellshape_global)
545
546 do n=1, mesh%LOCAL_MESH_NUM
547 lmesh2d => mesh%lcmesh_list(n)
548 elem => lmesh2d%refElem2D
549
550 allocate( r(elem%Np) )
551 !$omp parallel do private(r)
552 do ke2d=lmesh2d%NeS, lmesh2d%NeE
553 r(:) = rplanet / bell_r &
554 * acos( sin(bell_clat) * sin(lmesh2d%lat(:,ke2d)) &
555 + cos(bell_clat) * cos(lmesh2d%lat(:,ke2d)) * cos(lmesh2d%lon(:,ke2d) - bell_clon) )
556 topo%local(n)%val(:,ke2d) = bell_height / ( 1.0_rp + r(:)**2 )
557 end do
558 deallocate( r )
559 end do
560
561 return
562 end subroutine mktopo_bellshape_global
563
564 !-----------------------------------------------------------------------------
570 subroutine mktopo_schaer_type_global( mesh, topo )
571 implicit none
572
573 class(meshbase2d), intent(in), target :: mesh
574 class(meshfield2d), intent(inout) :: topo
575
576 ! bell-shaped mountain parameter
577 real(rp) :: schaer_clon
578 real(rp) :: schaer_clat
579 real(rp) :: schaer_r
580 real(rp) :: schaer_lambda
581 integer :: schaer_shape_id
582 real(rp) :: schaer_height
583 logical :: quasi_2d_flag
584
585 namelist / param_mktopo_schaer_global / &
586 schaer_clon, &
587 schaer_clat, &
588 schaer_r, &
589 schaer_lambda, &
590 schaer_shape_id, &
591 schaer_height, &
592 quasi_2d_flag
593
594 integer :: ierr
595 integer :: n
596 integer :: ke2d
597 type(localmesh2d), pointer :: lmesh2d
598 class(elementbase2d), pointer :: elem
599
600 real(rp), allocatable :: r(:)
601 real(rp), allocatable :: mwt_shape(:)
602
603 !--------------------------------
604
605 log_info("MKTOPO_Schaer_global",*) 'Setup'
606
607 schaer_clon = 1.0_rp / 4.0_rp * pi
608 schaer_clat = 0.0_rp
609 schaer_r = rplanet * 3.0_rp / 4.0_rp * pi
610 schaer_lambda = rplanet * pi / 16.0_rp
611 schaer_shape_id = 1
612 schaer_height = 2000.0_rp
613 quasi_2d_flag = .false.
614
615 !--- read namelist
616 rewind(io_fid_conf)
617 read(io_fid_conf,nml=param_mktopo_schaer_global,iostat=ierr)
618 if( ierr < 0 ) then !--- missing
619 log_info("MKTOPO_Schaer_global",*) 'Not found namelist. Default used.'
620 elseif( ierr > 0 ) then !--- fatal error
621 log_error("MKTOPO_Schaer_global",*) 'Not appropriate names in namelist PARAM_MKTOPO_SCHAER_GLOBAL. Check!'
622 call prc_abort
623 endif
624 log_nml(param_mktopo_schaer_global)
625
626 do n=1, mesh%LOCAL_MESH_NUM
627 lmesh2d => mesh%lcmesh_list(n)
628 elem => lmesh2d%refElem2D
629
630 allocate( r(elem%Np), mwt_shape(elem%Np) )
631 !$omp parallel do private(r, mwt_shape)
632 do ke2d=lmesh2d%NeS, lmesh2d%NeE
633 if ( quasi_2d_flag ) then
634! r(:) = RPlanet * ( min( abs(lmesh2D%lon(:,ke2d) - SCHAER_Clon), abs(lmesh2D%lon(:,ke2d) - (SCHAER_Clon + 2.0_RP * PI ) ) ) )
635 r(:) = rplanet * abs(lmesh2d%lon(:,ke2d) - schaer_clon)
636 else
637 r(:) = rplanet &
638 * acos( sin(schaer_clat) * sin(lmesh2d%lat(:,ke2d)) &
639 + cos(schaer_clat) * cos(lmesh2d%lat(:,ke2d)) * cos(lmesh2d%lon(:,ke2d) - schaer_clon) )
640 end if
641
642 if ( schaer_shape_id == 1 ) then
643 mwt_shape(:) = exp( - ( r(:) / schaer_r )**2 )
644 else if ( schaer_shape_id == 2 ) then
645 where ( r(:) < schaer_r )
646 mwt_shape(:) = 0.5_rp * ( 1.0_rp + cos( pi * r(:) / schaer_r ) )
647 elsewhere
648 mwt_shape(:) = 0.0_rp
649 end where
650 end if
651
652 if ( quasi_2d_flag ) then
653 mwt_shape(:) = mwt_shape(:) * cos(lmesh2d%lat(:,ke2d))
654 end if
655
656 topo%local(n)%val(:,ke2d) = schaer_height * mwt_shape(:) * cos( pi * r(:) / schaer_lambda )**2
657 end do
658 deallocate( r, mwt_shape )
659 end do
660
661 return
662 end subroutine mktopo_schaer_type_global
663
664 !-----------------------------------------------------------------------------
670!OCL SERIAL
671 subroutine mktopo_barocwave_global_jw2006( mesh, topo )
672 use mod_mktopo_util, only: &
674 implicit none
675
676 class(meshbase2d), intent(in), target :: mesh
677 class(meshfield2d), intent(inout) :: topo
678
679 ! parameters of baroclinic wave test case in JW2006
680 real(rp) :: eta0 = 0.252_rp
681 real(rp) :: u0 = 35.e0_rp
682
683 namelist / param_mktopo_barocwave_global_jw2006 / &
684 eta0, u0
685
686 integer :: ierr
687 integer :: n
688 integer :: ke2d
689 type(localmesh2d), pointer :: lmesh2d
690 class(elementbase2d), pointer :: elem
691 !--------------------------------
692
693 log_info("MKTOPO_barocwave_JW2006",*) 'Setup'
694
695 !--- read namelist
696 rewind(io_fid_conf)
697 read(io_fid_conf, nml=param_mktopo_barocwave_global_jw2006, iostat=ierr)
698 if( ierr < 0 ) then !--- missing
699 log_info("MKTOPO_barocwave_JW2006",*) 'Not found namelist. Default used.'
700 elseif( ierr > 0 ) then !--- fatal error
701 log_error("MKTOPO_barocwave_JW2006",*) 'Not appropriate names in namelist PARAM_MKTOPO_BAROCWAVE_GLOBAL_JW2006. Check!'
702 call prc_abort
703 endif
704 log_nml(param_mktopo_barocwave_global_jw2006)
705
706 do n=1, mesh%LOCAL_MESH_NUM
707 lmesh2d => mesh%lcmesh_list(n)
708 elem => lmesh2d%refElem2D
709
710 !$omp parallel do
711 do ke2d=lmesh2d%NeS, lmesh2d%NeE
712 call calc_topo( topo%local(n)%val(:,ke2d), & ! (out)
713 u0, eta0, lmesh2d%lat(:,ke2d), elem%Np ) ! (in)
714 end do
715 end do
716
717 return
718 end subroutine mktopo_barocwave_global_jw2006
719
720end module mod_mktopo
module ATMOSPHERE component
module Atmosphere / Mesh
module Utility for mktopo
subroutine, public mktopoutil_barocwave_global_jw2006_calc_topo(topo, u0, eta0, lat, np)
Calculate a topography for a global test case of idealized baroclinic wave in Jablonowski and William...
module INITIAL
integer, parameter, public i_input_file
integer, parameter, public i_barocwave_global_jw2006
integer, parameter, public i_bellshape_global
subroutine, public mktopo_setup
Setup.
subroutine, public mktopo(output, model_mesh, topography)
Driver.
integer, parameter, public i_ignore
subroutine, public mktopo_write(model_mesh, topography)
Output topography data.
integer, parameter, public i_scaer
integer, parameter, public i_flat
integer, public mktopo_type
integer, parameter, public i_scaher_global
integer, parameter, public i_bellshape
module FElib / Element / Base
module FElib / Element / hexahedron
module FElib / Mesh / Local 2D
module FElib / Mesh / Local 3D
module FElib / Mesh / Base 2D
integer, public meshbase2d_dimtypeid_xy
module FElib / Mesh / Base 3D
module FElib / Mesh / Cubed-sphere 2D domain
module FElib / Mesh / Rectangle 2D domain
module FElib / Mesh / Topography
module FElib / Data / base
FElib / model framework / variable manager.