ROMS
Loading...
Searching...
No Matches
set_2dfldr_mod Module Reference

Functions/Subroutines

subroutine set_2dfldr_tile (ng, tile, model, ifield, lbi, ubi, lbj, ubj, finp, fout, update, setbc)
 

Function/Subroutine Documentation

◆ set_2dfldr_tile()

subroutine set_2dfldr_mod::set_2dfldr_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model,
integer, intent(in) ifield,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
real(r8), dimension(lbi:,lbj:,:), intent(in) finp,
real(r8), dimension(lbi:,lbj:), intent(out) fout,
logical, intent(out) update,
logical, intent(in), optional setbc )

Definition at line 22 of file set_2dfldr.F.

26!***********************************************************************
27!
28 USE mod_param
29 USE mod_parallel
30 USE mod_iounits
31 USE mod_ncparam
32 USE mod_scalars
33!
35# ifdef DISTRIBUTE
37# endif
38!
39! Imported variable declarations.
40!
41 logical, intent(in), optional :: SetBC
42
43 logical, intent(out) :: update
44
45 integer, intent(in) :: ng, tile, model, ifield
46 integer, intent(in) :: LBi, UBi, LBj, UBj
47!
48# ifdef ASSUMED_SHAPE
49 real(r8), intent(in) :: Finp(LBi:,LBj:,:)
50 real(r8), intent(out) :: Fout(LBi:,LBj:)
51# else
52 real(r8), intent(in) :: Finp(LBi:UBi,LBj:UBj,2)
53 real(r8), intent(out) :: Fout(LBi:UBi,LBj:UBj)
54# endif
55!
56! Local variable declarations.
57!
58 logical :: LapplyBC, Lgrided, Lonerec
59
60 integer :: Tindex, gtype, i, it1, it2, j
61
62 real(dp) :: SecScale, fac, fac1, fac2
63 real(r8) :: Fval
64
65# include "set_bounds.h"
66!
67!----------------------------------------------------------------------
68! Set-up requested field for current tile.
69!----------------------------------------------------------------------
70!
71! Set switch to apply boundary conditions.
72!
73 IF (PRESENT(setbc)) THEN
74 lapplybc=setbc
75 ELSE
76 lapplybc=.true.
77 END IF
78!
79! Get requested field information from global storage.
80!
81 lgrided=linfo(1,ifield,ng)
82 lonerec=linfo(3,ifield,ng)
83 gtype =iinfo(1,ifield,ng)
84 tindex =iinfo(8,ifield,ng)
85 update=.true.
86!
87! Set linear, time interpolation factors. Fractional seconds are
88! rounded to the nearest milliseconds integer towards zero in the
89! time interpolation weights.
90!
91 secscale=1000.0_dp ! seconds to milliseconds
92 it1=3-tindex
93 it2=tindex
94 fac1=anint((time(ng)-tintrp(it2,ifield,ng))*secscale,dp)
95 fac2=anint((tintrp(it1,ifield,ng)-time(ng))*secscale,dp)
96!
97! Load time-invariant data. Time interpolation is not necessary.
98!
99 IF (lonerec) THEN
100 IF (lgrided) THEN
101 DO j=jstrr,jendr
102 DO i=istrr,iendr
103 fout(i,j)=finp(i,j,tindex)
104 END DO
105 END DO
106 ELSE
107 fval=fpoint(tindex,ifield,ng)
108 DO j=jstrr,jendr
109 DO i=istrr,iendr
110 fout(i,j)=fval
111 END DO
112 END DO
113 END IF
114!
115! Time-interpolate from gridded or point data.
116!
117 ELSE IF (((fac1*fac2).ge.0.0_dp).and. &
118 & ((fac1+fac2).gt.0.0_dp)) THEN
119 fac=1.0_dp/(fac1+fac2)
120 fac1=fac*fac1 ! nondimensional
121 fac2=fac*fac2 ! nondimensional
122 IF (lgrided) THEN
123 DO j=jstrr,jendr
124 DO i=istrr,iendr
125 fout(i,j)=fac1*finp(i,j,it1)+fac2*finp(i,j,it2)
126 END DO
127 END DO
128 ELSE
129 fval=fac1*fpoint(it1,ifield,ng)+fac2*fpoint(it2,ifield,ng)
130 DO j=jstrr,jendr
131 DO i=istrr,iendr
132 fout(i,j)=fval
133 END DO
134 END DO
135 END IF
136!
137! Activate synchronization flag if a new time record needs to be
138! read in at the next time step.
139!
140 IF ((time(ng)-dt(ng)).lt.tintrp(it2,ifield,ng)) THEN
141 IF (domain(ng)%SouthWest_Test(tile)) synchro_flag(ng)=.true.
142 END IF
143!
144! Unable to set-up requested field. Activate error flag to quit.
145!
146 ELSE
147 IF (domain(ng)%SouthWest_Test(tile)) THEN
148 IF (master) THEN
149 WRITE (stdout,10) trim(vname(1,ifield)), tdays(ng), &
150 & finfo(1,ifield,ng), finfo(2,ifield,ng), &
151 & finfo(3,ifield,ng), finfo(4,ifield,ng), &
152 & tintrp(it1,ifield,ng)*sec2day, &
153 & tintrp(it2,ifield,ng)*sec2day, &
154 & fac1*sec2day/secscale, &
155 & fac2*sec2day/secscale
156 END IF
157 10 FORMAT (/,' SET_2DFLDR - current model time', &
158 & ' exceeds ending value for variable: ',a, &
159 & /,14x,'TDAYS = ',f15.4, &
160 & /,14x,'Data Tmin = ',f15.4,2x,'Data Tmax = ',f15.4, &
161 & /,14x,'Data Tstr = ',f15.4,2x,'Data Tend = ',f15.4, &
162 & /,14x,'TINTRP1 = ',f15.4,2x,'TINTRP2 = ',f15.4, &
163 & /,14x,'FAC1 = ',f15.4,2x,'FAC2 = ',f15.4)
164 exit_flag=2
165 update=.false.
166 END IF
167 END IF
168!
169! Exchange boundary data.
170!
171 IF (update) THEN
172 IF (lapplybc.and.(ewperiodic(ng).or.nsperiodic(ng))) THEN
173 IF (gtype.eq.r2dvar) THEN
174 CALL exchange_r2d_tile (ng, tile, &
175 & lbi, ubi, lbj, ubj, &
176 & fout)
177 ELSE IF (gtype.eq.u2dvar) THEN
178 CALL exchange_u2d_tile (ng, tile, &
179 & lbi, ubi, lbj, ubj, &
180 & fout)
181 ELSE IF (gtype.eq.v2dvar) THEN
182 CALL exchange_v2d_tile (ng, tile, &
183 & lbi, ubi, lbj, ubj, &
184 & fout)
185 END IF
186 END IF
187
188# ifdef DISTRIBUTE
189 IF (.not.lapplybc) THEN
190 CALL mp_exchange2d (ng, tile, model, 1, &
191 & lbi, ubi, lbj, ubj, &
192 & nghostpoints, &
193 & .false., .false., &
194 & fout)
195 ELSE
196 CALL mp_exchange2d (ng, tile, model, 1, &
197 & lbi, ubi, lbj, ubj, &
198 & nghostpoints, &
199 & ewperiodic(ng), nsperiodic(ng), &
200 & fout)
201 END IF
202# endif
203 END IF
204
205 RETURN
subroutine exchange_r2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
subroutine exchange_u2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
subroutine exchange_v2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
integer stdout
logical, dimension(:,:,:), allocatable linfo
real(dp), dimension(:,:,:), allocatable tintrp
real(dp), dimension(:,:,:), allocatable fpoint
real(dp), dimension(:,:,:), allocatable finfo
character(len=maxlen), dimension(6, 0:nv) vname
integer, dimension(:,:,:), allocatable iinfo
logical master
integer nghostpoints
Definition mod_param.F:710
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, parameter u2dvar
Definition mod_param.F:718
integer, parameter r2dvar
Definition mod_param.F:717
integer, parameter v2dvar
Definition mod_param.F:719
real(dp), dimension(:), allocatable dt
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
logical, dimension(:), allocatable synchro_flag
real(dp), dimension(:), allocatable tdays
real(dp), parameter sec2day
integer exit_flag
real(dp), dimension(:), allocatable time
subroutine mp_exchange2d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, nghost, ew_periodic, ns_periodic, a, b, c, d)

References mod_param::domain, mod_scalars::dt, mod_scalars::ewperiodic, exchange_2d_mod::exchange_r2d_tile(), exchange_2d_mod::exchange_u2d_tile(), exchange_2d_mod::exchange_v2d_tile(), mod_scalars::exit_flag, mod_ncparam::finfo, mod_ncparam::fpoint, mod_ncparam::iinfo, mod_ncparam::linfo, mod_parallel::master, mp_exchange_mod::mp_exchange2d(), mod_param::nghostpoints, mod_scalars::nsperiodic, mod_param::r2dvar, mod_scalars::sec2day, mod_iounits::stdout, mod_scalars::synchro_flag, mod_scalars::tdays, mod_scalars::time, mod_ncparam::tintrp, mod_param::u2dvar, mod_param::v2dvar, and mod_ncparam::vname.

Referenced by ad_set_data_tile().

Here is the call graph for this function:
Here is the caller graph for this function: