source: palm/trunk/UTIL/chemistry/gasphase_preproc/kpp4palm/templates/kp4_compress_subroutines @ 3820

Last change on this file since 3820 was 3458, checked in by kanani, 5 years ago

Reintegrated fixes/changes from branch chemistry

File size: 4.6 KB
RevLine 
[3458]1
2CONTAINS
3! Public subroutines
4
5  SUBROUTINE kco_initialize (npoints_initial, f, IERR, Kacc, Krej)
6    implicit    none
7    integer,intent(IN)                                  :: npoints_initial
8    real(kind=ep),intent(IN),dimension(:,:)             :: f     !Only to get dimension for f_done
9    integer,intent(IN),dimension(:)                     :: IERR
10    integer,intent(IN),dimension(:)                     :: Kacc, Krej
11
12!-- local variables
13    integer                          :: jl
14
15    kpoints      = npoints_initial
16    kpoints_save = kpoints
17
18    allocate(index_org(kpoints))
19    allocate(done_check(kpoints))
20    allocate(index_step(kpoints))
21    allocate(cell_done(kpoints))
22   
23    allocate(f_done(size(f,1),size(f,2)))
24    allocate(Kacc_done(size(Kacc)))
25    allocate(Krej_done(size(Krej)))
26    allocate(ierr_done(size(ierr)))
27
28    cell_done     = .false.
29    done_check    = .false.
30
31    do jl=1,kpoints
32      index_org(jl) = jl
33    end do
34
35    compress_done = .FALSE.
36
37!    write(9,*) 'kco_initialize ',npoints_initial
38
39    return
40  end SUBROUTINE kco_initialize
41
42  SUBROUTINE kco_compress (n_points, T, H, Hnew, IERR, f, RCONST, FIX,     &
43                                   RejectLastH, RejectMoreH, Kacc, Krej)
44    implicit    none
45    integer,intent(INOUT)                       :: n_points
46    real(kind=ep),intent(INOUT),dimension(:)    :: T, H, Hnew
47    INTEGER,intent(INOUT),dimension(:)          :: IERR
48    real(kind=ep),intent(INOUT),dimension(:,:)  :: f
49    real(kind=ep),intent(INOUT),dimension(:,:)  :: RCONST, FIX
50    LOGICAL,intent(INOUT),dimension(:)          :: RejectLastH, RejectMoreH
51    integer,intent(INOUT),dimension(:)          :: Kacc, Krej
52
53!-- local variables
54    integer                          :: i,jl
55    integer                          :: ic
56
57!   Check, if converged variables exist
58
59!    write(9,*) 'kco_compress_1 ',cell_done(1:kpoints)
60
61    if(.not. any(cell_done(1:kpoints)))   return
62
63    compress_done = .TRUE.
64
65!   Move converged cells to original position
66
67    do i=1,size(f,2)
68!CDIR NODEP
69      do jl=1,kpoints
70        if(cell_done(jl))   then
71          f_done(index_org(jl),i) = f(jl,i)
72        end if
73      end do
74    end do
75   
76    do jl=1,kpoints
77      if(cell_done(jl))   then
78        ierr_done(index_org(jl)) = ierr(jl)
79        Kacc_done(index_org(jl)) = Kacc(jl)
80        Krej_done(index_org(jl)) = Krej(jl)
81      end if
82    end do
83
84!CDIR NODEP
85    do jl=1,kpoints
86      if(cell_done(jl))   then
87        done_check(index_org(jl)) = .true.
88      end if
89    end do
90
91!    write(9,*) 'done check ',kpoints,count(done_check)
92       
93!   Compute new index vector
94
95    ic = 0
96    do jl=1,kpoints
97      if(.not.cell_done(jl))   then
98        ic = ic+1
99        index_step(ic) = jl
100      end if
101    end do
102
103!    write(9,*) 'new number of points: ',ic,kpoints,count(done_check)
104
105!   Set new number of Points
106
107    kpoints  = ic
108    n_points = ic
109
110!   Compress arrays
111
112    index_org(1:ic) = index_org(index_step(1:ic))
113
114    do i=1,size(f,2)
115      f(1:ic,i) = f(index_step(1:ic),i)
116    end do
117
118!   compress argument list arrays
119
120    T(1:ic)           = T(index_step(1:ic))
121    IERR(1:ic)        = IERR(index_step(1:ic))
122    H(1:ic)           = H(index_step(1:ic))
123    Hnew(1:ic)        = Hnew(index_step(1:ic))
124
125    do i=1,size(RCONST,2)
126       RCONST(1:ic,i) = RCONST(index_step(1:ic),i)
127    end do
128
129    do i=1,size(FIX,2)
130       FIX(1:ic,i)    = FIX(index_step(1:ic),i)
131    end do
132
133    RejectLastH(1:ic) = RejectLastH(index_step(1:ic))
134    RejectMoreH(1:ic) = RejectMoreH(index_step(1:ic))
135    Kacc(1:ic)        = Kacc(index_step(1:ic))
136    Krej(1:ic)        = Krej(index_step(1:ic))
137    cell_done(1:ic)   = cell_done(index_step(1:ic))
138!    write(9,*) 'kco_compress_3 ',n_points
139
140    return
141  end SUBROUTINE kco_compress
142
143  SUBROUTINE kco_finalize ( f, IERR, Kacc, Krej)
144    implicit    none
145    real(kind=ep),intent(INOUT),dimension(:,:)   :: f
146    integer,      intent(INOUT),dimension(:)     :: IERR, Kacc, Krej
147
148!    write(9,*) 'kco_finalize ',kpoints,kpoints_save
149
150    if(compress_done)  then
151       kpoints = kpoints_save
152
153       f(1:kpoints,:)  = f_done(1:kpoints,:)   ! mz_pj_20070903 (1:kpoints) added
154       ierr(1:kpoints) = ierr_done(1:kpoints)  ! mz_pj_20070903 ...
155       Kacc(1:kpoints) = Kacc_done(1:kpoints)  ! mz_pj_20070903 ...
156       Krej(1:kpoints) = Krej_done(1:kpoints)  ! mz_pj_20070903 ...
157    end if
158
159    deallocate(index_org)
160    deallocate(index_step)
161    deallocate(cell_done)
162    deallocate(done_check)
163
164    deallocate(f_done)
165    deallocate(ierr_done)
166    deallocate(Kacc_done)
167    deallocate(Krej_done)
168
169    return
170  end SUBROUTINE kco_finalize
171
172end module kp4_compress
Note: See TracBrowser for help on using the repository browser.