| 1 |
C $Header: /u/gcmpack/MITgcm_contrib/dgoldberg/streamice/streamice_velmask_upd.F,v 1.4 2013/04/06 17:43:41 dgoldberg Exp $ |
| 2 |
C $Name: $ |
| 3 |
|
| 4 |
#include "STREAMICE_OPTIONS.h" |
| 5 |
|
| 6 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
| 7 |
|
| 8 |
CBOP |
| 9 |
SUBROUTINE STREAMICE_VELMASK_UPD ( myThid ) |
| 10 |
|
| 11 |
C /============================================================\ |
| 12 |
C | SUBROUTINE | |
| 13 |
C | o | |
| 14 |
C |============================================================| |
| 15 |
C | | |
| 16 |
C \============================================================/ |
| 17 |
IMPLICIT NONE |
| 18 |
|
| 19 |
C === Global variables === |
| 20 |
#include "SIZE.h" |
| 21 |
#include "GRID.h" |
| 22 |
#include "EEPARAMS.h" |
| 23 |
#include "PARAMS.h" |
| 24 |
#include "STREAMICE.h" |
| 25 |
#ifdef ALLOW_USE_MPI |
| 26 |
# include "EESUPPORT.h" |
| 27 |
#endif |
| 28 |
! #include "STREAMICE_ADV.h" |
| 29 |
|
| 30 |
INTEGER myThid |
| 31 |
|
| 32 |
#ifdef ALLOW_STREAMICE |
| 33 |
|
| 34 |
INTEGER i, j, bi, bj, ki, kj |
| 35 |
INTEGER maskFlag |
| 36 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
| 37 |
#ifdef ALLOW_USE_MPI |
| 38 |
integer mpiRC, mpiMyWid |
| 39 |
#endif |
| 40 |
#ifdef ALLOW_PETSC |
| 41 |
_RS DoFCount |
| 42 |
integer n_dofs_proc_loc (0:nPx*nPy-1) |
| 43 |
integer n_dofs_cum_sum (0:nPx*nPy-1) |
| 44 |
#endif |
| 45 |
|
| 46 |
_EXCH_XY_RL( H_streamice, myThid ) |
| 47 |
_EXCH_XY_RL( area_shelf_streamice, myThid ) |
| 48 |
_EXCH_XY_RL( STREAMICE_hmask, myThid ) |
| 49 |
|
| 50 |
DO bj=myByLo(myThid),myByHi(myThid) |
| 51 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
| 52 |
DO j=1-OLy,sNy+OLy |
| 53 |
DO i=1-OLx,sNx+OLx |
| 54 |
STREAMICE_umask(i,j,bi,bj) = 0. _d 0 |
| 55 |
STREAMICE_vmask(i,j,bi,bj) = 0. _d 0 |
| 56 |
STREAMICE_ufacemask(i,j,bi,bj) = 0. _d 0 |
| 57 |
STREAMICE_vfacemask(i,j,bi,bj) = 0. _d 0 |
| 58 |
ENDDO |
| 59 |
ENDDO |
| 60 |
ENDDO |
| 61 |
ENDDO |
| 62 |
|
| 63 |
DO bj=myByLo(myThid),myByHi(myThid) |
| 64 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
| 65 |
DO j=0,sNy+1 |
| 66 |
DO i=0,sNx+1 |
| 67 |
IF (STREAMICE_hmask(i,j,bi,bj) .eq. 1) THEN |
| 68 |
|
| 69 |
DO kj=0,1 |
| 70 |
DO ki=0,1 |
| 71 |
STREAMICE_umask (i+ki,j+kj,bi,bj) = 1.0 |
| 72 |
STREAMICE_vmask (i+ki,j+kj,bi,bj) = 1.0 |
| 73 |
ENDDO |
| 74 |
ENDDO |
| 75 |
|
| 76 |
DO ki=0,1 |
| 77 |
maskFlag=INT(STREAMICE_ufacemask_bdry(i+ki,j,bi,bj)) |
| 78 |
IF (maskFlag.EQ.3) THEN |
| 79 |
DO kj=0,1 |
| 80 |
STREAMICE_umask(i+ki,j+kj,bi,bj) = 3.0 |
| 81 |
STREAMICE_vmask(i+ki,j+kj,bi,bj) = 3.0 |
| 82 |
ENDDO |
| 83 |
STREAMICE_ufacemask(i+ki,j,bi,bj) = 3.0 |
| 84 |
ELSE IF (maskFlag.EQ.2) THEN |
| 85 |
!DO kj=0,1 |
| 86 |
STREAMICE_ufacemask(i+ki,j,bi,bj) = 2.0 |
| 87 |
!ENDDO |
| 88 |
ELSE IF (maskFlag.EQ.4) THEN |
| 89 |
DO kj=0,1 |
| 90 |
STREAMICE_umask(i+ki,j+kj,bi,bj) = 0.0 |
| 91 |
STREAMICE_vmask(i+ki,j+kj,bi,bj) = 0.0 |
| 92 |
ENDDO |
| 93 |
STREAMICE_ufacemask(i+ki,j,bi,bj) = 4.0 |
| 94 |
ELSE IF (maskFlag.EQ.0) THEN |
| 95 |
DO kj=0,1 |
| 96 |
STREAMICE_umask(i+ki,j+kj,bi,bj) = 0.0 |
| 97 |
STREAMICE_vmask(i+ki,j+kj,bi,bj) = 0.0 |
| 98 |
ENDDO |
| 99 |
STREAMICE_ufacemask(i+ki,j,bi,bj) = 0.0 |
| 100 |
ELSE IF (maskFlag.EQ.1) THEN |
| 101 |
DO kj=0,1 |
| 102 |
STREAMICE_umask(i+ki,j+kj,bi,bj) = 0.0 |
| 103 |
ENDDO |
| 104 |
END IF |
| 105 |
ENDDO |
| 106 |
|
| 107 |
DO kj=0,1 |
| 108 |
maskFlag=INT(STREAMICE_vfacemask_bdry(i,j+kj,bi,bj)) |
| 109 |
IF (maskFlag.EQ.3) THEN |
| 110 |
DO ki=0,1 |
| 111 |
STREAMICE_vmask(i+ki,j+kj,bi,bj) = 3.0 |
| 112 |
STREAMICE_umask(i+ki,j+kj,bi,bj) = 3.0 |
| 113 |
ENDDO |
| 114 |
STREAMICE_vfacemask(i,j+kj,bi,bj) = 3.0 |
| 115 |
ELSE IF (maskFlag.EQ.2) THEN |
| 116 |
!DO ki=0,1 |
| 117 |
STREAMICE_vfacemask(i,j+kj,bi,bj) = 2.0 |
| 118 |
!ENDDO |
| 119 |
ELSE IF (maskFlag.EQ.4) THEN |
| 120 |
DO ki=0,1 |
| 121 |
STREAMICE_umask(i+ki,j+kj,bi,bj) = 0.0 |
| 122 |
STREAMICE_vmask(i+ki,j+kj,bi,bj) = 0.0 |
| 123 |
ENDDO |
| 124 |
STREAMICE_vfacemask(i,j+kj,bi,bj) = 4.0 |
| 125 |
ELSE IF (maskFlag.EQ.0) THEN |
| 126 |
DO ki=0,1 |
| 127 |
STREAMICE_umask(i+ki,j+kj,bi,bj) = 0.0 |
| 128 |
STREAMICE_vmask(i+ki,j+kj,bi,bj) = 0.0 |
| 129 |
ENDDO |
| 130 |
STREAMICE_vfacemask(i+ki,j,bi,bj) = 0.0 |
| 131 |
ELSE IF (maskFlag.EQ.1) THEN |
| 132 |
DO ki=0,1 |
| 133 |
STREAMICE_vmask(i+ki,j+kj,bi,bj) = 0.0 |
| 134 |
ENDDO |
| 135 |
ENDIF |
| 136 |
ENDDO |
| 137 |
|
| 138 |
IF (i .lt. sNx+OLx) THEN |
| 139 |
IF ((STREAMICE_hmask(i+1,j,bi,bj) .eq. 0.0) .OR. |
| 140 |
& (STREAMICE_hmask(i+1,j,bi,bj) .eq. 2.0)) THEN |
| 141 |
!right boundary or adjacent to unfilled cell |
| 142 |
STREAMICE_ufacemask(i+1,j,bi,bj) = 2.0 |
| 143 |
ENDIF |
| 144 |
ENDIF |
| 145 |
|
| 146 |
IF (i .gt. 1-OLx) THEN |
| 147 |
IF ((STREAMICE_hmask(i-1,j,bi,bj) .eq. 0.0) .OR. |
| 148 |
& (STREAMICE_hmask(i-1,j,bi,bj) .eq. 2.0)) THEN |
| 149 |
!left boundary or adjacent to unfilled cell |
| 150 |
STREAMICE_ufacemask(i,j,bi,bj) = 2 |
| 151 |
ENDIF |
| 152 |
ENDIF |
| 153 |
|
| 154 |
IF (j .lt. sNy+OLy) THEN |
| 155 |
IF ((STREAMICE_hmask(i,j+1,bi,bj) .eq. 0.0) .OR. |
| 156 |
& (STREAMICE_hmask(i,j+1,bi,bj) .eq. 2.0)) THEN |
| 157 |
!top boundary or adjacent to unfilled cell |
| 158 |
STREAMICE_vfacemask(i,j+1,bi,bj) = 2 |
| 159 |
ENDIF |
| 160 |
ENDIF |
| 161 |
|
| 162 |
IF (j .gt. 1-OLy) THEN |
| 163 |
IF ((STREAMICE_hmask(i,j-1,bi,bj) .eq. 0.0) .OR. |
| 164 |
& (STREAMICE_hmask(i,j-1,bi,bj) .eq. 2.0)) THEN |
| 165 |
!bot boundary or adjacent to unfilled cell |
| 166 |
STREAMICE_vfacemask(i,j,bi,bj) = 2.0 |
| 167 |
ENDIF |
| 168 |
ENDIF |
| 169 |
|
| 170 |
ENDIF |
| 171 |
ENDDO |
| 172 |
ENDDO |
| 173 |
ENDDO |
| 174 |
ENDDO |
| 175 |
|
| 176 |
_EXCH_XY_RL( STREAMICE_ufacemask, myThid ) |
| 177 |
_EXCH_XY_RL( STREAMICE_vfacemask, myThid ) |
| 178 |
_EXCH_XY_RL( STREAMICE_umask, myThid ) |
| 179 |
_EXCH_XY_RL( STREAMICE_vmask, myThid ) |
| 180 |
|
| 181 |
! CALL WRITE_FULLARRAY_RL ("umask",STREAMICE_umask, |
| 182 |
! c 1,0,0,1,0,myThid) |
| 183 |
! CALL WRITE_FLD_XY_RL ("umask","",STREAMICE_umask,0,myThid) |
| 184 |
! CALL WRITE_FLD_XY_RL ("vmask","",STREAMICE_vmask,0,myThid) |
| 185 |
! CALL WRITE_FLD_XY_RL ("ufacemask","",STREAMICE_ufacemask,0,myThid) |
| 186 |
! CALL WRITE_FLD_XY_RL ("vfacemask","",STREAMICE_vfacemask,0,myThid) |
| 187 |
|
| 188 |
#ifdef ALLOW_PETSC |
| 189 |
|
| 190 |
DO bj = myByLo(myThid), myByHi(myThid) |
| 191 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
| 192 |
DO j=1,sNy |
| 193 |
DO i=1,sNx |
| 194 |
streamice_petsc_dofs_u (i,j,bi,bj) = -2.0 |
| 195 |
streamice_petsc_dofs_v (i,j,bi,bj) = -2.0 |
| 196 |
ENDDO |
| 197 |
ENDDO |
| 198 |
ENDDO |
| 199 |
ENDDO |
| 200 |
|
| 201 |
DoFCount = -1.0 |
| 202 |
DO bj = myByLo(myThid), myByHi(myThid) |
| 203 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
| 204 |
DO j=1,sNy |
| 205 |
DO i=1,sNx |
| 206 |
|
| 207 |
|
| 208 |
C DOFS ARE NUMBERED AS FOLLOWS ON PROCESSOR DOMAIN: |
| 209 |
C grid is stepped through in order bj, bi, j, i |
| 210 |
C 1) if umask(i,j,bi,bj)==1, the counter is updated by 1; |
| 211 |
C streamice_petsc_dofs_u is assigned the counter; |
| 212 |
C o/w streamice_petsc_dofs_u is assigned -1 |
| 213 |
C 2) if vmask(i,j,bi,bj)==1, the counter is updated by 1; |
| 214 |
C streamice_petsc_dofs_v is assigned the counter; |
| 215 |
C o/w streamice_petsc_dofs_v is assigned -1 |
| 216 |
C NOTE THESE NUMBERING ARRAYS ARE USED TO CONSTRUCT PETSC VECTORS AND MATRIX |
| 217 |
|
| 218 |
if (STREAMICE_umask (i,j,bi,bj).eq.1) THEN |
| 219 |
DoFCount = DoFCount + 1.0 |
| 220 |
streamice_petsc_dofs_u (i,j,bi,bj) = DoFCount |
| 221 |
else |
| 222 |
streamice_petsc_dofs_u (i,j,bi,bj) = -1.0 |
| 223 |
endif |
| 224 |
|
| 225 |
if (STREAMICE_vmask (i,j,bi,bj).eq.1) THEN |
| 226 |
DoFCount = DoFCount + 1.0 |
| 227 |
streamice_petsc_dofs_v (i,j,bi,bj) = DoFCount |
| 228 |
else |
| 229 |
streamice_petsc_dofs_v (i,j,bi,bj) = -1.0 |
| 230 |
endif |
| 231 |
|
| 232 |
ENDDO |
| 233 |
ENDDO |
| 234 |
ENDDO |
| 235 |
ENDDO |
| 236 |
|
| 237 |
#ifdef ALLOW_USE_MPI |
| 238 |
|
| 239 |
DO i=0,nPx*nPy-1 |
| 240 |
n_dofs_proc_loc (i) = 0 |
| 241 |
ENDDO |
| 242 |
|
| 243 |
CALL MPI_COMM_RANK( MPI_COMM_WORLD, mpiMyWId, mpiRC ) |
| 244 |
|
| 245 |
n_dofs_proc_loc (mpiMyWId) = INT(DoFCount)+1 |
| 246 |
|
| 247 |
CALL MPI_Allreduce(n_dofs_proc_loc,n_dofs_process,nPx*nPy, |
| 248 |
& MPI_INTEGER, MPI_SUM,MPI_COMM_MODEL,mpiRC) |
| 249 |
|
| 250 |
n_dofs_cum_sum(0) = 0 |
| 251 |
|
| 252 |
DO i=1,nPx*nPy-1 |
| 253 |
n_dofs_cum_sum(i) = n_dofs_cum_sum(i-1)+ |
| 254 |
& n_dofs_process(i-1) |
| 255 |
ENDDO |
| 256 |
|
| 257 |
#else /* ALLOW_USE_MPI */ |
| 258 |
|
| 259 |
n_dofs_process (0) = INT(DoFCount)+1 |
| 260 |
n_dofs_cum_sum (0) = INT(DoFCount)+1 |
| 261 |
|
| 262 |
#endif /* ALLOW_USE_MPI */ |
| 263 |
|
| 264 |
DO bj = myByLo(myThid), myByHi(myThid) |
| 265 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
| 266 |
DO j=1,sNy |
| 267 |
DO i=1,sNx |
| 268 |
IF (streamice_petsc_dofs_u(i,j,bi,bj).ge.0 ) THEN |
| 269 |
streamice_petsc_dofs_u(i,j,bi,bj) = |
| 270 |
& streamice_petsc_dofs_u(i,j,bi,bj) + |
| 271 |
& n_dofs_cum_sum(mpimywid) |
| 272 |
ENDIF |
| 273 |
IF (streamice_petsc_dofs_v(i,j,bi,bj).ge.0 ) THEN |
| 274 |
streamice_petsc_dofs_v(i,j,bi,bj) = |
| 275 |
& streamice_petsc_dofs_v(i,j,bi,bj) + |
| 276 |
& n_dofs_cum_sum(mpimywid) |
| 277 |
ENDIF |
| 278 |
ENDDO |
| 279 |
ENDDO |
| 280 |
ENDDO |
| 281 |
ENDDO |
| 282 |
|
| 283 |
_EXCH_XY_RS(streamice_petsc_dofs_u,myThid) |
| 284 |
_EXCH_XY_RS(streamice_petsc_dofs_v,myThid) |
| 285 |
|
| 286 |
|
| 287 |
#endif /* ALLOW_PETSC */ |
| 288 |
|
| 289 |
|
| 290 |
#endif |
| 291 |
RETURN |
| 292 |
END |