source: trunk/LMDZ.MARS/libf/phymars/inistats.F @ 1242

Last change on this file since 1242 was 1130, checked in by emillour, 11 years ago

Mars GCM:
Series of changes to enable running in parallel (using LMDZ.COMMON dynamics);
Current LMDZ.MARS can still notheless be compiled and run in serial mode
"as previously".
Summary of main changes:

  • Main programs (newstart, start2archive, xvik) that used to be in dyn3d have been moved to phymars.
  • dyn3d/control.h is now module control_mod.F90
  • rearanged input/outputs routines everywhere to handle serial/MPI cases. physdem.F => phyredem.F90 , phyetat0.F => phyetat0.F90 ; all read/write routines for startfi files are gathered in module iostart.F90
  • added parallelism related routines init_phys_lmdz.F90, comgeomphy.F90, dimphy.F90, iniphysiq.F90, mod_grid_phy_lmdz.F90, mod_phys_lmdz_mpi_data.F90, mod_phys_lmdz_mpi_transfert.F90, mod_phys_lmdz_omp_data.F90, mod_phys_lmdz_omp_transfert.F90, mod_phys_lmdz_para.F90, mod_phys_lmdz_transfert_para.F90 in phymars and mod_const_mpi.F90 in dyn3d (for compliance with parallel case)
  • created generic routines 'planetwide_maxval' and 'planetwide_minval', in module "planetwide_mod", that enable obtaining the min and max of a field over the whole planet.

EM

File size: 3.6 KB
RevLine 
[38]1      subroutine inistats(ierr)
2
[1130]3      use mod_phys_lmdz_para, only : is_master
[38]4      implicit none
5
6#include "dimensions.h"
7#include "paramet.h"
8#include "comgeom.h"
9#include "comvert.h"
10#include "comconst.h"
11#include "statto.h"
12#include "netcdf.inc"
13
14      integer,intent(out) :: ierr
15      integer :: nid
16      integer :: l,nsteppd
17      real, dimension(llm) ::  sig_s
18      integer :: idim_lat,idim_lon,idim_llm,idim_llmp1,idim_time
19      real, dimension(istime) :: lt
20      integer :: nvarid
21
22      write (*,*)
23      write (*,*) '                        || STATS ||'
24      write (*,*)
25      write (*,*) 'daysec',daysec
26      write (*,*) 'dtphys',dtphys
27      nsteppd=nint(daysec/dtphys)
28      write (*,*) 'nsteppd=',nsteppd
29      if (abs(float(nsteppd)-daysec/dtphys).gt.1.e-8*daysec)
30     &   stop'Dans Instat:  1jour .ne. n pas physiques'
31
32      if(mod(nsteppd,istime).ne.0)
33     &   stop'Dans Instat:  1jour .ne. n*istime pas physiques'
34
35      istats=nsteppd/istime
36      write (*,*) 'istats=',istats
37      write (*,*) 'Storing ',istime,'times per day'
38      write (*,*) 'thus every ',istats,'physical timestep '
39      write (*,*)
40
41      do l= 1, llm
42         sig_s(l)=((ap(l)+ap(l+1))/preff+bp(l)+bp(l+1))/2.
43         pseudoalt(l)=-10.*log(presnivs(l)/preff)   
44      enddo
45
[1130]46      if (is_master) then
47      ! only the master needs do this
48
[410]49      ierr = NF_CREATE("stats.nc",IOR(NF_CLOBBER,NF_64BIT_OFFSET),nid)
[38]50      if (ierr.ne.NF_NOERR) then
51         write (*,*) NF_STRERROR(ierr)
52         stop ""
53      endif
54
55      ierr = NF_DEF_DIM (nid, "latitude", jjp1, idim_lat)
56      ierr = NF_DEF_DIM (nid, "longitude", iip1, idim_lon)
57      ierr = NF_DEF_DIM (nid, "altitude", llm, idim_llm)
58      ierr = NF_DEF_DIM (nid, "llmp1", llm+1, idim_llmp1)
59      ierr = NF_DEF_DIM (nid, "Time", NF_UNLIMITED, idim_time)
60
61      ierr = NF_ENDDEF(nid)
62      call def_var_stats(nid,"Time","Time",
63     &            "days since 0000-00-0 00:00:00",1,
64     &            idim_time,nvarid,ierr)
65! Time is initialised later by mkstats subroutine
66
67      call def_var_stats(nid,"latitude","latitude",
68     &            "degrees_north",1,idim_lat,nvarid,ierr)
69#ifdef NC_DOUBLE
70      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,rlatu/pi*180)
71#else
72      ierr = NF_PUT_VAR_REAL (nid,nvarid,rlatu/pi*180)
73#endif
74      call def_var_stats(nid,"longitude","East longitude",
75     &            "degrees_east",1,idim_lon,nvarid,ierr)
76#ifdef NC_DOUBLE
77      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,rlonv/pi*180)
78#else
79      ierr = NF_PUT_VAR_REAL (nid,nvarid,rlonv/pi*180)
80#endif
81
82! Niveaux verticaux, aps et bps
83      ierr = NF_REDEF (nid)
84#ifdef NC_DOUBLE
85      ierr = NF_DEF_VAR (nid,"altitude", NF_DOUBLE, 1,idim_llm,nvarid)
86#else
87      ierr = NF_DEF_VAR (nid,"altitude", NF_FLOAT, 1,idim_llm,nvarid)
88#endif
89      ierr = NF_PUT_ATT_TEXT (nid,nvarid,"long_name",8,"altitude")
90      ierr = NF_PUT_ATT_TEXT (nid,nvarid,'units',2,"km")
91      ierr = NF_PUT_ATT_TEXT (nid,nvarid,'positive',2,"up")
92      ierr = NF_ENDDEF(nid)
93#ifdef NC_DOUBLE
94      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,pseudoalt)
95#else
96      ierr = NF_PUT_VAR_REAL (nid,nvarid,pseudoalt)
97#endif
[690]98      call def_var_stats(nid,"aps","hybrid pressure at midlayers"
99     & ," ",1,idim_llm,nvarid,ierr)
[38]100#ifdef NC_DOUBLE
101      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,aps)
102#else
103      ierr = NF_PUT_VAR_REAL (nid,nvarid,aps)
104#endif
105
[690]106      call def_var_stats(nid,"bps","hybrid sigma at midlayers"
107     & ," ",1,idim_llm,nvarid,ierr)
[38]108#ifdef NC_DOUBLE
109      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,bps)
110#else
111      ierr = NF_PUT_VAR_REAL (nid,nvarid,bps)
112#endif
113
114      ierr=NF_CLOSE(nid)
115
[1130]116      endif ! of if (is_master)
[38]117      end subroutine inistats
118
Note: See TracBrowser for help on using the repository browser.