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