source: trunk/LMDZ.PLUTO.old/libf/phypluto/testconserv.F90 @ 3436

Last change on this file since 3436 was 3175, checked in by emillour, 11 months ago

Pluto PCM:
Add the old Pluto LMDZ for reference (required prior step to making
an LMDZ.PLUTO using the same framework as the other physics packages).
TB+EM

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