summaryrefslogtreecommitdiffstats
path: root/ast/ast_tester/testframeset.f
blob: 0c9c52eaee6ad65237b935e753dd570d01af3ece (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
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
      program testframeset
      implicit none

      include 'AST_PAR'
      include 'AST_ERR'
      include 'SAE_PAR'

      integer status, pfrm, ffrm, p2fmap, fs, p2fmap2, result, orig
      double precision ina(2), inb(2), outa(2), outb(2), xout, yout
      character text*100

c      call ast_watchmemory(100)


      status = sai__ok
      call err_mark( status )
      call ast_begin( status )


      pfrm = ast_frame( 2, "Domain=PIXEL", status )
      ffrm = ast_frame( 2, "Domain=FPLANE", status )

      ina( 1 ) = 1.0
      ina( 2 ) = 1.0
      inb( 1 ) = 100.0
      inb( 2 ) = 200.0

      outa( 1 ) = -2.5
      outa( 2 ) = -1.0
      outb( 1 ) = 2.5
      outb( 2 ) = 1.0
      p2fmap = ast_winmap( 2, ina, inb, outa, outb, ' ', status )

      fs = ast_frameset( pfrm, ' ', status )
      call ast_addframe( fs, AST__CURRENT, p2fmap, ffrm, status )

      call ast_setc( fs, 'Base', 'Fplane', status )
      if( ast_geti( fs, 'Base', status ) .ne. 2 )
     :    call stopit( status, 'Error -3' )

      call ast_setc( fs, 'Base', 'pixel', status )
      if( ast_geti( fs, 'Base', status ) .ne. 1 )
     :    call stopit( status, 'Error -2' )

      call ast_setc( fs, 'Current', 'PIXEL', status )
      if( ast_geti( fs, 'Current', status ) .ne. 1 )
     :    call stopit( status, 'Error -1' )

      call ast_setc( fs, 'Current', 'fplane', status )
      if( ast_geti( fs, 'Current', status ) .ne. 2 )
     :    call stopit( status, 'Error 0' )

      text = ast_getc( fs, 'AllVariants', status )
      if( text .ne. 'FPLANE' ) call stopit( status, 'Error 1' )

      text = ast_getc( fs, 'Variant', status )
      if( text .ne. 'FPLANE' ) call stopit( status, 'Error 2' )

      if( ast_test( fs, 'Variant', status ) ) call stopit( status,
     :                                                     'Error 3' )


      call ast_addvariant( FS, ast__null, 'FP1', status )

      text = ast_getc( fs, 'AllVariants', status )
      if( text .ne. 'FP1' ) call stopit( status, 'Error 4' )

      text = ast_getc( fs, 'Variant', status )
      if( text .ne. 'FP1' ) call stopit( status, 'Error 5' )

      if( .not. ast_test( fs, 'Variant', status ) ) call stopit( status,
     :                                                       'Error 6' )

      call ast_clear( fs, 'Variant', status )

      text = ast_getc( fs, 'AllVariants', status )
      if( text .ne. 'FPLANE' ) call stopit( status, 'Error 7' )

      text = ast_getc( fs, 'Variant', status )
      if( text .ne. 'FPLANE' ) call stopit( status, 'Error 8' )

      if( ast_test( fs, 'Variant', status ) ) call stopit( status,
     :                                                     'Error 9' )

      call ast_addvariant( FS, ast__null, 'FP1', status )

      outa( 1 ) = 100.0
      outa( 2 ) = 100.0
      outb( 1 ) = 200.0
      outb( 2 ) = 200.0
      p2fmap2 = ast_winmap( 2, ina, inb, outa, outb, ' ', status )

      call ast_invert( p2fmap, status )
      call ast_addvariant( fs, ast_simplify(
     :                     ast_cmpmap( p2fmap, p2fmap2, 1, ' ',
     :                                 status ), status ),
     :                     'FP2', status )

      text = ast_getc( fs, 'AllVariants', status )
      if( text .ne. 'FP1 FP2' ) call stopit( status, 'Error 10' )

      text = ast_getc( fs, 'Variant', status )
      if( text .ne. 'FP2' ) call stopit( status, 'Error 11' )

      if( .not. ast_test( fs, 'Variant', status ) ) call stopit( status,
     :                                                      'Error 12' )
      call ast_tran2( fs, 1, 50.5D0, 100.5D0, .TRUE., xout, yout,
     :                status )
      if( abs( xout - 150.0D0 ) .gt. 1.0E-6 .OR.
     :    abs( yout - 150.0D0 ) .gt. 1.0E-6 )  call stopit( status,
     :                                                      'Error 13' )

      call ast_setc( fs, 'Variant', 'FP1', status )
      call ast_tran2( fs, 1, 50.5D0, 100.5D0, .TRUE., xout, yout,
     :                status )
      if( abs( xout - 0.0D0 ) .gt. 1.0E-6 .OR.
     :    abs( yout - 0.0D0 ) .gt. 1.0E-6 )  call stopit( status,
     :                                                      'Error 14' )

      outa( 1 ) = -100.0
      outa( 2 ) = -100.0
      outb( 1 ) = -200.0
      outb( 2 ) = -200.0
      p2fmap2 = ast_winmap( 2, ina, inb, outa, outb, ' ', status )

      p2fmap = ast_getmapping( fs, AST__CURRENT, AST__BASE, status )
      call ast_addvariant( fs, ast_simplify(
     :                     ast_cmpmap( p2fmap, p2fmap2, 1, ' ',
     :                                 status ), status ),
     :                     'FP3', status )

      text = ast_getc( fs, 'AllVariants', status )
      if( text .ne. 'FP1 FP2 FP3' ) call stopit( status, 'Error 15' )

      text = ast_getc( fs, 'Variant', status )
      if( text .ne. 'FP3' ) call stopit( status, 'Error 16' )

      if( .not. ast_test( fs, 'Variant', status ) ) call stopit( status,
     :                                                      'Error 17' )
      call ast_tran2( fs, 1, 50.5D0, 100.5D0, .TRUE., xout, yout,
     :                status )
      if( abs( xout + 150.0D0 ) .gt. 1.0E-6 .OR.
     :    abs( yout + 150.0D0 ) .gt. 1.0E-6 )  call stopit( status,
     :                                                      'Error 18' )

      call ast_setc( fs, 'Variant', 'FP2', status )
      call ast_tran2( fs, 1, 50.5D0, 100.5D0, .TRUE., xout, yout,
     :                status )
      if( abs( xout - 150.0D0 ) .gt. 1.0E-6 .OR.
     :    abs( yout - 150.0D0 ) .gt. 1.0E-6 )  call stopit( status,
     :                                                      'Error 19' )

      call checkdump( fs, result, status )

      text = ast_getc( result, 'AllVariants', status )
      if( text .ne. 'FP1 FP2 FP3' ) call stopit( status, 'Error 20' )

      text = ast_getc( result, 'Variant', status )
      if( text .ne. 'FP2' ) call stopit( status, 'Error 21' )

      call ast_tran2( result, 1, 50.5D0, 100.5D0, .TRUE., xout, yout,
     :                status )
      if( abs( xout - 150.0D0 ) .gt. 1.0E-6 .OR.
     :    abs( yout - 150.0D0 ) .gt. 1.0E-6 )  call stopit( status,
     :                                                      'Error 22' )



      orig = ast_geti( fs, 'Current', status )
      call ast_addframe( fs, AST__CURRENT, AST_UNITMAP( 2, '', status ),
     :                   AST_FRAME( 2, "Domain=DSB", status ), status )
      call ast_tran2( fs, 1, 50.5D0, 100.5D0, .TRUE., xout, yout,
     :                status )
      if( abs( xout - 150.0D0 ) .gt. 1.0E-6 .OR.
     :    abs( yout - 150.0D0 ) .gt. 1.0E-6 )  call stopit( status,
     :                                                      'Error 23' )

      if( status .eq. sai__ok ) then
         call ast_setc( fs, 'Variant', 'FP1', status )
         if( status .eq. ast__attin ) then
            call err_annul( status )
         else
            call err_flush( status )
            call stopit( status, 'Error 24' )
         end if
      end if

      text = ast_getc( fs, 'AllVariants', status )
      if( text .ne. 'DSB' ) call stopit( status, 'Error 25' )

      text = ast_getc( fs, 'Variant', status )
      if( text .ne. 'DSB' ) call stopit( status, 'Error 26' )

      if( ast_test( fs, 'Variant', status ) ) call stopit( status,
     :                                                     'Error 27' )

      call ast_mirrorvariants( fs, orig, status )

      text = ast_getc( fs, 'AllVariants', status )
      if( text .ne. 'FP1 FP2 FP3' ) call stopit( status, 'Error 28' )

      text = ast_getc( fs, 'Variant', status )
      if( text .ne. 'FP2' ) call stopit( status, 'Error 29' )

      if( .not. ast_test( fs, 'Variant', status ) ) call stopit( status,
     :                                                     'Error 30' )

      call ast_tran2( fs, 1, 50.5D0, 100.5D0, .TRUE., xout, yout,
     :                status )
      if( abs( xout - 150.0D0 ) .gt. 1.0E-6 .OR.
     :    abs( yout - 150.0D0 ) .gt. 1.0E-6 )  call stopit( status,
     :                                                      'Error 31' )

      call ast_set( fs, 'Variant=FP1', status )
      text = ast_getc( fs, 'Variant', status )
      if( text .ne. 'FP1' ) call stopit( status, 'Error 32' )

      call ast_tran2( fs, 1, 50.5D0, 100.5D0, .TRUE., xout, yout,
     :                status )
      if( abs( xout - 0.0D0 ) .gt. 1.0E-6 .OR.
     :    abs( yout - 0.0D0 ) .gt. 1.0E-6 )  call stopit( status,
     :                                                      'Error 33' )

      call checkdump( fs, result, status )

      text = ast_getc( result, 'AllVariants', status )
      if( text .ne. 'FP1 FP2 FP3' ) call stopit( status, 'Error 34' )

      text = ast_getc( result, 'Variant', status )
      if( text .ne. 'FP1' ) call stopit( status, 'Error 35' )

      call ast_tran2( result, 1, 50.5D0, 100.5D0, .TRUE., xout, yout,
     :                status )
      if( abs( xout - 0.0D0 ) .gt. 1.0E-6 .OR.
     :    abs( yout - 0.0D0 ) .gt. 1.0E-6 )  call stopit( status,
     :                                                      'Error 36' )

      result = ast_copy( fs, status )

      text = ast_getc( result, 'AllVariants', status )
      if( text .ne. 'FP1 FP2 FP3' ) call stopit( status, 'Error 37' )

      text = ast_getc( result, 'Variant', status )
      if( text .ne. 'FP1' ) call stopit( status, 'Error 38' )

      call ast_tran2( result, 1, 50.5D0, 100.5D0, .TRUE., xout, yout,
     :                status )
      if( abs( xout - 0.0D0 ) .gt. 1.0E-6 .OR.
     :    abs( yout - 0.0D0 ) .gt. 1.0E-6 )  call stopit( status,
     :                                                      'Error 39' )

      call ast_clear( fs, 'Variant', status )
      text = ast_getc( fs, 'Variant', status )
      if( text .ne. 'DSB' ) call stopit( status, 'Error 40' )






      call ast_end( status )
      call err_rlse( status )

      call ast_activememory( 'testframeset' )
      call ast_flushmemory( 1 );

      if( status .eq. sai__ok ) then
         write(*,*) 'All FrameSet tests passed'
      else
         write(*,*) 'FrameSet tests failed'
      end if

      end



      subroutine stopit( status, text )
      implicit none
      include 'SAE_PAR'
      integer status
      character text*(*)

      if( status .ne. sai__ok ) return
      status = sai__error
      write(*,*) text

      end

      subroutine checkdump( obj, result, status )
      implicit none
      include 'SAE_PAR'
      include 'AST_PAR'
      character key*30,txt1*50,txt2*50
      integer obj, status, next, end, ch, result, ll, overlap, size,
     :        i, type,obj1,obj2,l1,l2,nl,nrow,nrowold
      external mysource, mysink
      character buf*400000

      common /ss1/ buf
      common /ss2/ next, end, ll, nl

      if( status .ne. sai__ok ) return

      ch = ast_channel( mysource, mysink, ' ', status )

      nl = 0
      ll = 110
      next = 1
      if( ast_write( ch, obj, status ) .ne.1 ) then
         call stopit( status, 'Cannot write supplied object to '//
     :                'channel' )
      end if

      next = 1
      nl = 0
      result = ast_read( ch, status )

      if( result .eq. ast__null ) then
         call stopit( status, 'Cannot read object from channel' )
      end if

      if( .not. ast_isaframeset( result, status ) ) then
         call stopit( status, 'Object read from channel is not a '//
     :                'FrameSet')
      end if

      end

      subroutine mysource( status )
      implicit none
      include 'SAE_PAR'
      include 'AST_PAR'
      integer status, next, end, ll, nl
      character buf*400000

      common /ss1/ buf
      common /ss2/ next, end, ll,nl

      if( status .ne. sai__ok ) return

      if( next .ge. end ) then
         call ast_putline( buf, -1, status )
      else

c         write(*,*) buf( next : next + ll - 1 )
         call ast_putline( buf( next : ), ll, status )
         nl = nl + 1
      endif

      next = next + ll

      end

      subroutine mysink( status )
      implicit none
      include 'SAE_PAR'
      include 'AST_PAR'
      integer status, next, end, f, l, ll, nl
      character buf*400000
      character line*1000

      common /ss1/ buf
      common /ss2/ next, end, ll, nl

      if( status .ne. sai__ok ) return

      line = ' '
      call ast_getline( line, l, status )
      call chr_fandl( line( : l ), f, l )
      buf( next : ) = line( f : l )
      l = l - f + 1

      if( next + ll - 1 .ge. 400000 ) then
         write(*,*)
         call stopit( status, 'Buffer overflow in mysink!!' )
      else if( l .gt. ll ) then
         write(*,*)
         write(*,*) buf( next : next + l)
         write(*,*) 'Line length ',l
         call stopit( status, 'Line overflow in mysink!!' )
      else
         end = next + l
         buf( end : next + ll - 1 ) = ' '
         nl = nl + 1
      endif

      next = next + ll

      end