[3571] | 1 | MODULE get_timelen_PCM_mod |
---|
[2794] | 2 | |
---|
[3076] | 3 | use 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 |
---|
[2794] | 6 | |
---|
[3076] | 7 | implicit none |
---|
[2794] | 8 | |
---|
[3076] | 9 | character(256) :: msg, var, modname |
---|
| 10 | |
---|
[2794] | 11 | !======================================================================= |
---|
[3076] | 12 | contains |
---|
| 13 | !======================================================================= |
---|
| 14 | |
---|
[3571] | 15 | SUBROUTINE get_timelen_PCM(fichnom,timelen) |
---|
[3076] | 16 | |
---|
| 17 | implicit none |
---|
| 18 | |
---|
| 19 | !======================================================================= |
---|
[2794] | 20 | ! |
---|
[3149] | 21 | ! Purpose: Read in the data_PCM_Y*.nc the number of time steps |
---|
[2794] | 22 | ! |
---|
[2855] | 23 | ! Author: RV |
---|
[2794] | 24 | !======================================================================= |
---|
| 25 | |
---|
[3076] | 26 | include "dimensions.h" |
---|
[2794] | 27 | |
---|
[3076] | 28 | !======================================================================= |
---|
[2794] | 29 | ! Arguments: |
---|
[3076] | 30 | character(len = *), intent(in) :: fichnom !--- FILE NAME |
---|
| 31 | !======================================================================= |
---|
[3571] | 32 | ! Local variables |
---|
[3206] | 33 | integer :: fID, vID, ierr |
---|
[3076] | 34 | integer :: timelen ! number of times stored in the file |
---|
[2794] | 35 | !----------------------------------------------------------------------- |
---|
[3571] | 36 | modname = "get_timelen_PCM" |
---|
[2794] | 37 | |
---|
| 38 | ! Open initial state NetCDF file |
---|
[3076] | 39 | var = fichnom |
---|
| 40 | call error_msg(NF90_OPEN(var,NF90_NOWRITE,fID),"open",var) |
---|
[2794] | 41 | |
---|
[3076] | 42 | ierr = nf90_inq_varid (fID, "temps", vID) |
---|
| 43 | if (ierr /= nf90_noerr) then |
---|
[3096] | 44 | write(*,*)"read_data_PCM: Le champ <temps> est absent" |
---|
| 45 | write(*,*)"read_data_PCM: J essaie <time_counter>" |
---|
[3076] | 46 | ierr = nf90_inq_varid (fID, "time_counter", vID) |
---|
| 47 | if (ierr /= nf90_noerr) then |
---|
[3096] | 48 | write(*,*)"read_data_PCM: Le champ <time_counter> est absent" |
---|
| 49 | write(*,*)"read_data_PCM: J essaie <Time>" |
---|
[3591] | 50 | ierr = nf90_inq_varid (fID,"Time",vID) |
---|
[3076] | 51 | if (ierr /= nf90_noerr) then |
---|
[3096] | 52 | write(*,*)"read_data_PCM: Le champ <Time> est absent" |
---|
[2841] | 53 | write(*,*)trim(nf90_strerror(ierr)) |
---|
[3571] | 54 | call abort_gcm("get_timelen_PCM", "", 1) |
---|
[3076] | 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 |
---|
[2841] | 60 | ! Get the length of the "time_counter" dimension |
---|
| 61 | ierr = nf90_inq_dimid(fID,"time_counter",vID) |
---|
[3076] | 62 | ierr = nf90_inquire_dimension(fID,vID,len = timelen) |
---|
| 63 | endif |
---|
| 64 | else |
---|
| 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) |
---|
| 68 | endif |
---|
[2794] | 69 | |
---|
[3076] | 70 | call error_msg(NF90_CLOSE(fID),"close",fichnom) |
---|
[2794] | 71 | |
---|
[3076] | 72 | write(*,*) "The number of timestep of the PCM run data=", timelen |
---|
[2835] | 73 | |
---|
[3571] | 74 | END SUBROUTINE get_timelen_PCM |
---|
[2794] | 75 | |
---|
[3076] | 76 | !======================================================================= |
---|
[2794] | 77 | |
---|
[3076] | 78 | SUBROUTINE error_msg(ierr,typ,nam) |
---|
| 79 | |
---|
| 80 | implicit none |
---|
| 81 | |
---|
| 82 | integer, intent(in) :: ierr !--- NetCDF ERROR CODE |
---|
| 83 | character(len = *), intent(in) :: typ !--- TYPE OF OPERATION |
---|
| 84 | character(len = *), intent(in) :: nam !--- FIELD/FILE NAME |
---|
| 85 | |
---|
| 86 | if (ierr == nf90_noerr) return |
---|
| 87 | select 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 |
---|
| 95 | end select |
---|
| 96 | call abort_gcm(trim(modname),trim(msg),ierr) |
---|
| 97 | |
---|
| 98 | END SUBROUTINE error_msg |
---|
| 99 | |
---|
[3571] | 100 | END MODULE get_timelen_PCM_mod |
---|