source: trunk/LMDZ.COMMON/libf/evolution/get_timelen_PCM_mod.F90 @ 3578

Last change on this file since 3578 was 3571, checked in by jbclement, 10 days ago

PEM:

  • New way to manage the pressure: now the PEM manages only the average pressure and keeps the pressure deviation with the instantaneous pressure from the start to reconstruct the pressure at the end ('ps_avg = ps_start + ps_dev'). As a consequence, everything related to pressure in the PEM is modified accordingly.
  • Surface temperatures management is now simpler. It follows the strategy for the pressure (and soil temperature) described above.
  • Soil temperatures are now adapted to match the surface temperature changes occured during the PEM by modifying the soil temperature deviation at the end.
  • Few simplifications/optimizations: notably, the two PCM years are now read in one go in 'read_data_PCM_mod.F90' and only the needed variables are extracted.
  • Deletion of unused variables and unnecessary intermediate variables (memory saving and loop deletion in some cases).
  • Renaming of variables and subroutines to make everything clearer. In particular, the suffixes: '_avg' = average, '_start' = PCM start file, '_dev' = deviation, '_ini' or '0' = initial, '_dyn' = dynamical grid, '_timeseries' = daily average of last PCM year.
  • Cosmetic cleanings for readability.

JBC

File size: 3.3 KB
Line 
1MODULE get_timelen_PCM_mod
2
3use netcdf, only: nf90_open, NF90_NOWRITE, nf90_noerr, nf90_strerror, &
4                  nf90_get_var, nf90_inq_varid, nf90_inq_dimid,       &
5                  nf90_inquire_dimension, nf90_close
6
7implicit none
8
9character(256) :: msg, var, modname
10
11!=======================================================================
12contains
13!=======================================================================
14
15SUBROUTINE get_timelen_PCM(fichnom,timelen)
16
17implicit none
18
19!=======================================================================
20!
21! Purpose: Read in the data_PCM_Y*.nc the number of time steps
22!
23! Author: RV
24!=======================================================================
25
26include "dimensions.h"
27
28!=======================================================================
29! Arguments:
30character(len = *), intent(in) :: fichnom !--- FILE NAME
31!=======================================================================
32!   Local variables
33integer        :: fID, vID, ierr
34integer        :: timelen ! number of times stored in the file
35!-----------------------------------------------------------------------
36modname = "get_timelen_PCM"
37
38!  Open initial state NetCDF file
39var = fichnom
40call error_msg(NF90_OPEN(var,NF90_NOWRITE,fID),"open",var)
41
42ierr = nf90_inq_varid (fID, "temps", vID)
43if (ierr /= nf90_noerr) then
44    write(*,*)"read_data_PCM: Le champ <temps> est absent"
45    write(*,*)"read_data_PCM: J essaie <time_counter>"
46    ierr = nf90_inq_varid (fID, "time_counter", vID)
47    if (ierr /= nf90_noerr) then
48        write(*,*)"read_data_PCM: Le champ <time_counter> est absent"
49        write(*,*)"read_data_PCM: J essaie <Time>"
50        ierr = nf90_inq_varid (fID, "Time", vID)
51        if (ierr /= nf90_noerr) then
52            write(*,*)"read_data_PCM: Le champ <Time> est absent"
53            write(*,*)trim(nf90_strerror(ierr))
54            call abort_gcm("get_timelen_PCM", "", 1)
55        endif
56        ! Get the length of the "Time" dimension
57        ierr = nf90_inq_dimid(fID,"Time",vID)
58        ierr = nf90_inquire_dimension(fID,vID,len = timelen)
59    else
60        ! Get the length of the "time_counter" dimension
61        ierr = nf90_inq_dimid(fID,"time_counter",vID)
62        ierr = nf90_inquire_dimension(fID,vID,len = timelen)
63    endif
64else
65    ! Get the length of the "temps" dimension
66    ierr = nf90_inq_dimid(fID,"temps",vID)
67    ierr = nf90_inquire_dimension(fID,vID,len = timelen)
68endif
69
70call error_msg(NF90_CLOSE(fID),"close",fichnom)
71
72write(*,*) "The number of timestep of the PCM run data=", timelen
73
74END SUBROUTINE get_timelen_PCM
75
76!=======================================================================
77
78SUBROUTINE error_msg(ierr,typ,nam)
79
80implicit none
81
82integer,            intent(in) :: ierr !--- NetCDF ERROR CODE
83character(len = *), intent(in) :: typ  !--- TYPE OF OPERATION
84character(len = *), intent(in) :: nam  !--- FIELD/FILE NAME
85
86if (ierr == nf90_noerr) return
87select case(typ)
88    case('inq');   msg = "Field <"//trim(nam)//"> is missing"
89    case('get');   msg = "Reading failed for <"//trim(nam)//">"
90    case('open');  msg = "File opening failed for <"//trim(nam)//">"
91    case('close'); msg = "File closing failed for <"//trim(nam)//">"
92    case default
93        write(*,*) 'There is no message for this error.'
94        error stop
95end select
96call abort_gcm(trim(modname),trim(msg),ierr)
97
98END SUBROUTINE error_msg
99
100END MODULE get_timelen_PCM_mod
Note: See TracBrowser for help on using the repository browser.