| 1 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 2 | ! Interface pour ecrire en netcdf avec les routines d'enseignement |
|---|
| 3 | ! iotd de Frederic Hourdin |
|---|
| 4 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 5 | |
|---|
| 6 | subroutine iophys_ecrit(nom,lllm,titre,unite,px) |
|---|
| 7 | |
|---|
| 8 | USE mod_phys_lmdz_para, ONLY: klon_omp, is_mpi_root |
|---|
| 9 | USE mod_phys_lmdz_transfert_para, ONLY: gather |
|---|
| 10 | USE mod_grid_phy_lmdz, ONLY: klon_glo, nbp_lon, nbp_lat, grid1dto2d_glo |
|---|
| 11 | IMPLICIT NONE |
|---|
| 12 | |
|---|
| 13 | |
|---|
| 14 | |
|---|
| 15 | ! Ecriture de variables diagnostiques au choix dans la physique |
|---|
| 16 | ! dans un fichier NetCDF nomme 'diagfi'. Ces variables peuvent etre |
|---|
| 17 | ! 3d (ex : temperature), 2d (ex : temperature de surface), ou |
|---|
| 18 | ! 0d (pour un scalaire qui ne depend que du temps : ex : la longitude |
|---|
| 19 | ! solaire) |
|---|
| 20 | ! (ou encore 1d, dans le cas de testphys1d, pour sortir une colonne) |
|---|
| 21 | ! La periode d'ecriture est donnee par |
|---|
| 22 | ! "ecritphy " regle dans le fichier de controle de run : run.def |
|---|
| 23 | ! |
|---|
| 24 | ! writediagfi peut etre appele de n'importe quelle subroutine |
|---|
| 25 | ! de la physique, plusieurs fois. L'initialisation et la creation du |
|---|
| 26 | ! fichier se fait au tout premier appel. |
|---|
| 27 | ! |
|---|
| 28 | ! WARNING : les variables dynamique (u,v,t,q,ps) |
|---|
| 29 | ! sauvees par writediagfi avec une |
|---|
| 30 | ! date donnee sont legerement differentes que dans le fichier histoire car |
|---|
| 31 | ! on ne leur a pas encore ajoute de la dissipation et de la physique !!! |
|---|
| 32 | ! IL est RECOMMANDE d'ajouter les tendance physique a ces variables |
|---|
| 33 | ! avant l'ecriture dans diagfi (cf. physiq.F) |
|---|
| 34 | ! |
|---|
| 35 | ! Modifs: Aug.2010 Ehouarn: enforce outputs to be real*4 |
|---|
| 36 | ! |
|---|
| 37 | ! parametres (input) : |
|---|
| 38 | ! ---------- |
|---|
| 39 | ! unit : unite logique du fichier de sortie (toujours la meme) |
|---|
| 40 | ! nom : nom de la variable a sortir (chaine de caracteres) |
|---|
| 41 | ! titre: titre de la variable (chaine de caracteres) |
|---|
| 42 | ! unite : unite de la variable (chaine de caracteres) |
|---|
| 43 | ! px : variable a sortir (real 0, 1, 2, ou 3d) |
|---|
| 44 | ! |
|---|
| 45 | !================================================================= |
|---|
| 46 | |
|---|
| 47 | |
|---|
| 48 | ! Arguments on input: |
|---|
| 49 | integer lllm |
|---|
| 50 | character (len=*) :: nom,titre,unite |
|---|
| 51 | integer imjmax |
|---|
| 52 | parameter (imjmax=100000) |
|---|
| 53 | real px(klon_omp,lllm) |
|---|
| 54 | real xglo(klon_glo,lllm) |
|---|
| 55 | real zx(nbp_lon,nbp_lat,lllm) |
|---|
| 56 | |
|---|
| 57 | |
|---|
| 58 | CALL Gather(px,xglo) |
|---|
| 59 | !$OMP MASTER |
|---|
| 60 | IF (is_mpi_root) THEN |
|---|
| 61 | CALL Grid1Dto2D_glo(xglo,zx) |
|---|
| 62 | call iotd_ecrit(nom,lllm,titre,unite,zx) |
|---|
| 63 | ENDIF |
|---|
| 64 | !$OMP END MASTER |
|---|
| 65 | |
|---|
| 66 | return |
|---|
| 67 | end |
|---|
| 68 | |
|---|
| 69 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 70 | ! Version avec reindexation pour appeler depuis les routines internes |
|---|
| 71 | ! à la sous surface |
|---|
| 72 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 73 | |
|---|
| 74 | |
|---|
| 75 | subroutine iophys_ecrit_index(nom,lllm,titre,unite,knon,knindex,px) |
|---|
| 76 | |
|---|
| 77 | USE mod_phys_lmdz_para, ONLY: klon_omp |
|---|
| 78 | USE dimphy, ONLY : klon |
|---|
| 79 | USE mod_grid_phy_lmdz, ONLY: klon_glo |
|---|
| 80 | IMPLICIT NONE |
|---|
| 81 | |
|---|
| 82 | ! This subroutine returns the sea surface temperature already read from limit.nc |
|---|
| 83 | ! |
|---|
| 84 | |
|---|
| 85 | ! Arguments on input: |
|---|
| 86 | INTEGER lllm |
|---|
| 87 | CHARACTER (len=*) :: nom,titre,unite |
|---|
| 88 | REAL px(klon_omp,lllm) |
|---|
| 89 | INTEGER, INTENT(IN) :: knon ! nomber of points on compressed grid |
|---|
| 90 | INTEGER, DIMENSION(klon), INTENT(IN) :: knindex ! grid point number for compressed grid |
|---|
| 91 | REAL, DIMENSION(klon,lllm) :: varout |
|---|
| 92 | |
|---|
| 93 | INTEGER :: i,l |
|---|
| 94 | |
|---|
| 95 | IF (klon/=klon_omp) THEN |
|---|
| 96 | print*,'klon, klon_omp',klon,klon_omp |
|---|
| 97 | STOP'probleme de dimension parallele' |
|---|
| 98 | ENDIF |
|---|
| 99 | |
|---|
| 100 | varout(1:klon,1:lllm)=0. |
|---|
| 101 | DO l = 1, lllm |
|---|
| 102 | DO i = 1, knon |
|---|
| 103 | varout(knindex(i),l) = px(i,l) |
|---|
| 104 | END DO |
|---|
| 105 | END DO |
|---|
| 106 | CALL iophys_ecrit(nom,lllm,titre,unite,varout) |
|---|
| 107 | |
|---|
| 108 | END SUBROUTINE iophys_ecrit_index |
|---|
| 109 | |
|---|
| 110 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 111 | SUBROUTINE iophys_ini |
|---|
| 112 | USE mod_phys_lmdz_para, ONLY: is_mpi_root |
|---|
| 113 | USE vertical_layers_mod, ONLY: presnivs |
|---|
| 114 | USE regular_lonlat_mod, ONLY: lon_reg, lat_reg |
|---|
| 115 | USE dimphy, ONLY: klev |
|---|
| 116 | USE mod_grid_phy_lmdz, ONLY: klon_glo |
|---|
| 117 | |
|---|
| 118 | IMPLICIT NONE |
|---|
| 119 | |
|---|
| 120 | !======================================================================= |
|---|
| 121 | ! |
|---|
| 122 | ! Auteur: L. Fairhead , P. Le Van, Y. Wanherdrick, F. Forget |
|---|
| 123 | ! ------- |
|---|
| 124 | ! |
|---|
| 125 | ! Objet: |
|---|
| 126 | ! ------ |
|---|
| 127 | ! |
|---|
| 128 | ! 'Initialize' the diagfi.nc file: write down dimensions as well |
|---|
| 129 | ! as time-independent fields (e.g: geopotential, mesh area, ...) |
|---|
| 130 | ! |
|---|
| 131 | !======================================================================= |
|---|
| 132 | !----------------------------------------------------------------------- |
|---|
| 133 | ! Declarations: |
|---|
| 134 | ! ------------- |
|---|
| 135 | |
|---|
| 136 | real pi |
|---|
| 137 | INTEGER nlat_eff |
|---|
| 138 | |
|---|
| 139 | ! Arguments: |
|---|
| 140 | ! ---------- |
|---|
| 141 | |
|---|
| 142 | !$OMP MASTER |
|---|
| 143 | IF (is_mpi_root) THEN |
|---|
| 144 | |
|---|
| 145 | ! Bidouille pour gerer le fait que lat_reg contient deux latitudes |
|---|
| 146 | ! en version uni-dimensionnelle (chose qui pourrait être résolue |
|---|
| 147 | ! par ailleurs) |
|---|
| 148 | IF (klon_glo==1) THEN |
|---|
| 149 | nlat_eff=1 |
|---|
| 150 | ELSE |
|---|
| 151 | nlat_eff=size(lat_reg) |
|---|
| 152 | ENDIF |
|---|
| 153 | pi=2.*asin(1.) |
|---|
| 154 | call iotd_ini('phys.nc ', & |
|---|
| 155 | size(lon_reg),nlat_eff,klev,lon_reg(:)*180./pi,lat_reg*180./pi,presnivs) |
|---|
| 156 | ENDIF |
|---|
| 157 | !$OMP END MASTER |
|---|
| 158 | |
|---|
| 159 | END |
|---|
| 160 | |
|---|
| 161 | #ifdef und |
|---|
| 162 | SUBROUTINE gr_fi_ecrit(nfield,nlon,iim,jjmp1,fi,ecrit) |
|---|
| 163 | IMPLICIT none |
|---|
| 164 | |
|---|
| 165 | !======================================================================= |
|---|
| 166 | INTEGER nfield,nlon,iim,jjmp1, jjm |
|---|
| 167 | REAL fi(nlon,nfield), ecrit(iim*jjmp1,nfield) |
|---|
| 168 | |
|---|
| 169 | INTEGER i, n, ig |
|---|
| 170 | |
|---|
| 171 | jjm = jjmp1 - 1 |
|---|
| 172 | DO n = 1, nfield |
|---|
| 173 | DO i=1,iim |
|---|
| 174 | ecrit(i,n) = fi(1,n) |
|---|
| 175 | ecrit(i+jjm*iim,n) = fi(nlon,n) |
|---|
| 176 | ENDDO |
|---|
| 177 | DO ig = 1, nlon - 2 |
|---|
| 178 | ecrit(iim+ig,n) = fi(1+ig,n) |
|---|
| 179 | ENDDO |
|---|
| 180 | ENDDO |
|---|
| 181 | RETURN |
|---|
| 182 | END |
|---|
| 183 | |
|---|
| 184 | #endif |
|---|