/[MITgcm]/MITgcm/pkg/seaice/seaice_init.F
ViewVC logotype

Contents of /MITgcm/pkg/seaice/seaice_init.F

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


Revision 1.7 - (show annotations) (download)
Wed Apr 30 07:04:08 2003 UTC (21 years ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint50e_post, checkpoint50c_post, checkpoint50c_pre, checkpoint50d_pre, checkpoint50d_post, checkpoint50f_post, checkpoint50f_pre, checkpoint50e_pre
Changes since 1.6: +9 -9 lines
checkpoint50c_pre
Merging from release1_p13:
o bug fix for pkg/seaice dynamic solver
o Added SEAICE_initialHEFF to pkg/seaice

1 C $Header:
2
3 #include "SEAICE_OPTIONS.h"
4
5 CStartOfInterface
6 SUBROUTINE SEAICE_INIT( myThid )
7 C /==========================================================\
8 C | SUBROUTINE SEAICE_INIT |
9 C | o Initialization of sea ice model. |
10 C |==========================================================|
11 C \==========================================================/
12 IMPLICIT NONE
13
14 C === Global variables ===
15 #include "SIZE.h"
16 #include "EEPARAMS.h"
17 #include "PARAMS.h"
18 #include "GRID.h"
19 #include "SEAICE.h"
20 #include "SEAICE_GRID.h"
21 #include "SEAICE_DIAGS.h"
22 #include "SEAICE_PARAMS.h"
23 #include "SEAICE_EXTERNAL.h"
24
25 C === Routine arguments ===
26 C myThid - Thread no. that called this routine.
27 INTEGER myThid
28 CEndOfInterface
29
30 #ifdef ALLOW_SEAICE
31 C === Local variables ===
32 C i,j,k,bi,bj - Loop counters
33
34 INTEGER i, j, k, bi, bj
35 _RS mask_uice
36
37 #ifdef ALLOW_TIMEAVE
38 C Initialize averages to zero
39 DO bj = myByLo(myThid), myByHi(myThid)
40 DO bi = myBxLo(myThid), myBxHi(myThid)
41 CALL TIMEAVE_RESET(FUtave ,1,bi,bj,myThid)
42 CALL TIMEAVE_RESET(FVtave ,1,bi,bj,myThid)
43 CALL TIMEAVE_RESET(EmPmRtave,1,bi,bj,myThid)
44 CALL TIMEAVE_RESET(QNETtave ,1,bi,bj,myThid)
45 CALL TIMEAVE_RESET(QSWtave ,1,bi,bj,myThid)
46 CALL TIMEAVE_RESET(UICEtave ,1,bi,bj,myThid)
47 CALL TIMEAVE_RESET(VICEtave ,1,bi,bj,myThid)
48 CALL TIMEAVE_RESET(HEFFtave ,1,bi,bj,myThid)
49 CALL TIMEAVE_RESET(AREAtave ,1,bi,bj,myThid)
50 DO k=1,Nr
51 SEAICE_TimeAve(k,bi,bj)=ZERO
52 ENDDO
53 ENDDO
54 ENDDO
55 #endif /* ALLOW_TIMEAVE */
56
57 C--- initialize grid info
58 DO bj=myByLo(myThid),myByHi(myThid)
59 DO bi=myBxLo(myThid),myBxHi(myThid)
60 DO J=1,sNy
61 DO I=1,sNx
62 CSTICE(i,j,bi,bj) =cos(yC(I,J,bi,bj)*deg2rad)
63 CSUICE(i,j,bi,bj) =cos(yG(I,J,bi,bj)*deg2rad)
64 SINEICE(i,j,bi,bj)=sin(yG(I,J,bi,bj)*deg2rad)
65 TNGTICE(i,j,bi,bj)=tan(acos(CSTICE(i,j,bi,bj)))
66 TNGICE(i,j,bi,bj) =SINEICE(i,j,bi,bj)/CSUICE(i,j,bi,bj)
67 DXTICE(i,j,bi,bj)=dxF(i,j,bi,bj)/CSTICE(i,j,bi,bj)
68 DXUICE(i,j,bi,bj)=dxV(i,j,bi,bj)/CSUICE(i,j,bi,bj)
69 DYTICE(i,j,bi,bj)=dyF(i,j,bi,bj)
70 DYUICE(i,j,bi,bj)=dyU(i,j,bi,bj)
71 ENDDO
72 ENDDO
73 DO j=1-OLy,sNy+OLy
74 DO i=1-OLx,sNx+OLx
75 HEFFM(i,j,bi,bj)=ONE
76 IF (_hFacC(i,j,1,bi,bj).eq.0.) HEFFM(i,j,bi,bj)=ZERO
77 ENDDO
78 ENDDO
79 DO J=1,sNy
80 DO I=1,sNx
81 UVM(i,j,bi,bj)=ZERO
82 mask_uice=HEFFM(I,J, bi,bj)+HEFFM(I-1,J-1,bi,bj)
83 & +HEFFM(I,J-1,bi,bj)+HEFFM(I-1,J, bi,bj)
84 IF(mask_uice.GT.3.5) UVM(I,J,bi,bj)=ONE
85 ENDDO
86 ENDDO
87 DO j=1-OLy,sNy+OLy
88 DO i=1-OLx,sNx+OLx
89 TICE(I,J,bi,bj)=273.0 _d 0
90 UICEC(I,J,bi,bj)=ZERO
91 VICEC(I,J,bi,bj)=ZERO
92 AMASS(I,J,bi,bj)=1000.0 _d 0
93 ENDDO
94 ENDDO
95 ENDDO
96 ENDDO
97
98 C-- Update overlap regions
99 _EXCH_XY_R8(UVM, myThid)
100 _EXCH_XY_R8(TNGTICE, myThid)
101 _EXCH_XY_R8(TNGICE, myThid)
102 _EXCH_XY_R8(CSTICE, myThid)
103 _EXCH_XY_R8(CSUICE, myThid)
104 _EXCH_XY_R8(SINEICE, myThid)
105 _EXCH_XY_R8(DXTICE, myThid)
106 _EXCH_XY_R8(DXUICE, myThid)
107 _EXCH_XY_R8(DYTICE, myThid)
108 _EXCH_XY_R8(DYUICE, myThid)
109
110 C-- Set model variables to initial/restart conditions
111 IF ( nIter0 .NE. 0 ) THEN
112 CALL SEAICE_READ_PICKUP ( myThid )
113 ELSE
114 DO bj=myByLo(myThid),myByHi(myThid)
115 DO bi=myBxLo(myThid),myBxHi(myThid)
116 DO j=1-OLy,sNy+OLy
117 DO i=1-OLx,sNx+OLx
118 HSNOW(I,J,bi,bj)=0.2 _d 0
119 YNEG(I,J,bi,bj)=ZERO
120 TMIX(I,J,bi,bj)=TICE(I,J,bi,bj)
121 DO k=1,3
122 HEFF(I,J,k,bi,bj)=SEAICE_initialHEFF
123 AREA(I,J,k,bi,bj)=HEFFM(i,j,bi,bj)
124 UICE(I,J,k,bi,bj)=ZERO
125 VICE(I,J,k,bi,bj)=ZERO
126 ENDDO
127 ENDDO
128 ENDDO
129 ENDDO
130 ENDDO
131 ENDIF
132
133 C-- Read initial sea-ice thickness from file if available.
134 IF ( HeffFile .NE. ' ' ) THEN
135 _BEGIN_MASTER( myThid )
136 CALL READ_FLD_XY_RL( HeffFile, ' ', ZETA, 0, myThid )
137 _END_MASTER(myThid)
138 _EXCH_XY_R8(ZETA,myThid)
139 DO bj=myByLo(myThid),myByHi(myThid)
140 DO bi=myBxLo(myThid),myBxHi(myThid)
141 DO j=1-OLy,sNy+OLy
142 DO i=1-OLx,sNx+OLx
143 DO k=1,3
144 HEFF(I,J,k,bi,bj) = MAX(ZETA(i,j,bi,bj),ZERO)
145 IF ( ZETA(i,j,bi,bj).EQ.ZERO )
146 & AREA(I,J,k,bi,bj) = ZERO
147 ENDDO
148 ENDDO
149 ENDDO
150 ENDDO
151 ENDDO
152 ENDIF
153
154 C--- Complete initialization
155 DO bj=myByLo(myThid),myByHi(myThid)
156 DO bi=myBxLo(myThid),myBxHi(myThid)
157 DO j=1-OLy,sNy+OLy
158 DO i=1-OLx,sNx+OLx
159 ZETA(I,J,bi,bj)=HEFF(I,J,1,bi,bj)*(1.0 _d 11)
160 ETA(I,J,bi,bj)=ZETA(I,J,bi,bj)/4.0 _d 0
161 surfaceTendencyTice(i,j,bi,bj) = ZERO
162 ENDDO
163 ENDDO
164 ENDDO
165 ENDDO
166
167 #endif /* ALLOW_SEAICE */
168
169 RETURN
170 END

  ViewVC Help
Powered by ViewVC 1.1.22