18 use mom_io,
only : open_file
19 use mom_io,
only : append_file, ascii_file, multiple, single_file
20 use mom_time_manager,
only : time_type, get_time, get_date, set_date,
operator(-)
25 implicit none ;
private
27 #include <MOM_memory.h>
29 public write_u_accel, write_v_accel, pointaccel_init
33 character(len=200) :: u_trunc_file
35 character(len=200) :: v_trunc_file
39 integer :: cols_written
43 type(time_type),
pointer :: time => null()
49 real,
pointer,
dimension(:,:,:) :: &
56 u_accel_bt => null(), &
58 real,
pointer,
dimension(:,:,:) :: pbce => null()
69 subroutine write_u_accel(I, j, um, hin, ADp, CDp, dt, G, GV, US, CS, vel_rpt, str, a, hv)
70 integer,
intent(in) :: i
71 integer,
intent(in) :: j
75 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
77 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
83 real,
intent(in) :: dt
86 real,
intent(in) :: vel_rpt
87 real,
optional,
intent(in) :: str
89 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
90 optional,
intent(in) :: a
91 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
92 optional,
intent(in) :: hv
98 real :: inorm(szk_(g))
100 real :: h_scale, uh_scale
101 integer :: yr, mo, day, hr, minute, sec, yearday
104 logical :: do_k(szk_(g)+1)
105 logical :: prev_avail
108 angstrom = gv%Angstrom_H + gv%H_subroundoff
109 h_scale = gv%H_to_m ; uh_scale = gv%H_to_m
113 if (cs%cols_written < cs%max_writes)
then
114 cs%cols_written = cs%cols_written + 1
120 if (cs%u_file < 0)
then
121 if (len_trim(cs%u_trunc_file) < 1)
return
122 call open_file(cs%u_file, trim(cs%u_trunc_file), action=append_file, &
123 form=ascii_file, threading=multiple, fileset=single_file)
124 if (cs%u_file < 0)
then
125 call mom_error(note,
'Unable to open file '//trim(cs%u_trunc_file)//
'.')
131 prev_avail = (
associated(cs%u_prev) .and.
associated(cs%v_prev))
135 if (((max(cs%u_av(i,j,k),um(i,j,k)) >= vel_rpt) .or. &
136 (min(cs%u_av(i,j,k),um(i,j,k)) <= -vel_rpt)) .and. &
137 ((hin(i,j,k) + hin(i+1,j,k)) > 3.0*angstrom))
exit
141 if (((max(cs%u_av(i,j,k), um(i,j,k)) >= vel_rpt) .or. &
142 (min(cs%u_av(i,j,k), um(i,j,k)) <= -vel_rpt)) .and. &
143 ((hin(i,j,k) + hin(i+1,j,k)) > 3.0*angstrom))
exit
147 ks = 1; ke = nz;
write(file,
'("U: Unable to set ks & ke.")')
150 call get_date(cs%Time, yr, mo, day, hr, minute, sec)
151 call get_time((cs%Time - set_date(yr, 1, 1, 0, 0, 0)), sec, yearday)
152 write (file,
'(/,"--------------------------")')
153 write (file,
'(/,"Time ",i5,i4,F6.2," U-velocity violation at ",I4,": ",2(I3), &
154 & " (",F7.2," E "F7.2," N) Layers ",I3," to ",I3,". dt = ",1PG10.4)') &
155 yr, yearday, (real(sec)/3600.0), pe_here(), i, j, &
156 g%geoLonCu(i,j), g%geoLatCu(i,j), ks, ke, dt
158 if (ks <= gv%nk_rho_varies) ks = 1
160 if ((hin(i,j,k) + hin(i+1,j,k)) > 3.0*angstrom) do_k(k) = .true.
163 write(file,
'(/,"Layers:",$)')
164 do k=ks,ke ;
if (do_k(k))
write(file,
'(I10," ",$)') (k);
enddo
165 write(file,
'(/,"u(m): ",$)')
166 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (um(i,j,k));
enddo
168 write(file,
'(/,"u(mp): ",$)')
169 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (cs%u_prev(i,j,k));
enddo
171 write(file,
'(/,"u(3): ",$)')
172 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (cs%u_av(i,j,k));
enddo
174 write(file,
'(/,"CFL u: ",$)')
175 do k=ks,ke ;
if (do_k(k))
then
176 cfl = abs(um(i,j,k)) * dt * g%dy_Cu(i,j)
177 if (um(i,j,k) < 0.0)
then ; cfl = cfl * g%IareaT(i+1,j)
178 else ; cfl = cfl * g%IareaT(i,j) ;
endif
179 write(file,
'(ES10.3," ",$)') cfl
181 write(file,
'(/,"CFL0 u:",$)')
182 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
183 abs(um(i,j,k)) * dt * g%IdxCu(i,j) ;
enddo
186 write(file,
'(/,"du: ",$)')
187 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
188 ((um(i,j,k)-cs%u_prev(i,j,k)));
enddo
190 write(file,
'(/,"CAu: ",$)')
191 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (dt*adp%CAu(i,j,k));
enddo
192 write(file,
'(/,"PFu: ",$)')
193 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (dt*adp%PFu(i,j,k));
enddo
194 write(file,
'(/,"diffu: ",$)')
195 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (dt*us%s_to_T*adp%diffu(i,j,k));
enddo
197 if (
associated(adp%gradKEu))
then
198 write(file,
'(/,"KEu: ",$)')
199 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
200 (dt*adp%gradKEu(i,j,k));
enddo
202 if (
associated(adp%rv_x_v))
then
203 write(file,
'(/,"Coru: ",$)')
204 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
205 dt*(adp%CAu(i,j,k)-adp%rv_x_v(i,j,k));
enddo
207 if (
associated(adp%du_dt_visc))
then
208 write(file,
'(/,"ubv: ",$)')
209 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
210 (um(i,j,k)-dt*adp%du_dt_visc(i,j,k));
enddo
211 write(file,
'(/,"duv: ",$)')
212 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
213 (dt*adp%du_dt_visc(i,j,k));
enddo
215 if (
associated(adp%du_other))
then
216 write(file,
'(/,"du_other: ",$)')
217 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
218 (adp%du_other(i,j,k));
enddo
221 write(file,
'(/,"a: ",$)')
222 do k=ks,ke+1 ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') a(i,j,k)*us%Z_to_m*dt;
enddo
224 if (
present(hv))
then
225 write(file,
'(/,"hvel: ",$)')
226 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') hv(i,j,k);
enddo
228 write(file,
'(/,"Stress: ",ES10.3)') str
230 if (
associated(cs%u_accel_bt))
then
231 write(file,
'("dubt: ",$)')
232 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
233 (dt*cs%u_accel_bt(i,j,k)) ;
enddo
237 write(file,
'(/,"h--: ",$)')
238 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (h_scale*hin(i,j-1,k));
enddo
239 write(file,
'(/,"h+-: ",$)')
240 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (h_scale*hin(i+1,j-1,k));
enddo
241 write(file,
'(/,"h-0: ",$)')
242 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (h_scale*hin(i,j,k));
enddo
243 write(file,
'(/,"h+0: ",$)')
244 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (h_scale*hin(i+1,j,k));
enddo
245 write(file,
'(/,"h-+: ",$)')
246 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (h_scale*hin(i,j+1,k));
enddo
247 write(file,
'(/,"h++: ",$)')
248 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (h_scale*hin(i+1,j+1,k));
enddo
251 e(nz+1) = -us%Z_to_m*g%bathyT(i,j)
252 do k=nz,1,-1 ; e(k) = e(k+1) + h_scale*hin(i,j,k) ;
enddo
253 write(file,
'(/,"e-: ",$)')
254 write(file,
'(ES10.3," ",$)') e(ks)
255 do k=ks+1,ke+1 ;
if (do_k(k-1))
write(file,
'(ES10.3," ",$)') e(k);
enddo
257 e(nz+1) = -us%Z_to_m*g%bathyT(i+1,j)
258 do k=nz,1,-1 ; e(k) = e(k+1) + h_scale*hin(i+1,j,k) ;
enddo
259 write(file,
'(/,"e+: ",$)')
260 write(file,
'(ES10.3," ",$)') e(ks)
261 do k=ks+1,ke+1 ;
if (do_k(k-1))
write(file,
'(ES10.3," ",$)') e(k) ;
enddo
262 if (
associated(cs%T))
then
263 write(file,
'(/,"T-: ",$)')
264 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') cs%T(i,j,k);
enddo
265 write(file,
'(/,"T+: ",$)')
266 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') cs%T(i+1,j,k);
enddo
268 if (
associated(cs%S))
then
269 write(file,
'(/,"S-: ",$)')
270 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') cs%S(i,j,k);
enddo
271 write(file,
'(/,"S+: ",$)')
272 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') cs%S(i+1,j,k);
enddo
276 write(file,
'(/,"v--: ",$)')
277 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (cs%v_prev(i,j-1,k));
enddo
278 write(file,
'(/,"v-+: ",$)')
279 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (cs%v_prev(i,j,k));
enddo
280 write(file,
'(/,"v+-: ",$)')
281 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (cs%v_prev(i+1,j-1,k));
enddo
282 write(file,
'(/,"v++: ",$)')
283 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (cs%v_prev(i+1,j,k));
enddo
286 write(file,
'(/,"vh--: ",$)')
287 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
288 (uh_scale*cdp%vh(i,j-1,k)*g%IdxCv(i,j-1));
enddo
289 write(file,
'(/," vhC--:",$)')
290 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
291 (0.5*cs%v_av(i,j-1,k)*h_scale*(hin(i,j-1,k) + hin(i,j,k)));
enddo
293 write(file,
'(/," vhCp--:",$)')
294 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
295 (0.5*cs%v_prev(i,j-1,k)*h_scale*(hin(i,j-1,k) + hin(i,j,k)));
enddo
298 write(file,
'(/,"vh-+: ",$)')
299 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
300 (uh_scale*cdp%vh(i,j,k)*g%IdxCv(i,j));
enddo
301 write(file,
'(/," vhC-+:",$)')
302 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
303 (0.5*cs%v_av(i,j,k)*h_scale*(hin(i,j,k) + hin(i,j+1,k)));
enddo
305 write(file,
'(/," vhCp-+:",$)')
306 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
307 (0.5*cs%v_prev(i,j,k)*h_scale*(hin(i,j,k) + hin(i,j+1,k)));
enddo
310 write(file,
'(/,"vh+-: ",$)')
311 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
312 (uh_scale*cdp%vh(i+1,j-1,k)*g%IdxCv(i+1,j-1));
enddo
313 write(file,
'(/," vhC+-:",$)')
314 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
315 (0.5*cs%v_av(i+1,j-1,k)*h_scale*(hin(i+1,j-1,k) + hin(i+1,j,k)));
enddo
317 write(file,
'(/," vhCp+-:",$)')
318 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
319 (0.5*cs%v_prev(i+1,j-1,k)*h_scale*(hin(i+1,j-1,k) + hin(i+1,j,k)));
enddo
322 write(file,
'(/,"vh++: ",$)')
323 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
324 (uh_scale*cdp%vh(i+1,j,k)*g%IdxCv(i+1,j));
enddo
325 write(file,
'(/," vhC++:",$)')
326 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
327 (0.5*cs%v_av(i+1,j,k)*h_scale*(hin(i+1,j,k) + hin(i+1,j+1,k)));
enddo
329 write(file,
'(/," vhCp++:",$)')
330 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
331 (0.5*cs%v_av(i+1,j,k)*h_scale*(hin(i+1,j,k) + hin(i+1,j+1,k)));
enddo
334 write(file,
'(/,"D: ",2(ES10.3))') us%Z_to_m*g%bathyT(i,j),us%Z_to_m*g%bathyT(i+1,j)
339 du = um(i,j,k)-cs%u_prev(i,j,k)
340 if (abs(du) < 1.0e-6) du = 1.0e-6
344 write(file,
'(2/,"Norm: ",$)')
345 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') (1.0/inorm(k));
enddo
347 write(file,
'(/,"du: ",$)')
348 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
349 ((um(i,j,k)-cs%u_prev(i,j,k))*inorm(k));
enddo
351 write(file,
'(/,"CAu: ",$)')
352 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
353 (dt*adp%CAu(i,j,k)*inorm(k));
enddo
355 write(file,
'(/,"PFu: ",$)')
356 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
357 (dt*adp%PFu(i,j,k)*inorm(k));
enddo
359 write(file,
'(/,"diffu: ",$)')
360 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
361 (dt*us%s_to_T*adp%diffu(i,j,k)*inorm(k));
enddo
363 if (
associated(adp%gradKEu))
then
364 write(file,
'(/,"KEu: ",$)')
365 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
366 (dt*adp%gradKEu(i,j,k)*inorm(k));
enddo
368 if (
associated(adp%rv_x_v))
then
369 write(file,
'(/,"Coru: ",$)')
370 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
371 dt*(adp%CAu(i,j,k)-adp%rv_x_v(i,j,k))*inorm(k);
enddo
373 if (
associated(adp%du_dt_visc))
then
374 write(file,
'(/,"duv: ",$)')
375 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
376 (dt*adp%du_dt_visc(i,j,k))*inorm(k);
enddo
378 if (
associated(adp%du_other))
then
379 write(file,
'(/,"du_other: ",$)')
380 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
381 (adp%du_other(i,j,k))*inorm(k);
enddo
383 if (
associated(cs%u_accel_bt))
then
384 write(file,
'(/,"dubt: ",$)')
385 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
386 (dt*cs%u_accel_bt(i,j,k)*inorm(k)) ;
enddo
395 end subroutine write_u_accel
400 subroutine write_v_accel(i, J, vm, hin, ADp, CDp, dt, G, GV, US, CS, vel_rpt, str, a, hv)
401 integer,
intent(in) :: i
402 integer,
intent(in) :: j
406 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
408 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
414 real,
intent(in) :: dt
417 real,
intent(in) :: vel_rpt
418 real,
optional,
intent(in) :: str
420 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
421 optional,
intent(in) :: a
422 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
423 optional,
intent(in) :: hv
429 real :: inorm(szk_(g))
431 real :: h_scale, uh_scale
432 integer :: yr, mo, day, hr, minute, sec, yearday
435 logical :: do_k(szk_(g)+1)
436 logical :: prev_avail
439 angstrom = gv%Angstrom_H + gv%H_subroundoff
440 h_scale = gv%H_to_m ; uh_scale = gv%H_to_m
444 if (cs%cols_written < cs%max_writes)
then
445 cs%cols_written = cs%cols_written + 1
451 if (cs%v_file < 0)
then
452 if (len_trim(cs%v_trunc_file) < 1)
return
453 call open_file(cs%v_file, trim(cs%v_trunc_file), action=append_file, &
454 form=ascii_file, threading=multiple, fileset=single_file)
455 if (cs%v_file < 0)
then
456 call mom_error(note,
'Unable to open file '//trim(cs%v_trunc_file)//
'.')
462 prev_avail = (
associated(cs%u_prev) .and.
associated(cs%v_prev))
465 if (((max(cs%v_av(i,j,k), vm(i,j,k)) >= vel_rpt) .or. &
466 (min(cs%v_av(i,j,k), vm(i,j,k)) <= -vel_rpt)) .and. &
467 ((hin(i,j,k) + hin(i,j+1,k)) > 3.0*angstrom))
exit
471 if (((max(cs%v_av(i,j,k), vm(i,j,k)) >= vel_rpt) .or. &
472 (min(cs%v_av(i,j,k), vm(i,j,k)) <= -vel_rpt)) .and. &
473 ((hin(i,j,k) + hin(i,j+1,k)) > 3.0*angstrom))
exit
477 ks = 1; ke = nz;
write(file,
'("V: Unable to set ks & ke.")')
480 call get_date(cs%Time, yr, mo, day, hr, minute, sec)
481 call get_time((cs%Time - set_date(yr, 1, 1, 0, 0, 0)), sec, yearday)
482 write (file,
'(/,"--------------------------")')
483 write (file,
'(/,"Time ",i5,i4,F6.2," V-velocity violation at ",I4,": ",2(I3), &
484 & " (",F7.2," E ",F7.2," N) Layers ",I3," to ",I3,". dt = ",1PG10.4)') &
485 yr, yearday, (real(sec)/3600.0), pe_here(), i, j, &
486 g%geoLonCv(i,j), g%geoLatCv(i,j), ks, ke, dt
488 if (ks <= gv%nk_rho_varies) ks = 1
490 if ((hin(i,j,k) + hin(i,j+1,k)) > 3.0*angstrom) do_k(k) = .true.
493 write(file,
'(/,"Layers:",$)')
494 do k=ks,ke ;
if (do_k(k))
write(file,
'(I10," ",$)') (k);
enddo
495 write(file,
'(/,"v(m): ",$)')
496 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (vm(i,j,k));
enddo
499 write(file,
'(/,"v(mp): ",$)')
500 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (cs%v_prev(i,j,k));
enddo
503 write(file,
'(/,"v(3): ",$)')
504 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (cs%v_av(i,j,k));
enddo
505 write(file,
'(/,"CFL v: ",$)')
506 do k=ks,ke ;
if (do_k(k))
then
507 cfl = abs(vm(i,j,k)) * dt * g%dx_Cv(i,j)
508 if (vm(i,j,k) < 0.0)
then ; cfl = cfl * g%IareaT(i,j+1)
509 else ; cfl = cfl * g%IareaT(i,j) ;
endif
510 write(file,
'(ES10.3," ",$)') cfl
512 write(file,
'(/,"CFL0 v:",$)')
513 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
514 abs(vm(i,j,k)) * dt * g%IdyCv(i,j) ;
enddo
517 write(file,
'(/,"dv: ",$)')
518 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
519 ((vm(i,j,k)-cs%v_prev(i,j,k)));
enddo
522 write(file,
'(/,"CAv: ",$)')
523 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (dt*adp%CAv(i,j,k));
enddo
525 write(file,
'(/,"PFv: ",$)')
526 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (dt*adp%PFv(i,j,k));
enddo
528 write(file,
'(/,"diffv: ",$)')
529 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') (dt*us%s_to_T*adp%diffv(i,j,k));
enddo
531 if (
associated(adp%gradKEv))
then
532 write(file,
'(/,"KEv: ",$)')
533 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
534 (dt*adp%gradKEv(i,j,k));
enddo
536 if (
associated(adp%rv_x_u))
then
537 write(file,
'(/,"Corv: ",$)')
538 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
539 dt*(adp%CAv(i,j,k)-adp%rv_x_u(i,j,k));
enddo
541 if (
associated(adp%dv_dt_visc))
then
542 write(file,
'(/,"vbv: ",$)')
543 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
544 (vm(i,j,k)-dt*adp%dv_dt_visc(i,j,k));
enddo
546 write(file,
'(/,"dvv: ",$)')
547 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
548 (dt*adp%dv_dt_visc(i,j,k));
enddo
550 if (
associated(adp%dv_other))
then
551 write(file,
'(/,"dv_other: ",$)')
552 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
553 (adp%dv_other(i,j,k));
enddo
556 write(file,
'(/,"a: ",$)')
557 do k=ks,ke+1 ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') a(i,j,k)*us%Z_to_m*dt;
enddo
559 if (
present(hv))
then
560 write(file,
'(/,"hvel: ",$)')
561 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') hv(i,j,k);
enddo
563 write(file,
'(/,"Stress: ",ES10.3)') str
565 if (
associated(cs%v_accel_bt))
then
566 write(file,
'("dvbt: ",$)')
567 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
568 (dt*cs%v_accel_bt(i,j,k)) ;
enddo
572 write(file,
'("h--: ",$)')
573 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') h_scale*hin(i-1,j,k);
enddo
574 write(file,
'(/,"h0-: ",$)')
575 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') h_scale*hin(i,j,k);
enddo
576 write(file,
'(/,"h+-: ",$)')
577 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') h_scale*hin(i+1,j,k);
enddo
578 write(file,
'(/,"h-+: ",$)')
579 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') h_scale*hin(i-1,j+1,k);
enddo
580 write(file,
'(/,"h0+: ",$)')
581 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') h_scale*hin(i,j+1,k);
enddo
582 write(file,
'(/,"h++: ",$)')
583 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') h_scale*hin(i+1,j+1,k);
enddo
585 e(nz+1) = -us%Z_to_m*g%bathyT(i,j)
586 do k=nz,1,-1 ; e(k) = e(k+1) + h_scale*hin(i,j,k);
enddo
587 write(file,
'(/,"e-: ",$)')
588 write(file,
'(ES10.3," ",$)') e(ks)
589 do k=ks+1,ke+1 ;
if (do_k(k-1))
write(file,
'(ES10.3," ",$)') e(k);
enddo
591 e(nz+1) = -us%Z_to_m*g%bathyT(i,j+1)
592 do k=nz,1,-1 ; e(k) = e(k+1) + h_scale*hin(i,j+1,k) ;
enddo
593 write(file,
'(/,"e+: ",$)')
594 write(file,
'(ES10.3," ",$)') e(ks)
595 do k=ks+1,ke+1 ;
if (do_k(k-1))
write(file,
'(ES10.3," ",$)') e(k);
enddo
596 if (
associated(cs%T))
then
597 write(file,
'(/,"T-: ",$)')
598 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') cs%T(i,j,k);
enddo
599 write(file,
'(/,"T+: ",$)')
600 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') cs%T(i,j+1,k);
enddo
602 if (
associated(cs%S))
then
603 write(file,
'(/,"S-: ",$)')
604 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') cs%S(i,j,k);
enddo
605 write(file,
'(/,"S+: ",$)')
606 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') cs%S(i,j+1,k);
enddo
610 write(file,
'(/,"u--: ",$)')
611 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') cs%u_prev(i-1,j,k);
enddo
612 write(file,
'(/,"u-+: ",$)')
613 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') cs%u_prev(i-1,j+1,k);
enddo
614 write(file,
'(/,"u+-: ",$)')
615 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') cs%u_prev(i,j,k);
enddo
616 write(file,
'(/,"u++: ",$)')
617 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') cs%u_prev(i,j+1,k);
enddo
620 write(file,
'(/,"uh--: ",$)')
621 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
622 (uh_scale*cdp%uh(i-1,j,k)*g%IdyCu(i-1,j));
enddo
623 write(file,
'(/," uhC--: ",$)')
624 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
625 (cs%u_av(i-1,j,k) * h_scale*0.5*(hin(i-1,j,k) + hin(i,j,k)));
enddo
627 write(file,
'(/," uhCp--:",$)')
628 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
629 (cs%u_prev(i-1,j,k) * h_scale*0.5*(hin(i-1,j,k) + hin(i,j,k)));
enddo
632 write(file,
'(/,"uh-+: ",$)')
633 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
634 (uh_scale*cdp%uh(i-1,j+1,k)*g%IdyCu(i-1,j+1));
enddo
635 write(file,
'(/," uhC-+: ",$)')
636 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
637 (cs%u_av(i-1,j+1,k) * h_scale*0.5*(hin(i-1,j+1,k) + hin(i,j+1,k)));
enddo
639 write(file,
'(/," uhCp-+:",$)')
640 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
641 (cs%u_prev(i-1,j+1,k) * h_scale*0.5*(hin(i-1,j+1,k) + hin(i,j+1,k)));
enddo
644 write(file,
'(/,"uh+-: ",$)')
645 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
646 (uh_scale*cdp%uh(i,j,k)*g%IdyCu(i,j));
enddo
647 write(file,
'(/," uhC+-: ",$)')
648 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
649 (cs%u_av(i,j,k) * h_scale*0.5*(hin(i,j,k) + hin(i+1,j,k)));
enddo
651 write(file,
'(/," uhCp+-:",$)')
652 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
653 (cs%u_prev(i,j,k) * h_scale*0.5*(hin(i,j,k) + hin(i+1,j,k)));
enddo
656 write(file,
'(/,"uh++: ",$)')
657 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
658 (uh_scale*cdp%uh(i,j+1,k)*g%IdyCu(i,j+1));
enddo
659 write(file,
'(/," uhC++: ",$)')
660 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
661 (cs%u_av(i,j+1,k) * 0.5*h_scale*(hin(i,j+1,k) + hin(i+1,j+1,k)));
enddo
663 write(file,
'(/," uhCp++:",$)')
664 do k=ks,ke ;
if (do_k(k))
write(file,
'(ES10.3," ",$)') &
665 (cs%u_prev(i,j+1,k) * h_scale*0.5*(hin(i,j+1,k) + hin(i+1,j+1,k)));
enddo
668 write(file,
'(/,"D: ",2(ES10.3))') us%Z_to_m*g%bathyT(i,j),us%Z_to_m*g%bathyT(i,j+1)
673 dv = vm(i,j,k)-cs%v_prev(i,j,k)
674 if (abs(dv) < 1.0e-6) dv = 1.0e-6
678 write(file,
'(2/,"Norm: ",$)')
679 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') (1.0/inorm(k));
enddo
680 write(file,
'(/,"dv: ",$)')
681 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
682 ((vm(i,j,k)-cs%v_prev(i,j,k))*inorm(k));
enddo
683 write(file,
'(/,"CAv: ",$)')
684 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
685 (dt*adp%CAv(i,j,k)*inorm(k));
enddo
686 write(file,
'(/,"PFv: ",$)')
687 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
688 (dt*adp%PFv(i,j,k)*inorm(k));
enddo
689 write(file,
'(/,"diffv: ",$)')
690 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
691 (dt*us%s_to_T*adp%diffv(i,j,k)*inorm(k));
enddo
693 if (
associated(adp%gradKEu))
then
694 write(file,
'(/,"KEv: ",$)')
695 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
696 (dt*adp%gradKEv(i,j,k)*inorm(k));
enddo
698 if (
associated(adp%rv_x_u))
then
699 write(file,
'(/,"Corv: ",$)')
700 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
701 dt*(adp%CAv(i,j,k)-adp%rv_x_u(i,j,k))*inorm(k);
enddo
703 if (
associated(adp%dv_dt_visc))
then
704 write(file,
'(/,"dvv: ",$)')
705 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
706 (dt*adp%dv_dt_visc(i,j,k)*inorm(k));
enddo
708 if (
associated(adp%dv_other))
then
709 write(file,
'(/,"dv_other: ",$)')
710 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
711 (adp%dv_other(i,j,k)*inorm(k));
enddo
713 if (
associated(cs%v_accel_bt))
then
714 write(file,
'(/,"dvbt: ",$)')
715 do k=ks,ke ;
if (do_k(k))
write(file,
'(F10.6," ",$)') &
716 (dt*cs%v_accel_bt(i,j,k)*inorm(k)) ;
enddo
725 end subroutine write_v_accel
728 subroutine pointaccel_init(MIS, Time, G, param_file, diag, dirs, CS)
730 target,
intent(in) :: mis
733 type(time_type),
target,
intent(in) :: time
737 type(
diag_ctrl),
target,
intent(inout) :: diag
744 #include "version_variable.h"
745 character(len=40) :: mdl =
"MOM_PointAccel"
747 if (
associated(cs))
return
750 cs%diag => diag ; cs%Time => time
752 cs%T => mis%T ; cs%S => mis%S ; cs%pbce => mis%pbce
753 cs%u_accel_bt => mis%u_accel_bt ; cs%v_accel_bt => mis%v_accel_bt
754 cs%u_prev => mis%u_prev ; cs%v_prev => mis%v_prev
755 cs%u_av => mis%u_av;
if (.not.
associated(mis%u_av)) cs%u_av => mis%u(:,:,:)
756 cs%v_av => mis%v_av;
if (.not.
associated(mis%v_av)) cs%v_av => mis%v(:,:,:)
760 call get_param(param_file, mdl,
"U_TRUNC_FILE", cs%u_trunc_file, &
761 "The absolute path to the file where the accelerations "//&
762 "leading to zonal velocity truncations are written. \n"//&
763 "Leave this empty for efficiency if this diagnostic is "//&
764 "not needed.", default=
"", debuggingparam=.true.)
765 call get_param(param_file, mdl,
"V_TRUNC_FILE", cs%v_trunc_file, &
766 "The absolute path to the file where the accelerations "//&
767 "leading to meridional velocity truncations are written. \n"//&
768 "Leave this empty for efficiency if this diagnostic is "//&
769 "not needed.", default=
"", debuggingparam=.true.)
770 call get_param(param_file, mdl,
"MAX_TRUNC_FILE_SIZE_PER_PE", cs%max_writes, &
771 "The maximum number of colums of truncations that any PE "//&
772 "will write out during a run.", default=50, debuggingparam=.true.)
774 if (len_trim(dirs%output_directory) > 0)
then
775 if (len_trim(cs%u_trunc_file) > 0) &
776 cs%u_trunc_file = trim(dirs%output_directory)//trim(cs%u_trunc_file)
777 if (len_trim(cs%v_trunc_file) > 0) &
778 cs%v_trunc_file = trim(dirs%output_directory)//trim(cs%v_trunc_file)
779 call log_param(param_file, mdl,
"output_dir/U_TRUNC_FILE", cs%u_trunc_file)
780 call log_param(param_file, mdl,
"output_dir/V_TRUNC_FILE", cs%v_trunc_file)
782 cs%u_file = -1 ; cs%v_file = -1 ; cs%cols_written = 0
784 end subroutine pointaccel_init