Ignore:
Timestamp:
Feb 10, 2026, 9:46:10 AM (8 days ago)
Author:
emillour
Message:

Mars PCM:
Add a "output_diagfi" (default .true.) flag so that the user can decide if
a diagfi.nc file should be outputted.
Also do some cleaning around the few remaining calls to writediagfi
or reference to it througout the code; write_output() should be used.
EM

Location:
trunk/LMDZ.MARS/libf/aeronomars
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/libf/aeronomars/calchim_mod.F90

    r3726 r4058  
    4040      use iono_h, only: temp_elect
    4141      use wstats_mod, only: wstats
     42      use write_output_mod, only: write_output
    4243      use callkeys_mod, only: photochem
    4344
     
    223224!     for output:
    224225
    225       logical,save :: output        ! to issue calls to writediagfi and stats
     226      logical,save :: output        ! to issue calls to write_output and stats
    226227      real :: jo3_3d(ngrid,nlayer)  ! Photodissociation rate O3->O1D (s-1)
    227228      real :: jh2o_3d(ngrid,nlayer)  ! Photodissociation rate H2O->H+OH (s-1)
     
    956957      if (photochem .and. output) then
    957958         if (ngrid > 1) then
    958             call writediagfi(ngrid,'jo3','j o3->o1d',    &
    959                              's-1',3,jo3_3d(1,1))
    960             call writediagfi(ngrid,'jh2o','jh2o',    &
    961                              's-1',3,jh2o_3d(1,1))
    962             call writediagfi(ngrid,'iter','iterations',  &
    963                              ' ',3,iter_3d(1,1))
    964 
    965 !            if (.not. unichim) then
    966                call writediagfi(ngrid,'emission_no',        &
    967                     'NO nightglow emission rate','cm-3 s-1',3,emission_no)
    968                call writediagfi(ngrid,'emission_o2',        &
    969                     'O2 nightglow emission rate','cm-3 s-1',3,emission_o2)
    970 !            endif
     959           call write_output('jo3','j o3->o1d','s-1',jo3_3d)
     960           call write_output('jh2o','jh2o','s-1',jh2o_3d)
     961           call write_output('iter_chemistry','iterations',' ',iter_3d)
     962
     963           call write_output('emission_no', &
     964                             'NO nightglow emission rate','cm-3 s-1', &
     965                             emission_no)
     966           call write_output('emission_o2', &
     967                             'O2 nightglow emission rate','cm-3 s-1', &
     968                             emission_o2)
    971969           
    972                call wstats(ngrid,'jo3','j o3->o1d',       &
     970           call wstats(ngrid,'jo3','j o3->o1d',       &
    973971                           's-1',3,jo3_3d(1,1))
    974                call wstats(ngrid,'emission_no',           &
     972           call wstats(ngrid,'emission_no',           &
    975973                   'NO nightglow emission rate','cm-3 s-1',3,emission_no)
    976                call wstats(ngrid,'emission_o2',           &
     974           call wstats(ngrid,'emission_o2',           &
    977975                   'O2 nightglow emission rate','cm-3 s-1',3,emission_o2)
    978                call wstats(ngrid,'mmean','mean molecular mass',       &
     976           call wstats(ngrid,'mmean','mean molecular mass',       &
    979977                           'g.mole-1',3,mmean(1,1))
    980978         end if ! of if (ngrid.gt.1)
  • trunk/LMDZ.MARS/libf/aeronomars/surfacearea.F

    r3726 r4058  
    1616      use comcstfi_h, only: pi
    1717      use wstats_mod, only: wstats
     18      use write_output_mod, only: write_output
    1819      use callkeys_mod, only: microphys
    1920      implicit none
     
    126127     $            "micron2 cm-3",3,surfice*1.e6)
    127128
    128       call writediagfi(ngrid,"surfdust", "Dust surface area",
    129      $            "micron2 cm-3",3,surfdust*1.e6)
    130       call writediagfi(ngrid,"surfice", "Ice cloud surface area",
    131      $            "micron2 cm-3",3,surfice*1.e6)
     129      call write_output("surfdust","Dust surface area",
     130     &                  "micron2 cm-3",surfdust*1.e6)
     131      call write_output("surfice","Ice cloud surface area",
     132     &                  "micron2 cm-3",surfice*1.e6)
    132133
    133134      end subroutine surfacearea
Note: See TracChangeset for help on using the changeset viewer.