1 |
MODULE dynredem0_m |
2 |
|
3 |
IMPLICIT NONE |
4 |
|
5 |
CONTAINS |
6 |
|
7 |
SUBROUTINE dynredem0(fichnom, iday_end, phis) |
8 |
|
9 |
! From dyn3d/dynredem.F, version 1.2, 2004/06/22 11:45:30 |
10 |
! \'Ecriture du fichier de red\'emarrage au format NetCDF (initialisation) |
11 |
|
12 |
USE comconst, ONLY: cpp, daysec, dtvr, g, kappa, omeg, rad |
13 |
USE dimens_m, ONLY: iim, jjm, llm, nqmx |
14 |
USE disvert_m, ONLY: ap, bp, pa, preff, presnivs |
15 |
use dynetat0_m, only: day_ref, annee_ref, clat, clon, dzoomx, dzoomy, & |
16 |
grossismx, grossismy, taux, tauy, rlatu, rlatv, rlonu, rlonv, rlatu1, & |
17 |
rlatu2, yprimu1, yprimu2, xprimp025, xprimm025, xprimu, xprimv |
18 |
USE ener, ONLY: ang0, etot0, ptot0, stot0, ztot0 |
19 |
USE iniadvtrac_m, ONLY: tname, ttext |
20 |
USE ju2ymds_m, ONLY: ju2ymds |
21 |
USE netcdf, ONLY: nf90_clobber, nf90_float, nf90_global, nf90_unlimited |
22 |
USE netcdf95, ONLY: nf95_close, nf95_create, nf95_def_dim, nf95_def_var, & |
23 |
nf95_enddef, nf95_inq_varid, nf95_put_att, nf95_put_var |
24 |
USE paramet_m, ONLY: iip1, jjp1, llmp1 |
25 |
use ymds2ju_m, only: ymds2ju |
26 |
|
27 |
CHARACTER(len=*), INTENT(IN):: fichnom |
28 |
INTEGER, INTENT(IN):: iday_end |
29 |
REAL, INTENT(IN):: phis(:, :) |
30 |
|
31 |
! Local: |
32 |
|
33 |
INTEGER iq |
34 |
INTEGER, PARAMETER:: length = 100 |
35 |
REAL tab_cntrl(length) ! tableau des param\`etres du run |
36 |
|
37 |
! Pour NetCDF : |
38 |
INTEGER idim_index |
39 |
INTEGER idim_rlonu, idim_rlonv, idim_rlatu, idim_rlatv |
40 |
INTEGER idim_s, idim_sig |
41 |
INTEGER dimid_temps |
42 |
INTEGER ncid, varid |
43 |
integer varid_controle, varid_rlonu, varid_rlatu, varid_rlonv, varid_rlatv |
44 |
integer varid_xprimu, varid_xprimv, varid_xprimm025, varid_xprimp025 |
45 |
integer varid_rlatu1, varid_rlatu2, varid_yprimu1, varid_yprimu2, varid_ap |
46 |
integer varid_bp, varid_presnivs, varid_phisinit |
47 |
|
48 |
REAL zjulian, hours |
49 |
INTEGER yyears0, jjour0, mmois0 |
50 |
CHARACTER(len=30) unites |
51 |
|
52 |
!----------------------------------------------------------------------- |
53 |
|
54 |
PRINT *, 'Call sequence information: dynredem0' |
55 |
|
56 |
CALL ymds2ju(annee_ref, 1, iday_end, 0., zjulian) |
57 |
CALL ju2ymds(zjulian, yyears0, mmois0, jjour0, hours) |
58 |
|
59 |
tab_cntrl(1) = iim |
60 |
tab_cntrl(2) = jjm |
61 |
tab_cntrl(3) = llm |
62 |
tab_cntrl(4) = day_ref |
63 |
tab_cntrl(5) = annee_ref |
64 |
tab_cntrl(6) = rad |
65 |
tab_cntrl(7) = omeg |
66 |
tab_cntrl(8) = g |
67 |
tab_cntrl(9) = cpp |
68 |
tab_cntrl(10) = kappa |
69 |
tab_cntrl(11) = daysec |
70 |
tab_cntrl(12) = dtvr |
71 |
tab_cntrl(13) = etot0 |
72 |
tab_cntrl(14) = ptot0 |
73 |
tab_cntrl(15) = ztot0 |
74 |
tab_cntrl(16) = stot0 |
75 |
tab_cntrl(17) = ang0 |
76 |
tab_cntrl(18) = pa |
77 |
tab_cntrl(19) = preff |
78 |
|
79 |
! Param\`etres pour le zoom : |
80 |
tab_cntrl(20) = clon |
81 |
tab_cntrl(21) = clat |
82 |
tab_cntrl(22) = grossismx |
83 |
tab_cntrl(23) = grossismy |
84 |
tab_cntrl(24) = 1. |
85 |
tab_cntrl(25) = dzoomx |
86 |
tab_cntrl(26) = dzoomy |
87 |
tab_cntrl(27) = 0. |
88 |
tab_cntrl(28) = taux |
89 |
tab_cntrl(29) = tauy |
90 |
|
91 |
tab_cntrl(30) = iday_end |
92 |
tab_cntrl(31:) = 0. |
93 |
|
94 |
CALL nf95_create(fichnom, nf90_clobber, ncid) |
95 |
CALL nf95_put_att(ncid, nf90_global, 'title', & |
96 |
'start file for the dynamics code') |
97 |
|
98 |
! Definir les dimensions du fichiers: |
99 |
|
100 |
CALL nf95_def_dim(ncid, 'index', length, idim_index) |
101 |
CALL nf95_def_dim(ncid, 'rlonu', iip1, idim_rlonu) |
102 |
CALL nf95_def_dim(ncid, 'rlatu', jjp1, idim_rlatu) |
103 |
CALL nf95_def_dim(ncid, 'rlonv', iip1, idim_rlonv) |
104 |
CALL nf95_def_dim(ncid, 'rlatv', jjm, idim_rlatv) |
105 |
CALL nf95_def_dim(ncid, 'sigs', llm, idim_s) |
106 |
CALL nf95_def_dim(ncid, 'sig', llmp1, idim_sig) |
107 |
CALL nf95_def_dim(ncid, 'temps', nf90_unlimited, dimid_temps) |
108 |
|
109 |
! Definir et enregistrer certains champs invariants: |
110 |
|
111 |
CALL nf95_def_var(ncid, 'controle', nf90_float, idim_index, varid_controle) |
112 |
CALL nf95_put_att(ncid, varid_controle, 'title', 'Parametres de controle') |
113 |
|
114 |
CALL nf95_def_var(ncid, 'rlonu', nf90_float, idim_rlonu, varid_rlonu) |
115 |
CALL nf95_put_att(ncid, varid_rlonu, 'title', 'Longitudes des points U') |
116 |
|
117 |
CALL nf95_def_var(ncid, 'rlatu', nf90_float, idim_rlatu, varid_rlatu) |
118 |
CALL nf95_put_att(ncid, varid_rlatu, 'title', 'Latitudes des points U') |
119 |
|
120 |
CALL nf95_def_var(ncid, 'rlonv', nf90_float, idim_rlonv, varid_rlonv) |
121 |
CALL nf95_put_att(ncid, varid_rlonv, 'title', 'Longitudes des points V') |
122 |
|
123 |
CALL nf95_def_var(ncid, 'rlatv', nf90_float, idim_rlatv, varid_rlatv) |
124 |
CALL nf95_put_att(ncid, varid_rlatv, 'title', 'Latitudes des points V') |
125 |
|
126 |
CALL nf95_def_var(ncid, 'xprimu', nf90_float, idim_rlonu, varid_xprimu) |
127 |
CALL nf95_put_att(ncid, varid_xprimu, 'title', 'dx / dX aux points u') |
128 |
|
129 |
CALL nf95_def_var(ncid, 'xprimv', nf90_float, idim_rlonv, varid_xprimv) |
130 |
CALL nf95_put_att(ncid, varid_xprimv, 'title', 'dx / dX aux points v') |
131 |
|
132 |
CALL nf95_def_var(ncid, 'xprimm025', nf90_float, idim_rlonu, & |
133 |
varid_xprimm025) |
134 |
CALL nf95_def_var(ncid, 'xprimp025', nf90_float, idim_rlonu, & |
135 |
varid_xprimp025) |
136 |
|
137 |
CALL nf95_def_var(ncid, 'rlatu1', nf90_float, idim_rlatv, varid_rlatu1) |
138 |
CALL nf95_def_var(ncid, 'rlatu2', nf90_float, idim_rlatv, varid_rlatu2) |
139 |
CALL nf95_def_var(ncid, 'yprimu1', nf90_float, idim_rlatv, varid_yprimu1) |
140 |
CALL nf95_def_var(ncid, 'yprimu2', nf90_float, idim_rlatv, varid_yprimu2) |
141 |
|
142 |
CALL nf95_def_var(ncid, 'ap', nf90_float, idim_sig, varid_ap) |
143 |
CALL nf95_put_att(ncid, varid_ap, 'title', 'Coefficient A pour hybride') |
144 |
|
145 |
CALL nf95_def_var(ncid, 'bp', nf90_float, idim_sig, varid_bp) |
146 |
CALL nf95_put_att(ncid, varid_bp, 'title', 'Coefficient B pour hybride') |
147 |
|
148 |
CALL nf95_def_var(ncid, 'presnivs', nf90_float, idim_s, varid_presnivs) |
149 |
|
150 |
! Geopentiel au sol: |
151 |
|
152 |
CALL nf95_def_var(ncid, 'phisinit', nf90_float, & |
153 |
(/idim_rlonv, idim_rlatu/), varid_phisinit) |
154 |
CALL nf95_put_att(ncid, varid_phisinit, 'title', 'Geopotentiel au sol') |
155 |
|
156 |
! Definir les variables pour pouvoir les enregistrer plus tard: |
157 |
|
158 |
CALL nf95_def_var(ncid, 'temps', nf90_float, dimid_temps, varid) |
159 |
CALL nf95_put_att(ncid, varid, 'title', 'Temps de simulation') |
160 |
WRITE(unites, fmt = 200) yyears0, mmois0, jjour0 |
161 |
200 FORMAT ('days since ', I4, '-', I2.2, '-', I2.2, ' 00:00:00') |
162 |
CALL nf95_put_att(ncid, varid, 'units', unites) |
163 |
|
164 |
CALL nf95_def_var(ncid, 'ucov', nf90_float, & |
165 |
(/idim_rlonu, idim_rlatu, idim_s, dimid_temps/), varid) |
166 |
CALL nf95_put_att(ncid, varid, 'title', 'Vitesse U') |
167 |
|
168 |
CALL nf95_def_var(ncid, 'vcov', nf90_float, & |
169 |
(/idim_rlonv, idim_rlatv, idim_s, dimid_temps/), varid) |
170 |
CALL nf95_put_att(ncid, varid, 'title', 'Vitesse V') |
171 |
|
172 |
CALL nf95_def_var(ncid, 'teta', nf90_float, & |
173 |
(/idim_rlonv, idim_rlatu, idim_s, dimid_temps/), varid) |
174 |
CALL nf95_put_att(ncid, varid, 'title', 'Temperature') |
175 |
|
176 |
DO iq = 1, nqmx |
177 |
CALL nf95_def_var(ncid, tname(iq), nf90_float, & |
178 |
(/idim_rlonv, idim_rlatu, idim_s, dimid_temps/), varid) |
179 |
CALL nf95_put_att(ncid, varid, 'title', ttext(iq)) |
180 |
END DO |
181 |
|
182 |
CALL nf95_def_var(ncid, 'masse', nf90_float, & |
183 |
(/idim_rlonv, idim_rlatu, idim_s, dimid_temps/), varid) |
184 |
CALL nf95_put_att(ncid, varid, 'title', 'C est quoi ?') |
185 |
|
186 |
CALL nf95_def_var(ncid, 'ps', nf90_float, & |
187 |
(/idim_rlonv, idim_rlatu, dimid_temps/), varid) |
188 |
CALL nf95_put_att(ncid, varid, 'title', 'Pression au sol') |
189 |
|
190 |
CALL nf95_enddef(ncid) |
191 |
|
192 |
CALL nf95_put_var(ncid, varid_controle, tab_cntrl) |
193 |
CALL nf95_put_var(ncid, varid_rlonu, rlonu) |
194 |
CALL nf95_put_var(ncid, varid_rlatu, rlatu) |
195 |
CALL nf95_put_var(ncid, varid_rlonv, rlonv) |
196 |
CALL nf95_put_var(ncid, varid_rlatv, rlatv) |
197 |
CALL nf95_put_var(ncid, varid_xprimu, xprimu) |
198 |
CALL nf95_put_var(ncid, varid_xprimv, xprimv) |
199 |
CALL nf95_put_var(ncid, varid_xprimm025, xprimm025) |
200 |
CALL nf95_put_var(ncid, varid_xprimp025, xprimp025) |
201 |
call NF95_PUT_VAR(ncid, varid_rlatu1, rlatu1) |
202 |
call NF95_PUT_VAR(ncid, varid_rlatu2, rlatu2) |
203 |
CALL nf95_put_var(ncid, varid_yprimu1, yprimu1) |
204 |
CALL nf95_put_var(ncid, varid_yprimu2, yprimu2) |
205 |
CALL nf95_put_var(ncid, varid_ap, ap) |
206 |
CALL nf95_put_var(ncid, varid_bp, bp) |
207 |
CALL nf95_put_var(ncid, varid_presnivs, presnivs) |
208 |
CALL nf95_put_var(ncid, varid_phisinit, phis) |
209 |
|
210 |
CALL nf95_close(ncid) |
211 |
|
212 |
PRINT *, 'iim, jjm, llm, iday_end', iim, jjm, llm, iday_end |
213 |
PRINT *, 'rad, omeg, g, cpp, kappa', rad, omeg, g, cpp, kappa |
214 |
|
215 |
END SUBROUTINE dynredem0 |
216 |
|
217 |
END MODULE dynredem0_m |