1 | MODULE solver |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE solver *** |
---|
4 | !! Ocean solver : initialization of ocean solver |
---|
5 | !!===================================================================== |
---|
6 | |
---|
7 | !!---------------------------------------------------------------------- |
---|
8 | !! solver_init: solver initialization |
---|
9 | !!---------------------------------------------------------------------- |
---|
10 | !! * Modules used |
---|
11 | USE oce ! ocean dynamics and tracers variables |
---|
12 | USE dom_oce ! ocean space and time domain variables |
---|
13 | USE zdf_oce ! ocean vertical physics variables |
---|
14 | USE sol_oce ! solver variables |
---|
15 | USE solmat ! ??? |
---|
16 | USE solisl ! ??? |
---|
17 | USE obc_oce ! Lateral open boundary condition |
---|
18 | USE in_out_manager ! I/O manager |
---|
19 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
20 | USE lib_mpp |
---|
21 | USE dynspg_oce ! choice/control of key cpp for surface pressure gradient |
---|
22 | |
---|
23 | IMPLICIT NONE |
---|
24 | |
---|
25 | !!---------------------------------------------------------------------- |
---|
26 | !! OPA 9.0 , LOCEAN-IPSL (2005) |
---|
27 | !! $Header$ |
---|
28 | !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt |
---|
29 | !!---------------------------------------------------------------------- |
---|
30 | |
---|
31 | CONTAINS |
---|
32 | |
---|
33 | SUBROUTINE solver_init( kt ) |
---|
34 | !!---------------------------------------------------------------------- |
---|
35 | !! *** ROUTINE solver_init *** |
---|
36 | !! |
---|
37 | !! ** Purpose : Initialization for the solver of the elliptic equation: |
---|
38 | !! * default option: barotropic stream function system |
---|
39 | !! and islands initialization (if lk_isl=T) |
---|
40 | !! * lk_dynspg_flt = T : transport divergence system. No specific |
---|
41 | !! treatment of islands. |
---|
42 | !! |
---|
43 | !! ** Method : |
---|
44 | !! - Compute the local depth of the water column at u- and v-point |
---|
45 | !! (lk_dynspg_flt = T) or its inverse (lk_dynspg_rl = T). |
---|
46 | !! The local depth of the water column is computed by summing |
---|
47 | !! the vertical scale factors. For its inverse, the thickness of |
---|
48 | !! the first model level is imposed as lower bound. The inverse of |
---|
49 | !! this depth is THEN taken and masked, so that the inverse of the |
---|
50 | !! local depth is zero when the local depth is zero. |
---|
51 | !! - Construct the matrix of the elliptic system by a call to |
---|
52 | !! solmat.F routine. |
---|
53 | !! - island (if lk_isl=T) |
---|
54 | !! isl_dom: find islands from the bathymetry file |
---|
55 | !! isl_bsf: compute the island barotropic stream function |
---|
56 | !! isl_mat: compute the inverse island matrix |
---|
57 | !! set mbathy to the number of non-zero w-levels of a water |
---|
58 | !! column (the minimum value of mbathy is 2): |
---|
59 | !! mbathy = min( mbathy, 1 ) + 1 |
---|
60 | !! |
---|
61 | !! ** Action : - hur, hvr : masked inverse of the local depth at |
---|
62 | !! u- and v-point. (lk_dynspg_rl = T) |
---|
63 | !! - hu, hv : masked local depth at u- and v- points |
---|
64 | !! (lk_dynspg_flt = T) |
---|
65 | !! - c_solver_pt : nature of the gridpoint at which the |
---|
66 | !! solver is applied |
---|
67 | !! References : |
---|
68 | !! Jensen, 1986: adv. phys. oceanogr. num. mod.,ed. o brien,87-110. |
---|
69 | !! Madec & Marti, 1990: internal rep. LODYC, 90/03., 29pp. |
---|
70 | !! |
---|
71 | !! History : |
---|
72 | !! ! 90-10 (G. Madec) Original code |
---|
73 | !! ! 93-02 (O. Marti) |
---|
74 | !! ! 97-02 (G. Madec) local depth inverse computation |
---|
75 | !! ! 98-10 (G. Roullet, G. Madec) free surface |
---|
76 | !! 9.0 ! 03-07 (G. Madec) free form, F90 |
---|
77 | !! " ! 05-11 (V. Garnier) Surface pressure gradient organization |
---|
78 | !!---------------------------------------------------------------------- |
---|
79 | !! * Arguments |
---|
80 | INTEGER, INTENT(in) :: kt |
---|
81 | |
---|
82 | !! * Local declarations |
---|
83 | INTEGER :: ji, jj ! dummy loop indices |
---|
84 | CHARACTER(len=80) :: clname |
---|
85 | |
---|
86 | NAMELIST/namsol/ nsolv, nsol_arp, nmin, nmax, nmod, eps, resmax, sor, epsisl, nmisl, rnu |
---|
87 | !!---------------------------------------------------------------------- |
---|
88 | |
---|
89 | IF(lwp) WRITE(numout,*) |
---|
90 | IF(lwp) WRITE(numout,*) 'solver_init : solver to compute the surface pressure gradient' |
---|
91 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~~' |
---|
92 | |
---|
93 | ! open elliptic solver statistics file |
---|
94 | clname = 'solver.stat' |
---|
95 | CALL ctlopn( numsol, clname, 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL', & |
---|
96 | 1, numout, lwp, 1 ) |
---|
97 | |
---|
98 | |
---|
99 | ! 0. Define the solver parameters |
---|
100 | ! ---------------------------- |
---|
101 | ! Namelist namsol : elliptic solver / islands / free surface |
---|
102 | REWIND( numnam ) |
---|
103 | READ ( numnam, namsol ) |
---|
104 | |
---|
105 | #if defined key_feti |
---|
106 | ! FETI algorithm, we force nsolv at 3 |
---|
107 | nsolv = 3 |
---|
108 | #endif |
---|
109 | |
---|
110 | |
---|
111 | ! 0. Parameter control and print |
---|
112 | ! --------------------------- |
---|
113 | |
---|
114 | ! Control print |
---|
115 | IF(lwp) WRITE(numout,*) ' Namelist namsol : set solver parameters' |
---|
116 | |
---|
117 | IF(lwp) THEN |
---|
118 | WRITE(numout,*) ' type of elliptic solver nsolv = ', nsolv |
---|
119 | WRITE(numout,*) ' absolute/relative (0/1) precision nsol_arp = ', nsol_arp |
---|
120 | WRITE(numout,*) ' minimum iterations for solver nmin = ', nmin |
---|
121 | WRITE(numout,*) ' maximum iterations for solver nmax = ', nmax |
---|
122 | WRITE(numout,*) ' frequency for test nmod = ', nmod |
---|
123 | WRITE(numout,*) ' absolute precision of solver eps = ', eps |
---|
124 | WRITE(numout,*) ' absolute precision for SOR solver resmax = ', resmax |
---|
125 | WRITE(numout,*) ' optimal coefficient of sor sor = ', sor |
---|
126 | IF(lk_isl) WRITE(numout,*) ' absolute precision stream fct epsisl = ', epsisl |
---|
127 | IF(lk_isl) WRITE(numout,*) ' maximum pcg iterations island nmisl = ', nmisl |
---|
128 | WRITE(numout,*) ' free surface parameter rnu = ', rnu |
---|
129 | WRITE(numout,*) |
---|
130 | ENDIF |
---|
131 | |
---|
132 | IF( lk_dynspg_flt ) THEN |
---|
133 | IF(lwp) WRITE(numout,*) |
---|
134 | IF(lwp) WRITE(numout,*) ' free surface formulation' |
---|
135 | IF( lk_isl ) CALL ctl_stop( ' key_islands inconsistent with key_dynspg_flt' ) |
---|
136 | |
---|
137 | ELSEIF( lk_dynspg_rl ) THEN |
---|
138 | IF(lwp) WRITE(numout,*) |
---|
139 | IF(lwp) WRITE(numout,*) ' Rigid lid formulation' |
---|
140 | ELSE |
---|
141 | CALL ctl_stop( ' Choose only one surface pressure gradient calculation: filtered or rigid-lid', & |
---|
142 | & ' Should not call this routine if dynspg_exp or dynspg_ts has been chosen' ) |
---|
143 | ENDIF |
---|
144 | IF( lk_dynspg_flt .AND. lk_dynspg_rl ) & |
---|
145 | CALL ctl_stop( ' Chose between free surface or rigid-lid, not both' ) |
---|
146 | |
---|
147 | SELECT CASE ( nsolv ) |
---|
148 | |
---|
149 | CASE ( 1 ) ! preconditioned conjugate gradient solver |
---|
150 | IF(lwp) WRITE(numout,*) ' a preconditioned conjugate gradient solver is used' |
---|
151 | IF( jpr2di /= 0 .AND. jpr2dj /= 0 ) & |
---|
152 | CALL ctl_stop( ' jpr2di and jpr2dj should be equal to zero' ) |
---|
153 | |
---|
154 | CASE ( 2 ) ! successive-over-relaxation solver |
---|
155 | IF(lwp) WRITE(numout,*) ' a successive-over-relaxation solver is used' |
---|
156 | IF( jpr2di /= 0 .AND. jpr2dj /= 0 ) & |
---|
157 | CALL ctl_stop( ' jpr2di and jpr2dj should be equal to zero' ) |
---|
158 | |
---|
159 | CASE ( 3 ) ! FETI solver |
---|
160 | IF(lwp) WRITE(numout,*) ' the FETI solver is used' |
---|
161 | IF( jpr2di /= 0 .AND. jpr2dj /= 0 ) & |
---|
162 | CALL ctl_stop( ' jpr2di and jpr2dj should be equal to zero' ) |
---|
163 | |
---|
164 | IF( .NOT.lk_mpp ) THEN |
---|
165 | CALL ctl_stop( ' The FETI algorithm is used only with the key_mpp_... option' ) |
---|
166 | ELSE |
---|
167 | IF( jpnij == 1 ) THEN |
---|
168 | CALL ctl_stop( ' The FETI algorithm needs more than one processor' ) |
---|
169 | ENDIF |
---|
170 | ENDIF |
---|
171 | |
---|
172 | CASE ( 4 ) ! successive-over-relaxation solver with extra outer halo |
---|
173 | IF(lwp) WRITE(numout,*) ' a successive-over-relaxation solver with extra outer halo is used' |
---|
174 | IF(lwp) WRITE(numout,*) ' with jpr2di =', jpr2di, ' and jpr2dj =', jpr2dj |
---|
175 | IF( .NOT. lk_mpp .AND. jpr2di /= 0 .AND. jpr2dj /= 0 ) THEN |
---|
176 | CALL ctl_stop( ' jpr2di and jpr2dj are not equal to zero', & |
---|
177 | & ' In this case this algorithm should be used only with the key_mpp_... option' ) |
---|
178 | ELSE |
---|
179 | IF( ( ( jperio == 1 .OR. jperio == 4 .OR. jperio == 6 ) .OR. ( jpni /= 1 ) ) & |
---|
180 | & .AND. ( jpr2di /= jpr2dj ) ) CALL ctl_stop( ' jpr2di should be equal to jpr2dj' ) |
---|
181 | ENDIF |
---|
182 | |
---|
183 | CASE DEFAULT |
---|
184 | WRITE(ctmp1,*) ' bad flag value for nsolv = ', nsolv |
---|
185 | CALL ctl_stop( ctmp1 ) |
---|
186 | |
---|
187 | END SELECT |
---|
188 | |
---|
189 | ! Grid-point at which the solver is applied |
---|
190 | ! ----------------------------------------- |
---|
191 | |
---|
192 | IF( lk_dynspg_rl ) THEN ! rigid-lid |
---|
193 | IF( lk_mpp ) THEN |
---|
194 | c_solver_pt = 'G' ! G= F with special staff ??? which one? |
---|
195 | ELSE |
---|
196 | c_solver_pt = 'F' |
---|
197 | ENDIF |
---|
198 | ELSE ! free surface T-point |
---|
199 | IF( lk_mpp ) THEN |
---|
200 | c_solver_pt = 'S' ! S=T with special staff ??? which one? |
---|
201 | ELSE |
---|
202 | c_solver_pt = 'T' |
---|
203 | ENDIF |
---|
204 | ENDIF |
---|
205 | |
---|
206 | |
---|
207 | ! Construction of the elliptic system matrix |
---|
208 | ! ------------------------------------------ |
---|
209 | |
---|
210 | CALL sol_mat( kt ) |
---|
211 | |
---|
212 | |
---|
213 | IF( lk_isl ) THEN |
---|
214 | |
---|
215 | ! Islands in the domain |
---|
216 | ! --------------------- |
---|
217 | |
---|
218 | IF ( jpisl == 0 ) THEN |
---|
219 | WRITE(ctmp1,*) ' bad islands parameter jpisl =', jpisl |
---|
220 | CALL ctl_stop( ctmp1 ) |
---|
221 | ENDIF |
---|
222 | |
---|
223 | ! open Island streamfunction statistic file |
---|
224 | CALL ctlopn( numisp, 'islands.stat', 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL', & |
---|
225 | & 1 , numout , lwp , 1 ) |
---|
226 | |
---|
227 | CALL isl_dom ! Island identification |
---|
228 | |
---|
229 | CALL isl_bsf ! Island barotropic stream function |
---|
230 | |
---|
231 | CALL isl_mat ! Comput and invert the island matrix |
---|
232 | |
---|
233 | ! mbathy set to the number of w-level (minimum value 2) |
---|
234 | DO jj = 1, jpj |
---|
235 | DO ji = 1, jpi |
---|
236 | mbathy(ji,jj) = MAX( 1, mbathy(ji,jj) ) + 1 |
---|
237 | END DO |
---|
238 | END DO |
---|
239 | |
---|
240 | ENDIF |
---|
241 | |
---|
242 | END SUBROUTINE solver_init |
---|
243 | |
---|
244 | !!====================================================================== |
---|
245 | END MODULE solver |
---|