source: CONFIG/publications/ICOLMDZORINCA_CO2_Transport_GMD_2023/INCA/build/ppsrc/INCA_SRC/sulfint.f90 @ 6610

Last change on this file since 6610 was 6610, checked in by acosce, 10 months ago

INCA used for ICOLMDZORINCA_CO2_Transport_GMD_2023

File size: 10.4 KB
Line 
1
2
3
4
5
6
7
8
9
10
11
12!$Id: sulfint.F90 104 2008-12-23 10:28:51Z acosce $
13!! =========================================================================
14!! INCA - INteraction with Chemistry and Aerosols
15!!
16!! Copyright Laboratoire des Sciences du Climat et de l'Environnement (LSCE)
17!!           Unite mixte CEA-CNRS-UVSQ
18!!
19!! Contributors to this INCA subroutine:
20!!
21!! Didier Hauglustaine, LSCE, hauglustaine@cea.fr
22!!
23!! Anne Cozic, LSCE, anne.cozic@cea.fr
24!! Yann Meurdesoif, LSCE, yann.meurdesoif@cea.fr
25!!
26!! This software is a computer program whose purpose is to simulate the
27!! atmospheric gas phase and aerosol composition. The model is designed to be
28!! used within a transport model or a general circulation model. This version
29!! of INCA was designed to be coupled to the LMDz GCM. LMDz-INCA accounts
30!! for emissions, transport (resolved and sub-grid scale), photochemical
31!! transformations, and scavenging (dry deposition and washout) of chemical
32!! species and aerosols interactively in the GCM. Several versions of the INCA
33!! model are currently used depending on the envisaged applications with the
34!! chemistry-climate model.
35!!
36!! This software is governed by the CeCILL  license under French law and
37!! abiding by the rules of distribution of free software.  You can  use,
38!! modify and/ or redistribute the software under the terms of the CeCILL
39!! license as circulated by CEA, CNRS and INRIA at the following URL
40!! "http://www.cecill.info".
41!!
42!! As a counterpart to the access to the source code and  rights to copy,
43!! modify and redistribute granted by the license, users are provided only
44!! with a limited warranty  and the software's author,  the holder of the
45!! economic rights,  and the successive licensors  have only  limited
46!! liability.
47!!
48!! In this respect, the user's attention is drawn to the risks associated
49!! with loading,  using,  modifying and/or developing or reproducing the
50!! software by the user in light of its specific status of free software,
51!! that may mean  that it is complicated to manipulate,  and  that  also
52!! therefore means  that it is reserved for developers  and  experienced
53!! professionals having in-depth computer knowledge. Users are therefore
54!! encouraged to load and test the software's suitability as regards their
55!! requirements in conditions enabling the security of their systems and/or
56!! data to be ensured and,  more generally, to use and operate it in the
57!! same conditions as regards security.
58!!
59!! The fact that you are presently reading this means that you have had
60!! knowledge of the CeCILL license and that you accept its terms.
61!! =========================================================================
62
63SUBROUTINE SULF_INTI (filename)
64  !----------------------------------------------------------------------
65  !       ... SULFATE climatologies initialize
66  ! Didier Hauglustaine, IPSL, 2000.
67  !----------------------------------------------------------------------
68
69  USE AC_SULF 
70  USE INCA_DIM
71  USE MOD_INCA_PARA
72  IMPLICIT NONE
73
74  INCLUDE 'netcdf.inc'
75
76  !----------------------------------------------------------------------
77  !       ... Dummy args
78  !----------------------------------------------------------------------
79  CHARACTER(len=*), INTENT(in) :: filename
80
81  !----------------------------------------------------------------------
82  !       ... Local variables
83  !----------------------------------------------------------------------
84  INTEGER :: iret
85  INTEGER :: varid
86  INTEGER :: error
87
88  ! ... Open the file
89!$OMP MASTER
90  IF (is_mpi_root) THEN
91      iret = nf_open(filename, 0, ncid)
92      CALL check_err(iret, 'SULF_INTI', 'problem to open file')
93
94      ! ... Get lengths
95      iret = nf_inq_dimid(ncid, 'time', varid)
96      CALL check_err(iret, 'SULF_INTI', 'problem to check dimid time')
97      iret = nf_inq_dimlen(ncid, varid, ntimes)
98      CALL check_err(iret, 'SULF_INTI', 'problem to check dimlen time')
99      iret = nf_inq_dimid(ncid, 'lev', varid)
100      CALL check_err(iret, 'SULF_INTI', 'problem to check dimid lev')
101      iret = nf_inq_dimlen(ncid, varid, nlevs)
102      CALL check_err(iret, 'SULF_INTI', 'problem to check dimlen lev')
103      iret = nf_inq_dimid(ncid, 'vector', varid)
104      CALL check_err(iret, 'SULF_INTI', 'problem to check dimid vector')
105      iret = nf_inq_dimlen(ncid, varid, klons)
106      CALL check_err(iret, 'SULF_INTI', 'problem to check dimlen vector')
107  ENDIF
108!$OMP END MASTER
109
110  CALL bcast(ntimes)
111  CALL bcast(nlevs)
112  CALL bcast(klons)
113
114  IF (klons/=nbp_glo) THEN
115     write(lunout,*) 'sulf_inti : [klons] [nbp_glo]', klons, nbp_glo
116      call print_err(3, 'SULF_INTI','there a problem of vector size', 'check klons and nbp_glo', '')
117  ELSE
118      klons=PLON
119  ENDIF
120
121  ! ... Allocate variables
122  ALLOCATE(times(ntimes), STAT=error)
123  IF (error /= 0) call print_err(3,  'SULF_INTI', 'Space requested not possible for variable times ', '', '')
124  ALLOCATE(levs(nlevs), STAT=error)
125  IF (error /= 0) CALL print_err(3,  'SULF_INTI', 'Space requested not possible for variable levs', '', '')
126  ALLOCATE(so4bd(klons,nlevs,2), STAT=error)
127  IF (error /= 0) CALL print_err(3, 'SULF_INTI', 'Space requested not possible for variable so4bd ', '', '')
128
129  ! ... Get the coordinates
130!$OMP MASTER
131  IF (is_mpi_root) THEN
132      iret = nf_inq_varid(ncid, 'time', varid)
133      CALL check_err(iret, 'SULF_INTI', 'problem to check varid time')
134      iret = nf_get_var_double(ncid, varid, times)
135      CALL check_err(iret, 'SULF_INTI', 'problem to read time')
136
137
138      iret = nf_inq_varid(ncid, 'lev', varid)
139      CALL check_err(iret, 'SULF_INTI', 'problem to check varid lev')
140      iret = nf_get_var_double(ncid, varid, levs)
141      CALL check_err(iret, 'SULF_INTI', 'problem to read lev')
142
143
144
145      ! ... Variable ids
146      iret = nf_inq_varid(ncid, 'so4', so4_id)
147      CALL check_err(iret, 'SULF_INTI', 'problem to check varid so4')
148  ENDIF
149!$OMP END MASTER
150  CALL bcast(times)
151  CALL bcast(so4_id)
152  call bcast(levs)
153
154END SUBROUTINE SULF_INTI
155
156SUBROUTINE SULF_READ(calday)
157  !----------------------------------------------------------------------
158  !       ... Read SO4 distributions
159  ! Didier Hauglustaine, IPSL, 2000.
160  !----------------------------------------------------------------------
161
162  USE AC_SULF
163  USE INCA_DIM
164  USE MOD_INCA_PARA
165  USE PRINT_INCA
166  IMPLICIT NONE
167
168  INCLUDE 'netcdf.inc'
169
170  !----------------------------------------------------------------------
171  !       ... Dummy args
172  !----------------------------------------------------------------------
173  REAL, INTENT(IN) :: calday
174
175  !----------------------------------------------------------------------
176  !       ... Local variables
177  !----------------------------------------------------------------------
178  INTEGER :: iret
179  INTEGER :: varid
180  INTEGER :: oldlotime, oldhitime
181  INTEGER, DIMENSION(3) :: start3, count3
182  INTEGER, DIMENSION(2) :: start2, count2
183  REAL :: so4bd_glo(nbp_glo,nlevs,2)
184
185  ! ... Check to see if model time still bounded by data set
186  oldlotime = lotime
187  oldhitime = hitime
188  CALL FINDPLB (times, ntimes, calday, lotime)
189  IF ( cyclical ) THEN
190      hitime = MOD(lotime,ntimes) + 1
191  ELSE
192      hitime = lotime + 1
193  END IF
194
195  ! ... Read new data if needed
196  IF ( hitime /= oldhitime ) THEN
197!$OMP MASTER
198      IF (is_mpi_root) THEN
199          loind = hiind
200          hiind = MOD( loind, 2) + 1
201
202          count2(1) = nbp_glo   
203          count2(2) = 1 
204          start2(1) = 1
205          count3(1) = nbp_glo   
206          count3(2) = nlevs 
207          count3(3) = 1 
208          start3(1) = 1
209          start3(2) = 1
210
211          WRITE(lunout, *) 'SULF_READ : read new data for time ', times(hitime)
212          start3(3) = hitime
213          start2(2) = hitime
214          iret = nf_get_vara_double(ncid, so4_id, start3, count3, so4bd_glo(1,1,hiind))
215          CALL check_err(iret, 'SULF_READ', 'problem to read so4bd_glo hiind')
216
217      ENDIF
218!$OMP END MASTER
219      CALL scatter(so4bd_glo(:,:,hiind),so4bd(:,:,hiind))
220
221      IF ( lotime /= oldlotime ) THEN
222!$OMP MASTER
223          IF (is_mpi_root) THEN
224              WRITE(lunout, *) 'SULF_READ : read new data for time ', times(lotime)
225              start3(3) = lotime
226              start2(2) = lotime
227              iret = nf_get_vara_double(ncid, so4_id, start3, count3, so4bd_glo(1,1,loind))
228              CALL check_err(iret, 'SULF_READ', 'problem to read so4bd_glo loind')
229
230          ENDIF
231!$OMP END MASTER
232          CALL scatter(so4bd_glo(:,:,loind),so4bd(:,:,loind))
233
234      END IF
235
236  END IF
237
238END SUBROUTINE SULF_READ
239
240SUBROUTINE SULF_INTERP(calday, pmid, zmid)
241  !----------------------------------------------------------------------
242  !       ... Interpolate SO4 on model time
243  ! Didier Hauglustaine, IPSL, 2000.
244  !----------------------------------------------------------------------
245
246  USE AC_SULF
247  USE INCA_DIM
248  IMPLICIT NONE
249
250  !----------------------------------------------------------------------
251  !       ... Dummy args
252  !----------------------------------------------------------------------
253  REAL, INTENT(IN) :: calday
254  REAL, INTENT(IN) :: pmid(PLON,PLEV)
255  REAL, INTENT(IN) :: zmid(PLON,PLEV)
256
257  !----------------------------------------------------------------------
258  !       ... Local variables
259  !----------------------------------------------------------------------
260  REAL, PARAMETER :: dzs = 1. !specific to  Considine's input file
261  INTEGER :: k, i, m, j, l
262  INTEGER :: lonlev 
263  REAL :: dt, dt1, tint
264  REAL    :: wght
265  REAL, DIMENSION(PLON,PLEV) :: so4wrk
266
267  lonlev = PLON * PLEV
268
269  !----------------------------------------------------------------------
270  !     Note: 365 day year inconsistent with LMDz (360 days) !!!
271  !----------------------------------------------------------------------
272  ! ... interpolate linearly on current model time 
273
274  IF ( times(hitime) < times(lotime) ) THEN
275      dt = 365. - times(lotime) + times(hitime) 
276      IF (calday <= times(hitime) ) THEN
277          dt1 = 365. - times(lotime) + calday
278      ELSE
279          dt1 = calday - times(lotime)
280      END IF
281  ELSE
282      dt  = times(hitime) - times(lotime)
283      dt1 = calday - times(lotime)
284  END IF
285
286  tint = dt1/dt
287
288  ! ... Interpolate to local time
289
290  CALL linintp (& 
291     lonlev,           &
292     0.,               &
293     1.,               &
294     tint,             &
295     so4bd(1,1,loind), &
296     so4bd(1,1,hiind), &
297     so4wrk)             
298
299!--------------------------------------------------------
300!       ... Change units (from kg/m3 to g/cm3)
301!--------------------------------------------------------
302
303  DO i = 1, PLON
304    DO l = 1, PLEV
305     
306      so4(i,l) = so4wrk(i,l) *1.e-3
307     
308    END DO
309  END DO
310
311
312END SUBROUTINE SULF_INTERP
313
Note: See TracBrowser for help on using the repository browser.