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

Last change on this file since 3130 was 3130, checked in by jbclement, 14 months ago

PEM:
The perennial co2 ice is now taken into account with co2 frost (qsurf) to compute the tendency and to make the update + Rework of how co2 frost is converted to perennial co2 ice at the end of the PEM run + Correction of the value of 'threshold_co2_frost2perennial' to correspond to 10 m + Perennial co2 ice is now handled outside 'paleoclimate' in "phyetat0_mod.F90" of the Mars PCM + Some cleanings.

/!\ Commit for the PEM management of co2 ice before a rework of ice management in the PEM!
JBC

File size: 12.1 KB
Line 
1PROGRAM testphys1d
2
3use comsoil_h,           only: inertiedat, inertiesoil, nsoilmx, tsoil, nqsoil, qsoil
4use surfdat_h,           only: albedodat, perennial_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 ioipsl_getincom,     only: getin ! To use 'getin'
18use init_testphys1d_mod, only: init_testphys1d
19use writerestart1D_mod,  only: writerestart1D
20! Mostly for XIOS outputs:
21use mod_const_mpi,       only: init_const_mpi
22use parallel_lmdz,       only: init_parallel
23
24implicit none
25
26!=======================================================================
27!   subject:
28!   --------
29!   PROGRAM useful to run physical part of the martian GCM in a 1D column
30!
31! Can be compiled with a command like (e.g. for 25 layers)
32!   "makegcm -p mars -d 25 testphys1d"
33! It requires the files "testphys1d.def" "callphys.def"
34!   and a 'run.def' file (containing a "INCLUDEDEF=callphys.def" line)
35!   and a file describing the sigma layers (e.g. "z2sig.def")
36!
37!   author: Frederic Hourdin, R. Fournier, F. Forget
38!   -------
39!
40!   update: 12/06/2003, including chemistry (S. Lebonnois)
41!                       and water ice (F. Montmessin)
42!           27/09/2023, upgrade to F90 (JB Clément)
43!
44!=======================================================================
45
46include "dimensions.h"
47!#include "dimradmars.h"
48!#include "comgeomfi.h"
49!#include "surfdat.h"
50!#include "slope.h"
51!#include "comsoil.h"
52!#include "comdiurn.h"
53include "callkeys.h"
54!#include "comsaison.h"
55!#include "control.h"
56include "netcdf.inc"
57!#include "advtrac.h"
58
59!--------------------------------------------------------------
60! Declarations
61!--------------------------------------------------------------
62integer, parameter                  :: ngrid = 1     ! (2+(jjm-1)*iim - 1/jjm)
63integer, parameter                  :: nlayer = llm
64real, parameter                     :: odpref = 610. ! DOD reference pressure (Pa)
65integer                             :: unitstart     ! unite d'ecriture de "startfi"
66integer                             :: ndt, ilayer, ilevel, isoil, idt, iq
67logical                             :: firstcall, lastcall
68integer                             :: day0          ! initial (sol ; =0 at Ls=0)
69real                                :: day           ! date during the run
70real                                :: time          ! time (0<time<1 ; time=0.5 a midi)
71real                                :: dttestphys    ! testphys1d timestep
72real, dimension(nlayer)             :: play          ! Pressure at the middle of the layers (Pa)
73real, dimension(nlayer + 1)         :: plev          ! intermediate pressure levels (pa)
74real                                :: psurf         ! Surface pressure
75real, dimension(nlayer)             :: u, v          ! zonal, meridional wind
76real                                :: gru, grv      ! prescribed "geostrophic" background wind
77real, dimension(nlayer)             :: temp          ! temperature at the middle of the layers
78real, dimension(:,:,:), allocatable :: q             ! tracer mixing ratio (e.g. kg/kg)
79real, dimension(1)                  :: wstar = 0.    ! Thermals vertical velocity
80
81! Physical and dynamical tandencies (e.g. m.s-2, K/s, Pa/s)
82real, dimension(nlayer)             :: du, dv, dtemp, dudyn, dvdyn, dtempdyn
83real, dimension(1)                  :: dpsurf
84real, dimension(:,:,:), allocatable :: dq, dqdyn
85
86! Various intermediate variables
87integer                 :: aslun
88real                    :: zls, pks, ptif, qtotinit, qtot
89real, dimension(nlayer) :: phi, h, s, w
90integer                 :: nq = 1 ! number of tracers
91real, dimension(1)      :: latitude, longitude, cell_area
92character(2)            :: str2
93character(7)            :: str7
94character(44)           :: txt
95
96! RV & JBC: Use of starting files for 1D
97logical :: startfiles_1D, therestart1D, therestartfi, there
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
117startfiles_1D = .false.
118!------------------------------------------------------
119! Loading run parameters from "run.def" file
120!------------------------------------------------------
121! check if 'run.def' file is around. Otherwise reading parameters
122! from callphys.def via getin() routine won't work.
123inquire(file = 'run.def',exist = there)
124if (.not. there) then
125    write(*,*) 'Cannot find required file "run.def"'
126    write(*,*) '  (which should contain some input parameters along with the following line: INCLUDEDEF=callphys.def)'
127    write(*,*) ' ... might as well stop here ...'
128    error stop
129endif
130
131startfiles_1D = .false. ! Default value
132write(*,*) 'Do you want to use starting files and/or to write restarting files?'
133call getin('startfiles_1D',startfiles_1D)
134write(*,*) 'startfiles_1D =', startfiles_1D
135
136therestart1D = .false. ! Default value
137inquire(file = 'start1D.txt',exist = therestart1D)
138if (startfiles_1D .and. .not. therestart1D) then
139    write(*,*) 'There is no "start1D.txt" file!'
140    write(*,*) 'Initialization is done with default values.'
141endif
142therestartfi = .false. ! Default value
143inquire(file = 'startfi.nc',exist = therestartfi)
144if (.not. therestartfi) then
145    write(*,*) 'There is no "startfi.nc" file!'
146    write(*,*) 'Initialization is done with default values.'
147endif
148
149call init_testphys1d('start1D.txt','startfi.nc',startfiles_1D,therestart1D,therestartfi,ngrid,nlayer,odpref, &
150                     nq,q,time,psurf,u,v,temp,ndt,ptif,pks,dttestphys,zqsat,dq,dqdyn,day0,day,gru,grv,w,     &
151                     play,plev,latitude,longitude,cell_area,atm_wat_profile,atm_wat_tau)
152
153! Write a "startfi" file
154! ----------------------
155! This file will be read during the first call to "physiq".
156! It is needed to transfert physics variables to "physiq"...
157if (.not. therestartfi) then
158    call physdem0("startfi.nc",longitude,latitude,nsoilmx,ngrid, &
159                  llm,nq,dttestphys,float(day0),0.,cell_area,    &
160                  albedodat,inertiedat,def_slope,subslope_dist)
161    call physdem1("startfi.nc",nsoilmx,ngrid,llm,nq,nqsoil,dttestphys,time,      &
162                  tsurf,tsoil,inertiesoil,albedo,emis,q2,qsurf,qsoil,tauscaling, &
163                  totcloudfrac,wstar,watercap,perennial_co2ice)
164endif !(.not. therestartfi)
165
166!=======================================================================
167!  1D MODEL TIME STEPPING LOOP
168!=======================================================================
169firstcall = .true.
170lastcall = .false.
171do idt = 1,ndt
172    if (idt == ndt) lastcall = .true.
173!    if (idt == ndt - day_step - 1) then ! test
174!        lastcall = .true.
175!        call solarlong(day*1.0,zls)
176!        write(103,*) 'Ls=',zls*180./pi
177!        write(103,*) 'Lat=', latitude(1)*180./pi
178!        write(103,*) 'Tau=', tauvis/odpref*psurf
179!        write(103,*) 'RunEnd - Atmos. Temp. File'
180!        write(103,*) 'RunEnd - Atmos. Temp. File'
181!        write(104,*) 'Ls=',zls*180./pi
182!        write(104,*) 'Lat=', latitude(1)
183!        write(104,*) 'Tau=', tauvis/odpref*psurf
184!        write(104,*) 'RunEnd - Atmos. Temp. File'
185!    endif
186
187    ! Compute geopotential
188    ! ~~~~~~~~~~~~~~~~~~~~
189    s(:) = (aps(:)/psurf + bps(:))**rcp
190    h(:) = cpp*temp(:)/(pks*s(:))
191
192    phi(1) = pks*h(1)*(1. - s(1))
193    do ilayer = 2,nlayer
194        phi(ilayer) = phi(ilayer - 1) + pks*(h(ilayer - 1) + h(ilayer))*.5*(s(ilayer - 1) - s(ilayer))
195    enddo
196
197    ! Modify atmospheric water if asked
198    ! Added "atm_wat_profile" flag (JN + JBC)
199    ! ---------------------------------------
200    if (water) then
201        call watersat(nlayer,temp,play,zqsat)
202        if (0. <= atm_wat_profile .and. atm_wat_profile <= 1.) then
203        ! If atmospheric water is monitored
204            if (atm_wat_tau < 0.) then ! Prescribed atm_wat_profile: wet if >0, dry if =0
205                q(1,:,igcm_h2o_vap) = min(zqsat(:),atm_wat_profile*g/psurf)
206                q(1,:,igcm_h2o_ice) = 0. ! reset h2o ice
207            else ! Relaxation towards the value atm_wat_profile with relaxation time atm_wat_tau
208                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)
209                q(1,:,igcm_h2o_vap) = min(zqsat(:),q(1,:,igcm_h2o_vap))
210                q(1,:,igcm_h2o_ice) = 0. ! reset h2o ice
211            endif
212        endif
213    !EV profile
214!            IF(atm_wat_profile.eq.2) THEN
215!             DO ilayer=1,16
216!              q(1,ilayer,igcm_h2o_vap)=min(zqsat(ilayer),6e-5)
217!             ENDDO! ilayer=1,16
218!             DO ilayer=17,22
219!              q(1,ilayer,igcm_h2o_vap)=min(zqsat(ilayer),3e-5)
220!             ENDDO! ilayer=17,22
221!             DO ilayer=23,nlayer
222!              q(1,ilayer,igcm_h2o_vap)=0
223!             ENDDO! ilayer=23,nlayer
224!             endif
225     endif
226
227    ! Call physics
228    ! --------------------
229    call physiq(1,llm,nq,firstcall,lastcall,day,time,dttestphys,plev,play,phi,u,v,temp,q,w,du,dv,dtemp,dq,dpsurf)
230    !                                                                                      ^----- outputs -----^
231
232    ! Wind increment: specific for 1D
233    ! --------------------------------
234    ! The physics compute the tendencies on u and v,
235    ! here we just add Coriolos effect
236    !do ilayer = 1,nlayer
237    !    du(ilayer) = du(ilayer) + ptif*(v(ilayer) - grv)
238    !    dv(ilayer) = dv(ilayer) + ptif*(-u(ilayer) + gru)
239    !enddo
240    ! For some tests: No coriolis force at equator
241    !if (latitude(1) == 0.) then
242    du(:) = du(:) + (gru - u(:))/1.e4
243    dv(:) = dv(:) + (grv - v(:))/1.e4
244    !endif
245
246    ! Compute time for next time step
247    ! -------------------------------
248    firstcall = .false.
249    time = time + dttestphys/daysec
250    if (time > 1.) then
251        time = time - 1.
252        day = day + 1
253    endif
254
255    ! Compute winds and temperature for next time step
256    ! ------------------------------------------------
257    u(:) = u(:) + dttestphys*du(:)
258    v(:) = v(:) + dttestphys*dv(:)
259    temp(:) = temp(:) + dttestphys*dtemp(:)
260
261    ! Compute pressure for next time step
262    ! -----------------------------------
263    psurf = psurf + dttestphys*dpsurf(1) ! surface pressure change
264    plev(:) = ap(:) + psurf*bp(:)
265    play(:) = aps(:) + psurf*bps(:)
266
267    ! Increment tracers
268    q(1,:,:) = q(1,:,:) + dttestphys*dq(1,:,:)
269enddo ! End of time stepping loop (idt=1,ndt)
270
271! Writing the "restart1D.txt" file for the next run
272if (startfiles_1D) call writerestart1D('restart1D.txt',psurf,tsurf(1,:),nlayer,size(tsurf,2),temp,u,v,nq,noms,qsurf(1,:,:),q)
273
274write(*,*) "testphys1d: everything is cool."
275
276END PROGRAM testphys1d
277
278!***********************************************************************
279!***********************************************************************
280! Dummy subroutine used only in 3D, but required to
281! compile testphys1d (to cleanly use writediagfi)
282SUBROUTINE gr_fi_dyn(nfield,ngrid,im,jm,pfi,pdyn)
283
284implicit none
285
286integer                       :: im, jm, ngrid, nfield
287real, dimension(im,jm,nfield) :: pdyn
288real, dimension(ngrid,nfield) :: pfi
289
290if (ngrid /= 1) error stop 'gr_fi_dyn error: in 1D ngrid should be 1!!!'
291
292pdyn(1,1,1:nfield) = pfi(1,1:nfield)
293
294END SUBROUTINE gr_fi_dyn
295
296!***********************************************************************
297!***********************************************************************
Note: See TracBrowser for help on using the repository browser.