1 |
C $Header: $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "GGL90_OPTIONS.h" |
5 |
|
6 |
CBOP |
7 |
C !ROUTINE: GGL90_SOLVE |
8 |
C !INTERFACE: |
9 |
SUBROUTINE GGL90_SOLVE( bi, bj, iMin, iMax, jMin, jMax, |
10 |
I a, b, c, |
11 |
U gXNm1, |
12 |
I myThid ) |
13 |
C !DESCRIPTION: \bv |
14 |
C *==========================================================* |
15 |
C | S/R GGL90_SOLVE |
16 |
C | o Solve implicit diffusion equation for vertical |
17 |
C | diffusivity for turbulent kinetic energy. |
18 |
C | o Tridiagonal matrix solver. |
19 |
C *==========================================================* |
20 |
C \ev |
21 |
|
22 |
C !USES: |
23 |
IMPLICIT NONE |
24 |
C == Global data == |
25 |
#include "SIZE.h" |
26 |
#include "EEPARAMS.h" |
27 |
#include "PARAMS.h" |
28 |
#include "GRID.h" |
29 |
CML#ifdef ALLOW_AUTODIFF_TAMC |
30 |
CML#include "tamc_keys.h" |
31 |
CML#endif |
32 |
|
33 |
C !INPUT/OUTPUT PARAMETERS: |
34 |
C == Routine Arguments == |
35 |
INTEGER bi,bj,iMin,iMax,jMin,jMax |
36 |
INTEGER myThid |
37 |
_RL a(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
38 |
_RL b(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
39 |
_RL c(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
40 |
_RL gXnm1(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
41 |
|
42 |
#ifdef ALLOW_GGL90 |
43 |
C !LOCAL VARIABLES: |
44 |
C == Local variables == |
45 |
INTEGER i,j,k |
46 |
_RL gYnm1(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
47 |
_RL bet(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
48 |
_RL gam(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
49 |
CEOP |
50 |
|
51 |
cph( |
52 |
cph Not good for TAF: may create irreducible control flow graph |
53 |
cph IF (Nr.LE.1) RETURN |
54 |
cph) |
55 |
|
56 |
C-- Initialise |
57 |
DO k=1,Nr |
58 |
DO j=jMin,jMax |
59 |
DO i=iMin,iMax |
60 |
gYNm1(i,j,k) = 0. _d 0 |
61 |
ENDDO |
62 |
ENDDO |
63 |
ENDDO |
64 |
|
65 |
C-- Old and new gam, bet are the same |
66 |
DO k=1,Nr |
67 |
DO j=jMin,jMax |
68 |
DO i=iMin,iMax |
69 |
bet(i,j,k) = 0. _d 0 |
70 |
gam(i,j,k) = 0. _d 0 |
71 |
ENDDO |
72 |
ENDDO |
73 |
ENDDO |
74 |
|
75 |
C-- Only need do anything if Nr>1 |
76 |
IF (Nr.GT.1) THEN |
77 |
|
78 |
k = 1 |
79 |
C-- Beginning of forward sweep (top level) |
80 |
DO j=jMin,jMax |
81 |
DO i=iMin,iMax |
82 |
IF (b(i,j,1).NE.0.) bet(i,j,1) = 1. _d 0 / b(i,j,1) |
83 |
ENDDO |
84 |
ENDDO |
85 |
|
86 |
ENDIF |
87 |
|
88 |
C-- Middle of forward sweep |
89 |
IF (Nr.GE.2) THEN |
90 |
|
91 |
CADJ loop = sequential |
92 |
DO k=2,Nr |
93 |
|
94 |
DO j=jMin,jMax |
95 |
DO i=iMin,iMax |
96 |
gam(i,j,k) = c(i,j,k-1)*bet(i,j,k-1) |
97 |
IF ( ( b(i,j,k) - a(i,j,k)*gam(i,j,k) ) .NE. 0.) |
98 |
& bet(i,j,k) = 1. _d 0 / ( b(i,j,k) - a(i,j,k)*gam(i,j,k) ) |
99 |
ENDDO |
100 |
ENDDO |
101 |
|
102 |
ENDDO |
103 |
|
104 |
ENDIF |
105 |
|
106 |
|
107 |
DO j=jMin,jMax |
108 |
DO i=iMin,iMax |
109 |
gYNm1(i,j,1) = gXNm1(i,j,1)*bet(i,j,1) |
110 |
ENDDO |
111 |
ENDDO |
112 |
DO k=2,Nr |
113 |
DO j=jMin,jMax |
114 |
DO i=iMin,iMax |
115 |
gYnm1(i,j,k) = bet(i,j,k)* |
116 |
& (gXNm1(i,j,k) - a(i,j,k)*gYnm1(i,j,k-1)) |
117 |
ENDDO |
118 |
ENDDO |
119 |
ENDDO |
120 |
|
121 |
|
122 |
C-- Backward sweep |
123 |
CADJ loop = sequential |
124 |
DO k=Nr-1,1,-1 |
125 |
DO j=jMin,jMax |
126 |
DO i=iMin,iMax |
127 |
gYnm1(i,j,k)=gYnm1(i,j,k) - gam(i,j,k+1)*gYnm1(i,j,k+1) |
128 |
ENDDO |
129 |
ENDDO |
130 |
ENDDO |
131 |
|
132 |
DO k=1,Nr |
133 |
DO j=jMin,jMax |
134 |
DO i=iMin,iMax |
135 |
gXNm1(i,j,k)=gYnm1(i,j,k) |
136 |
ENDDO |
137 |
ENDDO |
138 |
ENDDO |
139 |
|
140 |
#endif /* ALLOW_GGL90 */ |
141 |
RETURN |
142 |
END |
143 |
|