ROMS
Loading...
Searching...
No Matches
ana_rain.h
Go to the documentation of this file.
1!!
2 SUBROUTINE ana_rain (ng, tile, model)
3!
4!! git $Id$
5!!======================================================================
6!! Copyright (c) 2002-2025 The ROMS Group !
7!! Licensed under a MIT/X style license !
8!! See License_ROMS.md !
9!=======================================================================
10! !
11! This routine sets precipitation rate (kg/m2/s) using an !
12! analytical expression. !
13! !
14!=======================================================================
15!
16 USE mod_param
17 USE mod_forces
18 USE mod_ncparam
19!
20! Imported variable declarations.
21!
22 integer, intent(in) :: ng, tile, model
23!
24! Local variable declarations.
25!
26 character (len=*), parameter :: MyFile = &
27 & __FILE__
28!
29#include "tile.h"
30!
31 CALL ana_rain_tile (ng, tile, model, &
32 & lbi, ubi, lbj, ubj, &
33 & imins, imaxs, jmins, jmaxs, &
34 & forces(ng) % rain)
35!
36! Set analytical header file name used.
37!
38#ifdef DISTRIBUTE
39 IF (lanafile) THEN
40#else
41 IF (lanafile.and.(tile.eq.0)) THEN
42#endif
43 ananame(21)=myfile
44 END IF
45!
46 RETURN
47 END SUBROUTINE ana_rain
48!
49!***********************************************************************
50 SUBROUTINE ana_rain_tile (ng, tile, model, &
51 & LBi, UBi, LBj, UBj, &
52 & IminS, ImaxS, JminS, JmaxS, &
53 & rain)
54!***********************************************************************
55!
56 USE mod_param
57 USE mod_ncparam
58 USE mod_scalars
59!
61#ifdef DISTRIBUTE
63#endif
64!
65! Imported variable declarations.
66!
67 integer, intent(in) :: ng, tile, model
68 integer, intent(in) :: LBi, UBi, LBj, UBj
69 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
70!
71#ifdef ASSUMED_SHAPE
72 real(r8), intent(out) :: rain(LBi:,LBj:)
73#else
74 real(r8), intent(out) :: rain(LBi:UBi,LBj:UBj)
75#endif
76!
77! Local variable declarations.
78!
79 integer :: i, j
80
81#include "set_bounds.h"
82!
83!-----------------------------------------------------------------------
84! Set analytical precipitation rate (kg/m2/s).
85!-----------------------------------------------------------------------
86!
87 DO j=jstrt,jendt
88 DO i=istrt,iendt
89 rain(i,j)=0.0_r8
90 END DO
91 END DO
92!
93! Exchange boundary data.
94!
95 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
96 CALL exchange_r2d_tile (ng, tile, &
97 & lbi, ubi, lbj, ubj, &
98 & rain)
99 END IF
100
101#ifdef DISTRIBUTE
102 CALL mp_exchange2d (ng, tile, model, 1, &
103 & lbi, ubi, lbj, ubj, &
104 & nghostpoints, &
105 & ewperiodic(ng), nsperiodic(ng), &
106 & rain)
107#endif
108!
109 RETURN
110 END SUBROUTINE ana_rain_tile
subroutine ana_rain_tile(ng, tile, model, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, rain)
Definition ana_rain.h:54
subroutine ana_rain(ng, tile, model)
Definition ana_rain.h:3
subroutine exchange_r2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
type(t_forces), dimension(:), allocatable forces
Definition mod_forces.F:554
logical lanafile
character(len=256), dimension(39) ananame
integer nghostpoints
Definition mod_param.F:710
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
subroutine mp_exchange2d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, nghost, ew_periodic, ns_periodic, a, b, c, d)