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 @ 16

Last change on this file since 16 was 16, checked in by opalod, 20 years ago

CT : UPDATE001 : First major NEMO update

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 9.0 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_rl       
22   USE dynspg_fsc     
23
24   IMPLICIT NONE
25
26   !!----------------------------------------------------------------------
27   !!   OPA 9.0 , LODYC-IPSL  (2003)
28   !!----------------------------------------------------------------------
29
30CONTAINS
31
32   SUBROUTINE solver_init
33      !!----------------------------------------------------------------------
34      !!                  ***  ROUTINE solver_init  ***
35      !!                   
36      !! ** Purpose :   Initialization for the solver of the elliptic equation:
37      !!       * default option: barotropic stream function system
38      !!         and islands initialization (if l_isl=T)
39      !!       * lk_dynspg_fsc = T : transport divergence system. No specific
40      !!         treatment of islands.
41      !!     
42      !! ** Method :
43      !!       - Compute the local depth of the water column at u- and v-point
44      !!      (lk_dynspg_fsc = T) or its inverse (lk_dynspg_rl = T).
45      !!      The local depth of the water column is computed by summing
46      !!      the vertical scale factors. For its inverse, the thickness of
47      !!      the first model level is imposed as lower bound. The inverse of
48      !!      this depth is THEN taken and masked, so that the inverse of the
49      !!      local depth is zero when the local depth is zero.
50      !!       - Construct the matrix of the elliptic system by a call to
51      !!      solmat.F routine.
52      !!       - island (if l_isl=T)
53      !!            isl_dom: find islands from the bathymetry file
54      !!            isl_bsf: compute the island barotropic stream function
55      !!            isl_mat: compute the inverse island matrix
56      !!            set mbathy to the number of non-zero w-levels of a water
57      !!            column (the minimum value of mbathy is 2):
58      !!                  mbathy = min( mbathy, 1 ) + 1
59      !!
60      !! ** Action : - hur, hvr : masked inverse of the local depth at
61      !!                                u- and v-point. (lk_dynspg_rl = T)
62      !!             - hu, hv   : masked local depth at u- and v- points
63      !!                                (lk_dynspg_fsc = T)
64      !!             - c_solver_pt : nature of the gridpoint at which the
65      !!                                solver is applied
66      !! References :
67      !!      Jensen, 1986: adv. phys. oceanogr. num. mod.,ed. o brien,87-110.
68      !!      Madec & Marti, 1990: internal rep. LODYC, 90/03., 29pp.
69      !!
70      !! History :
71      !!        !  90-10  (G. Madec)  Original code           
72      !!        !  93-02  (O. Marti)                         
73      !!        !  97-02  (G. Madec)  local depth inverse computation
74      !!        !  98-10  (G. Roullet, G. Madec)  free surface
75      !!   9.0  !  03-07  (G. Madec)  free form, F90
76      !!----------------------------------------------------------------------
77      !! * Local declarations
78      INTEGER :: ji, jj   ! dummy loop indices
79
80      NAMELIST/namsol/ nsolv, nmax, eps, sor, epsisl, nmisl, rnu
81      !!----------------------------------------------------------------------
82
83      IF(lwp) WRITE(numout,*)
84      IF(lwp) WRITE(numout,*) 'ini_sol : solver to compute the surface pressure gradient'
85      IF(lwp) WRITE(numout,*) '~~~~~~~'
86
87      ! open elliptic solver statistics file
88      CALL ctlopn( numsol, 'solver.stat', 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL',   &
89                   1, numout, lwp, 1 )
90
91
92      ! 0. Define the solver parameters
93      !    ----------------------------
94      ! Namelist namsol : elliptic solver / islands / free surface
95      REWIND( numnam )
96      READ  ( numnam, namsol )
97
98#if defined key_feti
99      ! FETI algorithm, we force nsolv at 3
100      nsolv = 3
101#endif
102
103
104      ! 0. Parameter control and print
105      !    ---------------------------
106
107      ! Control print
108      IF(lwp) WRITE(numout,*) '          Namelist namsol : set solver parameters'
109
110      IF(lwp) THEN
111         WRITE(numout,*)
112         WRITE(numout,*) '             type of elliptic solver        nsolv  = ', nsolv
113         WRITE(numout,*) '             maximum iterations for solver  nmax   = ', nmax
114         WRITE(numout,*) '             absolute precision of solver   eps    = ', eps
115         WRITE(numout,*) '             optimal coefficient of sor     sor    = ', sor
116         IF(l_isl) WRITE(numout,*) '             absolute precision stream fct  epsisl = ', epsisl
117         IF(l_isl) WRITE(numout,*) '             maximum pcg iterations island  nmisl  = ', nmisl
118         WRITE(numout,*) '             free surface parameter         rnu    = ', rnu
119         WRITE(numout,*)
120      ENDIF
121
122      IF( lk_dynspg_fsc ) THEN
123         IF(lwp) WRITE(numout,*)
124         IF(lwp) WRITE(numout,*) '          *** free surface formulation'
125         IF( l_isl ) THEN
126            IF(lwp) WRITE(numout,cform_err)
127            IF(lwp) WRITE(numout,*) ' key_islands inconsistent with key_dynspg_fsc'
128            nstop = nstop + 1
129         ENDIF
130      ELSEIF( lk_dynspg_rl ) THEN
131         IF(lwp) WRITE(numout,*)
132         IF(lwp) WRITE(numout,*) '          *** Rigid lid formulation'
133      ELSE
134         IF(lwp) WRITE(numout,cform_err)
135         IF(lwp) WRITE(numout,*) '          Chose at least one surface pressure gradient calculation: free surface or rigid-lid'
136         nstop = nstop + 1
137      ENDIF
138      IF( lk_dynspg_fsc .AND. lk_dynspg_rl ) THEN
139         IF(lwp) WRITE(numout,cform_err)
140         IF(lwp) WRITE(numout,*) '          Chose between free surface or rigid-lid, not both'
141         nstop = nstop + 1
142      ENDIF
143
144      SELECT CASE ( nsolv )
145
146      CASE ( 1 )                ! preconditioned conjugate gradient solver
147         IF(lwp) WRITE(numout,*) '          use a preconditioned conjugate gradient solver'
148
149      CASE ( 2 )                ! successive-over-relaxation solver
150         IF(lwp) WRITE(numout,*) '          use a successive-over-relaxation solver'
151
152      CASE ( 3 )                ! FETI solver
153         IF(lwp) WRITE(numout,*) '          use the FETI solver'
154         IF( .NOT.lk_mpp ) THEN
155            IF(lwp) WRITE(numout,*) ' The FETI algorithm is used only with the key_mpp_... option'
156            nstop = nstop + 1
157         ELSE
158            IF( jpnij == 1 ) THEN
159               IF(lwp) WRITE(numout,*) ' The FETI algorithm needs more than one processor'
160               nstop = nstop + 1
161            ENDIF
162         ENDIF
163         
164      CASE DEFAULT
165         IF(lwp) WRITE(numout,cform_err)
166         IF(lwp) WRITE(numout,*) '          bad flag value for nsolv = ', nsolv
167         nstop = nstop + 1
168         
169      END SELECT
170
171      ! Grid-point at which the solver is applied
172      ! -----------------------------------------
173
174      IF( lk_dynspg_rl ) THEN       ! rigid-lid
175         IF( lk_mpp ) THEN
176            c_solver_pt = 'G'   ! G= F with special staff ??? which one?
177         ELSE
178            c_solver_pt = 'F'
179         ENDIF
180      ELSE                          ! free surface T-point
181         IF( lk_mpp ) THEN
182            c_solver_pt = 'S'   ! S=T with special staff ??? which one?
183         ELSE
184            c_solver_pt = 'T'
185         ENDIF
186      ENDIF
187
188
189      ! Construction of the elliptic system matrix
190      ! ------------------------------------------
191
192      CALL sol_mat
193
194
195      IF( l_isl ) THEN
196     
197         ! Islands in the domain
198         ! ---------------------
199
200         IF ( jpisl == 0 ) THEN
201             IF(lwp)WRITE(numout,cform_err)
202             IF(lwp)WRITE(numout,*) ' bad islands parameter jpisl =', jpisl
203             nstop = nstop + 1
204         ENDIF
205
206         ! open Island streamfunction statistic file
207         CALL ctlopn( numisp, 'islands.stat', 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL',   &
208            &         1     , numout        , lwp      , 1                         )
209   
210         CALL isl_dom       ! Island identification
211
212         CALL isl_bsf       ! Island barotropic stream function
213
214         CALL isl_mat       ! Comput and invert the island matrix
215
216         ! mbathy set to the number of w-level (minimum value 2)
217         DO jj = 1, jpj
218            DO ji = 1, jpi
219               mbathy(ji,jj) = MAX( 1, mbathy(ji,jj) ) + 1
220            END DO
221         END DO
222
223      ENDIF
224
225   END SUBROUTINE solver_init
226
227   !!======================================================================
228END MODULE solver
Note: See TracBrowser for help on using the repository browser.