Back to home page

MITgcm

 
 

    


File indexing completed on 2023-08-04 05:10:46 UTC

view on githubraw file Latest commit 45315406 on 2023-08-03 16:50:12 UTC
809c36b928 Patr*0001 #include "SEAICE_OPTIONS.h"
fc004572cd Jean*0002 #ifdef ALLOW_OBCS
                0003 # include "OBCS_OPTIONS.h"
                0004 #endif
809c36b928 Patr*0005 
                0006 CStartOfInterface
7e4f7741f6 Patr*0007       SUBROUTINE SEAICE_INIT_VARIA( myThid )
e5b8ebbabe Jean*0008 C     *==========================================================*
7e4f7741f6 Patr*0009 C     | SUBROUTINE SEAICE_INIT_VARIA                             |
809c36b928 Patr*0010 C     | o Initialization of sea ice model.                       |
e5b8ebbabe Jean*0011 C     *==========================================================*
                0012 C     *==========================================================*
809c36b928 Patr*0013       IMPLICIT NONE
f4cb8ca7db Jean*0014 
809c36b928 Patr*0015 C     === Global variables ===
                0016 #include "SIZE.h"
                0017 #include "EEPARAMS.h"
                0018 #include "PARAMS.h"
                0019 #include "GRID.h"
36c84b38c7 Dimi*0020 #include "DYNVARS.h"
f3ce416a61 Jean*0021 #include "FFIELDS.h"
ccaa3c61f4 Patr*0022 #include "SEAICE_SIZE.h"
8c50aa8796 Jean*0023 #include "SEAICE_PARAMS.h"
                0024 #include "SEAICE.h"
ccaa3c61f4 Patr*0025 #include "SEAICE_TRACER.h"
8c50aa8796 Jean*0026 #include "SEAICE_TAVE.h"
fc004572cd Jean*0027 #ifdef OBCS_UVICE_OLD
a9eb030de2 Jean*0028 # include "OBCS_GRID.h"
ae1fb66b64 Dimi*0029 #endif
809c36b928 Patr*0030 
                0031 C     === Routine arguments ===
425e8efc36 Jean*0032 C     myThid :: Thread no. that called this routine.
809c36b928 Patr*0033       INTEGER myThid
                0034 CEndOfInterface
f4cb8ca7db Jean*0035 
809c36b928 Patr*0036 C     === Local variables ===
425e8efc36 Jean*0037 C     i,j,k,bi,bj :: Loop counters
809c36b928 Patr*0038 
edfdf5fa1d Jean*0039       INTEGER i, j, bi, bj
0320e25227 Mart*0040       INTEGER kSrf
8c50aa8796 Jean*0041       INTEGER k
f50f58ec54 Gael*0042 #ifdef ALLOW_SITRACER
                0043       INTEGER iTr, jTh
                0044 #endif
fc004572cd Jean*0045 #ifdef OBCS_UVICE_OLD
dc8f395c29 Dimi*0046       INTEGER I_obc, J_obc
                0047 #endif /* ALLOW_OBCS */
0dcda34642 Jean*0048 #ifdef SEAICE_CGRID
45315406aa Mart*0049       _RS mask_uice
2f5e8addfd Mart*0050       _RL recip_tensilDepth
0dcda34642 Jean*0051 #endif
486cba1d02 Mart*0052 
0320e25227 Mart*0053       IF ( usingPCoords ) THEN
                0054        kSrf = Nr
f4cb8ca7db Jean*0055       ELSE
0320e25227 Mart*0056        kSrf = 1
f4cb8ca7db Jean*0057       ENDIF
a6d2313c04 Dimi*0058 
f4cb8ca7db Jean*0059 C--   Initialise all variables in common blocks:
7109a141b2 Patr*0060       DO bj=myByLo(myThid),myByHi(myThid)
                0061        DO bi=myBxLo(myThid),myBxHi(myThid)
04016f2c47 Mart*0062         DO j=1-OLy,sNy+OLy
                0063          DO i=1-OLx,sNx+OLx
772590b63c Mart*0064           HEFF(i,j,bi,bj)=0. _d 0
                0065           AREA(i,j,bi,bj)=0. _d 0
45315406aa Mart*0066           HSNOW(i,j,bi,bj)   = 0. _d 0
86b84a92fc Patr*0067 #ifdef SEAICE_ITD
                0068           DO k=1,nITD
                0069            AREAITD(i,j,k,bi,bj) =0. _d 0
                0070            HEFFITD(i,j,k,bi,bj) =0. _d 0
45315406aa Mart*0071            HSNOWITD(i,j,k,bi,bj)=0. _d 0
86b84a92fc Patr*0072           ENDDO
                0073 #endif
772590b63c Mart*0074           UICE(i,j,bi,bj)=0. _d 0
                0075           VICE(i,j,bi,bj)=0. _d 0
                0076 C
5fe78992ba Mart*0077           uIceNm1(i,j,bi,bj) = 0. _d 0
                0078           vIceNm1(i,j,bi,bj) = 0. _d 0
45315406aa Mart*0079           DWATN  (i,j,bi,bj) = 0. _d 0
                0080 #ifdef SEAICE_CGRID
                0081           stressDivergenceX(i,j,bi,bj) = 0. _d 0
                0082           stressDivergenceY(i,j,bi,bj) = 0. _d 0
                0083 # ifdef SEAICE_ALLOW_EVP
                0084           seaice_sigma1 (i,j,bi,bj) = 0. _d 0
                0085           seaice_sigma2 (i,j,bi,bj) = 0. _d 0
                0086           seaice_sigma12(i,j,bi,bj) = 0. _d 0
                0087 # endif /* SEAICE_ALLOW_EVP */
                0088 #endif
                0089 #if ( defined SEAICE_CGRID || defined SEAICE_BGRID_DYNAMICS )
5fe78992ba Mart*0090           e11    (i,j,bi,bj) = 0. _d 0
                0091           e22    (i,j,bi,bj) = 0. _d 0
                0092           e12    (i,j,bi,bj) = 0. _d 0
                0093           deltaC (i,j,bi,bj) = 0. _d 0
45315406aa Mart*0094           PRESS  (i,j,bi,bj) = 0. _d 0
5fe78992ba Mart*0095           ETA    (i,j,bi,bj) = 0. _d 0
                0096           etaZ   (i,j,bi,bj) = 0. _d 0
                0097           ZETA   (i,j,bi,bj) = 0. _d 0
                0098           FORCEX (i,j,bi,bj) = 0. _d 0
                0099           FORCEY (i,j,bi,bj) = 0. _d 0
45315406aa Mart*0100           tensileStrFac(i,j,bi,bj) = 0. _d 0
                0101           PRESS0(i,j,bi,bj)  = 0. _d 0
                0102           FORCEX0(i,j,bi,bj) = 0. _d 0
                0103           FORCEY0(i,j,bi,bj) = 0. _d 0
                0104           SEAICE_zMax(i,j,bi,bj) = 0. _d 0
                0105           SEAICE_zMin(i,j,bi,bj) = 0. _d 0
                0106 #endif
f4cb8ca7db Jean*0107 #ifdef SEAICE_CGRID
                0108           seaiceMassC(i,j,bi,bj)=0. _d 0
                0109           seaiceMassU(i,j,bi,bj)=0. _d 0
                0110           seaiceMassV(i,j,bi,bj)=0. _d 0
45315406aa Mart*0111 # ifdef SEAICE_ALLOW_FREEDRIFT
                0112           uice_fd(i,j,bi,bj) = 0. _d 0
                0113           vice_fd(i,j,bi,bj) = 0. _d 0
                0114 # endif
                0115 # ifdef SEAICE_ALLOW_BOTTOMDRAG
                0116           CbotC(i,j,bi,bj)   = 0. _d 0
                0117 # endif /* SEAICE_ALLOW_BOTTOMDRAG */
                0118 #endif /* SEAICE_CGRID */
                0119 #ifdef SEAICE_BGRID_DYNAMICS
                0120           uIceB(i,j,bi,bj)   = 0. _d 0
                0121           vIceB(i,j,bi,bj)   = 0. _d 0
5fe78992ba Mart*0122           AMASS(i,j,bi,bj)   = 0. _d 0
                0123           DAIRN(i,j,bi,bj)   = 0. _d 0
                0124           WINDX(i,j,bi,bj)   = 0. _d 0
                0125           WINDY(i,j,bi,bj)   = 0. _d 0
                0126           GWATX(i,j,bi,bj)   = 0. _d 0
                0127           GWATY(i,j,bi,bj)   = 0. _d 0
45315406aa Mart*0128 #endif /* SEAICE_BGRID_DYNAMICS */
a98c4b8072 Ian *0129 #ifdef SEAICE_VARIABLE_SALINITY
f4cb8ca7db Jean*0130           HSALT(i,j,bi,bj)  = 0. _d 0
                0131 #endif
f50f58ec54 Gael*0132 #ifdef ALLOW_SITRACER
                0133           DO iTr = 1, SItrMaxNum
                0134            SItracer(i,j,bi,bj,iTr) = 0. _d 0
                0135            SItrBucket(i,j,bi,bj,iTr) = 0. _d 0
                0136 c "ice concentration" tracer that should remain .EQ.1.
                0137            if (SItrName(iTr).EQ.'one') SItracer(i,j,bi,bj,iTr)=1. _d 0
                0138           ENDDO
                0139           DO jTh = 1, 5
                0140            SItrHEFF (i,j,bi,bj,jTh) = 0. _d 0
                0141           ENDDO
bb24b8a3e6 Gael*0142           DO jTh = 1, 3
                0143            SItrAREA (i,j,bi,bj,jTh) = 0. _d 0
                0144           ENDDO
f50f58ec54 Gael*0145 #endif
f913c5a485 Mart*0146           DO k=1,nITD
f5282c5b03 Gael*0147             TICES(i,j,k,bi,bj)=0. _d 0
                0148           ENDDO
1cfef3b6ad Patr*0149           saltWtrIce(i,j,bi,bj) = 0. _d 0
                0150           frWtrIce(i,j,bi,bj)   = 0. _d 0
7109a141b2 Patr*0151          ENDDO
                0152         ENDDO
                0153        ENDDO
                0154       ENDDO
                0155 
c69ccc91fb antn*0156 #ifdef ALLOW_AUTODIFF
                0157 C- Note: To simplify dependency & avoid recomputations, when compiling
                0158 C        pkg/autodiff, we always call SEAICE_INIT_VARIA to initialise control
                0159 C        variables (as done above) without condition on useSEAICE.
                0160 C        Therefore, in this case, the "If useSEAICE" is added back here:
                0161       IF ( useSEAICE ) THEN
                0162 #endif
                0163 
8c50aa8796 Jean*0164 #ifdef ALLOW_TIMEAVE
                0165 C     Initialize averages to zero
                0166       DO bj = myByLo(myThid), myByHi(myThid)
                0167        DO bi = myBxLo(myThid), myBxHi(myThid)
                0168         CALL TIMEAVE_RESET( FUtave   , 1, bi, bj, myThid )
                0169         CALL TIMEAVE_RESET( FVtave   , 1, bi, bj, myThid )
                0170         CALL TIMEAVE_RESET( EmPmRtave, 1, bi, bj, myThid )
                0171         CALL TIMEAVE_RESET( QNETtave , 1, bi, bj, myThid )
                0172         CALL TIMEAVE_RESET( QSWtave  , 1, bi, bj, myThid )
                0173         CALL TIMEAVE_RESET( UICEtave , 1, bi, bj, myThid )
                0174         CALL TIMEAVE_RESET( VICEtave , 1, bi, bj, myThid )
                0175         CALL TIMEAVE_RESET( HEFFtave , 1, bi, bj, myThid )
                0176         CALL TIMEAVE_RESET( AREAtave , 1, bi, bj, myThid )
                0177         SEAICE_timeAve(bi,bj) = ZERO
                0178        ENDDO
                0179       ENDDO
                0180 #endif /* ALLOW_TIMEAVE */
                0181 
84cb676f82 Mart*0182 C--   Initialize (variable) grid info. As long as we allow masking of
                0183 C--   velocities outside of ice covered areas (in seaice_dynsolver)
                0184 C--   we need to re-initialize seaiceMaskU/V here for TAF/TAMC
809c36b928 Patr*0185       DO bj=myByLo(myThid),myByHi(myThid)
                0186        DO bi=myBxLo(myThid),myBxHi(myThid)
45315406aa Mart*0187 #ifdef SEAICE_CGRID
f4cb8ca7db Jean*0188         DO j=1-OLy+1,sNy+OLy
                0189          DO i=1-OLx+1,sNx+OLx
                0190           seaiceMaskU(i,j,bi,bj)=   0.0 _d 0
                0191           seaiceMaskV(i,j,bi,bj)=   0.0 _d 0
                0192           mask_uice=HEFFM(i,j,bi,bj)+HEFFM(i-1,j  ,bi,bj)
a808b72b4b Patr*0193           IF(mask_uice.GT.1.5 _d 0) seaiceMaskU(i,j,bi,bj)=1.0 _d 0
f4cb8ca7db Jean*0194           mask_uice=HEFFM(i,j,bi,bj)+HEFFM(i  ,j-1,bi,bj)
a808b72b4b Patr*0195           IF(mask_uice.GT.1.5 _d 0) seaiceMaskV(i,j,bi,bj)=1.0 _d 0
ab37184352 Mart*0196          ENDDO
                0197         ENDDO
45315406aa Mart*0198 # ifdef OBCS_UVICE_OLD
eb0337bb06 Dimi*0199         IF (useOBCS) THEN
b48f786726 Dimi*0200 C--   If OBCS is turned on, close southern and western boundaries
48474b8e08 Jean*0201          DO i=1-OLx,sNx+OLx
a6d2313c04 Dimi*0202 C Southern boundary
dc8f395c29 Dimi*0203           J_obc = OB_Js(i,bi,bj)
d6145066cc Jean*0204           IF ( J_obc.NE.OB_indexNone ) THEN
dc8f395c29 Dimi*0205            seaiceMaskU(i,J_obc,bi,bj)=   0.0 _d 0
                0206            seaiceMaskV(i,J_obc,bi,bj)=   0.0 _d 0
eb0337bb06 Dimi*0207           ENDIF
                0208          ENDDO
48474b8e08 Jean*0209          DO j=1-OLy,sNy+OLy
a6d2313c04 Dimi*0210 C Western boundary
d6145066cc Jean*0211           I_obc = OB_Iw(j,bi,bj)
                0212           IF ( I_obc.NE.OB_indexNone ) THEN
dc8f395c29 Dimi*0213            seaiceMaskU(I_obc,j,bi,bj)=   0.0 _d 0
                0214            seaiceMaskV(I_obc,j,bi,bj)=   0.0 _d 0
b48f786726 Dimi*0215           ENDIF
809c36b928 Patr*0216          ENDDO
b48f786726 Dimi*0217         ENDIF
45315406aa Mart*0218 # endif /* OBCS_UVICE_OLD */
a6d2313c04 Dimi*0219 #endif /* SEAICE_CGRID */
69c991a1ec Dimi*0220 
809c36b928 Patr*0221         DO j=1-OLy,sNy+OLy
                0222          DO i=1-OLx,sNx+OLx
f913c5a485 Mart*0223           DO k=1,nITD
f4cb8ca7db Jean*0224            TICES(i,j,k,bi,bj)=273.0 _d 0
09510da3bb Dimi*0225           ENDDO
45315406aa Mart*0226 #ifdef SEAICE_CGRID
f4cb8ca7db Jean*0227           seaiceMassC(i,j,bi,bj)=1000.0 _d 0
                0228           seaiceMassU(i,j,bi,bj)=1000.0 _d 0
                0229           seaiceMassV(i,j,bi,bj)=1000.0 _d 0
53092bcb42 Mart*0230 #endif
45315406aa Mart*0231 #ifdef SEAICE_BGRID_DYNAMICS
                0232           AMASS      (i,j,bi,bj)=1000.0 _d 0
                0233 #endif
809c36b928 Patr*0234          ENDDO
                0235         ENDDO
d5bc48d0a9 Dimi*0236 
809c36b928 Patr*0237        ENDDO
                0238       ENDDO
                0239 
                0240 C--   Update overlap regions
3ed8349f04 Mart*0241 #ifdef SEAICE_CGRID
b4f21203f7 Jean*0242       CALL EXCH_UV_XY_RL(seaiceMaskU,seaiceMaskV,.FALSE.,myThid)
45315406aa Mart*0243 #endif
                0244 #ifdef SEAICE_BGRID_DYNAMICS
48474b8e08 Jean*0245       _EXCH_XY_RS(UVM, myThid)
3ed8349f04 Mart*0246 #endif
809c36b928 Patr*0247 
a9a8db28c6 Alis*0248 C--   Now lets look at all these beasts
a24915ab1a Jean*0249       IF ( plotLevel.GE.debLevC ) THEN
3e0af136db Dimi*0250          CALL PLOT_FIELD_XYRL( HEFFM   , 'Current HEFFM   ' ,
23142459d0 Jean*0251      &        nIter0, myThid )
3ed8349f04 Mart*0252 #ifdef SEAICE_CGRID
                0253          CALL PLOT_FIELD_XYRL( seaiceMaskU, 'Current seaiceMaskU',
23142459d0 Jean*0254      &        nIter0, myThid )
3ed8349f04 Mart*0255          CALL PLOT_FIELD_XYRL( seaiceMaskV, 'Current seaiceMaskV',
23142459d0 Jean*0256      &        nIter0, myThid )
45315406aa Mart*0257 #endif
                0258 #ifdef SEAICE_BGRID_DYNAMICS
48474b8e08 Jean*0259          CALL PLOT_FIELD_XYRS( UVM     , 'Current UVM     ' ,
23142459d0 Jean*0260      &        nIter0, myThid )
3ed8349f04 Mart*0261 #endif
3e0af136db Dimi*0262       ENDIF
09510da3bb Dimi*0263 
809c36b928 Patr*0264 C--   Set model variables to initial/restart conditions
2a682dc462 Mart*0265       IF ( .NOT. ( startTime .EQ. baseTime .AND.  nIter0 .EQ. 0
                0266      &     .AND. pickupSuff .EQ. ' ') ) THEN
36f6faaf99 Dimi*0267 
6060ec2938 Dimi*0268          CALL SEAICE_READ_PICKUP ( myThid )
36f6faaf99 Dimi*0269 
6060ec2938 Dimi*0270       ELSE
36f6faaf99 Dimi*0271 
e7c33124ce Dimi*0272        DO bj=myByLo(myThid),myByHi(myThid)
                0273         DO bi=myBxLo(myThid),myBxHi(myThid)
                0274          DO j=1-OLy,sNy+OLy
                0275           DO i=1-OLx,sNx+OLx
772590b63c Mart*0276            HEFF(i,j,bi,bj)=SEAICE_initialHEFF*HEFFM(i,j,bi,bj)
                0277            UICE(i,j,bi,bj)=ZERO
                0278            VICE(i,j,bi,bj)=ZERO
e7c33124ce Dimi*0279           ENDDO
809c36b928 Patr*0280          ENDDO
e7c33124ce Dimi*0281         ENDDO
                0282        ENDDO
6060ec2938 Dimi*0283 
425e8efc36 Jean*0284 C--   Read initial sea-ice velocity from file (if available)
                0285        IF ( uIceFile .NE. ' ' )
                0286      &  CALL READ_FLD_XY_RL( uIceFile, ' ', uIce, 0, myThid )
                0287        IF ( vIceFile .NE. ' ' )
                0288      &  CALL READ_FLD_XY_RL( vIceFile, ' ', vIce, 0, myThid )
                0289        IF ( uIceFile .NE. ' ' .OR. vIceFile .NE. ' ' ) THEN
                0290 #ifdef SEAICE_CGRID
                0291         DO bj=myByLo(myThid),myByHi(myThid)
                0292          DO bi=myBxLo(myThid),myBxHi(myThid)
                0293           DO j=1-OLy,sNy+OLy
                0294            DO i=1-OLx,sNx+OLx
                0295             uIce(i,j,bi,bj) = uIce(i,j,bi,bj)*seaiceMaskU(i,j,bi,bj)
                0296             vIce(i,j,bi,bj) = vIce(i,j,bi,bj)*seaiceMaskV(i,j,bi,bj)
                0297            ENDDO
                0298           ENDDO
                0299          ENDDO
                0300         ENDDO
                0301 #endif /* SEAICE_CGRID */
                0302         CALL EXCH_UV_XY_RL( uIce, vIce, .TRUE., myThid )
                0303        ENDIF
                0304 
6060ec2938 Dimi*0305 C--   Read initial sea-ice thickness from file if available.
e7c33124ce Dimi*0306        IF ( HeffFile .NE. ' ' ) THEN
772590b63c Mart*0307         CALL READ_FLD_XY_RL( HeffFile, ' ', HEFF, 0, myThid )
                0308         _EXCH_XY_RL(HEFF,myThid)
e7c33124ce Dimi*0309         DO bj=myByLo(myThid),myByHi(myThid)
                0310          DO bi=myBxLo(myThid),myBxHi(myThid)
                0311           DO j=1-OLy,sNy+OLy
                0312            DO i=1-OLx,sNx+OLx
772590b63c Mart*0313             HEFF(i,j,bi,bj) = MAX(HEFF(i,j,bi,bj),ZERO)
e7c33124ce Dimi*0314            ENDDO
                0315           ENDDO
                0316          ENDDO
                0317         ENDDO
                0318        ENDIF
                0319 
                0320        DO bj=myByLo(myThid),myByHi(myThid)
                0321         DO bi=myBxLo(myThid),myBxHi(myThid)
                0322          DO j=1-OLy,sNy+OLy
                0323           DO i=1-OLx,sNx+OLx
772590b63c Mart*0324            IF(HEFF(i,j,bi,bj).GT.ZERO) AREA(i,j,bi,bj)=ONE
e7c33124ce Dimi*0325           ENDDO
                0326          ENDDO
                0327         ENDDO
                0328        ENDDO
f4cb8ca7db Jean*0329 
ff73878394 Dimi*0330 C--   Read initial sea-ice area from file if available.
09066b09cb Mart*0331        IF ( AreaFile .NE. ' ' ) THEN
772590b63c Mart*0332         CALL READ_FLD_XY_RL( AreaFile, ' ', AREA, 0, myThid )
                0333         _EXCH_XY_RL(AREA,myThid)
09066b09cb Mart*0334         DO bj=myByLo(myThid),myByHi(myThid)
                0335          DO bi=myBxLo(myThid),myBxHi(myThid)
                0336           DO j=1-OLy,sNy+OLy
                0337            DO i=1-OLx,sNx+OLx
772590b63c Mart*0338             AREA(i,j,bi,bj) = MAX(AREA(i,j,bi,bj),ZERO)
                0339             AREA(i,j,bi,bj) = MIN(AREA(i,j,bi,bj),ONE)
                0340             IF ( AREA(i,j,bi,bj) .LE. ZERO ) HEFF(i,j,bi,bj) = ZERO
                0341             IF ( HEFF(i,j,bi,bj) .LE. ZERO ) AREA(i,j,bi,bj) = ZERO
09066b09cb Mart*0342            ENDDO
                0343           ENDDO
                0344          ENDDO
                0345         ENDDO
                0346        ENDIF
                0347 
                0348        DO bj=myByLo(myThid),myByHi(myThid)
                0349         DO bi=myBxLo(myThid),myBxHi(myThid)
                0350          DO j=1-OLy,sNy+OLy
                0351           DO i=1-OLx,sNx+OLx
e5b8ebbabe Jean*0352            HSNOW(i,j,bi,bj) = 0.2 _d 0 * AREA(i,j,bi,bj)
09066b09cb Mart*0353           ENDDO
                0354          ENDDO
                0355         ENDDO
                0356        ENDDO
de31ea8481 Dimi*0357 
                0358 C--   Read initial snow thickness from file if available.
                0359        IF ( HsnowFile .NE. ' ' ) THEN
772590b63c Mart*0360         CALL READ_FLD_XY_RL( HsnowFile, ' ', HSNOW, 0, myThid )
                0361         _EXCH_XY_RL(HSNOW,myThid)
de31ea8481 Dimi*0362         DO bj=myByLo(myThid),myByHi(myThid)
                0363          DO bi=myBxLo(myThid),myBxHi(myThid)
                0364           DO j=1-OLy,sNy+OLy
                0365            DO i=1-OLx,sNx+OLx
772590b63c Mart*0366             HSNOW(i,j,bi,bj) = MAX(HSNOW(i,j,bi,bj),ZERO)
de31ea8481 Dimi*0367            ENDDO
                0368           ENDDO
                0369          ENDDO
                0370         ENDDO
                0371        ENDIF
fdfa8e151f Dimi*0372 
86b84a92fc Patr*0373 #ifdef SEAICE_ITD
                0374        DO bj=myByLo(myThid),myByHi(myThid)
                0375         DO bi=myBxLo(myThid),myBxHi(myThid)
                0376          DO j=1-OLy,sNy+OLy
                0377           DO i=1-OLx,sNx+OLx
5fe78992ba Mart*0378            AREAITD(i,j,1,bi,bj)   = AREA(i,j,bi,bj)
                0379            HEFFITD(i,j,1,bi,bj)   = HEFF(i,j,bi,bj)
                0380            HSNOWITD(i,j,1,bi,bj)  = HSNOW(i,j,bi,bj)
                0381            opnWtrFrac(i,j,bi,bj)  = 1. _d 0 - AREA(i,j,bi,bj)
                0382            fw2ObyRidge(i,j,bi,bj) = 0. _d 0
86b84a92fc Patr*0383           ENDDO
                0384          ENDDO
fd2abd71b8 Mart*0385          CALL SEAICE_ITD_REDIST(bi, bj, baseTime, nIter0, myThid)
86b84a92fc Patr*0386         ENDDO
                0387        ENDDO
                0388 #endif
                0389 
a98c4b8072 Ian *0390 #ifdef SEAICE_VARIABLE_SALINITY
3183247f03 Dimi*0391        DO bj=myByLo(myThid),myByHi(myThid)
                0392         DO bi=myBxLo(myThid),myBxHi(myThid)
f4cb8ca7db Jean*0393          DO j=1-OLy,sNy+OLy
                0394           DO i=1-OLx,sNx+OLx
0320e25227 Mart*0395            HSALT(i,j,bi,bj)=HEFF(i,j,bi,bj)*salt(i,j,kSrf,bi,bj)*
4eb4a54cba Jean*0396      &            SEAICE_rhoIce*SEAICE_saltFrac
                0397 cif     &            ICE2WATR*rhoConstFresh*SEAICE_saltFrac
a98c4b8072 Ian *0398 
3183247f03 Dimi*0399           ENDDO
                0400          ENDDO
                0401         ENDDO
                0402        ENDDO
                0403 
fdfa8e151f Dimi*0404 C--   Read initial sea ice salinity from file if available.
                0405        IF ( HsaltFile .NE. ' ' ) THEN
772590b63c Mart*0406         CALL READ_FLD_XY_RL( HsaltFile, ' ', HSALT, 0, myThid )
                0407         _EXCH_XY_RL(HSALT,myThid)
fdfa8e151f Dimi*0408        ENDIF
a98c4b8072 Ian *0409 #endif /* SEAICE_VARIABLE_SALINITY */
f4cb8ca7db Jean*0410 
e54fe3e1f9 Gael*0411 #ifdef ALLOW_SITRACER
f681b7f5d4 Dimi*0412 C--   Read initial sea ice age from file if available.
e54fe3e1f9 Gael*0413        DO iTr = 1, SItrMaxNum
                0414         IF ( SItrFile(iTr) .NE. ' ' ) THEN
                0415         CALL READ_FLD_XY_RL( siTrFile(iTr), ' ',
                0416      &   SItracer(1-OLx,1-OLy,1,1,iTr), 0, myThid )
                0417         _EXCH_XY_RL(SItracer(1-OLx,1-OLy,1,1,iTr),myThid)
ccaa3c61f4 Patr*0418         ENDIF
                0419        ENDDO
e54fe3e1f9 Gael*0420 #endif /* ALLOW_SITRACER */
f681b7f5d4 Dimi*0421 
809c36b928 Patr*0422       ENDIF
                0423 
4f6e464f38 Jean*0424 #ifdef ALLOW_OBCS
3352db77f4 Jean*0425 C--   In case we use scheme with a large stencil that extends into overlap:
                0426 C     no longer needed with the right masking in advection & diffusion S/R.
                0427 c     IF ( useOBCS ) THEN
                0428 c       DO bj=myByLo(myThid),myByHi(myThid)
                0429 c        DO bi=myBxLo(myThid),myBxHi(myThid)
48474b8e08 Jean*0430 c          CALL OBCS_COPY_TRACER( HEFF(1-OLx,1-OLy,bi,bj),
3352db77f4 Jean*0431 c    I                            1, bi, bj, myThid )
48474b8e08 Jean*0432 c          CALL OBCS_COPY_TRACER( AREA(1-OLx,1-OLy,bi,bj),
3352db77f4 Jean*0433 c    I                            1, bi, bj, myThid )
48474b8e08 Jean*0434 c          CALL OBCS_COPY_TRACER( HSNOW(1-OLx,1-OLy,bi,bj),
3352db77f4 Jean*0435 c    I                            1, bi, bj, myThid )
a98c4b8072 Ian *0436 #ifdef SEAICE_VARIABLE_SALINITY
48474b8e08 Jean*0437 c          CALL OBCS_COPY_TRACER( HSALT(1-OLx,1-OLy,bi,bj),
3352db77f4 Jean*0438 c    I                            1, bi, bj, myThid )
4f6e464f38 Jean*0439 #endif
3352db77f4 Jean*0440 c        ENDDO
                0441 c       ENDDO
                0442 c     ENDIF
4f6e464f38 Jean*0443 #endif /* ALLOW_OBCS */
                0444 
809c36b928 Patr*0445 C---  Complete initialization
                0446       DO bj=myByLo(myThid),myByHi(myThid)
                0447        DO bi=myBxLo(myThid),myBxHi(myThid)
45315406aa Mart*0448 #if ( defined SEAICE_CGRID || defined SEAICE_BGRID_DYNAMICS )
809c36b928 Patr*0449         DO j=1-OLy,sNy+OLy
                0450          DO i=1-OLx,sNx+OLx
8e32c48b8f Mart*0451           ZETA(i,j,bi,bj)        = HEFF(i,j,bi,bj)*(1.0 _d 11)
                0452           ETA(i,j,bi,bj)         = ZETA(i,j,bi,bj)/SEAICE_eccen**2
                0453           PRESS0(i,j,bi,bj)      = SEAICE_strength*HEFF(i,j,bi,bj)
ba6cfc5714 Mart*0454      &         *EXP(-SEAICE_cStar*(ONE-AREA(i,j,bi,bj)))
8e32c48b8f Mart*0455           SEAICE_zMax(i,j,bi,bj) = SEAICE_zetaMaxFac*PRESS0(i,j,bi,bj)
                0456           SEAICE_zMin(i,j,bi,bj) = SEAICE_zetaMin
                0457           PRESS0(i,j,bi,bj)      = PRESS0(i,j,bi,bj)*HEFFM(i,j,bi,bj)
809c36b928 Patr*0458          ENDDO
                0459         ENDDO
45315406aa Mart*0460 #endif
f3ce416a61 Jean*0461         IF ( useRealFreshWaterFlux .AND. .NOT.useThSIce ) THEN
                0462          DO j=1-OLy,sNy+OLy
                0463           DO i=1-OLx,sNx+OLx
772590b63c Mart*0464            sIceLoad(i,j,bi,bj) = HEFF(i,j,bi,bj)*SEAICE_rhoIce
13533b2b42 Mart*0465      &                         + HSNOW(i,j,bi,bj)*SEAICE_rhoSnow
f4cb8ca7db Jean*0466 
f3ce416a61 Jean*0467           ENDDO
                0468          ENDDO
                0469         ENDIF
809c36b928 Patr*0470        ENDDO
                0471       ENDDO
2f5e8addfd Mart*0472 #ifdef SEAICE_CGRID
                0473 C     compute tensile strength factor k: tensileStrength = k*PRESS
0dcda34642 Jean*0474 C     can be done in initialisation phase as long as it depends only
2f5e8addfd Mart*0475 C     on depth
                0476       IF ( SEAICE_tensilFac .NE. 0. _d 0 ) THEN
                0477        recip_tensilDepth = 0. _d 0
0dcda34642 Jean*0478        IF ( SEAICE_tensilDepth .GT. 0. _d 0 )
2f5e8addfd Mart*0479      &      recip_tensilDepth = 1. _d 0 / SEAICE_tensilDepth
                0480        DO bj=myByLo(myThid),myByHi(myThid)
                0481         DO bi=myBxLo(myThid),myBxHi(myThid)
                0482          DO j=1-OLy,sNy+OLy
                0483           DO i=1-OLx,sNx+OLx
                0484            tensileStrFac(i,j,bi,bj) = SEAICE_tensilFac*HEFFM(i,j,bi,bj)
                0485      &          *exp(-ABS(R_low(i,j,bi,bj))*recip_tensilDepth)
                0486           ENDDO
                0487          ENDDO
                0488         ENDDO
                0489        ENDDO
                0490       ENDIF
45315406aa Mart*0491 # if ( defined SEAICE_ALLOW_JFNK || defined SEAICE_ALLOW_KRYLOV )
                0492 C     Computing this metric cannot be done in S/R SEAICE_INIT_FIXED
                0493 C     where it belongs, because globalArea is only defined later after
                0494 C     S/R PACKAGES_INIT_FIXED, so we move this computation here.
                0495       CALL SEAICE_MAP_RS2VEC( nVec, rAw, rAs,
                0496      &     scalarProductMetric, .TRUE., myThid )
                0497       DO bj=myByLo(myThid),myByHi(myThid)
                0498        DO bi=myBxLo(myThid),myBxHi(myThid)
                0499         DO i=1,nVec
                0500          scalarProductMetric(i,1,bi,bj) =
                0501      &        scalarProductMetric(i,1,bi,bj)/globalArea
                0502         ENDDO
                0503        ENDDO
                0504       ENDDO
                0505 # endif /* SEAICE_ALLOW_JFNK or KRYLOV */
                0506 
2f5e8addfd Mart*0507 #endif /* SEAICE_CGRID */
809c36b928 Patr*0508 
c69ccc91fb antn*0509 #ifdef ALLOW_AUTODIFF
                0510 C-    end if useSEAICE block
                0511       ENDIF
                0512 #endif
                0513 
809c36b928 Patr*0514       RETURN
                0515       END