aboutsummaryrefslogtreecommitdiff
path: root/src/tov_lsoda.F
blob: b7a1e9f9c11cd50ae57d1bba5463d7af85c3e895 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
#ifdef HAVE_ODEPACK
c===========================================================
c     History: sode.f
c
c     Driver routine which integrates ODEs defining 
c     TOV star
c 
c    Based on bstar.f by Matthew Choptuik
c===========================================================
      subroutine     tov_lsoda ( rho, p, m, alf,
     &                          r,  nr, tol , n_succ, ssrc) 

      implicit       none


      integer       vsrc,  vsxynt
      real*8         rho2p,         p2rho

      integer        nr
      real*8         rho(nr),p(nr), m(nr), alf(nr)
      real*8         r(nr)
      real*8         tol
      integer        n_succ, ssrc

      character*12    codenm
      parameter    ( codenm = 'tov_lsoda' )

      integer        iargc,         indlnb,         i4arg
      real*8         r8arg

      real*8         r8_never
      parameter    ( r8_never = -1.0d-60 )

c-----------------------------------------------------------
c     Order of system.
c-----------------------------------------------------------
      integer        neq
      parameter    ( neq = 3 )

c-----------------------------------------------------------
c     Storage for solution at requested output radii.
c-----------------------------------------------------------
      integer        maxout
      parameter    ( maxout = 10 000 )

      real*8         y0(neq)
      real*8         vxout(maxout), vyout(maxout,neq),
     &               work(maxout)
      integer        nxout,         ixout,       nxout_succ

      integer        ieq

      logical        ltrace
      parameter    ( ltrace = .false. )

      integer        maxdump
      parameter    ( maxdump = 50 )

      logical        lsodatrace
      parameter    ( lsodatrace = .false. )


c-----------------------------------------------------------
c     LSODA Variables.
c-----------------------------------------------------------
      external       tov_fcn,        tov_jac
      real*8         y(neq)
      real*8         tbgn,       tend
      integer        itol
      real*8         rtol,       atol
      integer        itask,      istate,     iopt
      integer        lrw

c      parameter    ( lrw = 22 + neq * max(16 ,neq + 9) )
      parameter    ( lrw = 22 + neq * (16) )
      real*8         rwork(lrw)

      integer        liw
      parameter    ( liw = 20 + neq )
      integer        iwork(liw)
      integer        jt

      integer        i

c-----------------------------------------------------------
c     Common communication with routine 'tov_fcn' in 'tov_fcn.f' ...
c-----------------------------------------------------------
c   rho0:  Value of rho at r=0
c-----------------------------------------------------------
      real*8  rho0
      common / com_tov / rho0

c-----------------------------------------------------------
c     Parse arguments.  
c-----------------------------------------------------------
      rho0    = rho(1)

      nxout   = nr
      do ixout = 1, nxout
         vxout(ixout) = r(ixout)
      enddo

c-----------------------------------------------------------
c     Initialize those "boundary" conditions which are 
c     fixed by regularity, or which are arbitrary.
c-----------------------------------------------------------
      y0(1) = rho2p(rho0)
      y0(2) = 0.0d0
      y0(3) = 1.0


c-----------------------------------------------------------
c     Echo command line arguments if "local tracing" 
c     enabled ...
c-----------------------------------------------------------
      if( ltrace ) then
         write(0,*) 'rho0: ', rho0
         write(0,*) 'p(1): ', y0(1)
         write(0,*) 'm(1): ', y0(2)
         write(0,*) 'alf(1): ', y0(3)
      end if

c-----------------------------------------------------------
c     Get output radii 
c-----------------------------------------------------------
      if( ltrace ) then
         if( nxout .le. maxdump ) then
            call dvdump(vxout,nxout,codenm//': output radii',0)
         else 
            call dvdmp1(vxout,work,
     &                  int(1.0d0 * nxout / maxdump + 0.5d0),
     &                  nxout,codenm//': selected output radii',0)
         end if
         write(0,*)
         write(0,*) codenm//': Initial time: ', vxout(1)
         write(0,*)
      end if

c-----------------------------------------------------------
c     Set LSODA parameters ...
c-----------------------------------------------------------
      itol   = 1
      rtol   = tol
      atol   = tol
      itask  = 1 
      iopt   = 0
      jt     = 2

c-----------------------------------------------------------
c     Initialize the solution ...
c-----------------------------------------------------------
      do ieq = 1 , neq
         y(ieq)       = y0(ieq)
         vyout(1,ieq) = y0(ieq)
      end do

c-----------------------------------------------------------
c     Do the integration ...
c-----------------------------------------------------------
      ssrc = 0
      do ixout = 2 , nxout
         istate = 1
c-----------------------------------------------------------
c        Need these temporaries since lsoda overwrites 
c        tend ... 
c-----------------------------------------------------------
100      continue
         tbgn = vxout(ixout-1)
         tend = vxout(ixout)
         call lsoda(tov_fcn,neq,y,tbgn,tend,
     &           itol,rtol,atol,itask,    
     &           istate,iopt,rwork,lrw,iwork,liw,tov_jac,jt)
         if( lsodatrace ) then
            write(0,1000) codenm, ixout, vxout(ixout), 
     &      vxout(ixout+1)
1000        format(' ',a,': Step ',i4,'  t = ',1pe10.3,
     &      ' .. ',1pe10.3)
            write(0,*) codenm//': lsoda reurns ', istate
         end if

c         if (istate .eq. -2) then
c            write(0,*) 'istate = -2, trying again...'
c            itol   = 1
c            rtol   = tol
c            atol   = tol
c            itask  = 1
c            iopt   = 0
c            jt     = 2
c            istate = 3
c            goto 100
c         endif

         if( istate .lt. -1 ) then
            write(0,1500) codenm, istate, ixout, nxout, 
     &                    vxout(ixout-1), vxout(ixout)
1500        format(/' ',a,': Error return ',i2,
     &              ' from integrator LSODA.'/
     &              '       At output time ',i5,' of ',i5/
     &              '       Interval ',1pe11.3,' .. ',
     &             1pe11.3/)
             nxout_succ = ixout - 1
             ssrc = 1
            go to 500

         end if

         do ieq = 1 , neq
            vyout(ixout,ieq) = y(ieq)
         end do

      end do
      nxout_succ = nxout

500   continue

c-----------------------------------------------------------
c     OUTPUT
c-----------------------------------------------------------
      do ixout = 1 , nxout_succ
         p(ixout)    = vyout(ixout, 1)
         m(ixout)    = vyout(ixout, 2)
         alf(ixout)    = vyout(ixout, 3)
         rho(ixout)    = p2rho(p(ixout))
      end do

      n_succ = nxout_succ

      return

      end





c===========================================================
c     Some double precision vector utility routines. 
c
c     Originally from ~matt/vutil/dveclib.f
c===========================================================
 
c-----------------------------------------------------------
c     Dumps vector labelled with LABEL on UNIT.
c-----------------------------------------------------------
      subroutine dvdump(v,n,label,unit)
         implicit       none

         real*8         v(*)
         character*(*)  label
         integer        i, n, st, unit
 
         if( n .lt. 1 ) go to 130
            write(unit,100) label
 100        format(/' <<< ',A,' >>>'/)
            st = 1
 110        continue
               write(unit,120) ( v(i) , i = st , min(st+3,n))
 120           FORMAT(' ',4(1PE19.10))
               st = st + 4
            if( st .le. n ) go to 110
 
 130     continue
 
         return
      end
 
c-----------------------------------------------------------
c     Extension of dvdump which dumps every 'inc'th element.
c-----------------------------------------------------------
      subroutine dvdmp1(v,w,inc,n,label,unit)
         implicit       none

         real*8         v(*),    w(*)
         character*(*)  label
         integer        inc,     n,      unit
 
         call dvinj(v,w,inc,n)
         call dvdump(w,1+(n-1)/inc,label,unit)

         return

      end
 
c---------------------------------------------------------
c     Injects every incth element of v1 into v2
c---------------------------------------------------------
      subroutine dvinj(v1,v2,inc,n)
         implicit    none
   
         real*8      v1(*),     v2(*)
         integer     i,         inc,     j,      n
 
         j = 1
         do i = 1 , n , inc
            v2(j) = v1(i)
            j = j + 1
         end do

 
         return
      end

#else
      subroutine tov_lsoda()
      return 
      end
 
#endif