writehist_loc.f90 Source File


This file depends on

sourcefile~~writehist_loc.f90~~EfferentGraph sourcefile~writehist_loc.f90 writehist_loc.f90 sourcefile~com_io_dyn_mod.f90 com_io_dyn_mod.f90 sourcefile~writehist_loc.f90->sourcefile~com_io_dyn_mod.f90 sourcefile~parallel_lmdz.f90 parallel_lmdz.F90 sourcefile~writehist_loc.f90->sourcefile~parallel_lmdz.f90 sourcefile~comconst_mod.f90 comconst_mod.f90 sourcefile~writehist_loc.f90->sourcefile~comconst_mod.f90 sourcefile~paramet_mod_h.f90 paramet_mod_h.f90 sourcefile~writehist_loc.f90->sourcefile~paramet_mod_h.f90 sourcefile~misc_mod.f90 misc_mod.f90 sourcefile~writehist_loc.f90->sourcefile~misc_mod.f90 sourcefile~iniprint_mod_h.f90 iniprint_mod_h.f90 sourcefile~writehist_loc.f90->sourcefile~iniprint_mod_h.f90 sourcefile~temps_mod.f90 temps_mod.f90 sourcefile~writehist_loc.f90->sourcefile~temps_mod.f90 sourcefile~comgeom_mod_h.f90 comgeom_mod_h.f90 sourcefile~writehist_loc.f90->sourcefile~comgeom_mod_h.f90 sourcefile~infotrac.f90 infotrac.f90 sourcefile~writehist_loc.f90->sourcefile~infotrac.f90 sourcefile~parallel_lmdz.f90->sourcefile~paramet_mod_h.f90 sourcefile~parallel_lmdz.f90->sourcefile~iniprint_mod_h.f90 sourcefile~vampir.f90 vampir.F90 sourcefile~parallel_lmdz.f90->sourcefile~vampir.f90 sourcefile~lmdz_mpi.f90 lmdz_mpi.F90 sourcefile~parallel_lmdz.f90->sourcefile~lmdz_mpi.f90 sourcefile~mod_const_mpi.f90 mod_const_mpi.f90 sourcefile~parallel_lmdz.f90->sourcefile~mod_const_mpi.f90 sourcefile~control_mod.f90 control_mod.f90 sourcefile~parallel_lmdz.f90->sourcefile~control_mod.f90 sourcefile~wxios_mod.f90 wxios_mod.F90 sourcefile~parallel_lmdz.f90->sourcefile~wxios_mod.f90 sourcefile~comgeom_mod_h.f90->sourcefile~paramet_mod_h.f90 sourcefile~infotrac.f90->sourcefile~iniprint_mod_h.f90 sourcefile~lmdz_cppkeys_wrapper.f90 lmdz_cppkeys_wrapper.F90 sourcefile~infotrac.f90->sourcefile~lmdz_cppkeys_wrapper.f90 sourcefile~infotrac.f90->sourcefile~control_mod.f90 sourcefile~strings_mod.f90 strings_mod.f90 sourcefile~infotrac.f90->sourcefile~strings_mod.f90 sourcefile~readtracfiles_mod.f90 readTracFiles_mod.f90 sourcefile~infotrac.f90->sourcefile~readtracfiles_mod.f90 sourcefile~lmdz_reprobus_wrappers.f90 lmdz_reprobus_wrappers.F90 sourcefile~infotrac.f90->sourcefile~lmdz_reprobus_wrappers.f90 sourcefile~readtracfiles_mod.f90->sourcefile~strings_mod.f90 sourcefile~ioipsl_getin_p_mod.f90 ioipsl_getin_p_mod.f90 sourcefile~readtracfiles_mod.f90->sourcefile~ioipsl_getin_p_mod.f90 sourcefile~mod_grid_phy_lmdz.f90 mod_grid_phy_lmdz.f90 sourcefile~lmdz_reprobus_wrappers.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~wxios_mod.f90->sourcefile~iniprint_mod_h.f90 sourcefile~wxios_mod.f90->sourcefile~strings_mod.f90 sourcefile~dimphy.f90 dimphy.f90 sourcefile~wxios_mod.f90->sourcefile~dimphy.f90 sourcefile~mod_phys_lmdz_para.f90 mod_phys_lmdz_para.f90 sourcefile~wxios_mod.f90->sourcefile~mod_phys_lmdz_para.f90 sourcefile~lmdz_xios.f90 lmdz_xios.F90 sourcefile~wxios_mod.f90->sourcefile~lmdz_xios.f90 sourcefile~mod_phys_lmdz_mpi_data.f90 mod_phys_lmdz_mpi_data.f90 sourcefile~wxios_mod.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~geometry_mod.f90 geometry_mod.f90 sourcefile~wxios_mod.f90->sourcefile~geometry_mod.f90 sourcefile~infotrac_phy.f90 infotrac_phy.F90 sourcefile~wxios_mod.f90->sourcefile~infotrac_phy.f90 sourcefile~wxios_mod.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~print_control_mod.f90 print_control_mod.f90 sourcefile~wxios_mod.f90->sourcefile~print_control_mod.f90 sourcefile~mod_phys_lmdz_transfert_para.f90 mod_phys_lmdz_transfert_para.f90 sourcefile~wxios_mod.f90->sourcefile~mod_phys_lmdz_transfert_para.f90 sourcefile~nrtype.f90 nrtype.f90 sourcefile~wxios_mod.f90->sourcefile~nrtype.f90 sourcefile~wxios_mod.f90->sourcefile~ioipsl_getin_p_mod.f90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~print_control_mod.f90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~mod_phys_lmdz_transfert_para.f90 sourcefile~mod_phys_lmdz_omp_data.f90 mod_phys_lmdz_omp_data.F90 sourcefile~mod_phys_lmdz_para.f90->sourcefile~mod_phys_lmdz_omp_data.f90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~lmdz_mpi.f90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~lmdz_cppkeys_wrapper.f90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~print_control_mod.f90 sourcefile~geometry_mod.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~geometry_mod.f90->sourcefile~nrtype.f90 sourcefile~infotrac_phy.f90->sourcefile~iniprint_mod_h.f90 sourcefile~infotrac_phy.f90->sourcefile~lmdz_cppkeys_wrapper.f90 sourcefile~infotrac_phy.f90->sourcefile~strings_mod.f90 sourcefile~infotrac_phy.f90->sourcefile~readtracfiles_mod.f90 sourcefile~infotrac_phy.f90->sourcefile~lmdz_reprobus_wrappers.f90 sourcefile~infotrac_phy.f90->sourcefile~mod_phys_lmdz_para.f90 sourcefile~infotrac_phy.f90->sourcefile~ioipsl_getin_p_mod.f90 sourcefile~mod_phys_lmdz_transfert_para.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_omp_transfert.f90 mod_phys_lmdz_omp_transfert.f90 sourcefile~mod_phys_lmdz_transfert_para.f90->sourcefile~mod_phys_lmdz_omp_transfert.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90 mod_phys_lmdz_mpi_transfert.f90 sourcefile~mod_phys_lmdz_transfert_para.f90->sourcefile~mod_phys_lmdz_mpi_transfert.f90 sourcefile~ioipsl_getin_p_mod.f90->sourcefile~strings_mod.f90 sourcefile~ioipsl_getin_p_mod.f90->sourcefile~mod_phys_lmdz_para.f90 sourcefile~ioipsl_getin_p_mod.f90->sourcefile~mod_phys_lmdz_transfert_para.f90 sourcefile~mod_phys_lmdz_omp_data.f90->sourcefile~dimphy.f90 sourcefile~mod_phys_lmdz_omp_data.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_omp_data.f90->sourcefile~print_control_mod.f90 sourcefile~mod_phys_lmdz_omp_transfert.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_omp_transfert.f90->sourcefile~mod_phys_lmdz_omp_data.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90->sourcefile~lmdz_mpi.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90->sourcefile~mod_grid_phy_lmdz.f90

Contents

Source Code


Source Code

!
! $Id: writedynav_p.F 1279 2009-12-10 09:02:56Z fairhead $
!
subroutine writehist_loc( time, vcov, ucov,teta,ppk,phi,q, &
        masse,ps,phis)

  ! This routine needs IOIPSL
  USE iniprint_mod_h
  USE comgeom_mod_h
  USE ioipsl

  USE parallel_lmdz
  USE misc_mod
  USE infotrac, ONLY : nqtot
  use com_io_dyn_mod, only : histid,histvid,histuid
  USE comconst_mod, ONLY: cpp
  USE temps_mod, ONLY: itau_dyn

  USE dimensions_mod, ONLY: iim, jjm, llm, ndm
USE paramet_mod_h
implicit none

  !
  !   Ecriture du fichier histoire au format IOIPSL
  !
  !   Appels succesifs des routines: histwrite
  !
  !   Entree:
  !  histid: ID du fichier histoire
  !  time: temps de l'ecriture
  !  vcov: vents v covariants
  !  ucov: vents u covariants
  !  teta: temperature potentielle
  !  phi : geopotentiel instantane
  !  q   : traceurs
  !  masse: masse
  !  ps   :pression au sol
  !  phis : geopotentiel au sol
  !
  !
  !   Sortie:
  !  fileid: ID du fichier netcdf cree
  !
  !   L. Fairhead, LMD, 03/99
  !
  ! =====================================================================
  !
  !   Declarations


  !
  !   Arguments
  !

  REAL :: vcov(ijb_v:ije_v,llm),ucov(ijb_u:ije_u,llm)
  REAL :: teta(ijb_u:ije_u,llm),phi(ijb_u:ije_u,llm)
  REAL :: ppk(ijb_u:ije_u,llm)
  REAL :: ps(ijb_u:ije_u),masse(ijb_u:ije_u,llm)
  REAL :: phis(ijb_u:ije_u)
  REAL :: q(ijb_u:ije_u,llm,nqtot)
  integer :: time


  ! This routine needs IOIPSL
  !   Variables locales
  !
  INTEGER,SAVE,ALLOCATABLE :: ndex2d(:),ndexu(:),ndexv(:)
  INTEGER :: iq, ii, ll
  REAL,SAVE,ALLOCATABLE :: tm(:,:)
  REAL,SAVE,ALLOCATABLE :: vnat(:,:),unat(:,:)
  logical :: ok_sync
  integer :: itau_w
  integer :: ijb,ije,jjn
  LOGICAL,SAVE :: first=.TRUE.
!$OMP THREADPRIVATE(first)

  !
  !  Initialisations
  !
  if (adjust) return

  IF (first) THEN
!$OMP BARRIER
!$OMP MASTER
    ALLOCATE(unat(ijb_u:ije_u,llm))
    ALLOCATE(vnat(ijb_v:ije_v,llm))
    ALLOCATE(tm(ijb_u:ije_u,llm))
    ALLOCATE(ndex2d(ijnb_u*llm))
    ALLOCATE(ndexu(ijnb_u*llm))
    ALLOCATE(ndexv(ijnb_v*llm))
    ndex2d = 0
    ndexu = 0
    ndexv = 0
!$OMP END MASTER
!$OMP BARRIER
    first=.FALSE.
  ENDIF

  ok_sync = .TRUE.
  itau_w = itau_dyn + time

  ! Passage aux composantes naturelles du vent
  call covnat_loc(llm, ucov, vcov, unat, vnat)

  !
  !  Appels a histwrite pour l'ecriture des variables a sauvegarder
  !
  !  Vents U
  !

!$OMP BARRIER
!$OMP MASTER
  ijb=ij_begin
  ije=ij_end
  jjn=jj_nb

  call histwrite(histuid, 'u', itau_w, unat(ijb:ije,:), &
        iip1*jjn*llm, ndexu)
!$OMP END MASTER

  !
  !  Vents V
  !
  ije=ij_end
  if (pole_sud) jjn=jj_nb-1
  if (pole_sud) ije=ij_end-iip1
!$OMP BARRIER
!$OMP MASTER
  call histwrite(histvid, 'v', itau_w, vnat(ijb:ije,:), &
        iip1*jjn*llm, ndexv)
!$OMP END MASTER


  !
  !  Temperature potentielle
  !
  ijb=ij_begin
  ije=ij_end
  jjn=jj_nb
!$OMP MASTER
  call histwrite(histid, 'theta', itau_w, teta(ijb:ije,:), &
        iip1*jjn*llm, ndexu)
!$OMP END MASTER

  !
  !  Temperature
  !

!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
  do ll=1,llm
    do ii = ijb, ije
      tm(ii,ll) = teta(ii,ll) * ppk(ii,ll)/cpp
    enddo
  enddo
!$OMP ENDDO

!$OMP MASTER
  call histwrite(histid, 'temp', itau_w, tm(ijb:ije,:), &
        iip1*jjn*llm, ndexu)
!$OMP END MASTER


  !
  !  Geopotentiel
  !
!$OMP MASTER
  call histwrite(histid, 'phi', itau_w, phi(ijb:ije,:), &
        iip1*jjn*llm, ndexu)
!$OMP END MASTER


  !
  !  Traceurs
  !
  !!$OMP MASTER
  !    DO iq=1,nqtot
  !      call histwrite(histid, tracers(iq)%longName, itau_w,
  ! .                   q(ijb:ije,:,iq), iip1*jjn*llm, ndexu)
  !    enddo
  !!$OMP END MASTER


  !
  !  Masse
  !
!$OMP MASTER
   call histwrite(histid, 'masse', itau_w, masse(ijb:ije,:), &
         iip1*jjn*llm, ndexu)
!$OMP END MASTER


  !
  !  Pression au sol
  !
!$OMP MASTER
   call histwrite(histid, 'ps', itau_w, ps(ijb:ije), &
         iip1*jjn, ndex2d)
!$OMP END MASTER

  !
  !  Geopotentiel au sol
  !
!$OMP MASTER
    ! call histwrite(histid, 'phis', itau_w, phis(ijb:ije),
  ! .                 iip1*jjn, ndex2d)
!$OMP END MASTER

  !
  !  Fin
  !
!$OMP MASTER
  if (ok_sync) then
    call histsync(histid)
    call histsync(histvid)
    call histsync(histuid)
  endif
!$OMP END MASTER

end subroutine writehist_loc