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

Last change on this file since 1528 was 1528, checked in by rblod, 15 years ago

Suppress rigid-lid option, see ticket #486

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 8.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 obc_oce         ! Lateral open boundary condition
17   USE in_out_manager  ! I/O manager
18   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
19   USE lib_mpp
20   USE dynspg_oce      ! choice/control of key cpp for surface pressure gradient
21
22   IMPLICIT NONE
23
24   !!----------------------------------------------------------------------
25   !!   OPA 9.0 , LOCEAN-IPSL (2005)
26   !! $Id$
27   !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt
28   !!----------------------------------------------------------------------
29
30CONTAINS
31
32   SUBROUTINE solver_init( kt )
33      !!----------------------------------------------------------------------
34      !!                  ***  ROUTINE solver_init  ***
35      !!                   
36      !! ** Purpose :   Initialization for the solver of the elliptic equation:
37      !!       * lk_dynspg_flt = T : transport divergence system.
38      !!     
39      !! ** Method :
40      !!       - Compute the local depth of the water column at u- and v-point
41      !!      The local depth of the water column is computed by summing
42      !!      the vertical scale factors. For its inverse, the thickness of
43      !!      the first model level is imposed as lower bound. The inverse of
44      !!      this depth is THEN taken and masked, so that the inverse of the
45      !!      local depth is zero when the local depth is zero.
46      !!
47      !! ** Action : - hur, hvr : masked inverse of the local depth at
48      !!                                u- and v-point.
49      !!             - hu, hv   : masked local depth at u- and v- points
50      !!             - c_solver_pt : nature of the gridpoint at which the
51      !!                                solver is applied
52      !! References :
53      !!      Jensen, 1986: adv. phys. oceanogr. num. mod.,ed. o brien,87-110.
54      !!      Madec & Marti, 1990: internal rep. LODYC, 90/03., 29pp.
55      !!
56      !! History :
57      !!        !  90-10  (G. Madec)  Original code           
58      !!        !  93-02  (O. Marti)                         
59      !!        !  97-02  (G. Madec)  local depth inverse computation
60      !!        !  98-10  (G. Roullet, G. Madec)  free surface
61      !!   9.0  !  03-07  (G. Madec)  free form, F90
62      !!    "   !  05-11  (V. Garnier) Surface pressure gradient organization
63      !!----------------------------------------------------------------------
64      !! * Arguments
65      INTEGER, INTENT(in) :: kt
66
67      !! * Local declarations
68      CHARACTER(len=80) :: clname
69
70      NAMELIST/namsol/ nsolv, nsol_arp, nmin, nmax, nmod, eps, resmax, sor, rnu
71      !!----------------------------------------------------------------------
72
73      IF(lwp) WRITE(numout,*)
74      IF(lwp) WRITE(numout,*) 'solver_init : solver to compute the surface pressure gradient'
75      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
76
77      ! open elliptic solver statistics file
78      clname = 'solver.stat'
79      CALL ctlopn( numsol, clname, 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL',   &
80                   1, numout, lwp, 1 )
81
82
83      ! 0. Define the solver parameters
84      !    ----------------------------
85      ! Namelist namsol : elliptic solver / free surface
86      REWIND( numnam )
87      READ  ( numnam, namsol )
88
89#if defined key_feti
90      ! FETI algorithm, we force nsolv at 3
91      nsolv = 3
92#endif
93
94      ! 0. Parameter control and print
95      !    ---------------------------
96
97      ! Control print
98      IF(lwp) WRITE(numout,*) '          Namelist namsol : set solver parameters'
99
100      IF(lwp) THEN
101         WRITE(numout,*) '             type of elliptic solver            nsolv    = ', nsolv
102         WRITE(numout,*) '             absolute/relative (0/1) precision  nsol_arp = ', nsol_arp
103         WRITE(numout,*) '             minimum iterations for solver      nmin     = ', nmin
104         WRITE(numout,*) '             maximum iterations for solver      nmax     = ', nmax
105         WRITE(numout,*) '             frequency for test                 nmod     = ', nmod
106         WRITE(numout,*) '             absolute precision of solver       eps      = ', eps
107         WRITE(numout,*) '             absolute precision for SOR solver  resmax   = ', resmax
108         WRITE(numout,*) '             optimal coefficient of sor         sor      = ', sor
109         WRITE(numout,*) '             free surface parameter         rnu    = ', rnu
110         WRITE(numout,*)
111      ENDIF
112
113      IF( lk_dynspg_flt ) THEN
114         IF(lwp) WRITE(numout,*)
115         IF(lwp) WRITE(numout,*) '          free surface formulation'
116      ELSE
117         CALL ctl_stop( '          Choose only one surface pressure gradient calculation: filtered ',   &
118              &         '          Should not call this routine if dynspg_exp or dynspg_ts has been chosen' )
119      ENDIF
120
121      SELECT CASE ( nsolv )
122
123      CASE ( 1 )                ! preconditioned conjugate gradient solver
124         IF(lwp) WRITE(numout,*) '          a preconditioned conjugate gradient solver is used'
125         IF( jpr2di /= 0 .AND. jpr2dj /= 0 ) &
126            CALL ctl_stop( ' jpr2di and jpr2dj should be equal to zero' )
127
128      CASE ( 2 )                ! successive-over-relaxation solver
129         IF(lwp) WRITE(numout,*) '          a successive-over-relaxation solver with extra outer halo is used'
130         IF(lwp) WRITE(numout,*) '          with jpr2di =', jpr2di, ' and  jpr2dj =', jpr2dj
131         IF( .NOT. lk_mpp .AND. jpr2di /= 0 .AND. jpr2dj /= 0 ) THEN
132             CALL ctl_stop( ' jpr2di and jpr2dj are not equal to zero',   &
133             &              ' In this case this algorithm should be used only with the key_mpp_... option' )
134         ELSE
135            IF( ( ( jperio == 1 .OR. jperio == 4 .OR. jperio == 6 ) .OR. ( jpni /= 1 ) ) &
136              &  .AND. ( jpr2di /= jpr2dj ) ) CALL ctl_stop( '          jpr2di should be equal to jpr2dj' )
137         ENDIF
138
139      CASE ( 3 )                ! FETI solver
140         IF(lwp) WRITE(numout,*) '          the FETI solver is used'
141         IF( jpr2di /= 0 .AND. jpr2dj /= 0 ) &
142              CALL ctl_stop( ' jpr2di and jpr2dj should be equal to zero' )
143           
144         IF( .NOT.lk_mpp ) THEN
145            CALL ctl_stop( ' The FETI algorithm is used only with the key_mpp_... option' )
146         ELSE
147            IF( jpnij == 1 ) THEN
148               CALL ctl_stop( ' The FETI algorithm needs more than one processor' )
149            ENDIF
150         ENDIF
151         
152      CASE DEFAULT
153         WRITE(ctmp1,*) '          bad flag value for nsolv = ', nsolv
154         CALL ctl_stop( ctmp1 )
155         
156      END SELECT
157
158      IF( nbit_cmp == 1 ) THEN
159         IF( nsolv /= 2 ) THEN
160            CALL ctl_stop( ' Reproductibility tests (nbit_cmp=1) require the SOR solver: nsolv = 2' )
161         ELSE IF( MAX( jpr2di, jpr2dj ) > 0 ) THEN
162            CALL ctl_stop( ' Reproductibility tests (nbit_cmp=1) require jpr2di = jpr2dj = 0' )
163         END IF
164      END IF
165
166      ! Grid-point at which the solver is applied
167      ! -----------------------------------------
168
169      IF( lk_mpp ) THEN
170         c_solver_pt = 'S'   ! S=T with special staff ??? which one?
171      ELSE
172         c_solver_pt = 'T'
173      ENDIF
174
175      ! Construction of the elliptic system matrix
176      ! ------------------------------------------
177
178      CALL sol_mat( kt )
179
180   END SUBROUTINE solver_init
181
182   !!======================================================================
183END MODULE solver
Note: See TracBrowser for help on using the repository browser.