1 |
C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_get_gen.F,v 1.8 2005/01/12 20:33:13 heimbach Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CTRL_CPPOPTIONS.h" |
5 |
|
6 |
|
7 |
subroutine ctrl_get_gen( |
8 |
I xx_gen_file, xx_genstartdate, xx_genperiod, |
9 |
I genmask, genfld, xx_gen0, xx_gen1, xx_gen_dummy, |
10 |
I xx_gen_remo_intercept, xx_gen_remo_slope, |
11 |
I mytime, myiter, mythid |
12 |
& ) |
13 |
|
14 |
c ================================================================== |
15 |
c SUBROUTINE ctrl_get_gen |
16 |
c ================================================================== |
17 |
c |
18 |
c o new generic routine for reading time dependent control variables |
19 |
c heimbach@mit.edu 12-Jun-2003 |
20 |
c |
21 |
c ================================================================== |
22 |
c SUBROUTINE ctrl_get_gen |
23 |
c ================================================================== |
24 |
|
25 |
implicit none |
26 |
|
27 |
c == global variables == |
28 |
|
29 |
#include "EEPARAMS.h" |
30 |
#include "SIZE.h" |
31 |
#include "PARAMS.h" |
32 |
#include "GRID.h" |
33 |
|
34 |
#include "ctrl.h" |
35 |
#include "ctrl_dummy.h" |
36 |
#include "optim.h" |
37 |
#ifdef ALLOW_EXF |
38 |
# include "exf_fields.h" |
39 |
#endif |
40 |
|
41 |
c == routine arguments == |
42 |
|
43 |
character*(MAX_LEN_FNAM) xx_gen_file |
44 |
integer xx_genstartdate(4) |
45 |
_RL xx_genperiod |
46 |
_RL genmask(1-olx:snx+olx,1-oly:sny+oly,nr,nsx,nsy) |
47 |
_RL genfld(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy) |
48 |
_RL xx_gen0(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy) |
49 |
_RL xx_gen1(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy) |
50 |
_RL xx_gen_dummy |
51 |
_RL xx_gen_remo_intercept |
52 |
_RL xx_gen_remo_slope |
53 |
|
54 |
_RL mytime |
55 |
integer myiter |
56 |
integer mythid |
57 |
|
58 |
c == local variables == |
59 |
|
60 |
#ifdef ALLOW_EXF |
61 |
|
62 |
integer bi,bj |
63 |
integer i,j,k |
64 |
integer itlo,ithi |
65 |
integer jtlo,jthi |
66 |
integer jmin,jmax |
67 |
integer imin,imax |
68 |
integer ilgen |
69 |
|
70 |
_RL gensign |
71 |
_RL genfac |
72 |
logical doCtrlUpdate |
73 |
logical genfirst |
74 |
logical genchanged |
75 |
integer gencount0 |
76 |
integer gencount1 |
77 |
|
78 |
logical doglobalread |
79 |
logical ladinit |
80 |
|
81 |
character*(80) fnamegen |
82 |
|
83 |
c == external functions == |
84 |
|
85 |
integer ilnblnk |
86 |
external ilnblnk |
87 |
|
88 |
|
89 |
c == end of interface == |
90 |
|
91 |
jtlo = mybylo(mythid) |
92 |
jthi = mybyhi(mythid) |
93 |
itlo = mybxlo(mythid) |
94 |
ithi = mybxhi(mythid) |
95 |
jmin = 1-oly |
96 |
jmax = sny+oly |
97 |
imin = 1-olx |
98 |
imax = snx+olx |
99 |
|
100 |
c-- Now, read the control vector. |
101 |
doglobalread = .false. |
102 |
ladinit = .false. |
103 |
|
104 |
if (optimcycle .ge. 0) then |
105 |
ilgen=ilnblnk( xx_gen_file ) |
106 |
write(fnamegen(1:80),'(2a,i10.10)') |
107 |
& xx_gen_file(1:ilgen), '.', optimcycle |
108 |
endif |
109 |
|
110 |
c-- Get the counters, flags, and the interpolation factor. |
111 |
call ctrl_get_gen_rec( |
112 |
I xx_genstartdate, xx_genperiod, |
113 |
O genfac, genfirst, genchanged, |
114 |
O gencount0,gencount1, |
115 |
I mytime, myiter, mythid ) |
116 |
|
117 |
if ( genfirst ) then |
118 |
call active_read_xy_loc( fnamegen, xx_gen1, gencount0, |
119 |
& doglobalread, ladinit, optimcycle, |
120 |
& mythid, xx_gen_dummy ) |
121 |
#ifdef ALLOW_CTRL_SMOOTH |
122 |
if ( xx_gen_file .EQ. xx_tauu_file .OR. |
123 |
& xx_gen_file .EQ. xx_tauv_file ) |
124 |
& call ctrl_smooth(xx_gen1,genmask) |
125 |
#endif |
126 |
endif |
127 |
|
128 |
if (( genfirst ) .or. ( genchanged )) then |
129 |
call exf_SwapFFields( xx_gen0, xx_gen1, mythid ) |
130 |
|
131 |
call active_read_xy_loc( fnamegen, xx_gen1 , gencount1, |
132 |
& doglobalread, ladinit, optimcycle, |
133 |
& mythid, xx_gen_dummy ) |
134 |
#ifdef ALLOW_CTRL_SMOOTH |
135 |
if ( xx_gen_file .EQ. xx_tauu_file .OR. |
136 |
& xx_gen_file .EQ. xx_tauv_file ) |
137 |
& call ctrl_smooth(xx_gen1,genmask) |
138 |
#endif |
139 |
endif |
140 |
|
141 |
c-- Add control to model variable. |
142 |
cph( |
143 |
cph this flag ported from the SIO code |
144 |
cph Initial wind stress adjustments are too vigorous. |
145 |
if ( gencount0 .LE. 2 .AND. |
146 |
& ( xx_gen_file .EQ. xx_tauu_file .OR. |
147 |
& xx_gen_file .EQ. xx_tauv_file ) ) then |
148 |
doCtrlUpdate = .FALSE. |
149 |
else |
150 |
doCtrlUpdate = .TRUE. |
151 |
endif |
152 |
if ( xx_gen_file .EQ. xx_tauu_file .OR. |
153 |
& xx_gen_file .EQ. xx_tauv_file ) then |
154 |
gensign = -1. |
155 |
else |
156 |
gensign = 1. |
157 |
endif |
158 |
c |
159 |
cph since the above is ECCO specific, we undo it here: |
160 |
cph doCtrlUpdate = .TRUE. |
161 |
c |
162 |
if ( doCtrlUpdate ) then |
163 |
cph) |
164 |
do bj = jtlo,jthi |
165 |
do bi = itlo,ithi |
166 |
c-- Calculate mask for tracer cells (0 => land, 1 => water). |
167 |
k = 1 |
168 |
do j = 1,sny |
169 |
do i = 1,snx |
170 |
genfld(i,j,bi,bj) = genfld (i,j,bi,bj) |
171 |
& + gensign*genfac *xx_gen0(i,j,bi,bj) |
172 |
& + gensign*(1. _d 0 - genfac)*xx_gen1(i,j,bi,bj) |
173 |
genfld(i,j,bi,bj) = |
174 |
& genmask(i,j,k,bi,bj)*( genfld (i,j,bi,bj) - |
175 |
& ( xx_gen_remo_intercept + |
176 |
& xx_gen_remo_slope*(mytime-starttime) ) ) |
177 |
enddo |
178 |
enddo |
179 |
enddo |
180 |
enddo |
181 |
cph( |
182 |
endif |
183 |
cph) |
184 |
|
185 |
#endif /* ALLOW_EXF */ |
186 |
|
187 |
end |
188 |
|