source: palm/trunk/UTIL/chemistry/gasphase_preproc/kpp4palm/templates/messy_main_tools_kp4_compress.f90 @ 3298

Last change on this file since 3298 was 2696, checked in by kanani, 6 years ago

Merge of branch palm4u into trunk

File size: 5.5 KB
Line 
1module messy_main_tools_kp4_compress
2!
3! Module to compress data for ros_Integrator using variable time step
4! in vector mode. This module is not dependent on chemistry setup and
5! can be use for multiple chemistry setups in one job.
6!
7
8!     COPYRIGHT Klaus Ketelsen and MPI-CH              May 2007
9
10
11  implicit none
12  private
13  save
14
15! KPP DP - Double precision kind
16  INTEGER, PARAMETER :: dp = SELECTED_REAL_KIND(12,307)
17
18
19! variables
20  integer                                   :: kpoints
21  integer                                   :: kpoints_save
22  integer,dimension(:),allocatable,public   :: index_org
23  integer,dimension(:),allocatable          :: index_step
24  logical,dimension(:),allocatable,public   :: cell_done
25  logical,dimension(:),allocatable          :: done_check
26
27! Data arrays
28  real(kind=dp),dimension(:,:),allocatable  :: f_done
29  integer,      dimension(:),  allocatable  :: Kacc_done, Krej_done
30  integer,      dimension(:),  allocatable  :: ierr_done
31
32! interface block
33
34  interface kco_initialize
35    module procedure kco_initialize
36  end interface kco_initialize
37
38  interface kco_compress
39    module procedure kco_compress
40  end interface kco_compress
41
42  interface kco_finalize
43    module procedure kco_finalize
44  end interface kco_finalize
45
46  public kco_initialize,kco_compress
47  public kco_finalize
48
49 contains
50! Public subroutines
51
52  subroutine kco_initialize (npoints_initial, f, IERR, Kacc, Krej)
53    implicit    none
54    integer,intent(IN)                                  :: npoints_initial
55    real(kind=dp),intent(IN),dimension(:,:)             :: f     !Only to get dimension for f_done
56    integer,intent(IN),dimension(:)                     :: IERR
57    integer,intent(IN),dimension(:)                     :: Kacc, Krej
58
59!-- local variables
60    integer                          :: jl
61
62    kpoints      = npoints_initial
63    kpoints_save = kpoints
64
65    allocate(index_org(kpoints))
66    allocate(done_check(kpoints))
67    allocate(index_step(kpoints))
68    allocate(cell_done(kpoints))
69   
70    allocate(f_done(size(f,1),size(f,2)))
71    allocate(Kacc_done(size(Kacc)))
72    allocate(Krej_done(size(Krej)))
73    allocate(ierr_done(size(ierr)))
74
75    cell_done     = .false.
76    done_check    = .false.
77
78    do jl=1,kpoints
79      index_org(jl) = jl
80    end do
81
82    return
83  end subroutine kco_initialize
84
85  subroutine kco_compress (n_points, T, H, Hnew, IERR, f, RCONST, FIX,     &
86                                   RejectLastH, RejectMoreH, Kacc, Krej)
87    implicit    none
88    integer,intent(INOUT)                       :: n_points
89    real(kind=dp),intent(INOUT),dimension(:)    :: T, H, Hnew
90    INTEGER,intent(INOUT),dimension(:)          :: IERR
91    real(kind=dp),intent(INOUT),dimension(:,:)  :: f
92    real(kind=dp),intent(INOUT),dimension(:,:)  :: RCONST, FIX
93    LOGICAL,intent(INOUT),dimension(:)          :: RejectLastH, RejectMoreH
94    integer,intent(INOUT),dimension(:)          :: Kacc, Krej
95
96!-- local variables
97    integer                          :: i,jl
98    integer                          :: ic
99
100!   Check, if converged variables exist
101
102    if(.not. any(cell_done(1:kpoints)))   return
103
104!   Move converged cells to original position
105
106    do i=1,size(f,2)
107!CDIR NODEP
108      do jl=1,kpoints
109        if(cell_done(jl))   then
110          f_done(index_org(jl),i) = f(jl,i) 
111        end if
112      end do
113    end do
114   
115    do jl=1,kpoints
116      if(cell_done(jl))   then
117        ierr_done(index_org(jl)) = ierr(jl) 
118        Kacc_done(index_org(jl)) = Kacc(jl) 
119        Krej_done(index_org(jl)) = Krej(jl) 
120      end if
121    end do
122
123!CDIR NODEP
124    do jl=1,kpoints
125      if(cell_done(jl))   then
126        done_check(index_org(jl)) = .true.
127      end if
128    end do
129!kk    write(6,*) 'done check ',kpoints,count(done_check)
130       
131!   Compute new index vector
132
133    ic = 0
134    do jl=1,kpoints
135      if(.not.cell_done(jl))   then
136        ic = ic+1
137        index_step(ic) = jl
138      end if
139    end do
140
141!kk    write(6,*) 'new number of points: ',ic
142
143!   Set new number of Points
144
145    kpoints  = ic
146    n_points = ic
147
148!   Compress arrays
149
150    index_org(1:ic) = index_org(index_step(1:ic))
151
152    do i=1,size(f,2)
153      f(1:ic,i) = f(index_step(1:ic),i)
154    end do
155
156!   compress argument list arrays
157
158    T(1:ic)           = T(index_step(1:ic))
159    IERR(1:ic)        = IERR(index_step(1:ic))
160    H(1:ic)           = H(index_step(1:ic))
161    Hnew(1:ic)        = Hnew(index_step(1:ic))
162
163    do i=1,size(RCONST,2)
164       RCONST(1:ic,i) = RCONST(index_step(1:ic),i)
165    end do
166
167    do i=1,size(FIX,2)
168       FIX(1:ic,i)    = FIX(index_step(1:ic),i)
169    end do
170
171    RejectLastH(1:ic) = RejectLastH(index_step(1:ic))
172    RejectMoreH(1:ic) = RejectMoreH(index_step(1:ic))
173    Kacc(1:ic)        = Kacc(index_step(1:ic))
174    Krej(1:ic)        = Krej(index_step(1:ic))
175
176    return
177  end subroutine kco_compress
178
179  subroutine kco_finalize ( f, IERR, Kacc, Krej)
180    implicit    none
181    real(kind=dp),intent(INOUT),dimension(:,:)   :: f 
182    integer,      intent(INOUT),dimension(:)     :: IERR, Kacc, Krej
183
184    kpoints = kpoints_save
185
186    f(1:kpoints,:)  = f_done(1:kpoints,:)   ! mz_pj_20070903 (1:kpoints) added
187    ierr(1:kpoints) = ierr_done(1:kpoints)  ! mz_pj_20070903 ...
188    Kacc(1:kpoints) = Kacc_done(1:kpoints)  ! mz_pj_20070903 ...
189    Krej(1:kpoints) = Krej_done(1:kpoints)  ! mz_pj_20070903 ...
190
191    deallocate(index_org)
192    deallocate(index_step)
193    deallocate(cell_done)
194    deallocate(done_check)
195
196    deallocate(f_done)
197    deallocate(ierr_done)
198    deallocate(Kacc_done)
199    deallocate(Krej_done)
200
201    return
202  end subroutine kco_finalize
203
204end module messy_main_tools_kp4_compress
Note: See TracBrowser for help on using the repository browser.