ocnice_prep 1.14.0
Loading...
Searching...
No Matches
ocncalc_mod.F90
Go to the documentation of this file.
1
7module ocncalc_mod
8
9 use netcdf
10 use init_mod , only : nlevs, nxr, nyr
11 use init_mod, only : debug, logunit
12
13 use arrays_mod, only : hmin, maskspval, b3d, nbilin3d, rgb3d, eta
14 use utils_mod, only : getfield
15
16 implicit none
17
18 private
19
20 public calc_eta
21 public vfill
22
23contains
31 subroutine calc_eta(fname,dims,bathy)
32
33 character(len=*), intent(in) :: fname
34 integer, intent(in) :: dims(:)
35 real(kind=8), intent(in) :: bathy(:)
36
37 ! local variables
38 integer :: i,k
39 real(kind=8) :: denom
40 real(kind=8), allocatable, dimension(:) :: ssh,dilate
41 real(kind=8), allocatable, dimension(:,:) :: h
42 real(kind=8), allocatable, dimension(:,:) :: etmp
43 character(len=20) :: subname = 'calc_eta'
44 !----------------------------------------------------------------------------
45
46 if (debug)write(logunit,'(a)')'enter '//trim(subname)
47
48 allocate(ssh(dims(1)*dims(2))); ssh = 0.0
49 allocate(dilate(dims(1)*dims(2))); dilate = 0.0
50 allocate(h(dims(3),dims(1)*dims(2))); h = 0.0
51 ! note nlevs+1 for local array
52 allocate(etmp(dims(3)+1,dims(1)*dims(2))); etmp = 0.0
53
54 call getfield(trim(fname), 'sfc', (/dims(1),dims(2)/), ssh)
55 call getfield(trim(fname), 'h', (/dims(1),dims(2),dims(3)/), h)
56
57 ! ref: MOM_interface_heights.F90, SR find_eta_3d.F90, Boussinesq
58 etmp(dims(3)+1,:) = -bathy(:)
59 do k=dims(3),1,-1
60 etmp(k,:) = etmp(k+1,:) + h(k,:)
61 enddo
62
63 do i = 1,dims(1)*dims(2)
64 denom = etmp(1,i) + bathy(i)
65 if (denom .ne. 0.0) then
66 dilate(i) = (ssh(i) + bathy(i)) / (etmp(1,i) + bathy(i))
67 end if
68 end do
69
70 eta = 0.0
71 do k = 1,dims(3)
72 eta(k,:) = dilate(:)*(etmp(k,:) + bathy(:)) - bathy(:)
73 end do
74
75 if (debug)write(logunit,'(a)')'exit '//trim(subname)
76 end subroutine calc_eta
77
81 subroutine vfill()
82
83 ! local variables
84 integer :: n,i,k
85 integer :: idx1, klast
86 character(len=20) :: subname = 'vfill'
87 !----------------------------------------------------------------------------
88
89 if (debug)write(logunit,'(a)')'enter '//trim(subname)
90
91 idx1 = 0
92 do n = 1,nbilin3d
93 if (trim(b3d(n)%var_name) .eq. 'h') idx1 = n
94 end do
95
96 do i = 1,nxr*nyr
97 klast = nlevs
98 do k = 1,nlevs
99 if (rgb3d(idx1,k,i) .lt. maskspval)klast = k
100 end do
101 do n = 1,nbilin3d
102 do k = klast+1,nlevs
103 if (trim(b3d(n)%var_name) .eq. 'h') then
104 rgb3d(n,k,i) = hmin
105 else
106 rgb3d(n,k,i) = rgb3d(n,klast,i)
107 end if
108 end do
109 end do
110 end do
111
112 if (debug)write(logunit,'(a)')'exit '//trim(subname)
113 end subroutine vfill
114end module ocncalc_mod