LMDZ
readchlorophyll.F90
Go to the documentation of this file.
1 !
2 ! $Id$
3 !
4 
5 subroutine readchlorophyll(debut)
6 
7  use netcdf95, only: nf95_close, nf95_gw_var, nf95_inq_dimid, &
8  nf95_inq_varid, nf95_open
9  use netcdf, only: nf90_get_var, nf90_noerr, nf90_nowrite
10 
11  USE phys_cal_mod, ONLY : mth_cur
15  USE mod_phys_lmdz_para, ONLY: scatter
16  USE phys_state_var_mod, ONLY: chl_con
17 
18  implicit none
19 
20  include "YOMCST.h"
21 
22 ! Variable input
23  logical debut
24 
25 ! Variables locales
26  integer n_lat ! number of latitudes in the input data
27  integer n_lon ! number of longitudes in the input data
28  integer n_lev ! number of levels in the input data
29  integer n_month ! number of months in the input data
30  real, pointer:: latitude(:)
31  real, pointer:: longitude(:)
32  real, pointer:: time(:)
33  integer i, k
34  integer, save :: mth_pre
35 !$OMP THREADPRIVATE(mth_pre)
36 
37 ! Champs reconstitues
38  real, allocatable:: chlorocon(:, :, :)
39  real, allocatable:: chlorocon_mois(:, :)
40  real, allocatable:: chlorocon_mois_glo(:)
41 
42 ! For NetCDF:
43  integer ncid_in ! IDs for input files
44  integer varid, ncerr
45 
46 
47 !--------------------------------------------------------
48 
49 
50 !--only read file if beginning of run or start of new month
51  IF (debut.OR.mth_cur.NE.mth_pre) THEN
52 
53  IF (is_mpi_root) THEN
54 
55 
56  CALL nf95_open("chlorophyll.nc", nf90_nowrite, ncid_in)
57 
58  CALL nf95_inq_varid(ncid_in, "lon", varid)
59  CALL nf95_gw_var(ncid_in, varid, longitude)
60  n_lon = size(longitude)
61 ! print *, 'LON chlorophyll=', n_lon, longitude
62  IF (n_lon.NE.nbp_lon) THEN
63  print *,'Le nombre de lon n est pas egal a nbp_lon'
64  stop
65  ENDIF
66 
67 
68  CALL nf95_inq_varid(ncid_in, "lat", varid)
69  CALL nf95_gw_var(ncid_in, varid, latitude)
70  n_lat = size(latitude)
71 ! print *, 'LAT chlorophyll=', n_lat, latitude
72  IF (n_lat.NE.nbp_lat) THEN
73  print *,'Le nombre de lat n est pas egal a jnbp_lat'
74  stop
75  ENDIF
76 
77  CALL nf95_inq_varid(ncid_in, "time", varid)
78  CALL nf95_gw_var(ncid_in, varid, time)
79  n_month = size(time)
80 ! print *, 'TIME aerosol strato=', n_month, time
81  IF (n_month.NE.12) THEN
82  print *,'Le nombre de month n est pas egal a 12'
83  stop
84  ENDIF
85 
86  IF (.not.ALLOCATED(chlorocon)) ALLOCATE(chlorocon(n_lon, n_lat, n_month))
87  IF (.not.ALLOCATED(chlorocon_mois)) ALLOCATE(chlorocon_mois(n_lon, n_lat))
88  IF (.not.ALLOCATED(chlorocon_mois_glo)) ALLOCATE(chlorocon_mois_glo(klon_glo))
89 
90 !--reading stratospheric AOD at 550 nm
91  CALL nf95_inq_varid(ncid_in, "CHL", varid)
92  ncerr = nf90_get_var(ncid_in, varid, chlorocon)
93  print *,'code erreur readchlorophyll=', ncerr, varid
94 
95  CALL nf95_close(ncid_in)
96 
97 !---select the correct month
98  IF (mth_cur.LT.1.OR.mth_cur.GT.12) THEN
99  print *,'probleme avec le mois dans readchlorophyll =', mth_cur
100  ENDIF
101  chlorocon_mois(:,:) = chlorocon(:,:,mth_cur)
102 
103 !---reduce to a klon_glo grid
104  CALL grid2dto1d_glo(chlorocon_mois,chlorocon_mois_glo)
105 
106 
107  print*,"chrolophyll current month",mth_cur
108  do i=1,klon_glo
109 ! if(isnan(chlorocon_mois_glo(i)))then ! isnan() is not in the Fortran standard...
110 ! Another way to check for NaN:
111  if(chlorocon_mois_glo(i).ne.chlorocon_mois_glo(i)) then
112  chlorocon_mois_glo(i)=0.
113  endif
114  !print*,"high chl con",i,chlorocon_mois_glo(i)
115  enddo
116 
117 ! DEALLOCATE(chlorocon)
118 ! DEALLOCATE(chlorocon_mois)
119 ! DEALLOCATE(chlorocon_mois_glo)
120 
121  ENDIF !--is_mpi_root
122 
123 !--scatter on all proc
124  CALL scatter(chlorocon_mois_glo,chl_con)
125 
126 !--keep memory of previous month
127  mth_pre=mth_cur
128 
129  ENDIF !--debut ou nouveau mois
130 
131 end subroutine readchlorophyll
integer, save klon_glo
subroutine readchlorophyll(debut)
integer, save mth_cur
Definition: phys_cal_mod.F90:7
real, dimension(:), allocatable, save chl_con