Changeset 1610 for LMDZ5/trunk/libf/bibio
- Timestamp:
- Jan 23, 2012, 6:45:08 PM (13 years ago)
- File:
-
- 1 moved
Legend:
- Unmodified
- Added
- Removed
-
LMDZ5/trunk/libf/bibio/initdynav.F90
r1571 r1610 1 !2 1 ! $Id$ 3 ! 4 2 3 subroutine initdynav(day0,anne0,tstep,t_ops,t_wrt) 5 4 6 5 #ifdef CPP_IOIPSL 7 6 USE IOIPSL 8 7 #endif 9 10 use com_io_dyn_mod, only : histaveid,histvaveid,histuaveid,&11 &dynhistave_file,dynhistvave_file,dynhistuave_file12 8 USE infotrac, ONLY : nqtot, ttext 9 use com_io_dyn_mod, only : histaveid,histvaveid,histuaveid, & 10 dynhistave_file,dynhistvave_file,dynhistuave_file 11 implicit none 13 12 14 C15 C Routine d'initialisation des ecritures des fichiers histoires LMDZ16 C au format IOIPSL. Initialisation du fichier histoire moyenne.17 C18 C Appels succesifs des routines: histbeg19 C histhori20 C histver21 C histdef22 C histend23 C24 C Entree:25 C26 C infile: nom du fichier histoire a creer27 C day0,anne0: date de reference28 C tstep : frequence d'ecriture29 C t_ops: frequence de l'operation pour IOIPSL30 C t_wrt: frequence d'ecriture sur le fichier31 C32 C33 C L. Fairhead, LMD, 03/9934 C35 C =====================================================================36 C37 C Declarations38 #include "dimensions.h"39 #include "paramet.h"40 #include "comconst.h"41 #include "comvert.h"42 #include "comgeom.h"43 #include "temps.h"44 #include "ener.h"45 #include "logic.h"46 #include "description.h"47 #include "serre.h"48 #include "iniprint.h"49 13 50 C Arguments 51 C 52 integer day0, anne0 53 real tstep, t_ops, t_wrt 14 ! Routine d'initialisation des ecritures des fichiers histoires LMDZ 15 ! au format IOIPSL. Initialisation du fichier histoire moyenne. 16 17 ! Appels succesifs des routines: histbeg 18 ! histhori 19 ! histver 20 ! histdef 21 ! histend 22 23 ! Entree: 24 25 ! infile: nom du fichier histoire a creer 26 ! day0,anne0: date de reference 27 ! tstep : frequence d'ecriture 28 ! t_ops: frequence de l'operation pour IOIPSL 29 ! t_wrt: frequence d'ecriture sur le fichier 30 31 32 ! L. Fairhead, LMD, 03/99 33 34 include "dimensions.h" 35 include "paramet.h" 36 include "comconst.h" 37 include "comvert.h" 38 include "comgeom.h" 39 include "temps.h" 40 include "ener.h" 41 include "logic.h" 42 include "description.h" 43 include "serre.h" 44 include "iniprint.h" 45 46 ! Arguments 47 48 integer day0, anne0 49 real tstep, t_ops, t_wrt 54 50 55 51 #ifdef CPP_IOIPSL 56 ! This routine needs IOIPSL to work 57 C Variables locales 58 C 59 integer tau0 60 real zjulian 61 integer iq 62 real rlong(iip1,jjp1), rlat(iip1,jjp1) 63 integer uhoriid, vhoriid, thoriid, zvertiid 64 integer ii,jj 65 integer zan, dayref 66 C 67 C Initialisations 68 C 69 pi = 4. * atan (1.) 70 C 71 C Appel a histbeg: creation du fichier netcdf et initialisations diverses 72 C 52 ! This routine needs IOIPSL to work 53 ! Variables locales 73 54 74 zan = anne0 75 dayref = day0 76 CALL ymds2ju(zan, 1, dayref, 0.0, zjulian) 77 tau0 = itau_dyn 78 79 do jj = 1, jjp1 80 do ii = 1, iip1 81 rlong(ii,jj) = rlonv(ii) * 180. / pi 82 rlat(ii,jj) = rlatu(jj) * 180. / pi 83 enddo 84 enddo 85 86 ! Creation de 3 fichiers pour les differentes grilles horizontales 87 ! Restriction de IOIPSL: seulement 2 coordonnees dans le meme fichier 88 ! Grille Scalaire 89 call histbeg(dynhistave_file, iip1, rlong(:,1), jjp1, rlat(1,:), 90 . 1, iip1, 1, jjp1, 91 . tau0, zjulian, tstep, thoriid,histaveid) 55 integer tau0 56 real zjulian 57 integer iq 58 real rlong(iip1,jjp1), rlat(iip1,jjp1) 59 integer uhoriid, vhoriid, thoriid, zvertiid 60 integer ii,jj 61 integer zan, dayref 92 62 93 C Creation du fichier histoire pour les grilles en V et U (oblige pour l'instant, 94 C IOIPSL ne permet pas de grilles avec des nombres de point differents dans 95 C un meme fichier) 96 ! Grille V 97 do jj = 1, jjm 98 do ii = 1, iip1 99 rlong(ii,jj) = rlonv(ii) * 180. / pi 100 rlat(ii,jj) = rlatv(jj) * 180. / pi 101 enddo 102 enddo 63 !-------------------------------------------------------------------- 103 64 104 call histbeg(dynhistvave_file, iip1, rlong(:,1), jjm, rlat(1,:), 105 . 1, iip1, 1, jjm, 106 . tau0, zjulian, tstep, vhoriid,histvaveid) 107 ! Grille U 108 do jj = 1, jjp1 109 do ii = 1, iip1 110 rlong(ii,jj) = rlonu(ii) * 180. / pi 111 rlat(ii,jj) = rlatu(jj) * 180. / pi 112 enddo 113 enddo 65 ! Initialisations 114 66 115 call histbeg(dynhistuave_file, iip1, rlong(:,1),jjp1, rlat(1,:), 116 . 1, iip1, 1, jjp1, 117 . tau0, zjulian, tstep, uhoriid,histuaveid) 118 C 119 C Appel a histvert pour la grille verticale 120 C 121 call histvert(histaveid,'presnivs','Niveaux Pression 122 & approximatifs','mb',llm, presnivs/100., zvertiid,'down') 123 call histvert(histuaveid,'presnivs','Niveaux Pression 124 & approximatifs','mb',llm, presnivs/100., zvertiid,'down') 125 call histvert(histvaveid,'presnivs','Niveaux Pression 126 & approximatifs','mb',llm, presnivs/100., zvertiid,'down') 127 C 128 C Appels a histdef pour la definition des variables a sauvegarder 129 C 130 C Vents U 131 C 132 ! write(6,*)'inithistave',tstep 133 call histdef(histuaveid, 'u', 'vent u moyen ', 134 . 'm/s', iip1, jjp1, uhoriid, llm, 1, llm, zvertiid, 135 . 32, 'ave(X)', t_ops, t_wrt) 67 pi = 4. * atan (1.) 136 68 137 C Vents V 138 C 139 call histdef(histvaveid, 'v', 'vent v moyen', 140 . 'm/s', iip1, jjm, vhoriid, llm, 1, llm, zvertiid, 141 . 32, 'ave(X)', t_ops, t_wrt) 69 ! Appel a histbeg: creation du fichier netcdf et initialisations diverses 142 70 143 C 144 C Temperature 145 C 146 call histdef(histaveid, 'temp', 'temperature moyenne', 'K', 147 . iip1, jjp1, thoriid, llm, 1, llm, zvertiid, 148 . 32, 'ave(X)', t_ops, t_wrt) 149 C 150 C Temperature potentielle 151 C 152 call histdef(histaveid, 'theta', 'temperature potentielle', 'K', 153 . iip1, jjp1, thoriid, llm, 1, llm, zvertiid, 154 . 32, 'ave(X)', t_ops, t_wrt) 155 C 156 C Geopotentiel 157 C 158 call histdef(histaveid, 'phi', 'geopotentiel moyen', '-', 159 . iip1, jjp1, thoriid, llm, 1, llm, zvertiid, 160 . 32, 'ave(X)', t_ops, t_wrt) 161 C 162 C Traceurs 163 C 164 ! DO iq=1,nqtot 165 ! call histdef(histaveid, ttext(iq), ttext(iq), '-', 166 ! . iip1, jjp1, thoriid, llm, 1, llm, zvertiid, 167 ! . 32, 'ave(X)', t_ops, t_wrt) 168 ! enddo 169 C 170 C Masse 171 C 172 call histdef(histaveid, 'masse', 'masse', 'kg', 173 . iip1, jjp1, thoriid, llm, 1, llm, zvertiid, 174 . 32, 'ave(X)', t_ops, t_wrt) 175 C 176 C Pression au sol 177 C 178 call histdef(histaveid, 'ps', 'pression naturelle au sol', 'Pa', 179 . iip1, jjp1, thoriid, 1, 1, 1, -99, 180 . 32, 'ave(X)', t_ops, t_wrt) 181 C 182 C Geopotentiel au sol 183 C 184 ! call histdef(histaveid, 'phis', 'geopotentiel au sol', '-', 185 ! . iip1, jjp1, thoriid, 1, 1, 1, -99, 186 ! . 32, 'ave(X)', t_ops, t_wrt) 187 !C 188 C Fin 189 C 190 call histend(histaveid) 191 call histend(histuaveid) 192 call histend(histvaveid) 71 72 zan = anne0 73 dayref = day0 74 CALL ymds2ju(zan, 1, dayref, 0.0, zjulian) 75 tau0 = itau_dyn 76 77 do jj = 1, jjp1 78 do ii = 1, iip1 79 rlong(ii,jj) = rlonv(ii) * 180. / pi 80 rlat(ii,jj) = rlatu(jj) * 180. / pi 81 enddo 82 enddo 83 84 ! Creation de 3 fichiers pour les differentes grilles horizontales 85 ! Restriction de IOIPSL: seulement 2 coordonnees dans le meme fichier 86 ! Grille Scalaire 87 call histbeg(dynhistave_file, iip1, rlong(:,1), jjp1, rlat(1,:), & 88 1, iip1, 1, jjp1, & 89 tau0, zjulian, tstep, thoriid,histaveid) 90 91 ! Creation du fichier histoire pour les grilles en V et U (oblige 92 ! pour l'instant, IOIPSL ne permet pas de grilles avec des nombres 93 ! de point differents dans un meme fichier) 94 ! Grille V 95 do jj = 1, jjm 96 do ii = 1, iip1 97 rlong(ii,jj) = rlonv(ii) * 180. / pi 98 rlat(ii,jj) = rlatv(jj) * 180. / pi 99 enddo 100 enddo 101 102 call histbeg(dynhistvave_file, iip1, rlong(:,1), jjm, rlat(1,:), & 103 1, iip1, 1, jjm, & 104 tau0, zjulian, tstep, vhoriid,histvaveid) 105 ! Grille U 106 do jj = 1, jjp1 107 do ii = 1, iip1 108 rlong(ii,jj) = rlonu(ii) * 180. / pi 109 rlat(ii,jj) = rlatu(jj) * 180. / pi 110 enddo 111 enddo 112 113 call histbeg(dynhistuave_file, iip1, rlong(:,1),jjp1, rlat(1,:), & 114 1, iip1, 1, jjp1, & 115 tau0, zjulian, tstep, uhoriid,histuaveid) 116 117 ! Appel a histvert pour la grille verticale 118 119 call histvert(histaveid,'presnivs','Niveaux Pression & 120 approximatifs','mb',llm, presnivs/100., zvertiid,'down') 121 call histvert(histuaveid,'presnivs','Niveaux Pression & 122 approximatifs','mb',llm, presnivs/100., zvertiid,'down') 123 call histvert(histvaveid,'presnivs','Niveaux Pression & 124 approximatifs','mb',llm, presnivs/100., zvertiid,'down') 125 126 ! Appels a histdef pour la definition des variables a sauvegarder 127 128 ! Vents U 129 130 call histdef(histuaveid, 'u', 'vent u moyen ', & 131 'm/s', iip1, jjp1, uhoriid, llm, 1, llm, zvertiid, & 132 32, 'ave(X)', t_ops, t_wrt) 133 134 ! Vents V 135 136 call histdef(histvaveid, 'v', 'vent v moyen', & 137 'm/s', iip1, jjm, vhoriid, llm, 1, llm, zvertiid, & 138 32, 'ave(X)', t_ops, t_wrt) 139 140 141 ! Temperature 142 143 call histdef(histaveid, 'temp', 'temperature moyenne', 'K', & 144 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & 145 32, 'ave(X)', t_ops, t_wrt) 146 147 ! Temperature potentielle 148 149 call histdef(histaveid, 'theta', 'temperature potentielle', 'K', & 150 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & 151 32, 'ave(X)', t_ops, t_wrt) 152 153 ! Geopotentiel 154 155 call histdef(histaveid, 'phi', 'geopotentiel moyen', '-', & 156 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & 157 32, 'ave(X)', t_ops, t_wrt) 158 159 ! Traceurs 160 161 ! DO iq=1,nqtot 162 ! call histdef(histaveid, ttext(iq), ttext(iq), '-', & 163 ! iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & 164 ! 32, 'ave(X)', t_ops, t_wrt) 165 ! enddo 166 167 ! Masse 168 169 call histdef(histaveid, 'masse', 'masse', 'kg', & 170 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, & 171 32, 'ave(X)', t_ops, t_wrt) 172 173 ! Pression au sol 174 175 call histdef(histaveid, 'ps', 'pression naturelle au sol', 'Pa', & 176 iip1, jjp1, thoriid, 1, 1, 1, -99, & 177 32, 'ave(X)', t_ops, t_wrt) 178 179 ! Geopotentiel au sol 180 181 ! call histdef(histaveid, 'phis', 'geopotentiel au sol', '-', & 182 ! iip1, jjp1, thoriid, 1, 1, 1, -99, & 183 ! 32, 'ave(X)', t_ops, t_wrt) 184 185 call histend(histaveid) 186 call histend(histuaveid) 187 call histend(histvaveid) 193 188 #else 194 ! tell the user this routine should be run with ioipsl 195 write(lunout,*)"initdynav: Warning this routine should not be", 196 & " used without ioipsl" 189 write(lunout,*)"initdynav: Warning this routine should not be", & 190 " used without ioipsl" 197 191 #endif 198 ! of #ifdef CPP_IOIPSL199 return 200 end 192 ! of #ifdef CPP_IOIPSL 193 194 end subroutine initdynav
Note: See TracChangeset
for help on using the changeset viewer.