source: trunk/LMDZ.PLUTO/libf/phypluto/testconserv.F90

Last change on this file was 3412, checked in by afalco, 3 months ago

Pluto PCM: Imported glaciers & conserv mass routines from pluto.old.
AF

File size: 3.6 KB
Line 
1      subroutine testconserv(ngrid,nlayer,nq,igcm1,igcm2,ptimestep, &
2          pplev,zdq,zdqs,car1,car2)
3
4      use comgeomfi_h         
5      use comcstfi_mod, only: pi, g
6      use geometry_mod, only: latitude, longitude, cell_area
7      USE tracer_h, only: igcm_co_ice, igcm_ch4_ice
8      implicit none
9
10!==================================================================
11!     Purpose
12!     -------
13!     Test conservation of tracers
14
15!     Inputs
16!     ------
17!     ngrid                 Number of vertical columns
18!     nlayer                Number of layers
19!     pplev(ngrid,nlayer+1) Pressure levels
20!     zdq
21!     
22!     Outputs
23!     -------
24!     
25!     Both
26!     ----
27!
28!     Authors
29!     -------
30!     Tanguy Bertrand 2015
31!     
32!==================================================================
33
34#include "dimensions.h"
35
36!-----------------------------------------------------------------------
37!     Arguments
38
39      INTEGER ngrid, nlayer, nq
40      INTEGER igcm1, igcm2
41      REAL zdq(ngrid,nlayer,nq)
42      REAL zdqs(ngrid,nq)
43      REAL ptimestep
44      REAL pplev(ngrid,nlayer+1)
45      character(len=3) :: car1     
46      character(len=7) :: car2     
47
48      ! LOCAL VARIABLES
49      INTEGER l,ig,iq
50      REAL masse
51
52      ! OUTPUT
53      REAL dWtot   
54      REAL dWtots
55      REAL nconsMAX
56      INTEGER myig
57      REAL vdifcncons(ngrid)
58!-----------------------------------------------------------------------
59
60      write(*,*) 'TB igcm=',igcm1,igcm2
61      write(*,*) 'TB zdq1=',zdq(1,1,igcm1)
62      dWtot=0.0
63      dWtots=0.0
64      nconsMAX=0.0
65      do ig = 1, ngrid
66         vdifcncons(ig)=0.0
67         do l = 1, nlayer
68             masse = (pplev(ig,l) - pplev(ig,l+1))/g
69
70             iq    = igcm1
71             ! sum of atmospheric mass : kg
72             dWtot = dWtot + masse*zdq(ig,l,iq)*ptimestep*cell_area(ig)
73             ! for each column, total mass lost per sec : kg(tracer) / m2 /s
74             vdifcncons(ig)=vdifcncons(ig) + masse*zdq(ig,l,iq)
75 
76             ! if clouds :
77             IF ((igcm2.eq.igcm_co_ice).or.(igcm2.eq.igcm_ch4_ice)) THEN
78              iq    = igcm2
79              dWtot = dWtot + masse*zdq(ig,l,iq)*ptimestep*cell_area(ig)
80              vdifcncons(ig)=vdifcncons(ig) + masse*zdq(ig,l,iq)
81             ENDIF
82               
83         enddo
84
85             iq     = igcm1
86             dWtots = dWtots + zdqs(ig,iq)*ptimestep*cell_area(ig)
87             vdifcncons(ig)=vdifcncons(ig)+zdqs(ig,iq)
88
89             IF ((igcm2.eq.igcm_co_ice).or.(igcm2.eq.igcm_ch4_ice)) THEN
90              iq     = igcm2
91              dWtots = dWtots + zdqs(ig,iq)*ptimestep*cell_area(ig)
92              vdifcncons(ig)=vdifcncons(ig)+zdqs(ig,iq)
93             ENDIF
94             
95             ! vdifcncons is the total amount of material that appear or
96             ! disapear per second in the routine
97             ! it is the non conservative factor
98
99             if(vdifcncons(ig).gt.nconsMAX)then
100                nconsMAX=vdifcncons(ig)
101                myig=ig
102             endif
103
104      enddo
105      write(*,*) 'TB Dwtot=',dWtots,dWtot
106
107      dWtot  = dWtot/totarea
108      dWtots = dWtots/totarea
109      print*,'-------------------------------------------'
110      print*,'In ',car2,' atmospheric ',car1,' change=',dWtot,' kg m-2'
111      print*,'In ',car2,' surface ',car1,' change  =',dWtots,' kg m-2'
112      print*,'--> non-cons factor     =',dWtot+dWtots,' kg m-2'
113      print*,'--> MAX non-cons factor =',nconsMAX,' kg m-2 s-1'
114      IF (nconsMAX.gt.0.) then
115        print*,'--> obtained at lat/lon=',latitude(myig)*180./pi,longitude(myig)*180./pi
116      ENDIF
117      end subroutine testconserv
118
Note: See TracBrowser for help on using the repository browser.