source: trunk/LMDZ.PLUTO.old/libf/phypluto/testchange.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: 2.2 KB
Line 
1      subroutine testchange(ngrid,nlayer,nq,igcm1,igcm2,ptimestep, &
2          pplev,zdq,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 ptimestep
46      REAL pplev(ngrid,nlayer+1)
47      character(len=3) :: car1     
48      character(len=7) :: car2     
49
50      ! LOCAL VARIABLES
51      INTEGER l,ig,iq
52      REAL masse
53
54      ! OUTPUT
55      REAL dWtot
56      REAL nconsMAX
57      REAL vdifcncons(ngrid)
58!-----------------------------------------------------------------------
59
60
61      dWtot=0.0
62      nconsMAX=0.0
63      do ig = 1, ngrid
64         vdifcncons(ig)=0.0
65         do l = 1, nlayer
66               masse = (pplev(ig,l) - pplev(ig,l+1))/g
67
68               iq    = igcm1
69               dWtot = dWtot + masse*zdq(ig,l,iq)*ptimestep*area(ig)
70               vdifcncons(ig)=vdifcncons(ig) + masse*zdq(ig,l,iq)
71
72               iq    = igcm2
73               dWtot = dWtot + masse*zdq(ig,l,iq)*ptimestep*area(ig)
74               vdifcncons(ig)=vdifcncons(ig) + masse*zdq(ig,l,iq)
75               
76         enddo
77
78               if(vdifcncons(ig).gt.nconsMAX)then
79                  nconsMAX=vdifcncons(ig)
80               endif
81
82      enddo
83
84      dWtot  = dWtot/totarea
85      print*,'-------------------------------------------'
86      print*,'In ',car2,' atmospheric ',car1,' change=',dWtot,' kg m-2'
87      print*,'--> MAX non-cons factor =',nconsMAX,' kg m-2 s-1'
88 
89      end subroutine testchange
90
Note: See TracBrowser for help on using the repository browser.