/[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.14 - (show annotations) (download)
Thu Nov 13 06:35:15 2003 UTC (20 years, 7 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint52e_pre, checkpoint52e_post, checkpoint52d_pre, branch-netcdf, checkpoint52b_pre, checkpoint52d_post, checkpoint52a_post, checkpoint52b_post, checkpoint52f_post, checkpoint52c_post, checkpoint52f_pre, hrcube_1
Branch point for: netcdf-sm0
Changes since 1.13: +1 -3 lines
o modifications to make FREEZE flux visible to pkg/kpp
  - moved surfaceTendencyTice from pkg/seaice to main code
  - FREEZE moved to FORWARD_STEP
  - subroutine FREEZE now limits only surface temperature
    this means new output.txt for global_ocean.90x40x15,
    global_ocean.cs32x15, and global_with_exf, but note
    that results for these three experiments remain
    bit-identical to before if allowFreezing=.FALSE.
o added surface flux output variables to TIMEAVE_STATVARS
o time-averaged output for pkg/ptracers

1 C $Header: /usr/local/gcmpack/MITgcm/pkg/seaice/seaice_init.F,v 1.13 2003/11/06 22:13:00 heimbach Exp $
2 C $Name: $
3
4 #include "SEAICE_OPTIONS.h"
5
6 CStartOfInterface
7 SUBROUTINE SEAICE_INIT( myThid )
8 C /==========================================================\
9 C | SUBROUTINE SEAICE_INIT |
10 C | o Initialization of sea ice model. |
11 C |==========================================================|
12 C \==========================================================/
13 IMPLICIT NONE
14
15 C === Global variables ===
16 #include "SIZE.h"
17 #include "EEPARAMS.h"
18 #include "PARAMS.h"
19 #include "GRID.h"
20 #include "SEAICE.h"
21 #include "SEAICE_GRID.h"
22 #include "SEAICE_DIAGS.h"
23 #include "SEAICE_PARAMS.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 INTEGER myIter
37
38 #ifdef ALLOW_TIMEAVE
39 C Initialize averages to zero
40 DO bj = myByLo(myThid), myByHi(myThid)
41 DO bi = myBxLo(myThid), myBxHi(myThid)
42 CALL TIMEAVE_RESET(FUtave ,1,bi,bj,myThid)
43 CALL TIMEAVE_RESET(FVtave ,1,bi,bj,myThid)
44 CALL TIMEAVE_RESET(EmPmRtave,1,bi,bj,myThid)
45 CALL TIMEAVE_RESET(QNETtave ,1,bi,bj,myThid)
46 CALL TIMEAVE_RESET(QSWtave ,1,bi,bj,myThid)
47 CALL TIMEAVE_RESET(UICEtave ,1,bi,bj,myThid)
48 CALL TIMEAVE_RESET(VICEtave ,1,bi,bj,myThid)
49 CALL TIMEAVE_RESET(HEFFtave ,1,bi,bj,myThid)
50 CALL TIMEAVE_RESET(AREAtave ,1,bi,bj,myThid)
51 DO k=1,Nr
52 SEAICE_TimeAve(k,bi,bj)=ZERO
53 ENDDO
54 ENDDO
55 ENDDO
56 #endif /* ALLOW_TIMEAVE */
57
58 cph(
59 cph make sure TAF sees proper initialisation
60 cph to avoid partial recomputation issues
61 DO bj=myByLo(myThid),myByHi(myThid)
62 DO bi=myBxLo(myThid),myBxHi(myThid)
63 c
64 DO K=1,3
65 DO J=1,sNy
66 DO I=1,sNx
67 HEFF(I,J,k,bi,bj)=SEAICE_initialHEFF
68 AREA(I,J,k,bi,bj)=ZERO
69 UICE(I,J,k,bi,bj)=ZERO
70 VICE(I,J,k,bi,bj)=ZERO
71 ENDDO
72 ENDDO
73 ENDDO
74 c
75 DO J=1,sNy
76 DO I=1,sNx
77 HSNOW(I,J,bi,bj)=0.2 _d 0
78 ZETA(I,J,bi,bj)=ZERO
79 ENDDO
80 ENDDO
81 c
82 ENDDO
83 ENDDO
84 cph)
85
86
87 C--- initialize grid info
88 DO bj=myByLo(myThid),myByHi(myThid)
89 DO bi=myBxLo(myThid),myBxHi(myThid)
90 DO J=1,sNy
91 DO I=1,sNx
92 CSTICE(i,j,bi,bj) =cos(yC(I,J,bi,bj)*deg2rad)
93 CSUICE(i,j,bi,bj) =cos(yG(I,J,bi,bj)*deg2rad)
94 SINEICE(i,j,bi,bj)=sin(yC(I,J,bi,bj)*deg2rad)
95 TNGTICE(i,j,bi,bj)=SINEICE(i,j,bi,bj)/CSTICE(i,j,bi,bj)
96 SINEICE(i,j,bi,bj)=sin(yG(I,J,bi,bj)*deg2rad)
97 TNGICE(i,j,bi,bj) =SINEICE(i,j,bi,bj)/CSUICE(i,j,bi,bj)
98 DXTICE(i,j,bi,bj)=dxF(i,j,bi,bj)/CSTICE(i,j,bi,bj)
99 DXUICE(i,j,bi,bj)=dxV(i,j,bi,bj)/CSUICE(i,j,bi,bj)
100 DYTICE(i,j,bi,bj)=dyF(i,j,bi,bj)
101 DYUICE(i,j,bi,bj)=dyU(i,j,bi,bj)
102 ENDDO
103 ENDDO
104
105 DO j=1-OLy,sNy+OLy
106 DO i=1-OLx,sNx+OLx
107 HEFFM(i,j,bi,bj)=ONE
108 IF (_hFacC(i,j,1,bi,bj).eq.0.) HEFFM(i,j,bi,bj)=ZERO
109 ENDDO
110 ENDDO
111 DO J=1,sNy
112 DO I=1,sNx
113 UVM(i,j,bi,bj)=ZERO
114 mask_uice=HEFFM(I,J, bi,bj)+HEFFM(I-1,J-1,bi,bj)
115 & +HEFFM(I,J-1,bi,bj)+HEFFM(I-1,J, bi,bj)
116 IF(mask_uice.GT.3.5) UVM(I,J,bi,bj)=ONE
117 ENDDO
118 ENDDO
119
120 DO j=1-OLy,sNy+OLy
121 DO i=1-OLx,sNx+OLx
122 TICE(I,J,bi,bj)=273.0 _d 0
123 #ifdef SEAICE_MULTILEVEL
124 DO k=1,MULTDIM
125 TICES(I,J,k,bi,bj)=273.0 _d 0
126 ENDDO
127 #endif
128 UICEC(I,J,bi,bj)=ZERO
129 VICEC(I,J,bi,bj)=ZERO
130 AMASS(I,J,bi,bj)=1000.0 _d 0
131 ENDDO
132 ENDDO
133 ENDDO
134 ENDDO
135
136 C-- Update overlap regions
137 _EXCH_XY_R8(UVM, myThid)
138 _EXCH_XY_R8(TNGTICE, myThid)
139 _EXCH_XY_R8(TNGICE, myThid)
140 _EXCH_XY_R8(CSTICE, myThid)
141 _EXCH_XY_R8(CSUICE, myThid)
142 _EXCH_XY_R8(SINEICE, myThid)
143 _EXCH_XY_R8(DXTICE, myThid)
144 _EXCH_XY_R8(DXUICE, myThid)
145 _EXCH_XY_R8(DYTICE, myThid)
146 _EXCH_XY_R8(DYUICE, myThid)
147
148 myIter=0
149 CALL PLOT_FIELD_XYRL( CSTICE , 'Current CSTICE ' ,
150 & myIter, myThid )
151 CALL PLOT_FIELD_XYRL( CSUICE , 'Current CSUICE ' ,
152 & myIter, myThid )
153 CALL PLOT_FIELD_XYRL( TNGTICE , 'Current TNGTICE ' ,
154 & myIter, myThid )
155 CALL PLOT_FIELD_XYRL( TNGICE , 'Current TNGICE ' ,
156 & myIter, myThid )
157 CALL PLOT_FIELD_XYRL( SINEICE , 'Current SINEICE ' ,
158 & myIter, myThid )
159 CALL PLOT_FIELD_XYRL( DXTICE , 'Current DXTICE ' ,
160 & myIter, myThid )
161 CALL PLOT_FIELD_XYRL( DXUICE , 'Current DXUICE ' ,
162 & myIter, myThid )
163 CALL PLOT_FIELD_XYRL( DYTICE , 'Current DYTICE ' ,
164 & myIter, myThid )
165 CALL PLOT_FIELD_XYRL( DYUICE , 'Current DYUICE ' ,
166 & myIter, myThid )
167 CALL PLOT_FIELD_XYRL( HEFFM , 'Current HEFFM ' ,
168 & myIter, myThid )
169 CALL PLOT_FIELD_XYRL( UVM , 'Current UVM ' ,
170 & myIter, myThid )
171
172 C-- Set model variables to initial/restart conditions
173 IF ( nIter0 .NE. 0 ) THEN
174 CALL SEAICE_READ_PICKUP ( myThid )
175 ELSE
176 DO bj=myByLo(myThid),myByHi(myThid)
177 DO bi=myBxLo(myThid),myBxHi(myThid)
178 DO j=1-OLy,sNy+OLy
179 DO i=1-OLx,sNx+OLx
180 HSNOW(I,J,bi,bj)=0.2 _d 0
181 YNEG(I,J,bi,bj)=ZERO
182 TMIX(I,J,bi,bj)=TICE(I,J,bi,bj)
183 DO k=1,3
184 HEFF(I,J,k,bi,bj)=SEAICE_initialHEFF
185 AREA(I,J,k,bi,bj)=HEFFM(i,j,bi,bj)
186 UICE(I,J,k,bi,bj)=ZERO
187 VICE(I,J,k,bi,bj)=ZERO
188 ENDDO
189 ENDDO
190 ENDDO
191 ENDDO
192 ENDDO
193 ENDIF
194
195 C-- Read initial sea-ice thickness from file if available.
196 IF ( HeffFile .NE. ' ' ) THEN
197 _BEGIN_MASTER( myThid )
198 CALL READ_FLD_XY_RL( HeffFile, ' ', ZETA, 0, myThid )
199 _END_MASTER(myThid)
200 _EXCH_XY_R8(ZETA,myThid)
201 DO bj=myByLo(myThid),myByHi(myThid)
202 DO bi=myBxLo(myThid),myBxHi(myThid)
203 DO j=1-OLy,sNy+OLy
204 DO i=1-OLx,sNx+OLx
205 DO k=1,3
206 HEFF(I,J,k,bi,bj) = MAX(ZETA(i,j,bi,bj),ZERO)
207 IF ( ZETA(i,j,bi,bj).EQ.ZERO )
208 & AREA(I,J,k,bi,bj) = ZERO
209 ENDDO
210 ENDDO
211 ENDDO
212 ENDDO
213 ENDDO
214 ENDIF
215
216 C--- Complete initialization
217 DO bj=myByLo(myThid),myByHi(myThid)
218 DO bi=myBxLo(myThid),myBxHi(myThid)
219 DO j=1-OLy,sNy+OLy
220 DO i=1-OLx,sNx+OLx
221 ZETA(I,J,bi,bj)=HEFF(I,J,1,bi,bj)*(1.0 _d 11)
222 ETA(I,J,bi,bj)=ZETA(I,J,bi,bj)/4.0 _d 0
223 ENDDO
224 ENDDO
225 ENDDO
226 ENDDO
227
228 #endif /* ALLOW_SEAICE */
229
230 RETURN
231 END

  ViewVC Help
Powered by ViewVC 1.1.22