source: trunk/LMDZ.COMMON/libf/evolution/nb_time_step_PCM_mod.F90 @ 3170

Last change on this file since 3170 was 3149, checked in by jbclement, 15 months ago

PEM:

  • Simplification of the algorithm managing the stopping criteria;
  • Complete rework of the ice management in the PEM (H2O & CO2);

    Subroutines to evolve the H2O and CO2 ice are now in the same module "evol_ice_mod.F90".
    Tendencies are computed from the variation of "ice + frost" between the 2 PCM runs.
    Evolving ice in the PEM is now called 'h2o_ice' or 'co2_ice' (not anymore in 'qsurf' and free of 'water_reservoir').
    Default value 'ini_h2o_bigreservoir' (= 10 m) initializes the H2O ice of the first PEM run where there is 'watercap'. For the next PEM runs, initialization is done with the value kept in "startpem.nc". CO2 ice is taken from 'perennial_co2ice' of the PCM (paleoclimate flag must be true).
    Simplification of the condition to compute the surface ice cover needed for the stopping criteria.
    Frost ('qsurf') is not evolved by the PEM and given back to the PCM.
    New default threshold value 'inf_h2oice_threshold' (= 2 m) to decide at the end of the PEM run if the H2O ice should be 'watercap' or not for the next PCM runs. If H2O ice cannot be 'watercap', then the remaining H2O ice is transferred to the frost ('qsurf').

  • Renaming of variables/subroutines for clarity;
  • Some cleanings throughout the code;
  • Small updates in files of the deftank.

JBC

File size: 3.4 KB
Line 
1MODULE nb_time_step_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 nb_time_step_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        :: iq, fID, vID, idecal, ierr
34integer        :: timelen ! number of times stored in the file
35!-----------------------------------------------------------------------
36modname = "nb_time_step_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("nb_time_step_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 nb_time_step_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 nb_time_step_PCM_mod
Note: See TracBrowser for help on using the repository browser.