forked from altMITgcm/MITgcm
-
Notifications
You must be signed in to change notification settings - Fork 0
/
streamice_velmask_upd.F
234 lines (205 loc) · 7.38 KB
/
streamice_velmask_upd.F
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
C $Header: /u/gcmpack/MITgcm/pkg/streamice/streamice_velmask_upd.F,v 1.7 2014/04/30 12:18:32 dgoldberg Exp $
C $Name: $
#include "STREAMICE_OPTIONS.h"
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
CBOP
SUBROUTINE STREAMICE_VELMASK_UPD ( myThid )
C /============================================================\
C | SUBROUTINE |
C | o |
C |============================================================|
C | |
C \============================================================/
IMPLICIT NONE
C === Global variables ===
#include "SIZE.h"
#include "GRID.h"
#include "EEPARAMS.h"
#include "PARAMS.h"
#include "STREAMICE.h"
!#ifdef ALLOW_PETSC
!# ifdef ALLOW_USE_MPI
!# include "EESUPPORT.h"
!# endif
!#endif
! #include "STREAMICE_ADV.h"
#ifdef ALLOW_AUTODIFF_TAMC
# include "tamc.h"
#endif
INTEGER myThid
#ifdef ALLOW_STREAMICE
INTEGER i, j, bi, bj, ki, kj
INTEGER maskFlag, myThidCopy
CHARACTER*(MAX_LEN_MBUF) msgBuf
#ifdef ALLOW_USE_MPI
integer mpiRC, mpiMyWid
#endif
#ifdef ALLOW_PETSC
_RS DoFCount
integer n_dofs_proc_loc (0:nPx*nPy-1)
integer n_dofs_cum_sum (0:nPx*nPy-1)
#endif
_EXCH_XY_RL( H_streamice, myThid )
_EXCH_XY_RL( area_shelf_streamice, myThid )
_EXCH_XY_RS( STREAMICE_hmask, myThid )
DO bj=myByLo(myThid),myByHi(myThid)
DO bi=myBxLo(myThid),myBxHi(myThid)
DO j=1-OLy,sNy+OLy
DO i=1-OLx,sNx+OLx
STREAMICE_umask(i,j,bi,bj) = -1. _d 0
STREAMICE_vmask(i,j,bi,bj) = -1. _d 0
STREAMICE_ufacemask(i,j,bi,bj) = 0. _d 0
STREAMICE_vfacemask(i,j,bi,bj) = 0. _d 0
ENDDO
ENDDO
ENDDO
ENDDO
DO bj=myByLo(myThid),myByHi(myThid)
DO bi=myBxLo(myThid),myBxHi(myThid)
DO j=0,sNy+1
DO i=0,sNx+1
IF (STREAMICE_hmask(i,j,bi,bj) .eq. 1.0) THEN
DO kj=0,1
DO ki=0,1
if (STREAMICE_umask(i+ki,j+kj,bi,bj).eq.-1.0) then
STREAMICE_umask (i+ki,j+kj,bi,bj) = 1.0
endif
if (STREAMICE_vmask(i+ki,j+kj,bi,bj).eq.-1.0) then
STREAMICE_vmask (i+ki,j+kj,bi,bj) = 1.0
endif
ENDDO
ENDDO
DO ki=0,1
maskFlag=INT(STREAMICE_ufacemask_bdry(i+ki,j,bi,bj))
IF (maskFlag.EQ.3) THEN
DO kj=0,1
if (STREAMICE_umask(i+ki,j+kj,bi,bj).ne.0.0) then
STREAMICE_umask(i+ki,j+kj,bi,bj) = 3.0
endif
if(STREAMICE_vmask(i+ki,j+kj,bi,bj).ne.0.0) then
STREAMICE_vmask(i+ki,j+kj,bi,bj) = 3.0
endif
ENDDO
STREAMICE_ufacemask(i+ki,j,bi,bj) = 3.0
ELSE IF (maskFlag.EQ.2) THEN
!DO kj=0,1
STREAMICE_ufacemask(i+ki,j,bi,bj) = 2.0
!ENDDO
ELSE IF (maskFlag.EQ.4) THEN
DO kj=0,1
STREAMICE_umask(i+ki,j+kj,bi,bj) = 0.0
STREAMICE_vmask(i+ki,j+kj,bi,bj) = 0.0
ENDDO
STREAMICE_ufacemask(i+ki,j,bi,bj) = 4.0
ELSE IF (maskFlag.EQ.0) THEN
DO kj=0,1
STREAMICE_umask(i+ki,j+kj,bi,bj) = 0.0
STREAMICE_vmask(i+ki,j+kj,bi,bj) = 0.0
ENDDO
STREAMICE_ufacemask(i+ki,j,bi,bj) = 0.0
ELSE IF (maskFlag.EQ.1) THEN
DO kj=0,1
STREAMICE_umask(i+ki,j+kj,bi,bj) = 0.0
ENDDO
END IF
ENDDO
DO kj=0,1
maskFlag=INT(STREAMICE_vfacemask_bdry(i,j+kj,bi,bj))
IF (maskFlag.EQ.3) THEN
DO ki=0,1
if(STREAMICE_vmask(i+ki,j+kj,bi,bj).ne.0.0) then
STREAMICE_vmask(i+ki,j+kj,bi,bj) = 3.0
endif
if(STREAMICE_umask(i+ki,j+kj,bi,bj).ne.0.0) then
STREAMICE_umask(i+ki,j+kj,bi,bj) = 3.0
endif
ENDDO
STREAMICE_vfacemask(i,j+kj,bi,bj) = 3.0
ELSE IF (maskFlag.EQ.2) THEN
!DO ki=0,1
STREAMICE_vfacemask(i,j+kj,bi,bj) = 2.0
!ENDDO
ELSE IF (maskFlag.EQ.4) THEN
DO ki=0,1
STREAMICE_umask(i+ki,j+kj,bi,bj) = 0.0
STREAMICE_vmask(i+ki,j+kj,bi,bj) = 0.0
ENDDO
STREAMICE_vfacemask(i,j+kj,bi,bj) = 4.0
ELSE IF (maskFlag.EQ.0) THEN
DO ki=0,1
STREAMICE_umask(i+ki,j+kj,bi,bj) = 0.0
STREAMICE_vmask(i+ki,j+kj,bi,bj) = 0.0
ENDDO
STREAMICE_vfacemask(i+ki,j,bi,bj) = 0.0
ELSE IF (maskFlag.EQ.1) THEN
DO ki=0,1
STREAMICE_vmask(i+ki,j+kj,bi,bj) = 0.0
ENDDO
ENDIF
ENDDO
IF (i .lt. sNx+OLx) THEN
IF ((STREAMICE_hmask(i+1,j,bi,bj) .eq. 0.0) .OR.
& (STREAMICE_hmask(i+1,j,bi,bj) .eq. 2.0)) THEN
!right boundary or adjacent to unfilled cell
STREAMICE_ufacemask(i+1,j,bi,bj) = 2.0
ENDIF
ENDIF
IF (i .gt. 1-OLx) THEN
IF ((STREAMICE_hmask(i-1,j,bi,bj) .eq. 0.0) .OR.
& (STREAMICE_hmask(i-1,j,bi,bj) .eq. 2.0)) THEN
!left boundary or adjacent to unfilled cell
STREAMICE_ufacemask(i,j,bi,bj) = 2
ENDIF
ENDIF
IF (j .lt. sNy+OLy) THEN
IF ((STREAMICE_hmask(i,j+1,bi,bj) .eq. 0.0) .OR.
& (STREAMICE_hmask(i,j+1,bi,bj) .eq. 2.0)) THEN
!top boundary or adjacent to unfilled cell
STREAMICE_vfacemask(i,j+1,bi,bj) = 2
ENDIF
ENDIF
IF (j .gt. 1-OLy) THEN
IF ((STREAMICE_hmask(i,j-1,bi,bj) .eq. 0.0) .OR.
& (STREAMICE_hmask(i,j-1,bi,bj) .eq. 2.0)) THEN
!bot boundary or adjacent to unfilled cell
STREAMICE_vfacemask(i,j,bi,bj) = 2.0
ENDIF
ENDIF
ENDIF
ENDDO
ENDDO
ENDDO
ENDDO
!$TAF STORE streamice_umask = comlev1, key=ikey_dynamics
!$TAF STORE streamice_vmask = comlev1, key=ikey_dynamics
DO bj=myByLo(myThid),myByHi(myThid)
DO bi=myBxLo(myThid),myBxHi(myThid)
DO j=1,sNy
DO i=1,sNx
IF(streamice_umask(i,j,bi,bj).eq.-1.0) THEN
streamice_umask(i,j,bi,bj)=0.0
ENDIF
IF(streamice_vmask(i,j,bi,bj).eq.-1.0) THEN
streamice_vmask(i,j,bi,bj)=0.0
ENDIF
ENDDO
ENDDO
ENDDO
ENDDO
_EXCH_XY_RS( STREAMICE_ufacemask, myThid )
_EXCH_XY_RS( STREAMICE_vfacemask, myThid )
_EXCH_XY_RS( STREAMICE_umask, myThid )
_EXCH_XY_RS( STREAMICE_vmask, myThid )
! CALL WRITE_FULLARRAY_RL ("umask",STREAMICE_umask,
! c 1,0,0,1,0,myThid)
! CALL WRITE_FLD_XY_RL ("umask","",STREAMICE_umask,0,myThid)
! CALL WRITE_FLD_XY_RL ("vmask","",STREAMICE_vmask,0,myThid)
! CALL WRITE_FLD_XY_RL ("ufacemask","",STREAMICE_ufacemask,0,myThid)
! CALL WRITE_FLD_XY_RL ("vfacemask","",STREAMICE_vfacemask,0,myThid)
#ifdef ALLOW_PETSC
myThidCopy = myThid
call streamice_petsc_numerate (myThidCopy)
#endif
#endif
RETURN
END