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