source: trunk/LMDZ.PLUTO/libf/phypluto/testconservfast.F90 @ 3586

Last change on this file since 3586 was 3539, checked in by tbertrand, 8 weeks ago

LMDZ.PLUTO:
Update and Validation of the Volatile Transport Model (VTM, or "nogcm")
Merging missing elements with Pluto.old
TB

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