Line data Source code
1 : !----------------------------------------------------------------------
2 : ! this module computes the total advection tendencies of advected
3 : ! constituents for the finite volume dycore
4 : !----------------------------------------------------------------------
5 : module advect_tend
6 :
7 : use shr_kind_mod, only : r8 => shr_kind_r8
8 :
9 : save
10 : private
11 :
12 : public :: compute_adv_tends_xyz
13 :
14 : real(r8), allocatable :: adv_tendxyz(:,:,:,:,:)
15 :
16 : contains
17 :
18 : !----------------------------------------------------------------------
19 : ! computes the total advective tendencies
20 : ! called twice each time step:
21 : ! - first call sets the initial mixing ratios
22 : ! - second call computes and outputs the tendencies
23 : !----------------------------------------------------------------------
24 29184 : subroutine compute_adv_tends_xyz(elem,fvm,nets,nete,qn0,n0)
25 : use cam_history, only: outfld, hist_fld_active
26 : use time_manager, only: get_step_size
27 : use constituents, only: tottnam,pcnst
28 : use dimensions_mod, only: nc,np,nlev,use_cslam
29 : use element_mod, only: element_t
30 : use fvm_control_volume_mod, only: fvm_struct
31 : implicit none
32 :
33 : type (element_t), intent(in) :: elem(:)
34 : type(fvm_struct), intent(in) :: fvm(:)
35 : integer, intent(in) :: nets,nete,qn0,n0
36 : real(r8) :: dt,idt
37 : integer :: i,j,ic,nx,ie
38 : logical :: init
39 29184 : real(r8), allocatable, dimension(:,:) :: ftmp
40 :
41 29184 : if (use_cslam) then
42 : nx=nc
43 : else
44 0 : nx=np
45 : endif
46 87552 : allocate( ftmp(nx*nx,nlev) )
47 :
48 29184 : init = .false.
49 29184 : if ( .not. allocated( adv_tendxyz ) ) then
50 14592 : init = .true.
51 102144 : allocate( adv_tendxyz(nx,nx,nlev,pcnst,nets:nete) )
52 5090103192 : adv_tendxyz(:,:,:,:,:) = 0._r8
53 : endif
54 :
55 29184 : if (use_cslam) then
56 234384 : do ie=nets,nete
57 8647584 : do ic=1,pcnst
58 10180177200 : adv_tendxyz(:,:,:,ic,ie) = fvm(ie)%c(1:nc,1:nc,:,ic) - adv_tendxyz(:,:,:,ic,ie)
59 : end do
60 : end do
61 : else
62 0 : do ie=nets,nete
63 0 : do ic=1,pcnst
64 0 : adv_tendxyz(:,:,:,ic,ie) = elem(ie)%state%Qdp(:,:,:,ic,qn0)/elem(ie)%state%dp3d(:,:,:,n0) - adv_tendxyz(:,:,:,ic,ie)
65 : enddo
66 : end do
67 : end if
68 :
69 29184 : if ( .not. init ) then
70 14592 : dt = get_step_size()
71 14592 : idt = 1._r8/dt
72 :
73 117192 : do ie=nets,nete
74 4323792 : do ic = 1,pcnst
75 16826400 : do j=1,nx
76 54685800 : do i=1,nx
77 3571403400 : ftmp(i+(j-1)*nx,:) = adv_tendxyz(i,j,:,ic,ie)
78 : end do
79 : end do
80 4309200 : call outfld(tottnam(ic), ftmp,nx*nx, ie)
81 : end do
82 : end do
83 14592 : deallocate(adv_tendxyz)
84 : endif
85 29184 : deallocate(ftmp)
86 29184 : end subroutine compute_adv_tends_xyz
87 :
88 : end module advect_tend
|