source: trunk/LMDZ.MARS/libf/phymars/dyn1d/testphys1d.F90 @ 3074

Last change on this file since 3074 was 3074, checked in by jbclement, 13 months ago

Mars PCM:
Transformation of 'writerestart1D' into module.
JBC

File size: 10.3 KB
Line 
1PROGRAM testphys1d
2
3use comsoil_h,           only: inertiedat, inertiesoil, nsoilmx, tsoil
4use surfdat_h,           only: albedodat, perenial_co2ice, watercap, tsurf, emis, qsurf
5use comslope_mod,        only: def_slope, subslope_dist
6use phyredem,            only: physdem0, physdem1
7use watersat_mod,        only: watersat
8use tracer_mod,          only: igcm_h2o_vap, igcm_h2o_ice, noms
9use comcstfi_h,          only: pi, rad, omeg, g, mugaz, rcp, r, cpp
10use time_phylmdz_mod,    only: daysec, day_step
11use dimradmars_mod,      only: tauvis, totcloudfrac, albedo
12use dust_param_mod,      only: tauscaling
13use comvert_mod,         only: ap, bp, aps, bps, pa, preff, sig
14use physiq_mod,          only: physiq
15use turb_mod,            only: q2
16use write_output_mod,    only: write_output
17use init_testphys1d_mod, only: init_testphys1d
18use writerestart1D_mod,  only: writerestart1D
19! Mostly for XIOS outputs:
20use mod_const_mpi,       only: init_const_mpi
21use parallel_lmdz,       only: init_parallel
22
23implicit none
24
25!=======================================================================
26!   subject:
27!   --------
28!   PROGRAM useful to run physical part of the martian GCM in a 1D column
29!
30! Can be compiled with a command like (e.g. for 25 layers)
31!   "makegcm -p mars -d 25 testphys1d"
32! It requires the files "testphys1d.def" "callphys.def"
33!   and a 'run.def' file (containing a "INCLUDEDEF=callphys.def" line)
34!   and a file describing the sigma layers (e.g. "z2sig.def")
35!
36!   author: Frederic Hourdin, R.Fournier,F.Forget
37!   -------
38!
39!   update: 12/06/2003, including chemistry (S. Lebonnois)
40!                       and water ice (F. Montmessin)
41!           27/09/2023, upgrade to F90 (JB Clément)
42!
43!=======================================================================
44
45include "dimensions.h"
46!#include "dimradmars.h"
47!#include "comgeomfi.h"
48!#include "surfdat.h"
49!#include "slope.h"
50!#include "comsoil.h"
51!#include "comdiurn.h"
52include "callkeys.h"
53!#include "comsaison.h"
54!#include "control.h"
55include "netcdf.inc"
56!#include "advtrac.h"
57
58!--------------------------------------------------------------
59! Declarations
60!--------------------------------------------------------------
61integer, parameter                  :: ngrid = 1     ! (2+(jjm-1)*iim - 1/jjm)
62integer, parameter                  :: nlayer = llm
63real, parameter                     :: odpref = 610. ! DOD reference pressure (Pa)
64integer                             :: unitstart     ! unite d'ecriture de "startfi"
65integer                             :: ndt, ilayer, ilevel, isoil, idt, iq
66logical                             :: firstcall, lastcall
67integer                             :: day0          ! initial (sol ; =0 at Ls=0)
68real                                :: day           ! date during the run
69real                                :: time          ! time (0<time<1 ; time=0.5 a midi)
70real                                :: dttestphys    ! testphys1d timestep
71real, dimension(nlayer)             :: play          ! Pressure at the middle of the layers (Pa)
72real, dimension(nlayer + 1)         :: plev          ! intermediate pressure levels (pa)
73real                                :: psurf         ! Surface pressure
74real, dimension(nlayer)             :: u, v          ! zonal, meridional wind
75real                                :: gru, grv      ! prescribed "geostrophic" background wind
76real, dimension(nlayer)             :: temp          ! temperature at the middle of the layers
77real, dimension(:,:,:), allocatable :: q             ! tracer mixing ratio (e.g. kg/kg)
78real, dimension(1)                  :: wstar = 0.    ! Thermals vertical velocity
79
80! Physical and dynamical tandencies (e.g. m.s-2, K/s, Pa/s)
81real, dimension(nlayer)             :: du, dv, dtemp, dudyn, dvdyn, dtempdyn
82real, dimension(1)                  :: dpsurf
83real, dimension(:,:,:), allocatable :: dq, dqdyn
84
85! Various intermediate variables
86integer                 :: aslun
87real                    :: zls, pks, ptif, qtotinit, qtot
88real, dimension(nlayer) :: phi, h, s, w
89integer                 :: nq = 1 ! number of tracers
90real, dimension(1)      :: latitude, longitude, cell_area
91logical                 :: there
92character(len = 2)      :: str2
93character(len = 7)      :: str7
94character(len = 44)     :: txt
95
96! RV & JBC: Use of starting files for 1D
97logical :: startfiles_1D, therestart1D, therestartfi
98
99! JN & JBC: Force atmospheric water profiles
100real                            :: atm_wat_profile, atm_wat_tau
101real, dimension(:), allocatable :: zqsat ! useful for (atm_wat_profile=2)
102!=======================================================================
103
104!=======================================================================
105! INITIALISATION
106!=======================================================================
107#ifdef CPP_XIOS
108    call init_const_mpi
109    call init_parallel
110#endif
111
112! Initialize "serial/parallel" related stuff
113!call init_phys_lmdz(iim,jjp1,llm,1,(/(jjm-1)*iim+2/))
114!call init_phys_lmdz(1,1,llm,1,(/1/))
115!call initcomgeomphy
116
117call init_testphys1d(.false.,ngrid,nlayer,odpref,nq,q,time,psurf,u,v,temp,startfiles_1D,therestart1D, &
118                     therestartfi,ndt,ptif,pks,dttestphys,zqsat,dq,dqdyn,day0,day,gru,grv,w,          &
119                     play,plev,latitude,longitude,cell_area,atm_wat_profile,atm_wat_tau)
120
121! Write a "startfi" file
122! ----------------------
123! This file will be read during the first call to "physiq".
124! It is needed to transfert physics variables to "physiq"...
125if (.not. therestartfi) then
126    call physdem0("startfi.nc",longitude,latitude,nsoilmx,ngrid, &
127                  llm,nq,dttestphys,float(day0),0.,cell_area,    &
128                  albedodat,inertiedat,def_slope,subslope_dist)
129    call physdem1("startfi.nc",nsoilmx,ngrid,llm,nq,dttestphys,time,       &
130                  tsurf,tsoil,inertiesoil,albedo,emis,q2,qsurf,tauscaling, &
131                  totcloudfrac,wstar,watercap,perenial_co2ice)
132endif !(.not. therestartfi)
133
134!=======================================================================
135!  1D MODEL TIME STEPPING LOOP
136!=======================================================================
137firstcall = .true.
138lastcall = .false.
139do idt = 1,ndt
140    if (idt == ndt) lastcall = .true.
141!    if (idt == ndt - day_step - 1) then ! test
142!        lastcall = .true.
143!        call solarlong(day*1.0,zls)
144!        write(103,*) 'Ls=',zls*180./pi
145!        write(103,*) 'Lat=', latitude(1)*180./pi
146!        write(103,*) 'Tau=', tauvis/odpref*psurf
147!        write(103,*) 'RunEnd - Atmos. Temp. File'
148!        write(103,*) 'RunEnd - Atmos. Temp. File'
149!        write(104,*) 'Ls=',zls*180./pi
150!        write(104,*) 'Lat=', latitude(1)
151!        write(104,*) 'Tau=', tauvis/odpref*psurf
152!        write(104,*) 'RunEnd - Atmos. Temp. File'
153!    endif
154
155    ! Compute geopotential
156    ! ~~~~~~~~~~~~~~~~~~~~
157    s(:) = (aps(:)/psurf + bps(:))**rcp
158    h(:) = cpp*temp(:)/(pks*s(:))
159
160    phi(1) = pks*h(1)*(1. - s(1))
161    do ilayer = 2,nlayer
162        phi(ilayer) = phi(ilayer - 1) + pks*(h(ilayer - 1) + h(ilayer))*.5*(s(ilayer - 1) - s(ilayer))
163    enddo
164
165    ! Modify atmospheric water if asked
166    ! Added "atm_wat_profile" flag (JN + JBC)
167    ! ---------------------------------------
168    if (water) then
169        call watersat(nlayer,temp,play,zqsat)
170        if (0. <= atm_wat_profile .and. atm_wat_profile <= 1.) then
171        ! If atmospheric water is monitored
172            if (atm_wat_tau < 0.) then ! Prescribed atm_wat_profile: wet if >0, dry if =0
173                q(1,:,igcm_h2o_vap) = min(zqsat(:),atm_wat_profile*g/psurf)
174                q(1,:,igcm_h2o_ice) = 0. ! reset h2o ice
175            else ! Relaxation towards the value atm_wat_profile with relaxation time atm_wat_tau
176                q(1,:,igcm_h2o_vap) = atm_wat_profile*g/psurf + (q(1,:,igcm_h2o_vap) - atm_wat_profile*g/psurf)*dexp(-dttestphys/atm_wat_tau)
177                q(1,:,igcm_h2o_vap) = min(zqsat(:),q(1,:,igcm_h2o_vap))
178                q(1,:,igcm_h2o_ice) = 0. ! reset h2o ice
179            endif
180        endif
181    endif
182
183    ! Call physics
184    ! --------------------
185    call physiq(1,llm,nq,firstcall,lastcall,day,time,dttestphys,plev,play,phi,u,v,temp,q,w,du,dv,dtemp,dq,dpsurf)
186    !                                                                                      ^----- outputs -----^
187
188    ! Wind increment: specific for 1D
189    ! --------------------------------
190    ! The physics compute the tendencies on u and v,
191    ! here we just add Coriolos effect
192    !do ilayer = 1,nlayer
193    !    du(ilayer) = du(ilayer) + ptif*(v(ilayer) - grv)
194    !    dv(ilayer) = dv(ilayer) + ptif*(-u(ilayer) + gru)
195    !enddo
196    ! For some tests: No coriolis force at equator
197    !if (latitude(1) == 0.) then
198    du(:) = du(:) + (gru - u(:))/1.e4
199    dv(:) = dv(:) + (grv - v(:))/1.e4
200    !endif
201
202    ! Compute time for next time step
203    ! -------------------------------
204    firstcall = .false.
205    time = time + dttestphys/daysec
206    if (time > 1.) then
207        time = time - 1.
208        day = day + 1
209    endif
210
211    ! Compute winds and temperature for next time step
212    ! ------------------------------------------------
213    u(:) = u(:) + dttestphys*du(:)
214    v(:) = v(:) + dttestphys*dv(:)
215    temp(:) = temp(:) + dttestphys*dtemp(:)
216
217    ! Compute pressure for next time step
218    ! -----------------------------------
219    psurf = psurf + dttestphys*dpsurf(1) ! surface pressure change
220    plev(:) = ap(:) + psurf*bp(:)
221    play(:) = aps(:) + psurf*bps(:)
222
223    ! Increment tracers
224    q(1,:,:) = q(1,:,:) + dttestphys*dq(1,:,:)
225enddo ! End of time stepping loop (idt=1,ndt)
226
227! Writing the "restart1D.txt" file for the next run
228if (startfiles_1D) call writerestart1D('restart1D.txt',psurf,tsurf(1,:),nlayer,size(tsurf,2),temp,u,v,nq,noms,qsurf(1,:,:),q)
229
230write(*,*) "testphys1d: everything is cool."
231
232END PROGRAM testphys1d
233
234!***********************************************************************
235!***********************************************************************
236! Dummy subroutine used only in 3D, but required to
237! compile testphys1d (to cleanly use writediagfi)
238SUBROUTINE gr_fi_dyn(nfield,ngrid,im,jm,pfi,pdyn)
239
240implicit none
241
242integer                       :: im, jm, ngrid, nfield
243real, dimension(im,jm,nfield) :: pdyn
244real, dimension(ngrid,nfield) :: pfi
245
246if (ngrid /= 1) error stop 'gr_fi_dyn error: in 1D ngrid should be 1!!!'
247
248pdyn(1,1,1:nfield) = pfi(1,1:nfield)
249
250END SUBROUTINE gr_fi_dyn
251
252!***********************************************************************
253!***********************************************************************
Note: See TracBrowser for help on using the repository browser.