New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
solver.F90 in trunk/NEMO/OPA_SRC/SOL – NEMO

source: trunk/NEMO/OPA_SRC/SOL/solver.F90 @ 367

Last change on this file since 367 was 367, checked in by opalod, 18 years ago

nemo_v1_update_035 : CT : OBCs adapted to the new surface pressure gradient algorithms

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 11.5 KB
Line 
1MODULE 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
31CONTAINS
32
33   SUBROUTINE solver_init
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      !! * Local declarations
80      INTEGER :: ji, jj   ! dummy loop indices
81
82      NAMELIST/namsol/ nsolv, nsol_arp, nmin, nmax, nmod, eps, resmax, sor, epsisl, nmisl, rnu
83      !!----------------------------------------------------------------------
84
85      IF(lwp) WRITE(numout,*)
86      IF(lwp) WRITE(numout,*) 'solver_init : solver to compute the surface pressure gradient'
87      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
88
89      ! open elliptic solver statistics file
90      CALL ctlopn( numsol, 'solver.stat', 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL',   &
91                   1, numout, lwp, 1 )
92
93
94      ! 0. Define the solver parameters
95      !    ----------------------------
96      ! Namelist namsol : elliptic solver / islands / free surface
97      REWIND( numnam )
98      READ  ( numnam, namsol )
99
100#if defined key_feti
101      ! FETI algorithm, we force nsolv at 3
102      nsolv = 3
103#endif
104
105
106      ! 0. Parameter control and print
107      !    ---------------------------
108
109      ! Control print
110      IF(lwp) WRITE(numout,*) '          Namelist namsol : set solver parameters'
111
112      IF(lwp) THEN
113         WRITE(numout,*) '             type of elliptic solver            nsolv    = ', nsolv
114         WRITE(numout,*) '             absolute/relative (0/1) precision  nsol_arp = ', nsol_arp
115         WRITE(numout,*) '             minimum iterations for solver      nmin     = ', nmin
116         WRITE(numout,*) '             maximum iterations for solver      nmax     = ', nmax
117         WRITE(numout,*) '             frequency for test                 nmod     = ', nmod
118         WRITE(numout,*) '             absolute precision of solver       eps      = ', eps
119         WRITE(numout,*) '             absolute precision for SOR solver  resmax   = ', resmax
120         WRITE(numout,*) '             optimal coefficient of sor         sor      = ', sor
121         IF(lk_isl) WRITE(numout,*) '             absolute precision stream fct    epsisl   = ', epsisl
122         IF(lk_isl) WRITE(numout,*) '             maximum pcg iterations island    nmisl    = ', nmisl
123         WRITE(numout,*) '             free surface parameter         rnu    = ', rnu
124         WRITE(numout,*)
125      ENDIF
126
127      IF( lk_dynspg_flt ) THEN
128         IF(lwp) WRITE(numout,*)
129         IF(lwp) WRITE(numout,*) '          free surface formulation'
130         IF( lk_isl ) THEN
131            IF(lwp) WRITE(numout,cform_err)
132            IF(lwp) WRITE(numout,*) ' key_islands inconsistent with key_dynspg_flt'
133            nstop = nstop + 1
134         ENDIF
135      ELSEIF( lk_dynspg_rl ) THEN
136         IF(lwp) WRITE(numout,*)
137         IF(lwp) WRITE(numout,*) '          Rigid lid formulation'
138      ELSE
139         IF(lwp) WRITE(numout,cform_err)
140         IF(lwp) WRITE(numout,*) '          Choose only one surface pressure gradient calculation: filtered or rigid-lid'
141         IF(lwp) WRITE(numout,*) '          Should not call this routine if dynspg_exp or dynspg_ts has been chosen'
142         nstop = nstop + 1
143      ENDIF
144      IF( lk_dynspg_flt .AND. lk_dynspg_rl ) THEN
145         IF(lwp) WRITE(numout,cform_err)
146         IF(lwp) WRITE(numout,*) '          Chose between free surface or rigid-lid, not both'
147         nstop = nstop + 1
148      ENDIF
149
150      SELECT CASE ( nsolv )
151
152      CASE ( 1 )                ! preconditioned conjugate gradient solver
153         IF(lwp) WRITE(numout,*) '          a preconditioned conjugate gradient solver is used'
154         IF( jpr2di /= 0 .AND. jpr2dj /= 0 ) THEN
155            IF(lwp) WRITE(numout,cform_err)
156            IF(lwp) WRITE(numout,*) ' jpr2di and jpr2dj should be equal to zero'
157            nstop = nstop + 1
158         ENDIF
159
160      CASE ( 2 )                ! successive-over-relaxation solver
161         IF(lwp) WRITE(numout,*) '          a successive-over-relaxation solver is used'
162         IF( jpr2di /= 0 .AND. jpr2dj /= 0 ) THEN
163            IF(lwp) WRITE(numout,cform_err)
164            IF(lwp) WRITE(numout,*) ' jpr2di and jpr2dj should be equal to zero'
165            nstop = nstop + 1
166         ENDIF
167
168      CASE ( 3 )                ! FETI solver
169         IF(lwp) WRITE(numout,*) '          the FETI solver is used'
170         IF( jpr2di /= 0 .AND. jpr2dj /= 0 ) THEN
171            IF(lwp) WRITE(numout,cform_err)
172            IF(lwp) WRITE(numout,*) ' jpr2di and jpr2dj should be equal to zero'
173            nstop = nstop + 1
174         ENDIF
175         IF( .NOT.lk_mpp ) THEN
176            IF(lwp) WRITE(numout,cform_err)
177            IF(lwp) WRITE(numout,*) ' The FETI algorithm is used only with the key_mpp_... option'
178            nstop = nstop + 1
179         ELSE
180            IF( jpnij == 1 ) THEN
181               IF(lwp) WRITE(numout,cform_err)
182               IF(lwp) WRITE(numout,*) ' The FETI algorithm needs more than one processor'
183               nstop = nstop + 1
184            ENDIF
185         ENDIF
186         
187      CASE ( 4 )                ! successive-over-relaxation solver with extra outer halo
188         IF(lwp) WRITE(numout,*) '          a successive-over-relaxation solver with extra outer halo is used'
189         IF(lwp) WRITE(numout,*) '          with jpr2di =', jpr2di, ' and  jpr2dj =', jpr2dj
190         IF( .NOT. lk_mpp .AND. jpr2di /= 0 .AND. jpr2dj /= 0 ) THEN
191            IF(lwp) WRITE(numout,cform_err)
192            IF(lwp) WRITE(numout,*) ' jpr2di and jpr2dj are not equal to zero'
193            IF(lwp) WRITE(numout,*) ' In this case this algorithm should be used only with the key_mpp_... option'
194            nstop = nstop + 1
195         ELSE
196            IF( ( ( jperio == 1 .OR. jperio == 4 .OR. jperio == 6 ) .OR. ( jpni /= 1 ) ) &
197              &  .AND. ( jpr2di /= jpr2dj ) ) THEN 
198               IF(lwp) WRITE(numout,cform_err)
199               IF(lwp) WRITE(numout,*) '          jpr2di should be equal to jpr2dj'
200               nstop = nstop + 1
201            ENDIF
202         ENDIF
203
204      CASE DEFAULT
205         IF(lwp) WRITE(numout,cform_err)
206         IF(lwp) WRITE(numout,*) '          bad flag value for nsolv = ', nsolv
207         nstop = nstop + 1
208         
209      END SELECT
210
211      ! Grid-point at which the solver is applied
212      ! -----------------------------------------
213
214      IF( lk_dynspg_rl ) THEN       ! rigid-lid
215         IF( lk_mpp ) THEN
216            c_solver_pt = 'G'   ! G= F with special staff ??? which one?
217         ELSE
218            c_solver_pt = 'F'
219         ENDIF
220      ELSE                          ! free surface T-point
221         IF( lk_mpp ) THEN
222            c_solver_pt = 'S'   ! S=T with special staff ??? which one?
223         ELSE
224            c_solver_pt = 'T'
225         ENDIF
226      ENDIF
227
228
229      ! Construction of the elliptic system matrix
230      ! ------------------------------------------
231
232      CALL sol_mat
233
234
235      IF( lk_isl ) THEN
236     
237         ! Islands in the domain
238         ! ---------------------
239
240         IF ( jpisl == 0 ) THEN
241             IF(lwp)WRITE(numout,cform_err)
242             IF(lwp)WRITE(numout,*) ' bad islands parameter jpisl =', jpisl
243             nstop = nstop + 1
244         ENDIF
245
246         ! open Island streamfunction statistic file
247         CALL ctlopn( numisp, 'islands.stat', 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL',   &
248            &         1     , numout        , lwp      , 1                         )
249   
250         CALL isl_dom       ! Island identification
251
252         CALL isl_bsf       ! Island barotropic stream function
253
254         CALL isl_mat       ! Comput and invert the island matrix
255
256         ! mbathy set to the number of w-level (minimum value 2)
257         DO jj = 1, jpj
258            DO ji = 1, jpi
259               mbathy(ji,jj) = MAX( 1, mbathy(ji,jj) ) + 1
260            END DO
261         END DO
262
263      ENDIF
264
265   END SUBROUTINE solver_init
266
267   !!======================================================================
268END MODULE solver
Note: See TracBrowser for help on using the repository browser.