/[MITgcm]/MITgcm/pkg/ctrl/ctrl_getobcsw.F
ViewVC logotype

Annotation of /MITgcm/pkg/ctrl/ctrl_getobcsw.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.1.2.1 - (hide annotations) (download)
Tue Feb 5 20:23:58 2002 UTC (23 years, 5 months ago) by heimbach
Branch: ecco-branch
CVS Tags: ecco_c44_e19, ecco_c44_e18, ecco_c44_e17, ecco_c44_e16, ecco_c44_e22, ecco_c44_e23, ecco_c44_e20, ecco_c44_e21, ecco-branch-mod1, ecco-branch-mod2, ecco-branch-mod3, ecco-branch-mod4, ecco-branch-mod5
Changes since 1.1: +202 -0 lines
Starting from ecco-branch, replacing packages
cost, ctrl, ecco, obcs by ECCO packages.
Will create tag ecco-branch-mod1 after this modif.

1 heimbach 1.1.2.1
2     #include "CTRL_CPPOPTIONS.h"
3     #ifdef ALLOW_OBCS
4     # include "OBCS_OPTIONS.h"
5     #endif
6    
7    
8     subroutine ctrl_getobcsw(
9     I mytime,
10     I myiter,
11     I mythid
12     & )
13    
14     c ==================================================================
15     c SUBROUTINE ctrl_getobcsw
16     c ==================================================================
17     c
18     c o Get norhtern obc of the control vector and add it
19     c to dyn. fields
20     c
21     c started: heimbach@mit.edu, 29-Aug-2001
22     c
23     c ==================================================================
24     c SUBROUTINE ctrl_getobcsw
25     c ==================================================================
26    
27     implicit none
28    
29     #ifdef ALLOW_OBCSW_CONTROL
30    
31     c == global variables ==
32    
33     #include "EEPARAMS.h"
34     #include "SIZE.h"
35     #include "PARAMS.h"
36     #include "GRID.h"
37     #include "OBCS.h"
38    
39     #include "ctrl.h"
40     #include "ctrl_dummy.h"
41     #include "optim.h"
42    
43     c == routine arguments ==
44    
45     _RL mytime
46     integer myiter
47     integer mythid
48    
49     c == local variables ==
50    
51     integer bi,bj
52     integer i,j,k
53     integer itlo,ithi
54     integer jtlo,jthi
55     integer jmin,jmax
56     integer imin,imax
57     integer ilobcsw
58     integer iobcs
59    
60     _RL dummy
61     _RL obcswfac
62     logical obcswfirst
63     logical obcswchanged
64     integer obcswcount0
65     integer obcswcount1
66    
67     _RL maskyz (1-oly:sny+oly,nr,nsx,nsy)
68    
69     logical doglobalread
70     logical ladinit
71    
72     character*(80) fnameobcsw
73    
74    
75    
76     c == external functions ==
77    
78     integer ilnblnk
79     external ilnblnk
80    
81    
82     c == end of interface ==
83    
84     jtlo = mybylo(mythid)
85     jthi = mybyhi(mythid)
86     itlo = mybxlo(mythid)
87     ithi = mybxhi(mythid)
88     jmin = 1-oly
89     jmax = sny+oly
90     imin = 1-olx
91     imax = snx+olx
92    
93     c-- Now, read the control vector.
94     doglobalread = .false.
95     ladinit = .false.
96    
97     if (optimcycle .ge. 0) then
98     ilobcsw=ilnblnk( xx_obcsw_file )
99     write(fnameobcsw(1:80),'(2a,i10.10)')
100     & xx_obcsw_file(1:ilobcsw), '.', optimcycle
101     else
102     print*
103     print*,' ctrl_getobcsw: optimcycle not set correctly.'
104     print*,' ... stopped in ctrl_getobcsw.'
105     endif
106    
107     c-- Get the counters, flags, and the interpolation factor.
108     call ctrl_GetRec( 'xx_obcsw',
109     O obcswfac, obcswfirst, obcswchanged,
110     O obcswcount0,obcswcount1,
111     I mytime, myiter, mythid )
112    
113     do iobcs = 1,nobcs
114    
115     call active_read_yz( 'maskobcsw', maskyz,
116     & iobcs,
117     & doglobalread, ladinit, 0,
118     & mythid, dummy )
119    
120     call active_read_yz( fnameobcsw, tmpfldyz,
121     & (obcswcount0-1)*nobcs+iobcs,
122     & doglobalread, ladinit, optimcycle,
123     & mythid, xx_obcsw_dummy )
124    
125     do bj = jtlo,jthi
126     do bi = itlo,ithi
127     do k = 1,nr
128     do j = jmin,jmax
129     yz_obcs0(j,k,bi,bj) = tmpfldyz (j,k,bi,bj)
130     enddo
131     enddo
132     enddo
133     enddo
134    
135     call active_read_yz( fnameobcsw, tmpfldyz,
136     & (obcswcount1-1)*nobcs+iobcs,
137     & doglobalread, ladinit, optimcycle,
138     & mythid, xx_obcsw_dummy )
139    
140     do bj = jtlo,jthi
141     do bi = itlo,ithi
142     do k = 1,nr
143     do j = jmin,jmax
144     yz_obcs1 (j,k,bi,bj) = tmpfldyz (j,k,bi,bj)
145     enddo
146     enddo
147     enddo
148     enddo
149    
150     c-- Add control to model variable.
151     do bj = jtlo,jthi
152     do bi = itlo,ithi
153     c-- Calculate mask for tracer cells (0 => land, 1 => water).
154     do k = 1,nr
155     do j = 1,sny
156     if (iobcs .EQ. 1) then
157     OBWt(j,k,bi,bj) = OBWt (j,k,bi,bj)
158     & + obcswfac *yz_obcs0(j,k,bi,bj)
159     & + (1. _d 0 - obcswfac)*yz_obcs1(j,k,bi,bj)
160     OBWt(j,k,bi,bj) = OBWt(j,k,bi,bj)
161     & *maskyz(j,k,bi,bj)
162     else if (iobcs .EQ. 2) then
163     OBWs(j,k,bi,bj) = OBWs (j,k,bi,bj)
164     & + obcswfac *yz_obcs0(j,k,bi,bj)
165     & + (1. _d 0 - obcswfac)*yz_obcs1(j,k,bi,bj)
166     OBWs(j,k,bi,bj) = OBWs(j,k,bi,bj)
167     & *maskyz(j,k,bi,bj)
168     else if (iobcs .EQ. 3) then
169     OBWu(j,k,bi,bj) = OBWu (j,k,bi,bj)
170     & + obcswfac *yz_obcs0(j,k,bi,bj)
171     & + (1. _d 0 - obcswfac)*yz_obcs1(j,k,bi,bj)
172     OBWu(j,k,bi,bj) = OBWu(j,k,bi,bj)
173     & *maskyz(j,k,bi,bj)
174     else if (iobcs .EQ. 4) then
175     OBWv(j,k,bi,bj) = OBWv (j,k,bi,bj)
176     & + obcswfac *yz_obcs0(j,k,bi,bj)
177     & + (1. _d 0 - obcswfac)*yz_obcs1(j,k,bi,bj)
178     OBWv(j,k,bi,bj) = OBWv(j,k,bi,bj)
179     & *maskyz(j,k,bi,bj)
180     endif
181     enddo
182     enddo
183     enddo
184     enddo
185    
186     C-- End over iobcs loop
187     enddo
188    
189     #else /* ALLOW_OBCSW_CONTROL undefined */
190    
191     c == routine arguments ==
192    
193     _RL mytime
194     integer myiter
195     integer mythid
196    
197     c-- CPP flag ALLOW_OBCSW_CONTROL undefined.
198    
199     #endif /* ALLOW_OBCSW_CONTROL */
200    
201     end
202    

  ViewVC Help
Powered by ViewVC 1.1.22