!-------------------------------------------------------! ! liubo 2015-07-27 !-------------------------------------------------------! module restoring_mod use kinds_mod use domain use constants use broadcast use io use forcing_tools use time_management use prognostic use grid use exit_mod implicit none private save public :: init_restoring character (char_len) :: & forcing_infile_sst, & forcing_infile_spc real (r8), dimension(nx_block, ny_block,max_blocks_clinic), public :: & res_sst, res_spc contains subroutine init_restoring integer (int_kind) :: & i,j, k, n, iblock, &! loop indices nml_error, &! namelist error flag forcing_year_sst, & forcing_day_sst real(r8) :: forcing_time_sst, & forcing_time_spc, & forcing_data_inc character(char_len) :: & forcing_filename_sst, & forcing_filename_spc type (datafile) :: & res_data_file_sst, & res_data_file_spc type(io_field_desc) :: & res_datain_sst, & res_datain_spc type(io_dim) :: & i_dim, j_dim, t_dim !forcing_infile_sst is in the format 'root', such as !'/home/liub/sst' if the full name is !'/home/liub/sst.1880.001.12' namelist /res_nml/ & forcing_infile_sst, & forcing_infile_spc forcing_infile_sst='unknown-sst' forcing_infile_spc='unknown-spc' if (my_task == master_task) then open (nml_in, file=nml_filename, status='old', iostat=nml_error) if (nml_error /= 0) then nml_error = -1 else nml_error = 1 endif do while (nml_error > 0) read(nml_in, nml=res_nml,iostat=nml_error) end do if (nml_error == 0) close(nml_in) end if call broadcast_scalar(nml_error, master_task) if (nml_error /= 0)then call exit_POP(sigAbort,'ERROR: reading res_nml') endif call broadcast_scalar(forcing_infile_sst, master_task) call broadcast_scalar(forcing_infile_spc, master_task) ! forcing_time_sst = thour00 ! forcing_filename_sst is the output complete filename for ! sst restoring filename. ! forcing_infile_sst is expressed as '/../sst', that is the ! semi-finished sst restoring filename. ! forcing_time_sst is the accumulated hour since 01-01-0000 ! thour00 ! forcing_data_inc = 24 ! the following lines denoted that : ! for every timestep in the same day, using the same data file. forcing_year_sst = iyear forcing_time_sst = iyear*365*24 + (iday_of_year)*24 -12. ! if (my_task == master_task) then ! print*,"iyear: ", iyear ! print*, "iday_of_year: ", iday_of_year ! print*, "forcing_time_sst: ", forcing_time_sst ! print*, "thour00: ", thour00 ! endif !check if forcing_time_sst .eq. iday_of_year forcing_time_spc = forcing_time_sst forcing_data_inc = 24. !!!!!!!!!!!!!!!!!!!!!!!!!!!!read res_sst!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! call get_forcing_filename(forcing_filename_sst, & forcing_infile_sst, & forcing_time_sst, & forcing_data_inc) if (my_task == master_task) then print*,"forcing_filename_sst:",trim(forcing_filename_sst)//'.nc' endif res_data_file_sst = construct_file('nc', & full_name=trim(forcing_filename_sst)//'.nc' ) call data_set(res_data_file_sst, 'open_read') i_dim = construct_io_dim('i', nx_global) j_dim = construct_io_dim('j', ny_global) res_datain_sst = construct_io_field('sst', & dim1=i_dim, dim2=j_dim, & field_loc = field_loc_center, & field_type = field_type_scalar, & d2d_array = res_sst(:,:,:)) call data_set(res_data_file_sst, 'define', res_datain_sst) call data_set(res_data_file_sst, 'read', res_datain_sst) call data_set(res_data_file_sst, 'close') call destroy_io_field(res_datain_sst) call destroy_file(res_data_file_sst) !print*, "res_sst: ", res_sst !!!!!!!!!!!!!!!!!!!!!!!!!!!!read res_spc!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! call get_forcing_filename(forcing_filename_spc, & forcing_infile_spc, & forcing_time_spc, & forcing_data_inc) if (my_task == master_task) then print*,"forcing_filename_spc:",trim(forcing_filename_spc)//'.nc' endif res_data_file_spc = construct_file('nc', & full_name=trim(forcing_filename_spc)//'.nc' ) call data_set(res_data_file_spc, 'open_read') i_dim = construct_io_dim('i', nx_global) j_dim = construct_io_dim('j', ny_global) res_datain_spc = construct_io_field('X2', & dim1=i_dim, dim2=j_dim, & field_loc = field_loc_center, & field_type = field_type_scalar, & d2d_array = res_spc(:,:,:)) call data_set(res_data_file_spc, 'define', res_datain_spc) call data_set(res_data_file_spc, 'read', res_datain_spc) call data_set(res_data_file_spc, 'close') call destroy_io_field(res_datain_spc) call destroy_file(res_data_file_spc) end subroutine init_restoring end module restoring_mod