! module writehist_module ! ! Write a TIEGCM/TIMEGCM history. ! use readhist_module,only: history,handle_ncerr,field_3d,field_2d, | fminmaxspv,int_unread,input_file,real_unread implicit none #include "netcdf.inc" integer,private :: | id_time, ! dimension id for time (unlimited) | id_mtimedim, ! dimension id for model time (3) | id_lon,id_lat, ! horizontal dimension ids | id_lev,id_ilev, ! vertical dimension ids (midpoints,interfaces) | id_mlon,id_mlat, ! horizontal dimension ids for magnetic grid | id_mlev,id_imlev, ! vertical dimension ids for magnetic grid | id_latlon, ! dimension id for coords (2) | id_dtide, ! dimension id for dtide (2) | id_sdtide, ! dimension id for sdtide (10) | id_fake, ! fake dimension id (1) | id_datelen, ! length of date and time | id_filelen, ! max length of a file name | ids1(1),ids2(2),ids3(3),ids4(4),! vectors of dim id's | start_2d(2),count_2d(2), ! start loc and count for 2d vars | start_3d(3),count_3d(3), ! start loc and count for 3d vars | start_4d(4),count_4d(4), ! start loc and count for 4d vars ! ! Variable id's (idv_xxxx). Note lon, lat, lev, and time are ! coordinate variables (lon(lon), lat(lat), lev(lev), time(time)). ! Variables starting at mtime are "primary" variables. ! idv_f4d(nf4d) are id's for prognostic 4d fields. ! | idv_time,idv_mtime,idv_ut,idv_calendar_adv, | idv_lon,idv_lat,idv_lev,idv_ilev, | idv_lbc,idv_tlbc,idv_ulbc,idv_vlbc, | idv_tlbc_nm,idv_ulbc_nm,idv_vlbc_nm, | idv_mlon,idv_mlat,idv_mlev,idv_imlev,idv_p0_model, | idv_step,idv_iter,idv_year,idv_day,idv_p0,idv_grav, | idv_zg,idv_hpower,idv_ctpoten,idv_bximf,idv_byimf,idv_bzimf, | idv_colfac,idv_swden,idv_swvel,idv_al, | idv_gswm_mi_di_ncfile,idv_gswm_mi_sdi_ncfile, | idv_gswm_nm_di_ncfile,idv_gswm_nm_sdi_ncfile, | idv_see_ncfile,idv_gpi_ncfile,idv_ncep_ncfile, | idv_f107d,idv_f107a,idv_imf_ncfile,idv_mag,idv_dtide, | idv_sdtide,idv_e1,idv_e2,idv_h1,idv_h2,idv_alfac, | idv_ec,idv_alfad,idv_ed,idv_writedate integer,allocatable :: idv_sech(:) integer,parameter :: mxflds = 200 integer :: idv_f4d(mxflds),idv_f3d(mxflds) ! ! fieldlist: Single line of comma-separated field names written ! to the new history (useful for post-proc input) ! character(len=2048) :: fieldlist=' ' contains !----------------------------------------------------------------------- subroutine writehist(outfile,mtime,h,iprint,error) implicit none ! ! Args: ! History file to create: character(len=*),intent(in) :: outfile ! ! h: history structure to return (read from input file into h): type(history),intent(inout) :: h ! ! mtime(3): model time (day,hr,min) of history to be read: integer,intent(inout) :: mtime(3) integer,intent(in) :: iprint ! 0/1 flag for print to stdout integer,intent(out) :: error ! error return ! ! Local: integer :: istat,ncid,mtimes(3,1),ntimes,i istat = nf_create(outfile,NF_CLOBBER,ncid) write(6,"(/,'Created file ',a)") trim(outfile) mtimes(:,1) = mtime ntimes = 1 ! ! 10/12/08: Currently, this module writes only a single history ! to a single file, i.e., it will not write multiple histories ! to the same file. This capability should be added later. ! h%ihist = 1 ! ! Define history h on new nc file: h%history_file = trim(outfile) call nc_define(ncid,h,mtimes,ntimes) ! ! Write data to the file and close: call nc_wrhist(ncid,h,mtimes,ntimes) istat = nf_close(ncid) ! write(6,"(/,'4d fields (lon,lat,lev,time) written to file ',a, | ':',/)") trim(outfile) write(6,"(a)") trim(fieldlist) end subroutine writehist !----------------------------------------------------------------------- subroutine nc_define(ncid,h,mtimes,ntimes) ! ! Define dimensions, variables and attributes in a new history ! netcdf dataset (file has been opened, but no histories written yet). ! On input: ! ncid = file id for open netcdf dataset (in define mode). ! h = defined history structure (in hist_module) ! On output: ! Variables and dimensions for the netcdf file are defined. ! Values are given to the dimension variables, and the file ! is taken out of define mode (ready to receive histories) ! ! use hist_module,only: h,sh,nsource,nsecsource ! use cons_module,only: pi ! use input_module,only: potential_model,ncep_ncfile ! use params_module,only: tgcm_name,tgcm_version ! ! Args: integer,intent(in) :: ncid,ntimes integer,intent(in) :: mtimes(3,ntimes) type(history),intent(in) :: h ! ! Local: integer :: i,ii,istat,idum,ivar1(1),imo,ida,startmtime(4),nf4d character(len=120) :: char120 character(len=80) :: char80 character(len=24) :: create_date,create_time character(len=120) :: createdate character(len=1024) :: history_attrib real :: var1(1),rdays,rmins ! real,external :: mtime_to_datestr integer,parameter :: mxlen_filename=80 ! ! Define dimensions (field vars will be (time,lev,lat,lon) in ! CDL notation. Mag fields will be (time,[mlev or imlev],mlat,mlon): ! istat = nf_def_dim(ncid,"time",NF_UNLIMITED,id_time) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining time dimension') ! ! Geographic grid dimensions: istat = nf_def_dim(ncid,"lon",h%nlon,id_lon) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining longitude dimension') istat = nf_def_dim(ncid,"lat",h%nlat,id_lat) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining latitude dimension') ! ! Midpoint levels dimension "lev"=nlevp1: istat = nf_def_dim(ncid,"lev",h%nlev,id_lev) ! midpoint levels dimension if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining midpoint levels dimension') ! ! Interface levels dimension "ilev"=nlevp1: if (h%new) then istat = nf_def_dim(ncid,"ilev",h%nilev,id_ilev) ! interface levels dimension if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining interface levels dimension') else ! write(6,"('Did not define nilev on old history')") endif ! ! Magnetic grid dimensions: istat = nf_def_dim(ncid,"mlon",h%mlon,id_mlon) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining magnetic longitude dimension') istat = nf_def_dim(ncid,"mlat",h%mlat,id_mlat) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining magnetic latitude dimension') ! ! Mag midpoint levels dimension: istat = nf_def_dim(ncid,"mlev",h%mlev,id_mlev) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining magnetic midpoints levels dimension') ! ! Mag interface levels dimension: if (h%new) then istat = nf_def_dim(ncid,"imlev",h%imlev,id_imlev) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining magnetic interface levels dimension') endif ! istat = nf_def_dim(ncid,"mtimedim",3,id_mtimedim) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining mtime dimension') istat = nf_def_dim(ncid,"latlon",2,id_latlon) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining latlon dimension') istat = nf_def_dim(ncid,"dtidedim",2,id_dtide) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining dtide dimension') istat = nf_def_dim(ncid,"sdtidedim",10,id_sdtide) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining sdtide dimension') ! ! String length of date and time: if (h%old) then istat = nf_def_dim(ncid,"datelen",24,id_datelen) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining datelen=24 dimension') else istat = nf_def_dim(ncid,"datelen",h%datelen,id_datelen) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining datelen dimension') endif ! ! String length of file names (mxlen_filename is parameter in input.F): istat = nf_def_dim(ncid,"filelen",mxlen_filename,id_filelen) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining filelen dimension') ! ! Define dimension (coordinate) variables and attributes: ! ! Time (coordinate variable time(time)). This is days since ! the initial run's start time. The units string is: yyyy-m-d, ! where yyyy is the year, m is month, and d is day of the source ! start time. ! ids1(1) = id_time ! for 1d time-unlimited vars istat = nf_def_var(ncid,"time",NF_DOUBLE,1,ids1,idv_time) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining time dimension variable') ! rmins = mtime_to_datestr(h%initial_year,startmtime,imo,ida, ! | char80) ! real function mtime_to_datestr(iyear,mtime,imo,ida,datestr) rmins=mtime_to_datestr(h%year,mtimes(:,1),imo,ida,char80) istat = nf_put_att_text(ncid,idv_time,"long_name",4,"time") if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining long_name of time dimension variable') istat = nf_put_att_text(ncid,idv_time,"units",len_trim(char80), | trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining units of time dimension variable') istat = nf_put_att_int(ncid,idv_time,"initial_year",NF_INT,1, | h%year) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining initial_year attribute of time coord var') istat = nf_put_att_int(ncid,idv_time,"initial_mtime",NF_INT,3, | h%initial_mtime) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining initial_mtime attribute of time coord var') ! ! Geographic longitude (deg) (coordinate variable lon(lon)): ids1(1) = id_lon istat = nf_def_var(ncid,"lon",NF_DOUBLE,1,ids1,idv_lon) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining longitude dimension variable') write(char80,"('geographic longitude (-west, +east)')") istat = nf_put_att_text(ncid,idv_lon,"long_name", + len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining long_name of longitude dimension variable') istat = nf_put_att_text(ncid,idv_lon,"units",12,'degrees_east') if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining units of longitude dimension variable') ! ! Geographic latitude (deg) (coordinate variable lat(lat)): ids1(1) = id_lat istat = nf_def_var(ncid,"lat",NF_DOUBLE,1,ids1,idv_lat) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining latitude dimension variable') write(char80,"('geographic latitude (-south, +north)')") istat = nf_put_att_text(ncid,idv_lat,"long_name", + len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining long_name of latitude dimension variable') istat = nf_put_att_text(ncid,idv_lat,"units",13,'degrees_north') if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining units of latitude dimension variable') ! ! Midpoint levels coordinate variable lev(lev): ids1(1) = id_lev istat = nf_def_var(ncid,"lev",NF_DOUBLE,1,ids1,idv_lev) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining midpoint levels coordinate variable lev(lev)') ! long name of lev coord var: write(char80,"('midpoint levels')") istat = nf_put_att_text(ncid,idv_lev,"long_name", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, |'Error defining long_name of midpoint levels coordinate variable') ! short name = ln(p0/p) istat = nf_put_att_text(ncid,idv_lev,"short_name", | 8,"ln(p0/p)") if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining short_name of midpoint levels coord var') ! lev coord var is unitless: istat = nf_put_att_text(ncid,idv_lev,"units",0," ") if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining units of midpoint levels coord var') ! positive='up' istat = nf_put_att_text(ncid,idv_lev,"positive",2,'up') if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining positive attribute of midpoint levels coord var') ! standard name of midpoints lev coord var: write(char80,"('atmosphere_ln_pressure_coordinate')") istat = nf_put_att_text(ncid,idv_lev,"standard_name", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining standard_name of midpoint levels coord var') ! formula terms for midpoints lev coord var: write(char80,"('p0: p0 lev: lev')") istat = nf_put_att_text(ncid,idv_lev,"formula_terms", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining formula terms of midpoint levels coord var') ! formula to obtain pressure from midpoint lev: write(char80,"('p(k) = p0 * exp(-lev(k))')") istat = nf_put_att_text(ncid,idv_lev,"formula", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining formula for midpoint levels coord var') ! ! Interface levels coordinate array: if (h%new) then ids1(1) = id_ilev istat = nf_def_var(ncid,"ilev",NF_DOUBLE,1,ids1,idv_ilev) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining interface levels coord variable ilev(ilev)') ! long name of ilev coord var: write(char80,"('interface levels')") istat = nf_put_att_text(ncid,idv_ilev,"long_name", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining long_name of interface levels coordinate var') ! short name = ln(p0/p) istat = nf_put_att_text(ncid,idv_ilev,"short_name", | 8,"ln(p0/p)") if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining short_name of interface levels coord var') ! ilev coord var is unitless: istat = nf_put_att_text(ncid,idv_ilev,"units",0," ") if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining units of interface levels coord var') ! positive='up' istat = nf_put_att_text(ncid,idv_ilev,"positive",2,'up') if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining positive att of interface levels coord var') ! standard name of interface ilev coord var: write(char80,"('atmosphere_ln_pressure_coordinate')") istat = nf_put_att_text(ncid,idv_ilev,"standard_name", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining standard_name of interface levels coord var') ! formula terms for interface ilev coord var: write(char80,"('p0: p0 lev: ilev')") istat = nf_put_att_text(ncid,idv_ilev,"formula_terms", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining formula terms of interface levels coord var') ! formula to obtain pressure from interface ilev: write(char80,"('p(k) = p0 * exp(-ilev(k))')") istat = nf_put_att_text(ncid,idv_ilev,"formula", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining formula for interface levels coord var') endif ! h%new ! ! Magnetic longitude (deg) (coordinate variable mlon(mlon)): ids1(1) = id_mlon istat = nf_def_var(ncid,"mlon",NF_DOUBLE,1,ids1,idv_mlon) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining magnetic longitude dimension variable') write(char80,"('magnetic longitude (-west, +east)')") istat = nf_put_att_text(ncid,idv_mlon,"long_name", + len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining long_name of mag longitude dimension variable') istat = nf_put_att_text(ncid,idv_mlon,"units",12,'degrees_east') if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining units of magnetic longitude dimension variable') ! ! Magnetic latitude (deg) (coordinate variable mlat(mlat)): ids1(1) = id_mlat istat = nf_def_var(ncid,"mlat",NF_DOUBLE,1,ids1,idv_mlat) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining magnetic latitude dimension variable') write(char80,"('magnetic latitude (-south, +north)')") istat = nf_put_att_text(ncid,idv_mlat,"long_name", + len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining long_name of mag latitude dimension variable') istat = nf_put_att_text(ncid,idv_mlat,"units",13,'degrees_north') if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining units of magnetic latitude dimension variable') ! ! Magnetic midpoint coordinate variable mlev(mlev)): ids1(1) = id_mlev istat = nf_def_var(ncid,"mlev",NF_DOUBLE,1,ids1,idv_mlev) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining magnetic levels coord variable mlev') ! long name of mlev: write(char80,"('magnetic midpoint levels')") istat = nf_put_att_text(ncid,idv_mlev,"long_name", + len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining long_name of mlev coordinate variable') ! short name = ln(p0/p) istat = nf_put_att_text(ncid,idv_mlev,"short_name", | 8,"ln(p0/p)") if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining short_name of mlev levels coord var') ! units of mlev (unitless): istat = nf_put_att_text(ncid,idv_mlev,"units",0," ") if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining units of mlev coordinate variable') ! positive attribute of mlev: istat = nf_put_att_text(ncid,idv_mlev,"positive",2,'up') if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining positive attribute for mlev coord var') ! standard_name of mlev: write(char80,"('atmosphere_ln_pressure_coordinate')") istat = nf_put_att_text(ncid,idv_mlev,"standard_name", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining standard_name attribute for mlev coord var') ! formula terms for mlev: write(char80,"('p0: p0 lev: mlev')") istat = nf_put_att_text(ncid,idv_mlev,"formula_terms", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining formula terms of mlev coord var') ! formula to obtain pressure from mlev: write(char80,"('p(k) = p0 * exp(-mlev(k))')") istat = nf_put_att_text(ncid,idv_mlev,"formula", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining formula for mlev coord var') ! ! Magnetic interface coordinate variable imlev(imlev)): if (h%new) then ids1(1) = id_imlev istat = nf_def_var(ncid,"imlev",NF_DOUBLE,1,ids1,idv_imlev) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining magnetic levels coord variable imlev') ! long name of imlev: write(char80,"('magnetic interface levels')") istat = nf_put_att_text(ncid,idv_imlev,"long_name", + len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining long_name of imlev coordinate variable') ! short name = ln(p0/p) istat = nf_put_att_text(ncid,idv_imlev,"short_name", | 8,"ln(p0/p)") if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining short_name of imlev levels coord var') ! units of imlev (unitless): istat = nf_put_att_text(ncid,idv_imlev,"units",0," ") if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining units of imlev coordinate variable') ! positive attribute of imlev: istat = nf_put_att_text(ncid,idv_imlev,"positive",2,'up') if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining positive attribute for imlev coord var') ! standard_name of imlev: write(char80,"('atmosphere_ln_pressure_coordinate')") istat = nf_put_att_text(ncid,idv_imlev,"standard_name", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining standard_name attribute for imlev coord var') ! formula terms for imlev: write(char80,"('p0: p0 lev: imlev')") istat = nf_put_att_text(ncid,idv_imlev,"formula_terms", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining formula terms of imlev coord var') ! formula to obtain pressure from imlev: write(char80,"('p(k) = p0 * exp(-imlev(k))')") istat = nf_put_att_text(ncid,idv_imlev,"formula", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining formula for imlev coord var') endif ! h%new ! ! ids1(1) is id of time dimension for the following defines: ids1(1) = id_time ! for 1d time-unlimited vars ! ! Source date and time are attributes of time var ! (they are also global file attributes): ! istat = nf_put_att_int(ncid,idv_time,"initial_year",NF_INT,1, ! | sh%initial_year) ! istat = nf_put_att_int(ncid,idv_time,"initial_mtime",NF_INT,3, ! | sh%initial_mtime) ! istat = nf_put_att_int(ncid,idv_time,"initial_year",NF_INT,1, ! | h%year) ! istat = nf_put_att_int(ncid,idv_time,"initial_mtime",NF_INT,3, ! | h%modeltime(1:3)) ! ! Define "primary" (non-coordinate) time dependent variables: ! (these variables get a value for each history) ! ! Model time variable mtime(3,time), where 3 is the mtime dimension ! (id_mtimedim) (day,hour,minute), and time is the unlimited dimension ! (id_time): ! ids2(1) = id_mtimedim ids2(2) = id_time istat = nf_def_var(ncid,"mtime",NF_INT,2,ids2,idv_mtime) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining mtime variable') write(char80,"('model times (day, hour, minute)')") istat = nf_put_att_text(ncid,idv_mtime,"long_name", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining long_name attribute of mtime variable') istat = nf_put_att_text(ncid,idv_mtime,"units",17, | 'day, hour, minute') if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining units of mtime variable') ! ! Calendar year and day (this will be advanced if model is advancing ! in calendar time, otherwise will be constant from input) istat = nf_def_var(ncid,"year",NF_INT,1,ids1,idv_year) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining calendar year variable') istat = nf_put_att_text(ncid,idv_year,"long_name", | 13,'calendar year') if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining long_name of calendar year variable') ! istat = nf_def_var(ncid,"day",NF_INT,1,ids1,idv_day) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining calendar day variable') istat = nf_put_att_text(ncid,idv_day,"long_name", | 12,'calendar day') if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining long_name of calendar day variable') ! ! Calendar advance: istat = nf_def_var(ncid,"calendar_advance",NF_INT,1,ids1, | idv_calendar_adv) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining integer calendar advance variable') write(char80,"('calendar advance flag', | ' (1 if advancing calendar time)')") istat = nf_put_att_text(ncid,idv_calendar_adv,"long_name", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining long_name of calendar advance variable') ! ! Date each history was written: ids2(1) = id_datelen ids2(2) = id_time istat= nf_def_var(ncid,"write_date",NF_CHAR,2,ids2,idv_writedate) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining history write_date variable') write(char80,"('Date and time each history was written')") istat = nf_put_att_text(ncid,idv_writedate,"long_name", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining long_name of write_date') ! ! Number of time steps from model time 0,0,0 to current model time ! (iter(time)): istat = nf_def_var(ncid,"iter",NF_INT,1,ids1,idv_iter) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining iter variable') istat = nf_put_att_text(ncid,idv_iter,"long_name", | 42,'number of time steps from model time 0,0,0') if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining long_name of iter variable') ! ! Decimal ut(time) (ids1==id_time is unlimited time dimension): call defvar_time_dbl(ncid,"ut",ids1,idv_ut, | "universal time (from model time hour and minute)", | "hours") ! ! Time step (seconds) (integer scalar) (time-independent): istat = nf_def_var(ncid,"timestep",NF_INT,0,idum,idv_step) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining timestep variable') istat = nf_put_att_text(ncid,idv_step,"long_name",8,'timestep') if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining long_name of timestep variable') istat = nf_put_att_text(ncid,idv_step,"units",7,'seconds') if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining units of timestep variable') ! ! f107 daily and f107 average solar fluxes: ! 9/30/05: change units as suggested by John Caron (s.a. Udunits) call defvar_time_dbl(ncid,"f107d",ids1,idv_f107d, ! | "Daily 10.7 cm solar flux","W/m^2Hz*1.e-22") | "Daily 10.7 cm solar flux","1.e-22 kg.s-4") call defvar_time_dbl(ncid,"f107a",ids1,idv_f107a, ! | "81-day average 10.7 cm solar flux","W/m^2Hz*1.e-22") | "81-day average 10.7 cm solar flux","1.e-22 kg.s-4") ! ! Hemispheric power: call defvar_time_dbl(ncid,"hpower",ids1,idv_hpower, | "Hemispheric power","gw") ! ! Cross-tail potential: call defvar_time_dbl(ncid,"ctpoten",ids1,idv_ctpoten, | "Cross-tail potential","volts") ! ! BX,BY,BZ: call defvar_time_dbl(ncid,"bximf",ids1,idv_bximf, | "BX component of IMF","nT") call defvar_time_dbl(ncid,"byimf",ids1,idv_byimf, | "BY component of IMF","nT") call defvar_time_dbl(ncid,"bzimf",ids1,idv_bzimf, | "BZ component of IMF","nT") ! ! swvel,swden: call defvar_time_dbl(ncid,"swvel",ids1,idv_swvel, | "Solar wind velocity","km/s") call defvar_time_dbl(ncid,"swden",ids1,idv_swden, | "Solar wind density","cm3") ! ! al: call defvar_time_dbl(ncid,"al",ids1,idv_al, | "Lower magnetic auroral activity index","nT") ! ! Collision factor: call defvar_time_dbl(ncid,"colfac",ids1,idv_colfac, | "ion/neutral collision factor"," ") ! ! 1/31/01: Making tides time-dependent (dtide(time,2) and sdtide(time,10): ! Diurnal tide amp/phase: ids2(1) = id_dtide ids2(2) = id_time istat = nf_def_var(ncid,"dtide",NF_DOUBLE,2,ids2,idv_dtide) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining variable dtide') write(char80,"('amplitude and phase of diurnal tide mode (1,1)')") istat = nf_put_att_text(ncid,idv_dtide,"long_name", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining long_name attribute of variable dtide') ! ! Semi-diurnal tide amp/phase: ids2(1) = id_sdtide ids2(2) = id_time istat = nf_def_var(ncid,"sdtide",NF_DOUBLE,2,ids2,idv_sdtide) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining variable sdtide') write(char80,"('amplitudes and phases of semi-diurnal tide')") istat = nf_put_att_text(ncid,idv_sdtide,"long_name", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining long_name attribute of variable sdtide') ! ! Define ncep data file (used by timegcm only): call define_ncfile(ncid,"ncep_ncfile",idv_ncep_ncfile, | "long_name","ncep data file", | "data","10 mb lower boundary for TN and Z") ! ! Define gpi data file: call define_ncfile(ncid,"gpi_ncfile",idv_gpi_ncfile, | "long_name","GeoPhysical Indices data file", | "data","f107, f107a, ctpoten, power") ! ! Define TIMED SEE data file: call define_ncfile(ncid,"see_ncfile",idv_see_ncfile, | "long_name","TIMED SEE data file", | "data","Solar flux data") ! ! Define imf data file: call define_ncfile(ncid,"imf_ncfile",idv_imf_ncfile, | "long_name","IMF data file", | "data","Bx,By,Bz,Swvel,Swden,Kp,f107,f107a") ! ! Define gswm data files: call define_ncfile(ncid,"gswm_mi_di_ncfile", | idv_gswm_mi_di_ncfile, | "long_name","GSWM migrating diurnal tides data file", | "data","lower boundary perturbations for TN, UN, VN, and Z") call define_ncfile(ncid,"gswm_mi_sdi_ncfile", | idv_gswm_mi_sdi_ncfile, | "long_name","GSWM migrating semi-diurnal tides data file", | "data","lower boundary perturbations for TN, UN, VN, and Z") call define_ncfile(ncid,"gswm_nm_di_ncfile", | idv_gswm_nm_di_ncfile, | "long_name","GSWM non-migrating diurnal tides data file", | "data","lower boundary perturbations for TN, UN, VN, and Z") call define_ncfile(ncid,"gswm_nm_sdi_ncfile", | idv_gswm_nm_sdi_ncfile, | "long_name","GSWM non-migrating semi-diurnal tides data file", | "data","lower boundary perturbations for TN, UN, VN, and Z") ! ! Auroral parameters (alfa30, e30, alfad2, ed2 are timegcm only) ! ! alfa30: ! call defvar_time_dbl(ncid,"alfa30",ids1,idv_alfa30, ! | "Characteristic energy of high-energy auroral electrons", ! | "KeV") ! ! e30: ! call defvar_time_dbl(ncid,"e30",ids1,idv_e30, ! | "Energy flux of high-energy auroral electrons", ! | "ergs/cm2/s") ! ! alfad2: ! call defvar_time_dbl(ncid,"alfad2",ids1,idv_alfad2, ! | "Characteristic energy of solar protons in the polar cap", ! | "KeV") ! ! ed2: ! call defvar_time_dbl(ncid,"ed2",ids1,idv_ed2, ! | "Energy flux of solar protons in the polar cap", ! | "ergs/cm2/s") ! ! 1/24/08 btf: Auroral parameters: ! ! e1: call defvar_time_dbl(ncid,"e1",ids1,idv_e1, | "Peak energy flux in noon sector of the aurora", | "ergs/cm2/s") ! e2: call defvar_time_dbl(ncid,"e2",ids1,idv_e2, | "Peak energy flux in midnight sector of the aurora", | "ergs/cm2/s") ! h1: call defvar_time_dbl(ncid,"h1",ids1,idv_h1, | "Gaussian half-width of the noon auroral oval", | "degrees") ! h2: call defvar_time_dbl(ncid,"h2",ids1,idv_h2, | "Gaussian half-width of the midnight auroral oval", | "degrees") ! alfac: call defvar_time_dbl(ncid,"alfac",ids1,idv_alfac, | "Characteristic Maxwellian energy of polar cusp electrons", | "keV") ! ec: call defvar_time_dbl(ncid,"ec",ids1,idv_ec, | "Column energy input of polar cusp electrons", | "ergs/cm**2/s") ! alfad: call defvar_time_dbl(ncid,"alfad",ids1,idv_alfad, | "Characteristic Maxwellian energy of drizzle electrons", | "keV") ! ed: call defvar_time_dbl(ncid,"ed",ids1,idv_ed, | "Column energy input of drizzle electrons", | "ergs/cm**2/s") ! ! "New" histories do not write mag poles (apex code is always called) ! Magnetic pole coordinates (real(2,2)): ! ids2(1) = id_latlon ! ids2(2) = id_latlon ! istat = nf_def_var(ncid,"mag",NF_DOUBLE,2,ids2,idv_mag) ! if (istat /= NF_NOERR) call handle_ncerr(istat, ! + 'Error defining variable mag') ! write(char80,"('lat,lon coordinates of S,N magnetic poles')") ! istat = nf_put_att_text(ncid,idv_mag,"long_name", ! | len_trim(char80),trim(char80)) ! if (istat /= NF_NOERR) call handle_ncerr(istat, ! + 'Error defining long_name attribute of variable mag') ! istat = nf_put_att_text(ncid,idv_mag,"units",7, ! | 'degrees') ! if (istat /= NF_NOERR) call handle_ncerr(istat, ! + 'Error defining units of variable mag') ! ! Reference pressure used to convert ln(p0/p) to pressure (millibars (hPa)): istat = nf_def_var(ncid,"p0",NF_DOUBLE,0,idum,idv_p0) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining variable p0') istat = nf_put_att_text(ncid,idv_p0,"long_name", | 18,"Reference pressure") ! for conversion of vertical coordinates ! ! 4/26/06 btf: p0 units can be either hPa or millibars. Using the latter ! now to be consistent w/ microbar units of p0_model istat = nf_put_att_text(ncid,idv_p0,"units",9,"millibars") ! istat = nf_put_att_text(ncid,idv_p0,"units",3,"hPA") ! ! Reference pressure used by the model (dynes/cm2): ! p0_model units can be either dynes/cm2 or microbars. Using the latter ! to be consistent w/ 1981 paper by Dickinson, Ridley, and Roble ! istat = nf_def_var(ncid,"p0_model",NF_DOUBLE,0,idum,idv_p0_model) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining variable p0') char80 = ' ' write(char80,"('Reference pressure (as used by the model)')") istat = nf_put_att_text(ncid,idv_p0_model,"long_name", | len_trim(char80),trim(char80)) istat = nf_put_att_text(ncid,idv_p0_model,"units",9,"microbars") ! istat = nf_put_att_text(ncid,idv_p0_model,"units",9,"dynes/cm2") ! ! Gravity (constant used in the model, independent of height): istat = nf_def_var(ncid,"grav",NF_DOUBLE,0,idum,idv_grav) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error defining variable grav') istat = nf_put_att_text(ncid,idv_grav,"long_name", | 26,"gravitational acceleration") istat = nf_put_att_text(ncid,idv_grav,"units", | 4,"cm/s") write(char80,"('constant used in the model, independent of ', | 'height')") istat = nf_put_att_text(ncid,idv_grav,"info", | len_trim(char80),trim(char80)) ! ! Define 4d fields (3d plus time): if (h%nf3d_geo > 0) then do i=1,h%nf3d_geo call define_f4d(ncid,idv_f4d(i),h%f3d_geo(i)) enddo endif if (h%nf3d_mag > 0) then do i=1,h%nf3d_mag call define_f4d(ncid,idv_f4d(i),h%f3d_mag(i)) enddo endif ! ! Define 3d fields (2d plus time): if (h%nf2d_geo > 0) then do i=1,h%nf2d_geo call define_f3d(ncid,idv_f3d(i),h%f2d_geo(i)) enddo endif if (h%nf2d_mag > 0) then do i=1,h%nf2d_mag call define_f3d(ncid,idv_f3d(i),h%f2d_mag(i)) enddo endif ! ! Define LBC: bottom interface level (scalar coord of TLBC,ULBC,VLBC): if (h%new) then istat = nf_def_var(ncid,"LBC",NF_DOUBLE,0,idum,idv_lbc) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining var LBC') write(char80, | "('Interface level of t,u,v lower boundary condition')") istat = nf_put_att_text(ncid,idv_lbc,"long_name", | len_trim(char80),trim(char80)) if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error defining long_name attribute of variable LBC') endif ! ! Global file attributes (all history structures): ! istat = nf_put_att_text(ncid,NF_GLOBAL,"Conventions", | len_trim(h%conventions),h%conventions) istat = nf_put_att_text(ncid,NF_GLOBAL,"label", | len_trim(h%label),trim(h%label)) ! ! Echo create_date of original file (see history_attrib below ! for current creation date) istat=nf_put_att_text(ncid,NF_GLOBAL,"create_date", | len_trim(h%create_date),h%create_date) istat = nf_put_att_text(ncid,NF_GLOBAL,"logname", | len_trim(h%logname),trim(h%logname)) istat = nf_put_att_text(ncid,NF_GLOBAL,"host", | len_trim(h%host),trim(h%host)) istat = nf_put_att_text(ncid,NF_GLOBAL,"system", | len_trim(h%system),trim(h%system)) istat = nf_put_att_text(ncid,NF_GLOBAL,"model_name", | len_trim(h%model_name),trim(h%model_name)) istat = nf_put_att_text(ncid,NF_GLOBAL,"model_version", | len_trim(h%model_version),trim(h%model_version)) istat = nf_put_att_text(ncid,NF_GLOBAL,"output_file", | len_trim(h%output_file),trim(h%output_file)) istat = nf_put_att_text(ncid,NF_GLOBAL,"history_type", | len_trim(h%history_type),trim(h%history_type)) istat = nf_put_att_text(ncid,NF_GLOBAL,"run_type", | len_trim(h%run_type),trim(h%run_type)) istat = nf_put_att_text(ncid,NF_GLOBAL,"source_file", | len_trim(h%source_file),trim(h%source_file)) istat = nf_put_att_int (ncid,NF_GLOBAL,"source_mtime", | NF_INT,3,h%source_mtime) istat = nf_put_att_text(ncid,NF_GLOBAL,"initial_file", | len_trim(h%initial_file),trim(h%initial_file)) istat = nf_put_att_int (ncid,NF_GLOBAL,"initial_mtime", | NF_INT,3,h%initial_mtime) istat = nf_put_att_text(ncid,NF_GLOBAL,"lev_to_hPa_method1", | len_trim(h%lev_to_hPa_method1),trim(h%lev_to_hPa_method1)) istat = nf_put_att_text(ncid,NF_GLOBAL,"lev_to_hPa_method2", | len_trim(h%lev_to_hPa_method2),trim(h%lev_to_hPa_method2)) istat = nf_put_att_int (ncid,NF_GLOBAL,"nhist",NF_INT,1,h%nhist) istat = nf_put_att_int (ncid,NF_GLOBAL,"delhist_mins",NF_INT, | 1,h%delhist_mins) var1(1) = h%missing_value istat = nf_put_att_double(ncid,NF_GLOBAL,"missing_value", | NF_DOUBLE,1,var1) istat = nf_put_att_text(ncid,NF_GLOBAL,"potential_model", | len_trim(h%potential_model),trim(h%potential_model)) istat = nf_put_att_int (ncid,NF_GLOBAL,"tuv_lbc_intop",NF_INT,1, | h%tuv_lbc_intop) istat = nf_put_att_text(ncid,NF_GLOBAL,"contents", | len_trim(h%contents),trim(h%contents)) istat = nf_put_att_text(ncid,NF_GLOBAL,"contents_desc", | len_trim(h%contents_desc),trim(h%contents_desc)) ! ! Add history attribute to document that this file was made by interpolation ! of an original file. ! character(len=1024) :: history_attrib ! call datetime(create_date,create_time) ! create_date,create_time are output createdate = trim(create_date)//' '//trim(create_time) write(history_attrib,"('This file created on ',a, | ' by interpolation of original file ',a)") trim(createdate), | trim(input_file) istat = nf_put_att_text(ncid,NF_GLOBAL,"history", | len_trim(history_attrib),trim(history_attrib)) ! ! Exit define mode: istat = nf_enddef(ncid) if (istat /= NF_NOERR) call handle_ncerr(istat, + 'Error return from nf_enddef') ! ! Give values to coordinate variables: ! From readhist: ! real,pointer :: glat(:),glon(:) ! geographic lat lon grid ! real,pointer :: gmlat(:),gmlon(:) ! geomagnetic lat lon grid ! real,pointer :: lev(:),ilev(:) ! midpoint and interface geo levels ! real,pointer :: gmlev(:),gimlev(:) ! midpoint and interface mag levels ! ! Geographic horizontal: istat = nf_put_var_double(ncid,idv_lon,h%glon) if (istat /= NF_NOERR) | call handle_ncerr(istat,'Error giving values to glon coord var') istat = nf_put_var_double(ncid,idv_lat,h%glat) if (istat /= NF_NOERR) | call handle_ncerr(istat,'Error giving values to glat coord var') ! ! Vertical: istat = nf_put_var_double(ncid,idv_lev,h%lev) ! midpoint levels if (istat /= NF_NOERR) | call handle_ncerr(istat,'Error giving values to zpmid coord var') if (h%new) then istat = nf_put_var_double(ncid,idv_ilev,h%ilev) ! interface levels if (istat /= NF_NOERR) call handle_ncerr(istat, | 'Error giving values to zpint coord var') ! ! Bottom interface level LBC for t,u,v lbc: istat = nf_put_var_double(ncid,idv_lbc,h%ilev(1)) if (istat /= NF_NOERR) | call handle_ncerr(istat,'Error giving value to LBC scalar var') endif ! h%new ! ! Magnetic: istat = nf_put_var_double(ncid,idv_mlon,h%gmlon) istat = nf_put_var_double(ncid,idv_mlat,h%gmlat) ! ! 11/07 btf: mag midpoints are zpmag_mid(nmlev), and ! mag interfaces are zpmag_int(nimlev) ! (see init.F) ! istat = nf_put_var_double(ncid,idv_mlev,h%gmlev) if (istat /= NF_NOERR) then write(char120,"('Error return from nf_put_var_double ', + ' to assign values to gmlev')") call handle_ncerr(istat,char120) endif if (h%new) then istat = nf_put_var_double(ncid,idv_imlev,h%gimlev) if (istat /= NF_NOERR) then write(char120,"('Error return from nf_put_var_double ', + ' to assign values to gimlev')") call handle_ncerr(istat,char120) endif endif end subroutine nc_define !----------------------------------------------------------------------- subroutine define_ncfile(ncid,name,idv,att1name,att1val, | att2name,att2val) ! ! Args: integer,intent(in) :: ncid integer,intent(out) :: idv character(len=*),intent(in) :: name,att1name,att1val, | att2name,att2val ! ! Local: integer :: istat,ids2(2) character(len=120) :: char120 character(len=240) :: char240 ! if (len_trim(name) <= 0) then write(6,"('>>> define_ncfile: empty name! returning..')") return endif ! ids2(1) = id_filelen ids2(2) = id_time ! ! Define time-dependent char variable: istat= nf_def_var(ncid,name,NF_CHAR,2,ids2,idv) if (istat /= NF_NOERR) then write(char120,"('Error defining file var ',a)") trim(name) call handle_ncerr(istat,char120) endif ! ! First attribute: if (len_trim(att1name) > 0.and.len_trim(att1val) > 0) then istat = nf_put_att_text(ncid,idv,att1name,len_trim(att1val), | trim(att1val)) if (istat /= NF_NOERR) then write(char120,"('Error defining ',a,' attribute to var ', | a,': attribute value=',a)") trim(name),trim(att1name), | trim(att1val) call handle_ncerr(istat,char120) endif endif ! ! Second attribute: if (len_trim(att2name) > 0.and.len_trim(att2val) > 0) then istat = nf_put_att_text(ncid,idv,att2name,len_trim(att2val), | trim(att2val)) if (istat /= NF_NOERR) then ! write(char120,"('Error defining ',a,' attribute to var ', write(char240,"('Error defining ',a,' attribute to var ', | a,': attribute value=',a)") trim(name),trim(att2name), | trim(att2val) ! call handle_ncerr(istat,char120) call handle_ncerr(istat,char240) endif endif end subroutine define_ncfile !----------------------------------------------------------------------- subroutine defvar_time_dbl(ncid,name,idtime,idvar, | longname,units) ! ! Define time-dependent variable "name" on open netcdf file ncid. ! Variable has a single dimension in time (e.g., ut(time)). ! Also define long_name and units attributes. ! ! Args: integer,intent(in) :: ncid,idtime(1) integer,intent(out) :: idvar character(len=*),intent(in) :: name,longname,units ! ! Local: integer :: istat character(len=80) :: char80 ! ! Define the variable: istat = nf_def_var(ncid,name,NF_DOUBLE,1,idtime,idvar) if (istat /= NF_NOERR) then write(char80,"('Error defining double variable ',a)") name call handle_ncerr(istat,trim(char80)) endif ! ! Give long name attribute: if (len_trim(longname) > 0) then istat = nf_put_att_text(ncid,idvar,"long_name", | len_trim(longname),longname) if (istat /= NF_NOERR) then write(char80,"('Error defining long_name of double variable ', | a)") name call handle_ncerr(istat,trim(char80)) endif endif ! ! Give units attribute: if (len_trim(units) > 0) then istat = nf_put_att_text(ncid,idvar,"units",len_trim(units), | units) if (istat /= NF_NOERR) then write(char80,"('Error defining units of variable ',a)") name call handle_ncerr(istat,trim(char80)) endif endif end subroutine defvar_time_dbl !----------------------------------------------------------------------- subroutine define_f4d(ncid,idvar,f3d) ! ! e.g., call define_f4d(ncid,idv_f4d(i),h%f3d_geo(i)) ! ! Define a 4d (nlon,nlat,nlev,ntime) variable on current netcdf history file. ! ! type field_3d ! character(len=longname_len) :: long_name ! character(len=shortname_len) :: short_name ! character(len=units_len) :: units ! real :: missing_value ! real,pointer :: data(:,:,:) ! integer :: id1,id2,id3 ! usually nlon,nlat,nlev ! logical :: mag,geo ! character(len=8) :: levdimname ! lev, ilev, mlev, or imlev ! end type field_3d ! ! Args: integer,intent(in) :: ncid integer,intent(out) :: idvar type(field_3d),intent(in) :: f3d ! ! Local: integer :: istat,idimids(4) character(len=80) :: char80 character(len=160) :: char160 real :: var1(1) ! ! 4d fields are dimensioned (lon,lat,lev,time). ! The f3d data is (lon,lat,lev) for current time. ! if (f3d%geo) then idimids(1) = id_lon idimids(2) = id_lat idimids(3) = id_lev ! midpoints if (trim(f3d%levdimname)=='ilev') ! interfaces | idimids(3) = id_ilev else idimids(1) = id_mlon idimids(2) = id_mlat idimids(3) = id_mlev ! midpoints if (trim(f3d%levdimname)=='imlev') ! interfaces | idimids(3) = id_imlev endif idimids(4) = id_time ! ! This call returns the variable id, idvar. istat = nf_def_var(ncid,f3d%short_name,NF_DOUBLE,4,idimids,idvar) if (istat /= NF_NOERR) then write(char80,"('Error defining 4-d field ',a)") | f3d%short_name call handle_ncerr(istat,trim(char80)) endif ! ! Add long_name attribute: istat = nf_put_att_text(ncid,idvar,"long_name", | len_trim(f3d%long_name),trim(f3d%long_name)) ! ! Add units attribute: istat = nf_put_att_text(ncid,idvar,"units",len_trim(f3d%units), | f3d%units) if (istat /= NF_NOERR) then write(char80,"('Error defining units of 4-d', | ' variable ',a,': units = ',a)") trim(f3d%short_name), | trim(f3d%units) call handle_ncerr(istat,trim(char80)) endif ! ! 10/27/05 btf: Add former_name attribute to OMEGA to indicate it ! was previously known as W: if (f3d%short_name=='OMEGA') then istat = nf_put_att_text(ncid,idvar,"former_name",1,"W") endif ! ! Add missing_value attribute: if (f3d%missing_value /= real_unread) then var1(1) = f3d%missing_value istat = nf_put_att_double(ncid,idvar,"missing_value", | NF_DOUBLE,1,var1) if (istat /= NF_NOERR) then write(char160,"('define_f4d: Error return from ', | 'nf_put_att_double for missing_value var attribute')") call handle_ncerr(istat,char160) endif endif end subroutine define_f4d !----------------------------------------------------------------------- subroutine define_f3d(ncid,idvar,f2d) ! ! Define a 3d (nlon,nlat,ntime) variable on current netcdf history file. ! ! Args: integer,intent(in) :: ncid integer,intent(out) :: idvar type(field_2d),intent(in) :: f2d ! ! Local: integer :: istat,idimids(3) character(len=80) :: char80 character(len=160) :: char160 real :: var1(1) ! ! Var is dimensioned (lon,lat,time): if (f2d%geo) then idimids(1) = id_lon idimids(2) = id_lat else idimids(1) = id_mlon idimids(2) = id_mlat endif idimids(3) = id_time ! ! This call returns the variable id, idvar. istat = nf_def_var(ncid,f2d%short_name,NF_DOUBLE,3,idimids,idvar) if (istat /= NF_NOERR) then write(char80,"('Error defining 3-d field ',a)") f2d%short_name call handle_ncerr(istat,trim(char80)) endif ! ! Add long_name attribute: istat = nf_put_att_text(ncid,idvar,"long_name", | len_trim(f2d%long_name),trim(f2d%long_name)) ! ! Units: istat = nf_put_att_text(ncid,idvar,"units", | len_trim(f2d%units),trim(f2d%units)) ! ! 10/27/05 btf: Add former_name attribute to OMEGA to indicate it ! was previously known as W: if (f2d%short_name=='OMEGA') then istat = nf_put_att_text(ncid,idvar,"former_name",1,"W") endif ! ! Add missing_value attribute: if (f2d%missing_value /= real_unread) then var1(1) = f2d%missing_value istat = nf_put_att_double(ncid,idvar,"missing_value", | NF_DOUBLE,1,var1) if (istat /= NF_NOERR) then write(char160,"('define_f3d: Error return from ', | 'nf_put_att_double for missing_value var attribute')") call handle_ncerr(istat,char160) endif endif end subroutine define_f3d !----------------------------------------------------------------------- subroutine write_f4d(ncid,idvar,itime,f3d) ! ! call write_f4d(ncid,idv_f4d(i),h%f3d_geo(i)) ! implicit none ! ! Args: integer,intent(in) :: ncid,idvar,itime type(field_3d),intent(in) :: f3d ! ! Local: integer :: start_4d(4),count_4d(4),istat character(len=120) :: char120 real :: fmin,fmax start_4d(1:3) = 1 start_4d(4) = itime count_4d(1) = f3d%id1 count_4d(2) = f3d%id2 count_4d(3) = f3d%id3 count_4d(4) = 1 istat = nf_put_vara_double(ncid,idvar,start_4d,count_4d, | f3d%data(:,:,:)) if (istat /= NF_NOERR) then write(char120,"('Error return from nf_put_vara_double', | ' for field var ',a,' itime=',i2)") f3d%short_name,itime call handle_ncerr(istat,char120) else call fminmaxspv(f3d%data(:,:,:),f3d%id1*f3d%id2*f3d%id3, | fmin,fmax,f3d%missing_value) if (f3d%levdimname=='ilev') then write(6,"('Write 4d field ',a,' itime=',i3,' min,max=', | 2e12.4,' (interfaces)')") f3d%short_name,itime,fmin,fmax else write(6,"('Write 4d field ',a,' itime=',i3,' min,max=', | 2e12.4,' (midpoints)')") f3d%short_name,itime,fmin,fmax endif endif if (len_trim(fieldlist)==0) then write(fieldlist,"('''',a,'''')") trim(f3d%short_name) else write(fieldlist,"(a,',''',a,'''')") trim(fieldlist), | trim(f3d%short_name) endif end subroutine write_f4d !----------------------------------------------------------------------- subroutine write_f3d(ncid,idvar,itime,f2d) ! ! call write_f3d(ncid,idv_f3d(i),h%f2d_geo(i)) ! implicit none ! ! Args: integer,intent(in) :: ncid,idvar,itime type(field_2d),intent(in) :: f2d ! ! Local: integer :: start_3d(3),count_3d(3),istat character(len=120) :: char120 real :: fmin,fmax start_3d(1:2) = 1 start_3d(3) = itime count_3d(1) = f2d%id1 count_3d(2) = f2d%id2 count_3d(3) = 1 istat = nf_put_vara_double(ncid,idvar,start_3d,count_3d, | f2d%data(:,:)) if (istat /= NF_NOERR) then write(char120,"('Error return from nf_put_vara_double', | ' for field var ',a,' itime=',i2)") f2d%short_name,itime call handle_ncerr(istat,char120) else call fminmaxspv(f2d%data(:,:),f2d%id1*f2d%id2, | fmin,fmax,f2d%missing_value) write(6,"('Write 3d field ',a,' itime=',i3,' min,max=', | 2e12.4)") f2d%short_name,itime,fmin,fmax endif end subroutine write_f3d !----------------------------------------------------------------------- real function mtime_to_datestr(iyear,mtime,imo,ida,datestr) implicit none ! ! Given model time mtime (day,hr,min,sec) and iyear (4-digit year), ! return imo (2 digit month), ida (2 digit day of the month), and ! a date string "minutes since yyyy-m-d". Function value is real ! minutes into the day from modeltime(2:4) hour, min, and sec. ! ! Args: integer,intent(in) :: iyear,mtime(3) integer,intent(out) :: imo,ida character(len=*),intent(out) :: datestr ! ! Local: character(len=2) :: f_mon,f_day,f_hr,f_min character(len=120) :: format integer :: ndmon(12) = ! J F M A M J J A S O N D | (/31,28,31,30,31,30,31,31,30,31,30,31/) integer :: m,id,iday ! mtime_to_datestr = -1. datestr = ' ' ndmon(2) = 28 if (mod(iyear,4).eq.0) ndmon(2) = 29 ! iday = mtime(1) id = 0 do m=1,12 id = id+ndmon(m) if (id.eq.iday) then id = ndmon(m) goto 100 endif if (id.gt.iday) then id = ndmon(m)-id+iday goto 100 endif enddo write(6,"('>>> mtime_to_datestr: could not find date for mtime=', | 4i4,' iyear=',i4)") mtime,iyear imo = 0 ida = 0 return 100 continue imo = m ida = id ! ! Return real minutes: mtime_to_datestr = float(mtime(2))*60.+float(mtime(3)) ! ! Return date string: if (imo <= 9) then if (ida <= 9) then f_day = 'i1' f_mon = 'i1' else ! ida > 9 f_day = 'i2' f_mon = 'i1' endif else ! imo > 9 if (ida <= 9) then f_day = 'i1' f_mon = 'i2' else ! ida > 9 f_day = 'i2' f_mon = 'i2' endif endif if (mtime(2) <= 9) then if (mtime(3) <= 9) then f_hr = 'i1' f_min = 'i1' else f_hr = 'i1' f_min = 'i2' endif else if (mtime(3) <= 9) then f_hr = 'i2' f_min = 'i1' else f_hr = 'i2' f_min = 'i2' endif endif ! ! Example format: ! ('minutes since ',i4,'-',i1,'-',i2,' ',i1,':',i2,':0') ! Example date string: ! minutes since 1997-3-21 0:0:0 ! format = '(''minutes since '',i4,''-'','//f_mon//',''-'','//f_day format = trim(format)//','' '','//f_hr//','':'','//f_min//',' format = trim(format)//''':0'')' write(datestr,format) iyear,imo,ida,mtime(2),mtime(3) end function mtime_to_datestr !------------------------------------------------------------------- subroutine datetime(curdate,curtime) ! ! Return character*8 values for current date and time. ! (sub date_and_time is an f90 intrinsic) ! implicit none ! ! Args: character(len=*),intent(out) :: curdate,curtime ! ! Local: character(len=8) :: date character(len=10) :: time character(len=5) :: zone integer :: values(8) ! curdate = ' ' curtime = ' ' call date_and_time(date,time,zone,values) ! ! write(6,"('datetime: date=',a,' time=',a,' zone=',a)") ! | date,time,zone ! write(6,"('datetime: values=',8i8)") values ! curdate(1:2) = date(5:6) curdate(3:3) = '/' curdate(4:5) = date(7:8) curdate(6:6) = '/' curdate(7:8) = date(3:4) ! curtime(1:2) = time(1:2) curtime(3:3) = ':' curtime(4:5) = time(3:4) curtime(6:6) = ':' curtime(7:8) = time(5:6) ! end subroutine datetime !----------------------------------------------------------------------- real function mtime_delta(mtime0,mtime1) ! ! Return difference in time (minutes) from mtime0 to mtime1: ! implicit none integer,intent(in) :: mtime0(4),mtime1(4) real :: rmins0,rmins1 ! rmins0 = mtime0(1)*1440.+mtime0(2)*60.+mtime0(3)+mtime0(4)/60. rmins1 = mtime1(1)*1440.+mtime1(2)*60.+mtime1(3)+mtime1(4)/60. mtime_delta = rmins1-rmins0 ! if (mtime_delta < 0.) then ! write(6,"('>>> WARNING mtime_delta: negative delta: mtime0=', ! | 4i4,' mtime1=',4i4,' delta=',f10.2)") ! | mtime0,mtime1,mtime_delta ! endif end function mtime_delta !----------------------------------------------------------------------- subroutine nc_wrhist(ncid,h,mtimes,ntimes) ! ! Args: integer,intent(in) :: ncid,ntimes integer,intent(in) :: mtimes(3,ntimes) type(history),intent(in) :: h ! ! Local: ! ! Local: integer :: mins,istat,ivar1(1),i,j,iprog,imo,ida,mtimeinit(4), | ivar2(2),ii character(len=120) :: char120 character(len=240) :: char240 character(len=80) :: char80 integer :: ndims,nvars,ngatts,unlimdimid real :: var1(1),var22(2,2),var2(2),var10(10),rmins,rdays ! ! External: integer,external :: strloc ! ! Time: decimal days since start time of initial run: ! h%initial_mtime was set by sub hist_initype (if initial run), ! or was read from source history (if continuation run). ! mtimeinit(1:3) = h%initial_mtime(1:3) ; mtimeinit(4) = 0 rmins = mtime_delta(mtimeinit,h%modeltime) rdays = rmins/1440. if (rmins < 0.) then write(6,"(/,'>>> nc_wrhist: time from initial start must be', | ' >= 0: mtimeinit=',3i4,' h%modeltime=',3i4,' delta (mins)=', | f10.2)") mtimeinit(1:3),h%modeltime(1:3),rmins stop 'delta time from init run' else ! write(6,"('nc_wrhist: mtimeinit=',3i4,' h%modeltime=',3i4, ! | ' delta (mins)=',f10.2,' delta (days)=',f8.2)") ! | mtimeinit(1:3),h%modeltime(1:3),rmins,rdays endif if (idv_time <= 0) istat = nf_inq_varid(ncid,"time",idv_time) istat = nf_put_var1_double(ncid,idv_time,h%ihist,rmins) if (istat /= NF_NOERR) then write(char240,"('Error from nf_put_var1_double defining ', | 'time at h%modeltime=',4i4,' mtimeinit=',3i4, | ' delta (days)=',f8.2)") h%modeltime,rdays,mtimeinit(1:3) call handle_ncerr(istat,char240) endif ! ! Model time (day,hour,minute): start_2d(1) = 1 start_2d(2) = h%ihist count_2d(1) = 3 count_2d(2) = 1 if (idv_mtime <= 0) istat = nf_inq_varid(ncid,"mtime",idv_mtime) istat = nf_put_vara_int(ncid,idv_mtime,start_2d,count_2d, | h%modeltime(1:3)) if (istat /= NF_NOERR) then write(char120,"('Error defining mtime at h%ihist=',i2, + ' h%modeltime=',4i4,' istat=',i5)") h%ihist,h%modeltime,istat call handle_ncerr(istat,char120) endif ! ! Calendar year (initially from input, then updated if model is ! advancing calendar time): if (idv_year <= 0) istat = nf_inq_varid(ncid,"year",idv_year) istat = nf_put_var1_int(ncid,idv_year,h%ihist,h%year) if (istat /= NF_NOERR) then write(char120,"('Error from nf_put_var1_int defining ', | 'calendar year at h%ihist=',i3,': h%year=',i5)") | h%ihist,h%year call handle_ncerr(istat,char120) endif ! ! Calendar day (initially from input, then updated if model is ! advancing calendar time): if (idv_day <= 0) istat = nf_inq_varid(ncid,"day",idv_day) istat = nf_put_var1_int(ncid,idv_day,h%ihist,h%day) if (istat /= NF_NOERR) then write(char120,"('Error from nf_put_var1_int defining ', | 'calendar day at h%ihist=',i3,': h%day=',i5)") | h%ihist,h%day call handle_ncerr(istat,char120) endif ! ! Calendar advance: if (idv_calendar_adv <= 0) | istat = nf_inq_varid(ncid,"calendar_advance",idv_calendar_adv) istat = nf_put_var1_int(ncid,idv_calendar_adv,h%ihist, | h%calendar_advance) if (istat /= NF_NOERR) then write(char120,"('Error from nf_put_var1_int defining ', | 'calendar advance flag at h%ihist=',i3, | ': h%calendar_advance=',i5)") h%ihist,h%calendar_advance call handle_ncerr(istat,char120) endif ! ! Date and time history was written: if (idv_writedate <= 0) istat = nf_inq_varid(ncid,"write_date", | idv_writedate) start_2d(1) = 1 start_2d(2) = h%ihist count_2d(1) = len_trim(h%write_date) count_2d(2) = 1 istat = nf_put_vara_text(ncid,idv_writedate,start_2d,count_2d, | h%write_date) if (istat /= NF_NOERR) then write(char120,"('Error from nf_put_vara_text writing ', | 'h%write_date=',a)") trim(h%write_date) call handle_ncerr(istat,char120) endif ! ! Iteration number (number of steps from 0,0,0 to current mtime): if (idv_iter <= 0) istat = nf_inq_varid(ncid,"iter",idv_iter) istat = nf_put_var1_int(ncid,idv_iter,h%ihist,h%iter) if (istat /= NF_NOERR) then write(char120,"('Error from nf_put_var1_int defining ', + 'iter at h%ihist=',i3,': h%iter=',i8)") h%ihist,h%iter call handle_ncerr(istat,char120) endif ! ! Universal time (decimal hours): if (idv_ut <= 0) istat = nf_inq_varid(ncid,"ut",idv_ut) istat = nf_put_var1_double(ncid,idv_ut,h%ihist,h%ut) ! ! 1/31/01: tides are now time dependent 2d vars: ! Diurnal tide (dtide(time,2): start_2d(1) = 1 start_2d(2) = h%ihist count_2d(1) = 2 count_2d(2) = 1 if (idv_dtide <= 0) istat = nf_inq_varid(ncid,"dtide",idv_dtide) var2(:) = h%dtide(:) istat = nf_put_vara_double(ncid,idv_dtide,start_2d,count_2d,var2) ! ! Semi-diurnal tide (sdtide(time,10): start_2d(1) = 1 start_2d(2) = h%ihist count_2d(1) = 10 count_2d(2) = 1 if (idv_sdtide <= 0) istat = nf_inq_varid(ncid,"sdtide", | idv_sdtide) var10(:) = h%sdtide(:) istat = nf_put_vara_double(ncid,idv_sdtide,start_2d,count_2d, | var10) ! ! Update gpi_ncfile: call update_ncfile(ncid,"gpi_ncfile",h%gpi_ncfile,idv_gpi_ncfile, | h%ihist) ! ! Update ncep_ncfile: call update_ncfile(ncid,"ncep_ncfile",h%ncep_ncfile, | idv_ncep_ncfile,h%ihist) ! ! Update see_ncfile: call update_ncfile(ncid,"see_ncfile",h%see_ncfile,idv_see_ncfile, | h%ihist) ! ! Update imf_ncfile: call update_ncfile(ncid,"imf_ncfile",h%imf_ncfile,idv_imf_ncfile, | h%ihist) ! ! Update gswm files: call update_ncfile(ncid,"gswm_mi_di_ncfile",h%gswm_mi_di_ncfile, | idv_gswm_mi_di_ncfile,h%ihist) call update_ncfile(ncid,"gswm_mi_sdi_ncfile",h%gswm_mi_sdi_ncfile, | idv_gswm_mi_sdi_ncfile,h%ihist) call update_ncfile(ncid,"gswm_nm_di_ncfile",h%gswm_nm_di_ncfile, | idv_gswm_nm_di_ncfile,h%ihist) call update_ncfile(ncid,"gswm_nm_sdi_ncfile",h%gswm_nm_sdi_ncfile, | idv_gswm_nm_sdi_ncfile,h%ihist) ! ! f107d and f107a (daily and average solar fluxes): if (idv_f107d <= 0) istat = nf_inq_varid(ncid,"f107d", | idv_f107d) istat = nf_put_var1_double(ncid,idv_f107d,h%ihist,h%f107d) if (idv_f107a <= 0) istat = nf_inq_varid(ncid,"f107a", | idv_f107a) istat = nf_put_var1_double(ncid,idv_f107a,h%ihist,h%f107a) ! ! Hemispheric power: if (idv_hpower <= 0) istat = nf_inq_varid(ncid,"hpower", | idv_hpower) istat = nf_put_var1_double(ncid,idv_hpower,h%ihist,h%hpower) ! ! Cross-tail potential: if (idv_ctpoten <= 0) istat = nf_inq_varid(ncid,"ctpoten", | idv_ctpoten) istat = nf_put_var1_double(ncid,idv_ctpoten,h%ihist,h%ctpoten) ! ! BX,BY,BZ imf: if (idv_bximf <= 0) istat = nf_inq_varid(ncid,"bximf",idv_bximf) istat = nf_put_var1_double(ncid,idv_bximf,h%ihist,h%bximf) if (idv_byimf <= 0) istat = nf_inq_varid(ncid,"byimf",idv_byimf) istat = nf_put_var1_double(ncid,idv_byimf,h%ihist,h%byimf) if (idv_bzimf <= 0) istat = nf_inq_varid(ncid,"bzimf",idv_bzimf) istat = nf_put_var1_double(ncid,idv_bzimf,h%ihist,h%bzimf) ! ! swvel,swden: if (idv_swvel <= 0) istat = nf_inq_varid(ncid,"swvel",idv_swvel) istat = nf_put_var1_double(ncid,idv_swvel,h%ihist,h%swvel) if (idv_swden <= 0) istat = nf_inq_varid(ncid,"swden",idv_swden) istat = nf_put_var1_double(ncid,idv_swden,h%ihist,h%swden) ! ! al: if (idv_al <= 0) istat = nf_inq_varid(ncid,"swden",idv_al) istat = nf_put_var1_double(ncid,idv_al,h%ihist,h%al) ! ! Collision factor: if (idv_colfac <= 0) istat = nf_inq_varid(ncid,"colfac", | idv_colfac) istat = nf_put_var1_double(ncid,idv_colfac,h%ihist,h%colfac) ! ! e1,e2: if (idv_e1 <= 0) istat = nf_inq_varid(ncid,"e1",idv_e1) istat = nf_put_var1_double(ncid,idv_e1,h%ihist,h%e1) if (idv_e2 <= 0) istat = nf_inq_varid(ncid,"e2",idv_e2) istat = nf_put_var1_double(ncid,idv_e2,h%ihist,h%e2) ! h1,h2: if (idv_h1 <= 0) istat = nf_inq_varid(ncid,"h1",idv_h1) istat = nf_put_var1_double(ncid,idv_h1,h%ihist,h%h1) if (idv_h2 <= 0) istat = nf_inq_varid(ncid,"h2",idv_h2) istat = nf_put_var1_double(ncid,idv_h2,h%ihist,h%h2) ! alfac,ec: if (idv_alfac <= 0) istat = nf_inq_varid(ncid,"alfac",idv_alfac) istat = nf_put_var1_double(ncid,idv_alfac,h%ihist,h%alfac) if (idv_ec <= 0) istat = nf_inq_varid(ncid,"ec",idv_ec) istat = nf_put_var1_double(ncid,idv_ec,h%ihist,h%ec) ! alfad,ed: if (idv_alfad <= 0) istat = nf_inq_varid(ncid,"alfad",idv_alfad) istat = nf_put_var1_double(ncid,idv_alfad,h%ihist,h%alfad) if (idv_ed <= 0) istat = nf_inq_varid(ncid,"ed",idv_ed) istat = nf_put_var1_double(ncid,idv_ed,h%ihist,h%ed) ! ! Reference pressure used to convert ln(p0/p) to pressure (hPa): ! (h%p0 was set in sub define_hist, output.F) if (idv_p0 <= 0) istat = nf_inq_varid(ncid,"p0",idv_p0) istat = nf_put_var1_double(ncid,idv_p0,h%ihist,h%p0) ! ! Reference pressure used by the model (microbars): ! (h%p0_model was set in sub define_hist, output.F) if (idv_p0_model <= 0) istat = nf_inq_varid(ncid,"p0_model", | idv_p0_model) istat = nf_put_var1_double(ncid,idv_p0_model,h%ihist,h%p0_model) ! ! Gravitational constant used in the model: if (idv_grav <= 0) istat = nf_inq_varid(ncid,"grav",idv_grav) istat = nf_put_var1_double(ncid,idv_grav,h%ihist,h%grav) ! ! Timestep (seconds): if (idv_step <= 0) istat = nf_inq_varid(ncid,"timestep",idv_step) istat = nf_put_var1_int(ncid,idv_step,h%ihist,h%step) if (istat /= NF_NOERR) then write(char120,"('Error from nf_put_var1_int defining ', + 'ut at h%ihist=',i3,': h%step=',i4)") h%ihist,h%step call handle_ncerr(istat,char120) endif ! ! Update run_type (initial or continuation): istat = nf_redef(ncid) ! put dataset in define mode istat = nf_put_att_text(ncid,NF_GLOBAL,"run_type", | len_trim(h%run_type),trim(h%run_type)) if (istat /= NF_NOERR) then write(char120,"('nc_wrhist: Error return from nf_put_att_text ', | 'for run_type global attribute')") call handle_ncerr(istat,char120) endif ! ! Update source_file and source mtime: if (h%run_type(1:4)=='init') then ! initial source file istat = nf_put_att_text(ncid,NF_GLOBAL,"source_file", | len_trim(h%source_file)+10,trim(h%source_file)//' (initial)') elseif (h%run_type(1:4)=='cont') then ! continuation source file istat = nf_put_att_text(ncid,NF_GLOBAL,"source_file", | len_trim(h%source_file)+15, | trim(h%source_file)//' (continuation)') endif if (istat /= NF_NOERR) then write(char120,"('Error return from nf_put_att_text ', | 'for source_file global attribute')") call handle_ncerr(istat,char120) endif istat = nf_put_att_int(ncid,NF_GLOBAL,"source_mtime",NF_INT,3, | h%source_mtime) if (istat /= NF_NOERR) then write(char120,"('Error return from nf_put_att_text ', | 'for source_mtime global attribute')") call handle_ncerr(istat,char120) endif istat = nf_enddef(ncid) ! take out of define mode ! ! Update number of histories on the file (global attribute): ivar1(1) = h%ihist istat = nf_put_att_int(ncid,NF_GLOBAL,"nhist",NF_INT,1,ivar1) if (istat /= NF_NOERR) then write(char120,"('Error return from nf_put_att_int ', + 'for nhist global attribute')") call handle_ncerr(istat,char120) endif ! ! Write 4d fields (3d fields at current time): if (h%nf3d_geo > 0) then do i=1,h%nf3d_geo call write_f4d(ncid,idv_f4d(i),h%ihist,h%f3d_geo(i)) enddo endif if (h%nf2d_geo > 0) then do i=1,h%nf2d_geo call write_f3d(ncid,idv_f3d(i),h%ihist,h%f2d_geo(i)) enddo endif ! ! idv_tlbc = 0 ; idv_ulbc = 0 ; idv_vlbc = 0 ! forces wrlbc to inquire ! idv_tlbc_nm = 0 ; idv_ulbc_nm = 0 ; idv_vlbc_nm = 0 ! forces wrlbc to inquire ! if (h%hist_type(1:3)=='pri') then ! do i=1,nf4d_hist ! call wrf4d(ncid,f4d(i)%short_name,h%ihist,fakeflds,i) ! ! Write ZG after Z (zg was calculated by sub calczg in addiag.F, ! then gathered to fzg by mp_gather2root_f3d): ! 3/1/06 btf: do not write ZG to primary histories. ! if (trim(f4d(i)%short_name)=='Z') ! | call wrf3d(ncid,fzg,h%ihist,'ZG',idv_zg) ! enddo ! i=1,nf4d_hist end subroutine nc_wrhist !----------------------------------------------------------------------- subroutine update_ncfile(ncid,vname,ncfile,idv,ihist) ! ! Update file text var vname for history ihist. ! Value of the variable to write is ncfile. ! implicit none ! ! Args: integer,intent(in) :: ncid,ihist integer,intent(inout) :: idv character(len=*),intent(in) :: vname,ncfile ! ! Local: integer :: istat,start(2),count(2) character(len=len(ncfile)) :: file character(len=120) :: char120 ! if (idv <= 0) istat = nf_inq_varid(ncid,vname,idv) start(1) = 1 start(2) = ihist count(1) = len_trim(ncfile) count(2) = 1 istat = nf_put_vara_text(ncid,idv,start,count,trim(ncfile)) if (istat /= NF_NOERR) then write(char120,"('Error from nf_put_vara_text to update ncfile ', | a)") trim(file) call handle_ncerr(istat,char120) else ! write(6,"('Updated var ',a,' for ihist=',i3,' file=',a)") ! | trim(vname),ihist,trim(file) endif end subroutine update_ncfile !----------------------------------------------------------------------- end module writehist_module