~ubuntu-branches/ubuntu/warty/petsc/warty

« back to all changes in this revision

Viewing changes to src/ksp/pc/impls/samg/samgtools1.F

  • Committer: Bazaar Package Importer
  • Author(s): Adam C. Powell, IV
  • Date: 2004-06-07 13:41:43 UTC
  • mfrom: (1.1.1 upstream)
  • Revision ID: james.westby@ubuntu.com-20040607134143-92p586zrauvie0le
Tags: 2.2.0-2
* Upstream patch level 2.
* New PETSC_BOPT_EXTRA option for different BOPT and lib names, with _c++
  symlinks only for plain and single (closes: #249617).
* New DEBIAN_DIST=contrib option to link with hypre, parmetis (closes:
  #249619).
* Combined petsc-c and petsc-fortran substvars into petsc-compilers.
* Extra quote in -dev prerm eliminates "too many arguments" problem.

Show diffs side-by-side

added added

removed removed

Lines of Context:
 
1
      subroutine samgpetsc_apply_shift(ia, nnu, ia_shift, ja, nna,       &
 
2
     &                                 ja_shift)
 
3
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
 
4
! This routine applies a ia_shift to all elements in the ia array,! 
 
5
! and a ja_shift to all elements in the ja array.                 ! 
 
6
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
 
7
 
 
8
      implicit none 
 
9
      integer ia(*), nnu, ia_shift, ja(*), nna, ja_shift
 
10
      integer i 
 
11
 
 
12
      do i=1,nnu+1 
 
13
         ia(i) = ia(i) - ia_shift 
 
14
      enddo
 
15
 
 
16
      do i=1,nna
 
17
         ja(i) = ja(i) - ja_shift          
 
18
      enddo 
 
19
 
 
20
      end subroutine samgpetsc_apply_shift
 
21
 
 
22
      subroutine samgpetsc_get_levels(levelscp)
 
23
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
 
24
! Routine to get the numbers of levels created by SAMG            ! 
 
25
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
 
26
     
 
27
      implicit none 
 
28
    
 
29
      integer levelscp, levels 
 
30
 
 
31
      common /samg_mlev/ levels
 
32
 
 
33
      levelscp = levels 
 
34
 
 
35
      return 
 
36
      end subroutine samgpetsc_get_levels
 
37
 
 
38
      subroutine samgpetsc_get_dim_operator(k, nnu, nna)
 
39
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
 
40
! Routine do get number of unknowns and nonzeros of coarse grid   ! 
 
41
! matrix on level k                                               !
 
42
! input:  k:        coarse grid level                             !
 
43
! output: nnu, nna: number of unknowns and nonzeros               !
 
44
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
 
45
 
 
46
! The number of unknowns is determined from the imin and imax arrays. 
 
47
! These arrays are stored in common blocks 
 
48
! The number of nonzeros is determined from the ia array. 
 
49
! This array is stored in the u_wkspace work array. 
 
50
! For information on this latter array, see module samg_wkspace_status
 
51
! in amg_mod.f 
 
52
  
 
53
      use u_wkspace
 
54
 
 
55
      implicit none 
 
56
      integer              imin(25),imax(25),ipmin(25),ipmax(25)
 
57
      integer              levels
 
58
      common /samg_minx/   imin(25),imax(25),ipmin(25),ipmax(25)
 
59
      common /samg_mlev/   levels
 
60
      integer              k, ilo, ihi, n1, n2, nnu, nna
 
61
 
 
62
!     check level input parameter 
 
63
      if (k.lt.2.or.k.gt.levels) then 
 
64
          write(*,*) 'ERROR in samggetdimmat: k out of range'
 
65
          write(*,*) 'Specified value for k on input = ',k
 
66
          write(*,*) 'The current number of levels = ',levels
 
67
          stop 'k should satisfy: 2 <= k <= levels!'
 
68
      endif  
 
69
 
 
70
!     determine number of unknowns 
 
71
      ilo = imin(k)       
 
72
      ihi = imax(k)   
 
73
      nnu = ihi-ilo+1
 
74
!     determine number of nonzeros 
 
75
      n1  = ia(ilo)       
 
76
      n2  = ia(ihi+1)-1
 
77
      nna = n2-n1+1 
 
78
 
 
79
      write(*,*) 'The current level k = ', k 
 
80
      write(*,*) 'The number of levels = ', levels 
 
81
      write(*,*) 'De dimensie van de matrix = ', nnu
 
82
      write(*,*) 'Het aantal nullen van de matrix = ', nna 
 
83
      
 
84
      return
 
85
      end subroutine samgpetsc_get_dim_operator 
 
86
 
 
87
      subroutine samgpetsc_get_operator(k, aout, iaout, jaout)
 
88
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
 
89
! Routine do get coarse grid matrix on level k                    !
 
90
! input:        k:              coarse grid level                 !
 
91
! input/output: aut, ia, jaout: coarse grid matrix in skyline     ! 
 
92
!                               format                            !
 
93
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
 
94
 
 
95
! WARNING: This routine assumes memory allocation for aout, iaout, 
 
96
! and jaout OUTSIDE routine
 
97
 
 
98
! The coarse grid matrices are stored in (a, ia, ja). The array ia is 
 
99
! stored in the u_wkspace work array, while the arrays a and ja are stored 
 
100
! in the a_wkspace work array. See module samg_wkspace_status
 
101
! in amg_mod.f for details. 
 
102
 
 
103
      use u_wkspace
 
104
      use a_wkspace  
 
105
 
 
106
      implicit none 
 
107
      integer               imin(25),imax(25),ipmin(25),ipmax(25)
 
108
      integer               levels
 
109
      common /samg_minx/    imin(25),imax(25),ipmin(25),ipmax(25)
 
110
      common /samg_mlev/    levels
 
111
      double precision      aout(*)
 
112
      integer               iaout(*), jaout(*) 
 
113
      integer               k, ilo, ihi, n1, n2, nnu, nna 
 
114
 
 
115
!     check level input parameter 
 
116
      if (k.lt.2.or.k.gt.levels) then 
 
117
          write(*,*) 'ERROR in samggetmat: k out of range'
 
118
          write(*,*) 'Specified value for k on input = ',k
 
119
          write(*,*) 'The current number of levels = ',levels
 
120
          stop 'k should satisfy: 2 <= k <= levels!'
 
121
      endif  
 
122
 
 
123
!     determine number of unknowns 
 
124
      ilo = imin(k)       
 
125
      ihi = imax(k)   
 
126
      nnu = ihi-ilo+1
 
127
!     determine number of nonzeros 
 
128
      n1  = ia(ilo)       
 
129
      n2  = ia(ihi+1)-1
 
130
      nna = n2-n1+1 
 
131
 
 
132
!     get matrix values 
 
133
      iaout(1:nnu+1)   = ia(ilo:ihi+1)
 
134
      aout(1:nna)      = a(n1:n2)
 
135
      jaout(1:nna)     = ja(n1:n2)  
 
136
 
 
137
      return 
 
138
      end subroutine samgpetsc_get_operator 
 
139
 
 
140
      subroutine samgpetsc_get_dim_interpol(k, nnu, nna)
 
141
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
 
142
! Routine to get size and number of nonzeros of interpolation operator  ! 
 
143
! from grid k+1 (coarse grid) to k (finer grid).                        !
 
144
! input:  k:           fine grid level                                  ! 
 
145
! output: nnu and nna: size and number of nonzeros of interpolation     ! 
 
146
!                      operator                                         ! 
 
147
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
 
148
 
 
149
! The the size is determined from the imin and imax arrays. 
 
150
! These arrays are stored in common blocks 
 
151
! The number of nonzeros is determined from the iw array. 
 
152
! This array is stored in the u_wkspace work array. 
 
153
! For information on this latter array, see module samg_wkspace_status
 
154
! in amg_mod.f 
 
155
 
 
156
      use u_wkspace   
 
157
 
 
158
      implicit none 
 
159
      integer              imin(25),imax(25),ipmin(25),ipmax(25)
 
160
      integer              levels
 
161
      common /samg_minx/   imin(25),imax(25),ipmin(25),ipmax(25)
 
162
      common /samg_mlev/   levels
 
163
      integer              k, nnu, nna, ilo, ihi, n1, n2 
 
164
 
 
165
!     check level input parameter 
 
166
      if (k.lt.1.or.k.gt.levels-1) then 
 
167
          write(*,*) 'ERROR in samggetdimint: k out of range'
 
168
          write(*,*) 'Specified value for k on input = ',k
 
169
          write(*,*) 'The current number of levels = ',levels
 
170
          stop 'k should satisfy: 1 <= k <= levels-1!'
 
171
      endif  
 
172
 
 
173
!     determine size  
 
174
      ilo = imin(k)       
 
175
      ihi = imax(k)   
 
176
      nnu = ihi-ilo+1
 
177
!     determine number of nonzeros 
 
178
      n1  = iw(ilo)       
 
179
      n2  = iw(ihi+1)-1
 
180
      nna = n2-n1+1 
 
181
 
 
182
      write(*,*) 'The current level k = ', k 
 
183
      write(*,*) 'The number of levels = ', levels 
 
184
      write(*,*) 'De dimensie van de interp. op. k -> k+1= ', nnu
 
185
      write(*,*) 'Het aantal nullen van de interp. op = ', nna 
 
186
      
 
187
      return
 
188
      end subroutine samgpetsc_get_dim_interpol
 
189
 
 
190
      subroutine samgpetsc_get_interpol(k, wout, iwout, jwout)
 
191
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
 
192
! Routine do get interpolation operator from grid k+1 (coarse grid)     !
 
193
! to k (finer grid).                                                    !
 
194
! input:  k:                  fine grid level                           ! 
 
195
! output: wout, iwout, jwout: interpolation in skyline format           ! 
 
196
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
 
197
 
 
198
! WARNING: This routine assumes memory allocation for wout, iwout, 
 
199
! and jwout OUTSIDE routine
 
200
 
 
201
! The interpolation opreators are stored in (w, iw, jw). The array iw is 
 
202
! stored in the u_wkspace work array, while the arrays w and jw are stored 
 
203
! in the w_wkspace work array. See module samg_wkspace_status
 
204
! in amg_mod.f for details. 
 
205
 
 
206
      use u_wkspace   
 
207
      use w_wkspace   
 
208
 
 
209
      implicit none 
 
210
      integer              imin(25),imax(25),ipmin(25),ipmax(25)
 
211
      integer              levels
 
212
      common /samg_minx/   imin(25),imax(25),ipmin(25),ipmax(25)
 
213
      common /samg_mlev/   levels
 
214
      double precision     wout(*)
 
215
      integer              iwout(*), jwout(*) 
 
216
      integer              k, ilo, ihi, n1, n2, nnu, nna 
 
217
 
 
218
!     determine number of unknowns 
 
219
      ilo = imin(k)       
 
220
      ihi = imax(k)   
 
221
      nnu = ihi-ilo+1
 
222
!     determine number of nonzeros 
 
223
      n1  = iw(ilo)       
 
224
      n2  = iw(ihi+1)-1
 
225
      nna = n2-n1+1 
 
226
 
 
227
!     get interpolation values 
 
228
      iwout(1:nnu+1)   = iw(ilo:ihi+1)
 
229
      wout(1:nna)      = w(n1:n2)
 
230
      jwout(1:nna)     = jw(n1:n2)  
 
231
 
 
232
      end subroutine samgpetsc_get_interpol