psmile_enddef_metadata.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2006-2010, SGI Germany, Munich, Germany.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !BOP
00006 ! !ROUTINE: PSMILE_Enddef_Metadata
00007 ! !INTERFACE:
00008       subroutine psmile_enddef_metadata ( ierror )
00009 
00010 ! !USES:
00011       use PSMILe, dummy_interface => PSMILE_Enddef_Metadata
00012       use psmile_io_utils,only:indexi
00013 
00014       implicit none
00015       include 'prism.inc' 
00016 ! !OUTPUT PARAMETERS:
00017       integer, Intent (Out)               :: ierror
00018 
00019 !
00020 ! !LOCAL VARIABLES:
00021 !
00022       integer                             :: ierrp(2)
00023       integer                         :: ii,no_of_io_fields
00024       integer                         :: itemp_id,ino_of_glob_grids
00025       integer                             :: jj,kk,iid
00026       integer                             :: ino_of_items
00027       integer                             :: ilp,il_tskid
00028       integer                             :: icltn
00029       integer,allocatable                 :: itemp(:)
00030       integer,allocatable                 :: itemp1(:)
00031       integer,allocatable                 :: ivar_ids(:)
00032       integer,allocatable                 :: iglobal_grid_ids(:)
00033       integer,allocatable                 :: indices(:)
00034       integer,allocatable                 :: ino_grids(:)
00035       character(len=max_name)              :: cl_temp_name
00036 !
00037 ! !DESCRIPTION:
00038 !  
00039 ! 0th) Creates a table of related id's per var_id in case of a GME grid.
00040 ! 1st) Defines domain partitioning informations.
00041 ! 2nd) Open files.
00042 ! 3rd) Writes metadata informations to files.
00043 !
00044 ! !REVISION HISTORY:
00045 ! 
00046 !   Date      Programmer    Description
00047 ! ----------  -----------   -----------
00048 ! 31.10.2003  R. Vogelsang  created
00049 ! 09.12.2003  R. Vogelsang  call to psmile_write_meta inserted
00050 ! 21.06.2004  R. Vogelsang  Correction of the construction of the table of
00051 !                           related IDs.  
00052 ! ---------------------------------------------------------------------
00053 !EOP
00054       Character(len=len_cvs_string),save :: def_meta_cvs= 
00055 '$Id: psmile_enddef_metadata.F90 2687 2010-10-28 15:15:52Z coquart $'
00056 
00057       ierror = 0
00058 #ifdef __PSMILE_WITH_IO
00059 #ifdef VERBOSE
00060       print*,trim(ch_id),' : psmile_enddef_metadata: start'
00061       call psmile_flushstd
00062 #endif
00063 !
00064 !RV   Loop over var_ids find all associated grid_ids belonging to a certain
00065 !RV   global_grid_id
00066 !
00067       allocate(ivar_ids(Number_of_fields_allocated),stat=ierror)
00068       allocate(iglobal_grid_ids(Number_of_fields_allocated),stat=ierror)
00069 
00070       no_of_io_fields=0
00071       do ii=1,Number_of_fields_allocated
00072         if (Fields(ii)%status .eq.  PSMILe_status_defined) then
00073           if (associated(Fields(ii)%io_chan_infos) ) then
00074             if(Fields(ii)%io_chan_infos(1)%status.eq.PSMILe_status_defined) then
00075 !RV   At this point a Field is subject to I/O.
00076               no_of_io_fields=no_of_io_fields+1 
00077               ivar_ids(no_of_io_fields)=ii
00078               iglobal_grid_ids(no_of_io_fields)= &
00079                       Grids(Methods(Fields(ii)%method_id)%grid_id)%global_grid_id
00080             endif
00081           endif
00082         endif
00083       enddo
00084        
00085 
00086 #ifdef VERBOSE
00087       print*,trim(ch_id),' : psmile_enddef_metadata:0: ',no_of_io_fields,iglobal_grid_ids(1:no_of_io_fields)
00088       call psmile_flushstd
00089 #endif
00090 
00091       if(no_of_io_fields.gt.0) then
00092 !
00093 !RV   Sort arrays according to their global_grid_id
00094 !
00095       allocate(indices(no_of_io_fields),stat=ierror)
00096       allocate(itemp(no_of_io_fields),stat=ierror)
00097       allocate(itemp1(Number_of_Grids_allocated),stat=ierror)
00098       allocate(ino_grids(no_of_io_fields),stat=ierror)
00099 
00100       call indexi(no_of_io_fields,iglobal_grid_ids,indices)
00101 
00102       itemp(1:no_of_io_fields)=iglobal_grid_ids(indices(1:no_of_io_fields))
00103       iglobal_grid_ids(1:no_of_io_fields)=itemp(1:no_of_io_fields)
00104 
00105       itemp(1:no_of_io_fields)=ivar_ids(indices(1:no_of_io_fields))
00106       ivar_ids(1:no_of_io_fields)=itemp(1:no_of_io_fields)
00107 !
00108 !RV   Find chunks with identical global_grid_id
00109 !
00110 
00111       ino_of_glob_grids =1
00112       iid=iglobal_grid_ids(1)
00113       ino_grids(ino_of_glob_grids)=0
00114       do ii=1,no_of_io_fields
00115         if(iid.eq.iglobal_grid_ids(ii)) then
00116           ino_grids(ino_of_glob_grids)=ino_grids(ino_of_glob_grids)+1
00117         else
00118           iid=iglobal_grid_ids(ii)
00119           ino_of_glob_grids =ino_of_glob_grids +1
00120           ino_grids(ino_of_glob_grids)=1
00121         endif
00122       enddo
00123 #ifdef VERBOSE
00124       print*,trim(ch_id),' : psmile_enddef_metadata:1:' &
00125                        ,'ino_of_glob_grids,ino_grids(1:ino_of_glob_grids)' &
00126                        ,ino_of_glob_grids,ino_grids(1:ino_of_glob_grids)
00127       call psmile_flushstd
00128 #endif
00129 !
00130 !RV   Build the relationship for var_ids on a certain global grid
00131 !RV   having the same local_name. If several items ocure with the same
00132 !RV   name we have presumably a block structured grid where a certain name
00133 !RV   lives on several blocks! (Please refer to the parallelization
00134 !RV   of the GME grid! The parallelization is done for each logical
00135 !RV   rectangular block rather than distributing blocks over tasks!
00136 !
00137 
00138       iid=1
00139       do ii=1,ino_of_glob_grids
00140 
00141         call indexi(ino_grids(ii),ivar_ids(iid),indices)
00142         itemp(1:ino_grids(ii))=ivar_ids(indices(1:ino_grids(ii))+iid-1)
00143         
00144         icltn=1
00145         cl_temp_name=Fields(itemp(icltn))%local_name
00146 
00147         do while(icltn.le.ino_grids(ii))
00148 
00149 #ifdef VERBOSE
00150           print*,trim(ch_id),' : psmile_enddef_metadata:2:' &
00151                         ,'ii,iid,itemp(1:ino_grids(ii)),trim(cl_temp_name)' &
00152                         ,ii,iid,itemp(1:ino_grids(ii)),trim(cl_temp_name)
00153           call psmile_flushstd
00154 #endif
00155           ino_of_items=0
00156           do jj=icltn,ino_grids(ii)
00157 
00158             if(trim(cl_temp_name).eq.trim(Fields(itemp(jj))%local_name) .and. &
00159               jj.lt.ino_grids(ii)) then
00160 
00161               ino_of_items=ino_of_items+1
00162               itemp1(ino_of_items)=itemp(jj)
00163 
00164             else
00165               if(trim(cl_temp_name).eq.trim(Fields(itemp(jj))%local_name)) then
00166                 ino_of_items=ino_of_items+1
00167                 itemp1(ino_of_items)=itemp(jj)
00168               endif
00169 
00170               do kk=1,ino_of_items
00171 !
00172 !             This is an loop over multiple I/O tasks per varid
00173 !
00174                 do il_tskid=1,size(Fields(itemp1(kk))%io_task_lookup)
00175               
00176                   ilp=Fields(itemp1(kk))%io_task_lookup(il_tskid)
00177 !
00178 !             Let io_infos point to the real location of the I/O infos
00179 !
00180                   if(ilp .gt. 0 ) then
00181                     Fields(itemp1(kk))%io_infos =>  &
00182                      Fields(itemp1(kk))%io_chan_infos(ilp)
00183 !
00184                   if(.not.associated(Fields(itemp1(kk))%io_infos%related_ids)) &
00185                      allocate(Fields(itemp1(kk)) &
00186                               %io_infos%related_ids(ino_of_items))
00187 !
00188 !RV   The block_id in order to find the data at the correct record
00189 !RV   of the NetCDF file
00190 !
00191                   Fields(itemp1(kk)) &
00192                           %io_infos%block_id=kk
00193 !
00194 !RV   The relation ship between var_ids of the same local name.
00195 !
00196                   Fields(itemp1(kk)) &
00197                           %io_infos%related_ids(1:ino_of_items) &
00198                           =itemp1(1:ino_of_items)
00199 !
00200 !RV   opened implements a status register. If one of the entries
00201 !RV   is true a block member of the Field already opened the file.
00202 !
00203                   Fields(itemp1(kk)) &
00204                           %io_infos%opened=.false.
00205 !
00206 !RV   If done is all true all data which are  local to a process were 
00207 !RV   assembled in a 3d field and can be written to file. 
00208 !RV   Likewise on read a process reads its pile of blocks and sets
00209 !RV   done of all related var_ids to true.  
00210 !
00211                   Fields(itemp1(kk)) &
00212                           %io_infos%done=.false.
00213 #ifdef VERBOSE
00214                   print*,trim(ch_id),' : psmile_enddef_metadata: <'
00215                   print*,trim(ch_id),' : varid,local_name,blockid,relation' 
00216                   print*,trim(ch_id),' : ',itemp1(kk),trim(cl_temp_name) &
00217                                      ,Fields(itemp1(kk))%io_infos%block_id &
00218                         ,Fields(itemp1(kk))%io_infos%related_ids(1:ino_of_items)
00219                   print*,trim(ch_id),' : psmile_enddef_metadata: >'
00220       
00221                   call psmile_flushstd
00222 
00223 #endif
00224                 endif
00225                 enddo
00226               enddo
00227               icltn=icltn+ino_of_items
00228               ino_of_items=0
00229               cl_temp_name=Fields(itemp(jj))%local_name
00230             endif
00231 !
00232 !           Leave the search loop for a next iteration.
00233 !
00234             if(ino_of_items.eq.0) exit
00235           enddo
00236         enddo ! End of Do While
00237         iid=ino_grids(ii)+1
00238       enddo
00239 
00240       deallocate(indices,stat=ierror)
00241       deallocate(itemp,stat=ierror)
00242       deallocate(itemp1,stat=ierror)
00243       deallocate(ino_grids,stat=ierror)
00244 
00245       endif
00246     
00247       
00248       
00249       call psmile_def_domains(ierror)
00250       call psmile_open_files(ierror)
00251       call psmile_write_meta(ierror)
00252 
00253 
00254       deallocate(ivar_ids,stat=ierror)
00255       deallocate(iglobal_grid_ids,stat=ierror)
00256       
00257 #ifdef VERBOSE
00258       print*,trim(ch_id),' : psmile_enddef_metadata: end'
00259       call psmile_flushstd
00260       
00261 #endif      
00262 #endif
00263 
00264       end subroutine PSMILe_Enddef_Metadata

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1