UPP  001
 All Data Structures Files Functions Pages
CALMCVG.f
Go to the documentation of this file.
1 
33  SUBROUTINE calmcvg(Q1D,U1D,V1D,QCNVG)
34 
35 !
36 !
37 !
38  use masks, only: dx, dy, hbm2
39  use params_mod, only: d00, d25
40  use ctlblk_mod, only: jsta_2l, jend_2u, spval, jsta_m, jend_m, &
41  jsta_m2, jend_m2, im, jm
42  use gridspec_mod, only: gridtype
43 !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
44  implicit none
45 !
46 ! DECLARE VARIABLES.
47 !
48  REAL,dimension(IM,jsta_2l:jend_2u),intent(in) :: q1d, u1d, v1d
49  REAL,dimension(IM,jsta_2l:jend_2u),intent(inout) :: qcnvg
50 
51  REAL r2dy, r2dx
52  REAL, dimension(im,jsta_2l:jend_2u) :: uwnd, vwnd, qv
53  INTEGER ihe(jm),ihw(jm),ive(jm),ivw(jm)
54  integer i,j,ista,iend
55  real qvdy,qudx
56 !
57 !***************************************************************************
58 ! START CALMCVG HERE.
59 !
60 !
61 ! INITIALIZE MOISTURE CONVERGENCE ARRAY. LOAD TEMPORARY WIND ARRAYS.
62 !
63 !$omp parallel do private(i,j)
64  DO j=jsta_2l,jend_2u
65  DO i=1,im
66  IF(u1d(i,j)<spval)THEN
67  qcnvg(i,j) = 0.
68  ELSE
69  qcnvg(i,j) = spval
70  ENDIF
71  uwnd(i,j) = u1d(i,j)
72  vwnd(i,j) = v1d(i,j)
73  IF (uwnd(i,j) == spval) uwnd(i,j) = d00
74  IF (vwnd(i,j) == spval) vwnd(i,j) = d00
75  ENDDO
76  ENDDO
77 
78  CALL exch_f(q1d)
79  CALL exch_f(vwnd)
80 !
81  IF(gridtype == 'A')THEN
82 !$omp parallel do private(i,j,qudx,qvdy,r2dx,r2dy)
83  DO j=jsta_m,jend_m
84  DO i=2,im-1
85  IF(q1d(i,j+1)<spval.AND.q1d(i,j-1)<spval.AND. &
86  q1d(i+1,j)<spval.AND.q1d(i-1,j)<spval.AND. &
87  q1d(i,j)<spval) THEN
88  r2dx = 1./(2.*dx(i,j)) !MEB DX?
89  r2dy = 1./(2.*dy(i,j)) !MEB DY?
90  qudx = (q1d(i+1,j)*uwnd(i+1,j)-q1d(i-1,j)*uwnd(i-1,j))*r2dx
91  qvdy = (q1d(i,j+1)*vwnd(i,j+1)-q1d(i,j-1)*vwnd(i,j-1))*r2dy
92  qcnvg(i,j) = -(qudx + qvdy)
93  ELSE
94  qcnvg(i,j) = spval
95  ENDIF
96  ENDDO
97  ENDDO
98  ELSE IF(gridtype == 'E')THEN
99 
100  DO j=jsta_m,jend_m
101  ihe(j) = mod(j+1,2)
102  ihw(j) = ihe(j)-1
103  ive(j) = mod(j,2)
104  ivw(j) = ive(j)-1
105  END DO
106 
107 !$omp parallel do private(i,j,ista,iend)
108  DO j=jsta_m,jend_m
109  ista = 1+mod(j+1,2)
110  iend = im-mod(j,2)
111  DO i=ista,iend
112  IF(q1d(i,j-1)<spval.AND.q1d(i+ivw(j),j)<spval.AND.&
113  q1d(i+ive(j),j)<spval.AND.q1d(i,j+1)<spval) THEN
114  qv(i,j) = d25*(q1d(i,j-1)+q1d(i+ivw(j),j) &
115  +q1d(i+ive(j),j)+q1d(i,j+1))
116  ELSE
117  qv(i,j) = spval
118  ENDIF
119  END DO
120  END DO
121 
122  CALL exch_f(qv)
123 ! CALL EXCH_F(VWND)
124 
125 !
126 !$omp parallel do private(i,j,iend,qudx,qvdy,r2dx,r2dy)
127  DO j=jsta_m2,jend_m2
128  iend = im-1-mod(j,2)
129  DO i=2,iend
130  IF(qv(i+ihe(j),j)<spval.AND.uwnd(i+ihe(j),j)<spval.AND.&
131  qv(i+ihw(j),j)<spval.AND.uwnd(i+ihw(j),j)<spval.AND.&
132  qv(i,j)<spval.AND.qv(i,j-1)<spval.AND.qv(i,j+1)<spval) THEN
133  r2dx = 1./(2.*dx(i,j))
134  r2dy = 1./(2.*dy(i,j))
135  qudx = (qv(i+ihe(j),j)*uwnd(i+ihe(j),j) &
136  -qv(i+ihw(j),j)*uwnd(i+ihw(j),j))*r2dx
137  qvdy = (qv(i,j+1)*vwnd(i,j+1)-qv(i,j-1)*vwnd(i,j-1))*r2dy
138 
139  qcnvg(i,j) = -(qudx + qvdy) * hbm2(i,j)
140  ELSE
141  qcnvg(i,j) = spval
142  ENDIF
143  ENDDO
144  ENDDO
145  ELSE IF(gridtype=='B')THEN
146 
147  CALL exch_f(uwnd)
148 !
149 !$omp parallel do private(i,j,qudx,qvdy,r2dx,r2dy)
150  DO j=jsta_m,jend_m
151  DO i=2,im-1
152  IF(uwnd(i,j)<spval.AND.uwnd(i,j-1)<spval.AND.&
153  uwnd(i-1,j)<spval.AND.uwnd(i-1,j-1)<spval.AND.&
154  q1d(i,j)<spval.AND.q1d(i+1,j)<spval.AND.q1d(i-1,j)<spval.AND.&
155  vwnd(i,j)<spval.AND.vwnd(i-1,j)<spval.AND.&
156  vwnd(i,j-1)<spval.AND.vwnd(i-1,j-1)<spval.AND.&
157  q1d(i,j+1)<spval.AND.q1d(i,j-1)<spval) THEN
158  r2dx = 1./dx(i,j)
159  r2dy = 1./dy(i,j)
160  qudx=(0.5*(uwnd(i,j)+uwnd(i,j-1))*0.5*(q1d(i,j)+q1d(i+1,j)) &
161  -0.5*(uwnd(i-1,j)+uwnd(i-1,j-1))*0.5*(q1d(i,j)+q1d(i-1,j)))*r2dx
162  qvdy=(0.5*(vwnd(i,j)+vwnd(i-1,j))*0.5*(q1d(i,j)+q1d(i,j+1)) &
163  -0.5*(vwnd(i,j-1)+vwnd(i-1,j-1))*0.5*(q1d(i,j)+q1d(i,j-1)))*r2dy
164 
165  qcnvg(i,j) = -(qudx + qvdy)
166  ELSE
167  qcnvg(i,j) = spval
168  ENDIF
169 ! print*,'mcvg=',i,j,r2dx,r2dy,QCNVG(I,J)
170  ENDDO
171  ENDDO
172  ENDIF
173 !meb not sure about the indexing for the c-grid
174 !
175 ! END OF ROUTINE.
176 !
177  RETURN
178  END
179 
Definition: MASKS_mod.f:1
subroutine exch_f(a)
Definition: EXCH.f:65