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

Last change on this file since 3574 was 3574, checked in by jbclement, 32 hours ago

COMMON:
The compilation generates a Fortran subroutine to track compilation and version (SVN or Git) details through the executable file. Put the argument 'version' as an option when executing the code to display these information and create a file "version_diff.txt" containing the diff result.
It can work with every program but it has been implemented only for: 'gcm', 'parallel gcm', 'pem', '1D Mars PCM', 'Mars newstart', 'Mars start2archive' and '1D Generic PCM'.
JBC

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