-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathftst_find_angq.F90
350 lines (314 loc) · 15.7 KB
/
ftst_find_angq.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
!> @file
!! @brief unit test for angles
!! @author [email protected]
!!
!! Given a 5x5 block of coords, calculate the angle, anglet and angchk
!! and compare against values from MOM6 and CICE6 history output. The
!! 3 5x5 blocks are for the 1deg tripole grid. The first two blocks contain
!! the two polar points. The final block is located in the lower right quad
!! of the tripole region (i=300:304,j=306:310). Because the angle calculations
!! involve points outside of the 5x5 blocks (eg, i+1), not all 5x5 points are
!! checked.
!!
!! @author [email protected]
program ftst_find_angq
use gengrid_kinds, only : dbl_kind, real_kind
use angles , only : find_angq, find_ang, find_angchk
implicit none
integer , parameter :: ni = 5, nj = 5, nblocks=3
! verification values from model history files, mx100
! blocks 1 and 2 center on the left and right tripole
! block 3 is in lower right quad of the arctic at (300:304,306:310)
real(dbl_kind), dimension(ni,nj,nblocks) :: cice6_angle
real(dbl_kind), dimension(ni,nj,nblocks) :: cice6_anglet
real(dbl_kind), dimension(ni,nj,nblocks) :: mom6_sinrot
! block grid point values
real(dbl_kind), dimension(ni,nj,nblocks) :: lonct, lonbu, latbu
! test values
real(dbl_kind), dimension(ni,nj,nblocks) :: anglet, angle, angchk
real(dbl_kind), dimension(ni,nblocks) :: xangCt
real(dbl_kind) :: maxval, diff
character(len=1) :: ck
integer :: i,j,i1,i2,j1,j2,k
integer :: ipole(nblocks)
logical :: debug = .true.
!CICE6
data ((cice6_angle(i,j,1), i=1,ni), j=1,nj) / &
0.55923, 0.29966, 0.00000, -0.29966, -0.55923, &
0.70064, 0.38859, 0.00000, -0.38859, -0.70064, &
0.90943, 0.53371, 0.00000, -0.53371, -0.90943, &
1.20822, 0.75580, 0.00000, -0.75580, -1.20822, &
1.38117, 0.89328, 0.00000, -0.89328, -1.38117/
data ((cice6_angle(i,j,2), i=1,ni), j=1,nj) / &
0.55923, 0.29966, 0.00000, -0.29966, -0.55923, &
0.70064, 0.38859, 0.00000, -0.38859, -0.70064, &
0.90943, 0.53371, 0.00000, -0.53371, -0.90943, &
1.20822, 0.75580, 0.00000, -0.75580, -1.20822, &
1.38117, 0.89328, 0.00000, -0.89328, -1.38117/
data ((cice6_angle(i,j,3), i=1,ni), j=1,nj) / &
-1.14949, -1.16044, -1.17076, -1.18049, -1.18968, &
-1.18188, -1.19217, -1.20185, -1.21097, -1.21957, &
-1.21425, -1.22383, -1.23284, -1.24131, -1.24930, &
-1.24649, -1.25534, -1.26364, -1.27145, -1.27880, &
-1.27854, -1.28661, -1.29419, -1.30130, -1.30799 /
data ((cice6_anglet(i,j,1), i=1,ni), j=1,nj) / &
0.60483, 0.38996, 0.13521, -0.13521, -0.38996, &
0.73204, 0.48692, 0.17198, -0.17198, -0.48692, &
0.90492, 0.63285, 0.23026, -0.23026, -0.63285, &
1.13597, 0.85126, 0.32134, -0.32134, -0.85126, &
1.34026, 1.05945, 0.41175, -0.41175, -1.05945/
data ((cice6_anglet(i,j,2), i=1,ni), j=1,nj) / &
0.60483, 0.38996, 0.13521, -0.13521, -0.38996, &
0.73204, 0.48692, 0.17198, -0.17198, -0.48692, &
0.90492, 0.63285, 0.23026, -0.23026, -0.63285, &
1.13597, 0.85126, 0.32134, -0.32134, -0.85126, &
1.34026, 1.05945, 0.41175, -0.41175, -1.05945/
data ((cice6_anglet(i,j,3), i=1,ni), j=1,nj) / &
-1.12735, -1.13895, -1.14989, -1.16020, -1.16995, &
-1.16004, -1.17100, -1.18130, -1.19102, -1.20018, &
-1.19278, -1.20303, -1.21267, -1.22174, -1.23029, &
-1.22546, -1.23498, -1.24391, -1.25231, -1.26021, &
-1.25800, -1.26675, -1.27495, -1.28264, -1.28988 /
! MOM6 sin_rot
data ((mom6_sinrot(i,j,1), i=1,ni), j=1,nj) / &
-0.57231, -0.38355, -0.13608, 0.13608, 0.38355, &
-0.67404, -0.47460, -0.17342, 0.17342, 0.47460, &
-0.79354, -0.60678, -0.23239, 0.23239, 0.60678, &
-0.91250, -0.79512, -0.32692, 0.32692, 0.79512, &
-0.98949, -0.97272, -0.43486, 0.43486, 0.97272/
data ((mom6_sinrot(i,j,2), i=1,ni), j=1,nj) / &
-0.57231, -0.38355, -0.13608, 0.13608, 0.38355, &
-0.67404, -0.47460, -0.17342, 0.17342, 0.47460, &
-0.79354, -0.60678, -0.23239, 0.23239, 0.60678, &
-0.91250, -0.79512, -0.32692, 0.32692, 0.79512, &
-0.98949, -0.97272, -0.43486, 0.43486, 0.97272/
data ((mom6_sinrot(i,j,3), i=1,ni), j=1,nj) / &
0.90334, 0.90826, 0.91278, 0.91694, 0.92079, &
0.91689, 0.92120, 0.92516, 0.92881, 0.93216, &
0.92946, 0.93320, 0.93662, 0.93976, 0.94264, &
0.94103, 0.94420, 0.94711, 0.94977, 0.95221, &
0.95154, 0.95419, 0.95661, 0.95882, 0.96085 /
! lon Ct,Bu and latBu
data ((lonct(i,j,1), i=1,ni), j=1,nj) / &
-245.26937, -232.98519, -218.04609, -201.95391, -187.01481, &
-252.74561, -239.00064, -220.46665, -199.53335, -180.99936, &
-262.76769, -248.28139, -224.73626, -195.26374, -171.71861, &
-275.86250, -263.23478, -234.03925, -185.96075, -156.76522, &
-291.65329, -286.25383, -263.72015, -156.27985, -133.74617 /
data ((lonbu(i,j,1), i=1,ni), j=1,nj) / &
-242.68153, -227.77943, -210.00000, -192.22057, -177.31847, &
-251.01234, -233.49458, -210.00000, -186.50542, -168.98766, &
-262.99894, -243.55611, -210.00000, -176.44389, -157.00106, &
-279.68531, -263.47519, -210.00000, -156.52481, -140.31469, &
-300.00000, -300.00000, -210.00000, -120.00000, -120.00000 /
data ((latbu(i,j,1), i=1,ni), j=1,nj) / &
88.99001, 89.10706, 89.14964, 89.10706, 88.99001, &
89.16921, 89.31628, 89.37292, 89.31628, 89.16921, &
89.31750, 89.50699, 89.58912, 89.50699, 89.31750, &
89.41886, 89.66093, 89.79818, 89.66093, 89.41886, &
89.45503, 89.72754, 90.00000, 89.72754, 89.45503 /
data ((lonct(i,j,2), i=1,ni), j=1,nj) / &
-65.26937, -52.98519, -38.04609, -21.95391, -7.01481, &
-72.74561, -59.00064, -40.46665, -19.53335, -0.99936, &
-82.76769, -68.28139, -44.73626, -15.26374, 8.28139, &
-95.86250, -83.23478, -54.03925, -5.96075, 23.23478, &
-111.65329, -106.25383, -83.72015, 23.72015, 46.25383 /
data ((lonbu(i,j,2), i=1,ni), j=1,nj) / &
-62.68153, -47.77943, -30.00000, -12.22057, 2.68153, &
-71.01234, -53.49458, -30.00000, -6.50542, 11.01234, &
-82.99894, -63.55611, -30.00000, 3.55611, 22.99894, &
-99.68531, -83.47519, -30.00000, 23.47519, 39.68531, &
-120.00000, -120.00000, -30.00000, 60.00000, 60.00000 /
data ((latbu(i,j,2), i=1,ni), j=1,nj) / &
88.99001, 89.10706, 89.14964, 89.10706, 88.99001, &
89.16921, 89.31628, 89.37292, 89.31628, 89.16921, &
89.31750, 89.50699, 89.58912, 89.50699, 89.31750, &
89.41886, 89.66093, 89.79818, 89.66093, 89.41886, &
89.45503, 89.72754, 90.00000, 89.72754, 89.45503 /
data ((lonct(i,j,3), i=1,ni), j=1,nj) / &
38.08144, 38.86965, 39.62025, 40.33605, 41.01961, &
39.67340, 40.41426, 41.11877, 41.78975, 42.42974, &
41.27339, 41.96446, 42.62075, 43.24505, 43.83985, &
42.87655, 43.51566, 44.12188, 44.69789, 45.24611, &
44.47789, 45.06321, 45.61778, 46.14417, 46.64470 /
data ((lonbu(i,j,3), i=1,ni), j=1,nj) / &
39.26336, 40.00933, 40.71974, 41.39725, 42.04431, &
40.83545, 41.53354, 42.19746, 42.82988, 43.43321, &
42.41197, 43.06002, 43.67560, 44.26131, 44.81949, &
43.98813, 44.58428, 45.14993, 45.68757, 46.19946, &
45.55910, 46.10182, 46.61624, 47.10473, 47.56942 /
data ((latbu(i,j,3), i=1,ni), j=1,nj) / &
80.97195, 80.70286, 80.43134, 80.15742, 79.88114, &
81.07710, 80.80479, 80.53019, 80.25334, 79.97425, &
81.17218, 80.89690, 80.61947, 80.33992, 80.05826, &
81.25741, 80.97942, 80.69942, 80.41741, 80.13341, &
81.33309, 81.05265, 80.77032, 80.48610, 80.19999 /
print *,"Starting test of cpld_gridgen routine find_angq"
angle = 0.0
anglet = 0.0
angchk = 0.0
ipole = 0
j = nj
do i = 1,ni
if(latBu(i,j,1) .eq. 90.00)ipole(1) = i
if(latBu(i,j,2) .eq. 90.00)ipole(2) = i
enddo
!------------------------------------------
! find anglet and test against mom6 values
!------------------------------------------
call find_ang((/2,5/),(/2,5/),lonBu(:,:,1),latBu(:,:,1),lonCt(:,:,1),anglet(:,:,1))
call find_ang((/2,5/),(/2,5/),lonBu(:,:,2),latBu(:,:,2),lonCt(:,:,2),anglet(:,:,2))
call find_ang((/2,5/),(/2,5/),lonBu(:,:,3),latBu(:,:,3),lonCt(:,:,3),anglet(:,:,3))
!i1 = 2; j1 = 2
!i2 = 5; j2 = 5
i1 = 1; i2 = ni
j1 = 1; j2 = nj
do k = 1,nblocks
write(ck,'(i1.1)')k
maxval = 1.0e-30
do j = j1,j2
do i = i1,i2
if (abs(anglet(i,j,k)) .gt. 0.0)then
diff = abs(anglet(i,j,k) - asin(mom6_sinrot(i,j,k)))
maxval = max(diff,maxval)
end if
end do
end do
call passfail(maxval, 1.0e-4, 'MOM6 anglet, block '//trim(ck))
end do
if (debug) then
print *,'anglet '
print *,'left'
do j = 1,nj
print '(5f15.6)',(anglet(i,j,1), i = 1,ni)
end do
print *,'right'
do j = 1,nj
print '(5f15.6)',(anglet(i,j,2), i = 1,ni)
end do
print *,'quad'
do j = 1,nj
print '(5f15.6)',(anglet(i,j,3), i = 1,ni)
end do
end if
!--------------------------------------------------------
!
!---------------------------------------------------------
xangCt = 0.0
do i = 2,ni
i2 = ipole(2)+(ipole(1)-i)+1
xangCt(i,1) = -anglet(i2,nj,1) ! angle changes sign across seam
xangCt(i,2) = -anglet(i2,nj,2)
xangCt(i,3) = anglet(i,nj,3)
print *,i,i2,xangCt(i,:)
end do
!-------------------------------------------
! find angle and test against cice6 values
!-------------------------------------------
call find_angq((/2,5/),(/2,5/),xangCt(:,1),anglet(:,:,1),angle(:,:,1))
call find_angq((/2,5/),(/2,5/),xangCt(:,2),anglet(:,:,2),angle(:,:,2))
call find_angq((/2,5/),(/2,4/),xangCt(:,3),anglet(:,:,3),angle(:,:,3))
angle(:,:,:) = -angle(:,:,:)
!i1 = 2; j1 = 2
!i2 = 4; j2 = 5
i1 = 1; i2 = ni
j1 = 1; j2 = nj
do k = 1,nblocks
write(ck,'(i1.1)')k
maxval = 1.0e-30
! note reduced check, xangCt(:,3) is zero
!if (k .eq. nblocks) j2 = 4
do j = j1,j2
do i = i1,i2
if (abs(angle(i,j,k)) .gt. 0.0) then
diff = abs(angle(i,j,k) - cice6_angle(i,j,k))
maxval = max(diff,maxval)
end if
end do
end do
call passfail(maxval, 1.0e-4, 'CICE6 angle, block '//trim(ck))
end do
if (debug) then
print *,'angle'
print *,'left'
do j = 1,nj
print '(5f15.6)',(angle(i,j,1), i = 1,ni)
end do
print *,'right'
do j = 1,nj
print '(5f15.6)',(angle(i,j,2), i = 1,ni)
end do
print *,'quad'
do j = 1,nj
print '(5f15.6)',(angle(i,j,3), i = 1,ni)
end do
end if
!-----------------------------------------------------------------------
! find anglet calculated by CICE and test against cice6 and mom6 values
!-----------------------------------------------------------------------
call find_angchk((/3,4/),(/3,5/),angle(:,:,1),angchk(:,:,1))
call find_angchk((/3,4/),(/3,5/),angle(:,:,2),angchk(:,:,2))
call find_angchk((/3,4/),(/3,4/),angle(:,:,3),angchk(:,:,3))
i1 = 1; i2 = ni
j1 = 1; j2 = nj
!i1 = 3; j1 = 3
!i2 = 4; j2 = 5
do k = 1,nblocks
write(ck,'(i1.1)')k
maxval = 1.0e-30
!if (k .eq. nblocks) j2 = 4
do j = j1,j2
do i = i1,i2
if (abs(angchk(i,j,k)) .gt. 0.0) then
diff = abs(angchk(i,j,k) - cice6_anglet(i,j,k))
maxval = max(diff,maxval)
end if
end do
end do
call passfail(maxval, 1.0e-4, 'angchk vs CICE6 anglet, block '//trim(ck))
end do
!i1 = 3; j1 = 3
!i2 = 4; j2 = 4
do k = 1,nblocks
write(ck,'(i1.1)')k
maxval = 1.0e-30
!if (k .eq. nblocks) j2 = 4
do j = j1,j2
do i = i1,i2
if (abs(angchk(i,j,k)) .gt. 0.0) then
diff = abs(-angchk(i,j,k) - asin(mom6_sinrot(i,j,k)))
maxval = max(diff,maxval)
end if
end do
end do
call passfail(maxval, 1.2e-2, 'angchk vs MOM6 anglet, block '//trim(ck))
end do
if (debug) then
print *,'angchk'
print *,'left'
do j = 1,nj
print '(5f15.6)',(angchk(i,j,1), i = 1,ni)
end do
print *,'right'
do j = 1,nj
print '(5f15.6)',(angchk(i,j,2), i = 1,ni)
end do
print *,'quad'
do j = 1,nj
print '(5f15.6)',(angchk(i,j,3), i = 1,ni)
end do
end if
end program
subroutine passfail(maxval,tolerance,msg)
use gengrid_kinds, only : dbl_kind
implicit none
real(dbl_kind), intent(in) :: maxval
real(dbl_kind), intent(in) :: tolerance
character(len=*), intent(in) :: msg
if (maxval .le. tolerance) then
print *,'SUCCESS! '//trim(msg),' ',maxval,tolerance
else
print *,'FAIL! '//trim(msg),' ',maxval,tolerance
!stop 1
endif
end subroutine passfail