aboutsummaryrefslogtreecommitdiff
path: root/src/xoper.f
blob: 3ffcfdba33c2967098b9b8d633daf3d9403ec659 (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
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
1767
1768
1769
1770
1771
1772
1773
1774
1775
1776
1777
1778
1779
1780
1781
1782
1783
1784
1785
1786
1787
1788
1789
1790
1791
1792
1793
1794
1795
1796
1797
1798
1799
1800
1801
1802
1803
1804
1805
1806
1807
1808
1809
1810
1811
1812
1813
1814
1815
1816
1817
1818
1819
1820
1821
1822
1823
1824
1825
1826
1827
1828
1829
1830
1831
1832
1833
1834
1835
1836
1837
1838
1839
1840
1841
1842
1843
1844
1845
1846
1847
1848
1849
1850
1851
1852
1853
1854
1855
1856
1857
1858
1859
1860
1861
1862
1863
1864
1865
1866
1867
1868
1869
1870
1871
1872
1873
1874
1875
1876
1877
1878
1879
1880
1881
1882
1883
1884
1885
1886
1887
1888
1889
1890
1891
1892
1893
1894
1895
1896
1897
1898
1899
1900
1901
1902
1903
1904
1905
1906
1907
1908
1909
1910
1911
1912
1913
1914
1915
1916
1917
1918
1919
1920
1921
1922
1923
1924
1925
1926
1927
1928
1929
1930
1931
1932
1933
1934
1935
1936
1937
1938
1939
1940
1941
1942
1943
1944
1945
1946
1947
1948
1949
1950
1951
1952
1953
1954
1955
1956
1957
1958
1959
1960
1961
1962
1963
1964
1965
1966
1967
1968
1969
1970
1971
1972
1973
1974
1975
1976
1977
1978
1979
1980
1981
1982
1983
1984
1985
1986
1987
1988
1989
1990
1991
1992
1993
1994
1995
1996
1997
1998
1999
2000
2001
2002
2003
2004
2005
2006
2007
2008
2009
2010
2011
2012
2013
2014
2015
2016
2017
2018
2019
2020
2021
2022
2023
2024
2025
2026
2027
2028
2029
2030
2031
2032
2033
2034
2035
2036
2037
2038
2039
2040
2041
2042
2043
2044
2045
2046
2047
2048
2049
2050
2051
2052
2053
2054
2055
2056
2057
2058
2059
2060
2061
2062
2063
2064
2065
2066
2067
2068
2069
2070
2071
2072
2073
2074
2075
2076
2077
2078
2079
2080
2081
2082
2083
2084
2085
2086
2087
2088
2089
2090
2091
2092
2093
2094
2095
2096
2097
2098
2099
2100
2101
2102
2103
2104
2105
2106
2107
2108
2109
2110
2111
2112
2113
2114
2115
2116
2117
2118
2119
2120
2121
2122
2123
2124
2125
2126
2127
2128
2129
2130
2131
2132
2133
2134
2135
2136
2137
2138
2139
2140
2141
2142
2143
2144
2145
2146
2147
2148
2149
2150
2151
2152
2153
2154
2155
2156
2157
2158
2159
2160
2161
2162
2163
2164
2165
2166
2167
2168
2169
2170
2171
2172
2173
2174
2175
2176
2177
2178
2179
2180
2181
2182
2183
2184
2185
2186
2187
2188
2189
2190
2191
2192
2193
2194
2195
2196
2197
2198
2199
2200
2201
2202
2203
2204
2205
2206
2207
2208
2209
2210
2211
2212
2213
2214
2215
2216
2217
2218
2219
2220
2221
2222
2223
2224
2225
2226
2227
2228
2229
2230
2231
2232
2233
2234
2235
2236
2237
2238
2239
2240
2241
2242
2243
2244
2245
2246
2247
2248
2249
2250
2251
2252
2253
2254
2255
2256
2257
2258
2259
2260
2261
2262
2263
2264
2265
2266
2267
2268
2269
2270
2271
2272
2273
2274
2275
2276
2277
2278
2279
2280
2281
2282
2283
2284
2285
2286
2287
2288
2289
2290
2291
2292
2293
2294
2295
2296
2297
2298
2299
2300
2301
2302
2303
2304
2305
2306
2307
2308
2309
2310
2311
2312
2313
2314
2315
2316
2317
2318
2319
2320
2321
2322
2323
2324
2325
2326
2327
2328
2329
2330
2331
2332
2333
2334
2335
2336
2337
2338
2339
2340
2341
2342
2343
2344
2345
2346
2347
2348
2349
2350
2351
2352
2353
2354
2355
2356
2357
2358
2359
2360
2361
2362
2363
2364
2365
2366
2367
2368
2369
2370
2371
2372
2373
2374
2375
2376
2377
2378
2379
2380
2381
2382
2383
2384
2385
2386
2387
2388
2389
2390
2391
2392
2393
2394
2395
2396
2397
2398
2399
2400
2401
2402
2403
2404
2405
2406
2407
2408
2409
2410
2411
2412
2413
2414
2415
2416
2417
2418
2419
2420
2421
2422
2423
2424
2425
2426
2427
2428
2429
2430
2431
2432
2433
2434
2435
2436
2437
2438
2439
2440
2441
2442
2443
2444
2445
2446
2447
2448
2449
2450
2451
2452
2453
2454
2455
2456
2457
2458
2459
2460
2461
2462
2463
2464
2465
2466
2467
2468
2469
2470
2471
2472
2473
2474
2475
2476
2477
2478
2479
2480
2481
2482
2483
2484
2485
2486
2487
2488
2489
2490
2491
2492
2493
2494
2495
2496
2497
2498
2499
2500
2501
2502
2503
2504
2505
2506
2507
2508
2509
2510
2511
2512
2513
2514
2515
2516
2517
2518
2519
2520
2521
2522
2523
2524
2525
2526
2527
2528
2529
2530
2531
2532
2533
2534
2535
2536
2537
2538
2539
2540
2541
2542
2543
2544
2545
2546
2547
2548
2549
2550
2551
2552
2553
2554
2555
2556
2557
2558
2559
2560
2561
2562
2563
2564
2565
2566
2567
2568
2569
2570
2571
2572
2573
2574
2575
2576
2577
2578
2579
2580
2581
2582
2583
2584
2585
2586
2587
2588
2589
2590
2591
2592
2593
2594
2595
2596
2597
2598
2599
2600
2601
2602
2603
2604
2605
2606
2607
2608
2609
2610
2611
2612
2613
2614
2615
2616
2617
2618
2619
2620
2621
2622
2623
2624
2625
2626
2627
2628
2629
2630
2631
2632
2633
2634
2635
2636
2637
2638
2639
2640
2641
2642
2643
2644
2645
2646
2647
2648
2649
2650
2651
2652
2653
2654
2655
2656
2657
2658
2659
2660
2661
2662
2663
2664
2665
2666
2667
2668
2669
2670
2671
2672
2673
2674
2675
2676
2677
2678
2679
2680
2681
2682
2683
2684
2685
2686
2687
2688
2689
2690
2691
2692
2693
2694
2695
2696
2697
2698
2699
2700
2701
2702
2703
2704
2705
2706
2707
2708
2709
2710
2711
2712
2713
2714
2715
2716
2717
2718
2719
2720
2721
2722
2723
2724
2725
2726
2727
2728
2729
2730
2731
2732
2733
2734
2735
2736
2737
2738
2739
2740
2741
2742
2743
2744
2745
2746
2747
2748
2749
2750
2751
2752
2753
2754
2755
2756
2757
2758
2759
2760
2761
2762
2763
2764
2765
2766
2767
2768
2769
2770
2771
2772
2773
2774
2775
2776
2777
2778
2779
2780
C***********************************************************************
C    Module:  xoper.f
C 
C    Copyright (C) 2000 Mark Drela 
C 
C    This program is free software; you can redistribute it and/or modify
C    it under the terms of the GNU General Public License as published by
C    the Free Software Foundation; either version 2 of the License, or
C    (at your option) any later version.
C
C    This program is distributed in the hope that it will be useful,
C    but WITHOUT ANY WARRANTY; without even the implied warranty of
C    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
C    GNU General Public License for more details.
C
C    You should have received a copy of the GNU General Public License
C    along with this program; if not, write to the Free Software
C    Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
C***********************************************************************
C
      SUBROUTINE OPER
      INCLUDE 'XFOIL.INC'
      CHARACTER*1 ANS
      CHARACTER*4 COMAND, COMOLD
      LOGICAL LRECALC, LCPX
C
      CHARACTER*128 COMARG, ARGOLD, LINE
C
      PARAMETER (NPRX = 101)
      DIMENSION XPR(NPRX), YPR(NPRX), FPR(NPRX)
C
      DIMENSION NBLP(NPX)
      DIMENSION IPPAI(NPX), NAPOLT(NPX)
C
      DIMENSION IINPUT(20)
      DIMENSION RINPUT(20)
      LOGICAL ERROR
C
C---- retain last-command info if OPER is exited and then re-entered
      SAVE COMOLD, ARGOLD
C
C---- logical units for  polar save file,  polar dump file
      LUPLR = 9
      LUPLX = 11
C
      COMAND = '****'
      COMARG = ' '
      LRECALC = .FALSE.
      LCPX = .FALSE.
C
      IF(N.EQ.0) THEN
       WRITE(*,*)
       WRITE(*,*) '***  No airfoil available  ***'
       RETURN
      ENDIF
C
      IF(IPACT.NE.0) THEN
       WRITE(*,5000) IPACT
 5000  FORMAT(/'  Polar', I3,'  is active')
      ENDIF
C
ccc 500  CONTINUE
      COMOLD = COMAND
      ARGOLD = COMARG
C
C====================================================
C---- start of menu loop
 500  CONTINUE
C
      IF(LVISC) THEN
        IF(LPACC) THEN
         CALL ASKC('.OPERva^',COMAND,COMARG)
        ELSE
         CALL ASKC('.OPERv^',COMAND,COMARG)
        ENDIF
      ELSE
        IF(LPACC) THEN
         CALL ASKC('.OPERia^',COMAND,COMARG)
        ELSE
         CALL ASKC('.OPERi^',COMAND,COMARG)
        ENDIF
      ENDIF
C
C---- process previous command ?
      IF(COMAND(1:1).EQ.'!') THEN
        IF(COMOLD.EQ.'****') THEN
          WRITE(*,*) 'Previous .OPER command not valid'
          GO TO 500
        ELSE
          COMAND = COMOLD
          COMARG = ARGOLD
          LRECALC = .TRUE.
        ENDIF
      ELSE
        LRECALC = .FALSE.
      ENDIF
C
      IF(COMAND.EQ.'    ') THEN
C----- just <return> was typed... clean up plotting and exit OPER
       IF(LPLOT) CALL PLEND
       LPLOT = .FALSE.
       CALL CLRZOOM
       RETURN
      ENDIF
C
C---- extract command line numeric arguments
      DO I=1, 20
        IINPUT(I) = 0
        RINPUT(I) = 0.0
      ENDDO
      NINPUT = 20
      CALL GETFLT(COMARG,RINPUT,NINPUT,ERROR)
C
C---- don't try to read integers, since might get integer overflow
      DO I=1, NINPUT
        IF(ABS(RINPUT(I)) .GT. 2.1E9) THEN
         IINPUT(I) = 2**31
        ELSE
         IINPUT(I) = INT(RINPUT(I))
        ENDIF
      ENDDO
C
ccc      NINPUT = 20
ccc      CALL GETINT(COMARG,IINPUT,NINPUT,ERROR)
C
C--------------------------------------------------------
      IF(COMAND.EQ.'?   ') THEN
       WRITE(*,1050)
 1050  FORMAT(
     & /'   <cr>     Return to Top Level'
     & /'   !        Redo last ALFA,CLI,CL,ASEQ,CSEQ,VELS'
     &//'   Visc r   Toggle Inviscid/Viscous mode'
     & /'  .VPAR     Change BL parameter(s)'
     & /'   Re   r   Change Reynolds number'
     & /'   Mach r   Change Mach number'
     & /'   Type i   Change type of Mach,Re variation with CL'
     & /'   ITER     Change viscous-solution iteration limit'
     & /'   INIT     Toggle BL initialization flag'
     &//'   Alfa r   Prescribe alpha'
     & /'   CLI  r   Prescribe inviscid CL'
     & /'   Cl   r   Prescribe CL'
     & /'   ASeq rrr Prescribe a sequence of alphas'
     & /'   CSeq rrr Prescribe a sequence of CLs'
     &//'   SEQP     Toggle polar/Cp(x) sequence plot display'
     & /'   CINC     Toggle  minimum Cp  inclusion in polar'
     & /'   HINC     Toggle hinge moment inclusion in polar'
     & /'   Pacc i   Toggle auto point accumulation to active polar'
     & /'   PGET f   Read new polar from save file'
     & /'   PWRT i   Write polar to save file'
     & /'   PSUM     Show summary of stored polars'
     & /'   PLIS i   List stored polar(s)'
     & /'   PDEL i   Delete stored polar'
     & /'   PSOR i   Sort stored polar'
     & /'   PPlo ii. Plot stored polar(s)'
     & /'   APlo ii. Plot stored airfoil(s) for each polar'
     & /'   ASET i   Copy stored airfoil into current airfoil'
     & /'   PREM ir. Remove point(s) from stored polar'
     & /'   PNAM i   Change airfoil name of stored polar'
     & /'   PPAX     Change polar plot axis limits'
     &//'   RGET f   Read new reference polar from file'
     & /'   RDEL i   Delete stored reference polar'
     &//'   GRID     Toggle Cp vs x grid overlay'
     & /'   CREF     Toggle reference Cp data overlay'
     & /'   FREF     Toggle reference CL,CD.. data display'
     &//'   CPx      Plot Cp vs x'
     & /'   CPV      Plot airfoil with pressure vectors (gee wiz)'
     & /'  .VPlo     BL variable plots'
     & /'  .ANNO     Annotate current plot'
     & /'   HARD     Hardcopy current plot'
     & /'   SIZE r   Change plot-object size'
     & /'   CPMI r   Change minimum Cp axis annotation'
     &//'   BL   i   Plot boundary layer velocity profiles'
     & /'   BLC      Plot boundary layer velocity profiles at cursor'
     & /'   BLWT r   Change velocity profile scale weight'
     &//'   FMOM     Calculate flap hinge moment and forces'
     & /'   FNEW rr  Set new flap hinge point'
     & /'   VELS rr  Calculate velocity components at a point'
     & /'   DUMP f   Output Ue,Dstar,Theta,Cf vs s,x,y to file'
     & /'   CPWR f   Output x vs Cp to file'
     & /'   CPMN     Report minimum surface Cp'
     & /'   NAME s   Specify new airfoil name'
     & /'   NINC     Increment name version number')
c     &//'   IMAG    Toggle image-airfoil'
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'Z   ') THEN
       CALL USETZOOM(.TRUE.,.TRUE.)
       CALL REPLOT(IDEV)
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'U   ') THEN
       CALL CLRZOOM
       CALL REPLOT(IDEV)
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'VISC' .OR.
     &       COMAND.EQ.'V   '      ) THEN
       IF(LPACC) THEN
        WRITE(*,2100)
        GO TO 500
       ENDIF
C
       LVISC = .NOT. LVISC
C
       IF(LVISC) THEN
         IF(NINPUT.GE.1) THEN
           REINF1 = RINPUT(1)
         ELSE IF(REINF1 .EQ. 0.0) THEN
           CALL ASKR('Enter Reynolds number^',REINF1)
         ENDIF
C
         CALL MRSHOW(.TRUE.,.TRUE.)
       ENDIF
       LVCONV = .FALSE.
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'HARD') THEN
       IF(LPLOT) CALL PLEND
       LPLOT = .FALSE.
       CALL REPLOT(IDEVRP)
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'SIZE') THEN
       IF(NINPUT.GE.1) THEN
        SIZE = RINPUT(1)
       ELSE
        WRITE(*,*) 'Current plot-object size =', SIZE
        CALL ASKR('Enter new plot-object size^',SIZE)
       ENDIF
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'CPMI') THEN
       IF(NINPUT.GE.1) THEN
        CPMIN = RINPUT(1)
       ELSE
        WRITE(*,*) 'Current CPmin =', CPMIN
        CALL ASKR('Enter new CPmin^',CPMIN)
       ENDIF
C
       PFAC = PLOTAR/(CPMAX-CPMIN)
       CPDEL = -0.5
       IF(CPMIN .LT. -4.01) CPDEL = -1.0
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'UEMA') THEN
       IF(NINPUT.GE.1) THEN
        UEMAX = RINPUT(1)
       ELSE
        WRITE(*,*) 'Current UEmax =', UEMAX
        CALL ASKR('Enter new UEMAX^',UEMAX)
       ENDIF
C
       UFAC = PLOTAR/(UEMAX-UEMIN)
       UEDEL = 0.2
       IF((UEMAX-UEMIN) .GT. 2.51) UEDEL = 0.5
       IF((UEMAX-UEMIN) .GT. 5.01) UEDEL = 1.0
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'UEMI') THEN
       IF(NINPUT.GE.1) THEN
        UEMIN = RINPUT(1)
       ELSE
        WRITE(*,*) 'Current UEmin =', UEMIN
        CALL ASKR('Enter new UEMIN^',UEMIN)
       ENDIF
C
       UFAC = PLOTAR/(UEMAX-UEMIN)
       UEDEL = 0.2
       IF((UEMAX-UEMIN) .GT. 2.51) UEDEL = 0.5
       IF((UEMAX-UEMIN) .GT. 5.01) UEDEL = 1.0
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'VPAR') THEN
       CALL VPAR
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'RE  ' .OR.
     &       COMAND.EQ.'R   '      ) THEN
       IF(LPACC .AND. LVISC) THEN
        WRITE(*,2100)
        GO TO 500
       ENDIF
C
       IF(NINPUT.GE.1) THEN
        REINF1 = RINPUT(1)
       ELSE
        WRITE(*,*)
        WRITE(*,*) 'Currently...'
        CALL MRSHOW(.FALSE.,.TRUE.)
        CALL ASKR('Enter new Reynolds number^',REINF1)
       ENDIF
C
ccc    CALL MRSHOW(.FALSE.,.TRUE.)
       CALL MRCL(1.0,MINF_CL,REINF_CL)
       LVCONV = .FALSE.
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'MACH' .OR.
     &       COMAND.EQ.'M   '      ) THEN
       IF(LPACC) THEN
        WRITE(*,2100)
        GO TO 500
       ENDIF
C
 15    CONTINUE
       IF(NINPUT.GE.1) THEN
        MINF1 = RINPUT(1)
       ELSE
        WRITE(*,*)
        WRITE(*,*) 'Currently...'
        CALL MRSHOW(.TRUE.,.FALSE.)
        CALL ASKR('Enter Mach number^',MINF1)
       ENDIF
C
       IF(MINF1.GE.1.0) THEN
        WRITE(*,*) 'Supersonic freestream not allowed'
        NINPUT = 0
        GO TO 15
       ENDIF
ccc    CALL MRSHOW(.TRUE.,.FALSE.)
       CALL MRCL(1.0,MINF_CL,REINF_CL)
       CALL COMSET
C
       IF(MINF.GT.0.0) WRITE(*,1300) CPSTAR, QSTAR/QINF
 1300  FORMAT(/' Sonic Cp =', F10.2, '      Sonic Q/Qinf =', F10.3/)
C
       CALL CPCALC(N,QINV,QINF,MINF,CPI)
       IF(LVISC) CALL CPCALC(N+NW,QVIS,QINF,MINF,CPV)
       CALL CLCALC(N,X,Y,GAM,GAM_A,ALFA,MINF,QINF, XCMREF,YCMREF,
     &             CL,CM,CDP, CL_ALF,CL_MSQ)
       CALL CDCALC
       LVCONV = .FALSE.
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'TYPE' .OR.
     &       COMAND.EQ.'T'         ) THEN
       IF(LPACC) THEN
        WRITE(*,2100)
        GO TO 500
       ENDIF
C
   17  CONTINUE
       IF(NINPUT.GE.1) THEN
        ITYP = IINPUT(1)
       ELSE
        WRITE(*,1105)
 1105   FORMAT(
     & /' Type   parameters held constant       varying      fixed   '
     & /' ----   ------------------------       -------   -----------'
     & /'   1    M          , Re            ..   lift     chord, vel.'
     & /'   2    M sqrt(CL) , Re sqrt(CL)   ..   vel.     chord, lift'
     & /'   3    M          , Re CL         ..   chord    lift , vel.')
        CALL ASKI('Enter type of Mach,Re variation with CL^',ITYP)
       ENDIF
C
       IF(ITYP.EQ.1) THEN
        MATYP = 1
        RETYP = 1
       ELSE IF(ITYP.EQ.2) THEN
        MATYP = 2
        RETYP = 2
       ELSE IF(ITYP.EQ.3) THEN
        MATYP = 1
        RETYP = 3
       ENDIF
C
       IF(ITYP.LT.1 .OR. ITYP.GT.3) THEN
        NINPUT = 0
        GO TO 17
       ENDIF
C
       CALL MRSHOW(.TRUE.,.TRUE.)
       CALL MRCL(1.0,MINF_CL,REINF_CL)
       CALL COMSET
       LVCONV = .FALSE.
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'ITER') THEN
 18    CONTINUE
       IF(NINPUT.GE.1) THEN
        ITMAX = IINPUT(1)
       ELSE
        WRITE(*,*) 'Current iteration limit:', ITMAX
        CALL ASKI('Enter new iteration limit^',ITMAX)
       ENDIF
C
       IF(ITMAX.LT.1) THEN
        NINPUT = 0
        GO TO 18
       ENDIF
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'INIT') THEN
       LBLINI = .NOT.LBLINI
       IF(LBLINI) THEN
        WRITE(*,*) 'BLs are assumed to be initialized'
       ELSE
        WRITE(*,*) 'BLs will be initialized on next point'
        LIPAN = .FALSE.
       ENDIF
C
C--------------------------------------------------------
c      ELSEIF(COMAND.EQ.'IMAG') THEN
c       LIMAGE = .NOT.LIMAGE
c       IF(LIMAGE) THEN
c        CALL ASKR('Enter y-position of image plane^',YIMAGE)
c        CALL ASKI('Specify image type (1=wall -1=free jet)^',KIMAGE)
c       ELSE
c        WRITE(*,*) 'Image airfoil removed'
c       ENDIF
c       LGAMU = .FALSE.
c       LQAIJ = .FALSE.
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'ALFA' .OR.
     &       COMAND.EQ.'A   '     ) THEN
       IF(.NOT.LRECALC) THEN
C------- set inviscid solution only if point is not being recalculated
         IF(NINPUT.GE.1) THEN
          ADEG = RINPUT(1)
         ELSE
          ADEG = ALFA/DTOR          
          CALL ASKR('Enter angle of attack (deg)^',ADEG)
         ENDIF
         LALFA = .TRUE.
         ALFA = DTOR*ADEG
         QINF = 1.0
         CALL SPECAL
         IF(ABS(ALFA-AWAKE) .GT. 1.0E-5) LWAKE  = .FALSE.
         IF(ABS(ALFA-AVISC) .GT. 1.0E-5) LVCONV = .FALSE.
         IF(ABS(MINF-MVISC) .GT. 1.0E-5) LVCONV = .FALSE. 
       ENDIF
C
       IF(LVISC) CALL VISCAL(ITMAX)
       CALL CPX
       CALL FCPMIN
C
ccc    IF( LVISC .AND. LPACC .AND. LVCONV ) THEN
       IF( LPACC .AND. (LVCONV .OR. .NOT.LVISC)) THEN
        CALL PLRADD(LUPLR,IPACT)
        CALL PLXADD(LUPLX,IPACT)
       ENDIF
C
       IF(LVISC .AND. .NOT.LPACC .AND. .NOT.LVCONV) THEN
        WRITE(*,*) 'Type "!" to continue iterating'
       ENDIF
C
C      WRITE(*,*) 'N NW =', N, NW
C      call aski('Enter i^',ioff)
C      call askr('Enter dmass^',dms)
Cc
C      do 43 is=1, 2
C        do 430 ibl=2, nbl(is)
C          i = ipan(ibl,is)
C          mass(ibl,is) = 0.
C          if(i.eq.ioff) mass(ibl,is) = dms
C  430   continue
C   43 continue
Cc
C      call ueset
C      call qvfue
C      call gamqv
C      call cpcalc(N+NW,QVIS,QINF,MINF,CPV)
C      call cdcalc
c      CALL CLCALC(N,X,Y,GAM,GAM_A,ALFA,MINF,QINF,XCMREF,YCMREF,
c     &            CL,CM,CDP, CL_ALF,CL_MSQ)
C      call cpx
Cc
C
       COMOLD = COMAND
       ARGOLD = COMARG
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'CLI ') THEN
       IF(.NOT.LRECALC) THEN
        IF(NINPUT.GE.1) THEN
         CLSPEC = RINPUT(1)
        ELSE
         CALL ASKR('Enter inviscid lift coefficient^',CLSPEC)
        ENDIF
        LALFA = .TRUE.
        ALFA = 0.0
        QINF = 1.0
        CALL SPECCL
        ADEG = ALFA/DTOR
        IF(ABS(ALFA-AWAKE) .GT. 1.0E-5) LWAKE  = .FALSE.
        IF(ABS(ALFA-AVISC) .GT. 1.0E-5) LVCONV = .FALSE.
        IF(ABS(MINF-MVISC) .GT. 1.0E-5) LVCONV = .FALSE.
       ENDIF
C
       IF(LVISC) CALL VISCAL(ITMAX)
       CALL CPX
       CALL FCPMIN
C
ccc    IF( LVISC .AND. LPACC .AND. LVCONV ) THEN
       IF( LPACC .AND. (LVCONV .OR. .NOT.LVISC)) THEN
        CALL PLRADD(LUPLR,IPACT)
        CALL PLXADD(LUPLX,IPACT)
       ENDIF
C
       IF(LVISC .AND. .NOT.LPACC .AND. .NOT.LVCONV) THEN
        WRITE(*,*) 'Type "!" to continue iterating'
       ENDIF
C
       COMOLD = COMAND
       ARGOLD = COMARG
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'CL  ' .OR.
     &       COMAND.EQ.'C   '     ) THEN
       IF(.NOT.LRECALC) THEN
        IF(NINPUT.GE.1) THEN
         CLSPEC = RINPUT(1)
        ELSE
         CALL ASKR('Enter lift coefficient^',CLSPEC)
        ENDIF
        LALFA = .FALSE.
        ALFA = 0.0
        QINF = 1.0
        CALL SPECCL
        ADEG = ALFA/DTOR
        IF(ABS(ALFA-AWAKE) .GT. 1.0E-5) LWAKE  = .FALSE.
        IF(ABS(ALFA-AVISC) .GT. 1.0E-5) LVCONV = .FALSE.
        IF(ABS(MINF-MVISC) .GT. 1.0E-5) LVCONV = .FALSE.
       ENDIF
       IF(LVISC) CALL VISCAL(ITMAX)
       CALL FCPMIN
C
       CALL CPX
ccc    IF( LVISC .AND. LPACC .AND. LVCONV ) THEN
       IF( LPACC .AND. (LVCONV .OR. .NOT.LVISC)) THEN
        CALL PLRADD(LUPLR,IPACT)
        CALL PLXADD(LUPLX,IPACT)
       ENDIF
C
       IF(LVISC .AND. .NOT.LPACC .AND. .NOT.LVCONV) THEN
        WRITE(*,*) 'Type "!" to continue iterating'
       ENDIF
C
       COMOLD = COMAND
       ARGOLD = COMARG
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'ASEQ' .OR.
     &       COMAND.EQ.'AS  ' .OR.
     &       COMAND.EQ.'CSEQ' .OR.
     &       COMAND.EQ.'CS  '      ) THEN
       LALFA = COMAND.EQ.'ASEQ' .OR.
     &         COMAND.EQ.'AS  '
C
       IF(LALFA) THEN
         IF    (NINPUT.GE.3) THEN
          AA1 = RINPUT(1)
          AA2 = RINPUT(2)
          DAA = RINPUT(3)
         ELSEIF(NINPUT.GE.2) THEN
          AA1 = RINPUT(1)
          AA2 = RINPUT(2)
          DAA = DAA/DTOR
          CALL ASKR('Enter alfa increment   (deg)^',DAA)
         ELSEIF(NINPUT.GE.1) THEN
          AA1 = RINPUT(1)
          AA2 = AA2/DTOR
          CALL ASKR('Enter last  alfa value (deg)^',AA2)
          DAA = DAA/DTOR
          CALL ASKR('Enter alfa increment   (deg)^',DAA)
         ELSE
          AA1 = AA1/DTOR
          CALL ASKR('Enter first alfa value (deg)^',AA1)
          AA2 = AA2/DTOR
          CALL ASKR('Enter last  alfa value (deg)^',AA2)
          DAA = DAA/DTOR
          CALL ASKR('Enter alfa increment   (deg)^',DAA)
         ENDIF
         IF(AA2.LT.AA1) THEN
          DAA = -ABS(DAA)
         ELSE
          DAA =  ABS(DAA)
         ENDIF
         AA1 = AA1*DTOR
         AA2 = AA2*DTOR
         DAA = DAA*DTOR
         NPOINT = 1
         IF(DAA .NE. 0.0) NPOINT = INT((AA2-AA1)/DAA + 0.5) + 1
C
       ELSE
         IF    (NINPUT.GE.3) THEN
          CL1 = RINPUT(1)
          CL2 = RINPUT(2)
          DCL = RINPUT(3)
         ELSEIF(NINPUT.GE.2) THEN
          CL1 = RINPUT(1)
          CL2 = RINPUT(2)
          CALL ASKR('Enter CL increment  ^',DCL)
         ELSEIF(NINPUT.GE.1) THEN
          CL1 = RINPUT(1)
          CALL ASKR('Enter last  CL value^',CL2)
          CALL ASKR('Enter CL increment  ^',DCL)
         ELSE
          CALL ASKR('Enter first CL value^',CL1)
          CALL ASKR('Enter last  CL value^',CL2)
          CALL ASKR('Enter CL increment  ^',DCL)
         ENDIF
         IF(CL2.LT.CL1) THEN
          DCL = -ABS(DCL)
         ELSE
          DCL =  ABS(DCL)
         ENDIF
         NPOINT = 1
         IF(DCL .NE. 0.0) NPOINT = INT((CL2-CL1)/DCL + 0.5) + 1
       ENDIF
C
C- - - - - - - - - - - - - - - - - - 
C
C----- initialize plot
       CALL PLTINI
C
       IF(LPPSHO) THEN
C------ set up for polar plot
C
       ELSE
C------ set up for Cp(x) plot
C
C------ Cp scaling factor
        PFAC = PLOTAR/(CPMAX-CPMIN)
C
C------ determine airfoil box size and location
        CALL AIRLIM(N,X,Y,XMIN,XMAX,YMIN,YMAX)
C
C------ y-offset for airfoil in  Cp vs x  plot
        FACA = FACAIR/(XMAX-XMIN)
        XOFA = XOFAIR*(XMAX-XMIN) - XMIN
        YOFA = YOFAIR*(XMAX-XMIN) - YMAX - CPMAX*PFAC/FACA
C
C------ re-origin for  Cp vs x  plot
        CALL PLOT(0.09 , 0.04 + CPMAX*PFAC + (YMAX-YMIN)*FACA, -3)
C
C------ draw axes and airfoil picture for Cp vs x plot
        CALL CPAXES(LCPGRD,
     &              N,X,Y,XOFA,YOFA,FACA,
     &              CPMIN,CPMAX,CPDEL,PFAC,CH,
     &              'XFOIL',VERSION)
C
C------ set initial x,y-positions of sequence plot label top
        XL = 0.65
        IF(LVISC) XL = 0.48
        YL = -CPMIN*PFAC
C
C------ draw sequence plot label
        CALL SEQLAB(XL,YL,XL1,XL2,XL3,XL4,XL5,XL6,CHSEQ,1,LVISC)
C
        CALL PLFLUSH
C
C------ set label y position
        YL = YL - 0.2*CH
       ENDIF
C
C----- initialize unconverged-point counter
       ISEQEX = 0
       ALAST = ADEG
       CLAST = CL
C
C----- calculate each point, add Cp distribution to plot, and save to polar
       DO 115 IPOINT=1, NPOINT
C
C------- set proper alpha for this point
         IF(LALFA) THEN
          ALFA = AA1 + DAA*FLOAT(IPOINT-1)
         ELSE
          CLSPEC = CL1 + DCL*FLOAT(IPOINT-1)
          CALL SPECCL
         ENDIF
C
         IF(ABS(ALFA-AWAKE) .GT. 1.0E-5) LWAKE  = .FALSE.
         IF(ABS(ALFA-AVISC) .GT. 1.0E-5) LVCONV = .FALSE.
         IF(ABS(MINF-MVISC) .GT. 1.0E-5) LVCONV = .FALSE.
         CALL SPECAL
         ITMAXS = ITMAX + 5
         IF(LVISC) CALL VISCAL(ITMAXS)
C
         ADEG = ALFA/DTOR
C
         CALL FCPMIN
C
C------- add point to buffer polar and/or disk files
ccc      IF( LVISC .AND. LPACC .AND. LVCONV ) THEN
         IF( LPACC .AND. (LVCONV .OR. .NOT.LVISC)) THEN
          CALL PLRADD(LUPLR,IPACT)
          CALL PLXADD(LUPLX,IPACT)
         ENDIF
C
         IF(LPPSHO) THEN
          CALL PLTINI
ccc          CALL PLOTABS(0.5,0.5,-3)
          PSIZE = 1.0*SIZE
          CALL NEWFACTOR(PSIZE)
          CALL PLOT(5.0*CH,7.0*CH,-3)
C
          CH1 = CH*0.90
          CH2 = CH*0.75
          CLEXP = 1.0
C
          DO IP=1, NPOL
           NBLP(IP) = 1
          ENDDO
C
          CALL POLPLT(NAX,NPOL,NAPOL,CPOL,
     &            REYNP1,MACHP1,ACRITP,PTRATP,ETAPP,
     &            NAMEPOL,ICOLP,ILINP,
     &            NFX,NPOLREF,NDREF,CPOLREF,NAMEREF,ICOLR,ISYMR,
     &            ISX,NBLP,CPOLSD ,IMATYP,IRETYP,
     &            ' ','XFOIL',VERSION,
     &            PLOTAR,XCDWID,XALWID,XOCWID,CH1,CH2,CLEXP,
     &            LPGRID,LPCDW,LPLIST,LPLEGN,LAECEN,LPCDH,LPCMDOT,
     &            CPOLPLF,' ',0)
         ELSE
C-------- add alpha, CL, etc. to plot
          CALL SEQPLT(YL,XL1,XL2,XL3,XL4,XL5,XL6,CHSEQ,ADEG,CL,CM,LVISC)
C
C-------- add sonic Cp dashed line if within plot
          IF(CPSTAR.GE.CPMIN) CALL DASH(0.0,XL-CH,-CPSTAR*PFAC)
C
          CALL NEWPEN(2)
          IF(LVISC) THEN
C--------- Plot viscous -Cp distribution on airfoil
           CALL XYLINE(N+NW,X,CPV,-XOFA,FACA,0.0,-PFAC,1)
          ELSE
C--------- Plot inviscid -Cp distribution on airfoil
           CALL XYLINE(N,X,CPI,-XOFA,FACA,0.0,-PFAC,1)
          ENDIF
         ENDIF
C
         CALL PLFLUSH
c###
ccc    call dcpout
C
         IF(LVISC .AND. .NOT.LVCONV) THEN
C-------- increment unconverged-point counter
          ISEQEX = ISEQEX + 1
          IF(ISEQEX .GE. NSEQEX) THEN
           WRITE(*,1150) ISEQEX, ALAST, CLAST
 1150      FORMAT(
     & /' Sequence halted since previous',I3,' points did not converge'
     & /' Last-converged  alpha =', F8.3, '    CL =', F10.5)
           GO TO 116
          ENDIF
         ELSE
C-------- converged OK... reset unconverged-point counter
          ISEQEX = 0
          ALAST = ADEG
          CLAST = CL
         ENDIF
C
  115  CONTINUE
  116  CONTINUE
ccc      CALL ASKC('hit <cr>^',DUMMY,COMARG)
C
       COMOLD = COMAND
       ARGOLD = COMARG
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'SEQP') THEN
       LPPSHO = .NOT.LPPSHO
       IF(LPPSHO) THEN
        WRITE(*,*) 'Polar will be plotted during point sequence'
       ELSE
        WRITE(*,*) 'Cp(x) will be plotted during point sequence'
       ENDIF
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'PACC' .OR.
     &       COMAND.EQ.'P   '      ) THEN
       LPACC = .NOT.LPACC
C
       IF(LPACC) THEN
        IF(NINPUT.GE.1) THEN
C------- slot into which accumulated polar will go
         IP = MIN( MAX( IINPUT(1) , 0 ) , NPOL+1 )
        ELSE
C------- no command argument was given... just use next available slot
         IP = NPOL+1
         PFNAME(IP) = ' '
         PFNAMX(IP) = ' '
        ENDIF
C
        IF(IP.GT.NPOL) THEN
         IF(NPOL.EQ.NPX) THEN
          WRITE(*,*)
          WRITE(*,*) 'Number of polars is at array limit'
          WRITE(*,*) 'New polar will not be stored'
          IPACT = 0
         ELSE
          IPACT = NPOL + 1
          PFNAME(IPACT) = ' '
          PFNAMX(IPACT) = ' '
         ENDIF
C
        ELSE
         IPACT = IP
C
        ENDIF
C
C------ set up for appending to new or existing polar (if IPACT > 0)
        CALL PLRSET(IPACT)
C
C------ jump out if decision was made to abort polar accumulation
        IF(IPACT.LE.0) THEN
         LPACC = .FALSE.
         GO TO 500
        ENDIF
C
        CALL PLRINI(LUPLR,IPACT)
        CALL PLXINI(LUPLX,IPACT)
        WRITE(*,*)
        WRITE(*,*) 'Polar accumulation enabled'
C
       ELSE
        WRITE(*,*)
        WRITE(*,*) 'Polar accumulation disabled'
        IPACT = 0
C
       ENDIF
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'PGET') THEN
        IF(NPOL.GE.NPX) THEN
         WRITE(*,*)
         WRITE(*,*) 'Number of polars is at array limit'
         WRITE(*,*) 'Delete with PDEL if necessary'
         GO TO 500
        ENDIF
C
        IP = NPOL+1
C
        IF(COMARG.EQ.' ') THEN
         CALL ASKS('Enter polar filename^',FNAME)
        ELSE
         FNAME = COMARG
        ENDIF
C
        LU = 17
        CALL POLREAD(LU,FNAME,ERROR,
     &             NAX,NAPOL(IP),CPOL(1,1,IP), 
     &             REYNP1(IP),MACHP1(IP),ACRITP(IP),XSTRIPP(1,IP),
     &             PTRATP(IP),ETAPP(IP),
     &             NAMEPOL(IP),IRETYP(IP),IMATYP(IP),
     &             ISX,NBLP(IP),CPOLSD(1,1,1,IP),
     &             CODEPOL(IP),VERSPOL(IP) )
        IF(ERROR) THEN
         WRITE(*,*) 'Polar file READ error'
        ELSE
         NPOL = IP
         NXYPOL(IP) = 0
         CALL STRIP(NAMEPOL(IP),NNAMEP)
         NEL = 1
         CALL POLWRIT(6,' ',ERROR, .TRUE.,
     &             NAX, 1,NAPOL(IP), CPOL(1,1,IP), IPOL,NIPOL,
     &             REYNP1(IP),MACHP1(IP),ACRITP(IP),XSTRIPP(1,IP),
     &             PTRATP(IP),ETAPP(IP),
     &             NAMEPOL(IP),IRETYP(IP),IMATYP(IP),
     &             ISX,NEL,CPOLSD(1,1,1,IP), JPOL,NJPOL,
     &             CODEPOL(IP),VERSPOL(IP), .FALSE. )
         PFNAME(IP) = FNAME
         WRITE(*,5500) IP
 5500    FORMAT(/' Stored as  Polar', I4)
        ENDIF
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'PWRT') THEN
 75     CONTINUE
        IF(NPOL.EQ.1) THEN
         IP = 1
        ELSEIF(NINPUT.EQ.0) THEN
         CALL PLRSUM(1,NPOL,IPACT)
         CALL ASKI(
     &      'Enter index of polar to write (0=all, -1=abort)^',IP)
         IF(IP.EQ.-1) GO TO 500
        ELSE
         IP = IINPUT(1)
        ENDIF
C
        IF(IP.EQ.0) THEN
         IP1 = 1
         IP2 = NPOL
        ELSEIF(IP.GE.1 .AND. IP.LE.NPOL) THEN
         IP1 = IP
         IP2 = IP
        ELSE
         NINPUT = 0
         GO TO 75
        ENDIF
C
        NEL = 1
        DO IP = IP1, IP2
          LU = 19
          CALL PLRSUM(IP,IP,IPACT)
          CALL STRIP(PFNAME(IP),NPF)
          IF(NPF.EQ.0) THEN
           LINE = 'Enter polar output filename^' 
          ELSE
           LINE = 'Enter polar output filename [' 
     &            // PFNAME(IP)(1:NPF) // ']^'
          ENDIF
          CALL ASKS(LINE,FNAME)
          IF(NPF.NE.0 .AND. FNAME.EQ.' ') FNAME = PFNAME(IP)
C
          NIPOL = NIPOL0
          IF(LCMINP) THEN
           NIPOL = NIPOL + 1
           IPOL(IMC) = NIPOL
          ENDIF
          IF(LHMOMP) THEN
           NIPOL = NIPOL + 1
           IPOL(ICH) = NIPOL
          ENDIF
C
          CALL POLWRIT(LU,FNAME,ERROR, .TRUE.,
     &      NAX, 1,NAPOL(IP),CPOL(1,1,IP), IPOL,NIPOL,
     &      REYNP1(IP),MACHP1(IP),ACRITP(IP),XSTRIPP(1,IP),
     &      PTRATP(IP),ETAPP(IP),
     &      NAMEPOL(IP),IRETYP(IP),IMATYP(IP),
     &      ISX,NEL,CPOLSD(1,1,1,IP), JPOL,NJPOL,
     &      'XFOIL',VERSION, .TRUE. )
          IF(ERROR) THEN
           WRITE(*,1075) IP
 1075      FORMAT(' Polar', I3,'  not written')
          ELSE
           PFNAME(IP) = FNAME
          ENDIF
        ENDDO
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'RGET') THEN
        IF(NPOLREF.GE.NPX) THEN
         WRITE(*,*)
         WRITE(*,*) 'Number of reference polars is at array limit'
         WRITE(*,*) 'Delete with RDEL if necessary'
         GO TO 500
        ENDIF
C
        IR = NPOLREF+1
C
        IF(COMARG.EQ.' ') THEN
         CALL ASKS('Enter reference polar filename^',FNAME)
        ELSE
         FNAME = COMARG
        ENDIF
C
        LU = 9
        OPEN(LU,FILE=FNAME,STATUS='OLD',ERR=27)
        CALL POLREF(LU, FNAME, ERROR,
     &              NFX, NDREF(1,IR), CPOLREF(1,1,1,IR), NAMEREF(IR))
        CLOSE(LU)
        IF(ERROR) GO TO 27
C
        NPOLREF = IR
C
        CALL STRIP(NAMEREF(IR),NNREF)
        IF(NNREF.EQ.0) THEN
          CALL ASKS('Enter label for reference polar^',NAMEREF(IR))
          CALL STRIP(NAMEREF(IR),NNREF)
        ELSE
          WRITE(*,*)
          WRITE(*,*) NAMEREF(IR)
        ENDIF
C
ccc     ICOLR(IR) = NCOLOR - IR + 1
        ICOLR(IR) = 2 + IR
        ISYMR(IR) = MOD(IR,10)
 25   CONTINUE
      GO TO 500
C
 27   CONTINUE
      WRITE(*,*) 'File OPEN error'
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'RDEL') THEN
        IF(NPOLREF.EQ.0) THEN
         WRITE(*,*) 'No reference polars are stored'
         GO TO 500
        ENDIF
C
        IF(NINPUT.GE.1) THEN
         IR = IINPUT(1)
        ELSE
         IR = NPOLREF+1
        ENDIF
C
 35     CONTINUE
C
        IF(IR.EQ.0) THEN
C------- delete all polars
         NPOLREF = 0
C
        ELSEIF(IR.EQ.-1) THEN
C------- abort
         GO TO 500
C
        ELSEIF(IR.LT.-1 .OR. IR.GT.NPOLREF) THEN
         CALL PRFSUM(1,NPOLREF)
         CALL ASKI(
     &      'Specify ref. polar to delete (0 = all, -1 = abort)^',IR)
         GO TO 35
C
        ELSE
C------- delete ref. polar IR
         DO JR = IR+1, NPOLREF
           CALL PRFCOP(JR,JR-1)
           WRITE(*,1310) JR, JR-1
 1410      FORMAT(' Ref.polar',I3,'  moved into ref.polar',I3)
         ENDDO
         NPOLREF = NPOLREF-1
        ENDIF
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'PSUM') THEN
        IF(NPOL.EQ.0) THEN
         WRITE(*,*)
         WRITE(*,*) 'No polars are stored'
         GO TO 500
        ENDIF
C
        CALL PLRSUM(1,NPOL,IPACT)
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'PLIS') THEN
        IF(NPOL.EQ.0) THEN
         WRITE(*,*)
         WRITE(*,*) 'No polars are stored'
         GO TO 500
        ENDIF
C
        IF(NINPUT.EQ.0) THEN
         IP1 = 1
         IP2 = NPOL
        ELSE
         IP = IINPUT(1)
         IF(IP.EQ.0) THEN
          IP1 = 1
          IP2 = NPOL
         ELSEIF(IP.GE.1 .AND. IP.LE.NPOL) THEN
          IP1 = IP
          IP2 = IP
         ELSE
          WRITE(*,*)
          WRITE(*,*) 'Specified stored polar does not exist'
          GO TO 500
         ENDIF
        ENDIF
C
        NIPOL = NIPOL0
        IF(LCMINP) THEN
         NIPOL = NIPOL + 1
         IPOL(IMC) = NIPOL
        ENDIF
        IF(LHMOMP) THEN
         NIPOL = NIPOL + 1
         IPOL(ICH) = NIPOL
        ENDIF
C
        NEL = 1
        DO IP = IP1, IP2
          WRITE(*,3100) IP
 3100     FORMAT(
     &/' =============================================================='
     &/' Polar', I3)
          IA1 = 1
          IA2 = NAPOL(IP)
          CALL POLWRIT(6,' ',ERROR, .TRUE.,
     &                 NAX, IA1,IA2, CPOL(1,1,IP), IPOL,NIPOL,
     &                 REYNP1(IP),MACHP1(IP),ACRITP(IP),XSTRIPP(1,IP),
     &                 PTRATP(IP),ETAPP(IP),
     &                 NAMEPOL(IP), IRETYP(IP),IMATYP(IP),
     &                 ISX,NEL,CPOLSD(1,1,1,IP), JPOL,NJPOL,
     &                 'XFOIL',VERSION, .FALSE.)
        ENDDO
        NIPOL = NIPOL0
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'PDEL') THEN
        IF(NPOL.EQ.0) THEN
         WRITE(*,*) 'No polars are stored'
         GO TO 500
        ENDIF
C
        IF(NINPUT.GE.1) THEN
C------- use command argument
         IP = IINPUT(1)
        ELSE
C------- no argument given... set up for user query test below
         IP = NPOL+1
        ENDIF
C
 40     CONTINUE
        IF(IP.EQ.0) THEN
C------- delete all polars
         NPOL = 0
         IPACT = 0
         LPACC = .FALSE.
C
        ELSEIF(IP.EQ.-1) THEN
C------- abort
         GO TO 500
C
        ELSEIF(IP.LT.-1 .OR. IP.GT.NPOL) THEN
         CALL PLRSUM(1,NPOL,IPACT)
         CALL ASKI(
     &      'Specify polar to delete (0 = all, -1 = abort)^',IP)
         GO TO 40
C
        ELSE
C------- delete polar IP
         IF(IPACT.EQ.IP) THEN
          WRITE(*,*) 'Active polar deleted.  Accumulation turned off'
          IPACT = 0
          LPACC = .FALSE.
         ENDIF
C
         DO JP = IP+1, NPOL
           CALL PLRCOP(JP,JP-1)
           WRITE(*,1310) JP, JP-1
 1310      FORMAT(' Polar',I3,'  moved into polar',I3)
           IF(IPACT.EQ.JP) THEN
            IPACT = JP-1
           ENDIF
         ENDDO
         NPOL = NPOL-1
C
        ENDIF
C
        IF(IPACT.GT.0) THEN
         WRITE(*,1320) IPACT
 1320    FORMAT(' Polar',I3,'  is now active')
        ENDIF
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'PSOR') THEN
        IF(NPOL.EQ.0) THEN
         WRITE(*,*) 'No polars are stored'
         GO TO 500
        ENDIF
C
        IF(NINPUT.GE.1) THEN
C------- use command argument
         IP = IINPUT(1)
        ELSE
C------- no argument given... set up for user query test below
         IP = NPOL+1
        ENDIF
C
C------ sort polars in increasing alpha
        IDSORT = IAL
C
 42     CONTINUE
        IF    (IP.EQ.-1) THEN
C------- abort
         GO TO 500
C
        ELSEIF(IP.LT.-1 .OR. IP.GT.NPOL) THEN
         CALL PLRSUM(1,NPOL,IPACT)
         CALL ASKI(
     &      'Specify polar to sort (0 = all, -1 = abort)^',IP)
         GO TO 42
C
        ELSE
C------- sort polar(s)
         IF(IP.EQ.0) THEN
          IP1 = 1
          IP2 = NPOL
         ELSE
          IP1 = IP
          IP2 = IP
         ENDIF
         DO JP = IP1, IP2
           CALL PLRSRT(JP,IDSORT)
         ENDDO
        ENDIF
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'PPLO' .OR.
     &       COMAND.EQ.'PP  '      ) THEN
C------ set temporary polar-size array to plot only selected polars
        IF(NINPUT.EQ.0) THEN
C------- no polars specified... plot all of them
         DO IP=1, NPOL
           NAPOLT(IP) = NAPOL(IP)
         ENDDO
        ELSE
C------- set up to plot only specified polars
         DO IP=1, NPOL
           NAPOLT(IP) = 0
         ENDDO
         DO K=1, NINPUT
           IP = IINPUT(K)
           IF(IP.GE.1 .AND. IP.LE.NPOL) NAPOLT(IP) = NAPOL(IP)
         ENDDO
        ENDIF
C
        CALL PLTINI
ccc        CALL PLOTABS(0.5,0.5,-3)
        PSIZE = 1.0*SIZE
        CALL NEWFACTOR(PSIZE)
        CALL PLOT(5.0*CH,7.0*CH,-3)
C
        CH1 = CH*0.90
        CH2 = CH*0.75
        CLEXP = 1.0
        DO IP=1, NPOL
          NBLP(IP) = 1
        ENDDO
C
        CALL POLPLT(NAX,NPOL,NAPOLT,CPOL,
     &          REYNP1,MACHP1,ACRITP,PTRATP,ETAPP,
     &          NAMEPOL,ICOLP,ILINP,
     &          NFX,NPOLREF,NDREF,CPOLREF,NAMEREF,ICOLR,ISYMR,
     &          ISX,NBLP,CPOLSD ,IMATYP,IRETYP,
     &          ' ','XFOIL',VERSION,
     &          PLOTAR,XCDWID,XALWID,XOCWID,CH1,CH2,CLEXP,
     &          LPGRID,LPCDW,LPLIST,LPLEGN,LAECEN,LPCDH,LPCMDOT,
     &          CPOLPLF,' ',0)
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'APLO' .OR.
     &       COMAND.EQ.'AP  '      ) THEN
        IF(NPOL.EQ.0) THEN
         WRITE(*,*)
         WRITE(*,*) 'No polars are stored'
         GO TO 500
        ENDIF
C
        IF(NINPUT.EQ.0) THEN
         NPPAI = NPOL
         DO K=1, NPPAI
           IPPAI(K) = K
         ENDDO
        ELSE
         NPPAI = MIN( NINPUT , NPX )
         DO K=1, NPPAI
           IINP = IINPUT(K)
           IF(IINP.GE.1 .AND. IINP.LE.NPOL) THEN
            IPPAI(K) = IINP
           ELSE
            IPPAI(K) = 0
           ENDIF
         ENDDO
        ENDIF
C
        CALL PPAPLT(NPPAI,IPPAI)
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'ASET') THEN
       IF(NPOL.EQ.0) THEN
        WRITE(*,*)
        WRITE(*,*) 'No polar airfoils are stored'
        GO TO 500
       ENDIF
C
 50    CONTINUE
       IF(NINPUT.EQ.0) THEN
        IF(NPOL.EQ.1) THEN
         IP = 1
        ELSE
         CALL PLRSUM(1,NPOL,IPACT)
         CALL ASKI('Enter index of polar airfoil to set^',IP)
        ENDIF
       ELSE
        IP = IINPUT(1)
       ENDIF
C
       IF(IP.EQ.0) THEN
        GO TO 500
       ELSEIF(IP.LT.1 .OR. IP.GT.NPOL) THEN
        WRITE(*,*)
        WRITE(*,*) 'Specified polar airfoil does not exist'
        NINPUT = 0
        GO TO 50
       ENDIF
C
       WRITE(*,*)
       WRITE(*,*) 'Current airfoil will be overwritten.  Proceed?  Y'
       READ(*,1000) ANS
 1000  FORMAT(A)
C
       IF(INDEX('Nn',ANS) .NE. 0) THEN
        WRITE(*,*) 'No action taken'
        GO TO 500
       ELSE
        CALL APCOPY(IP)
       ENDIF
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'PREM') THEN
       IF(NPOL.EQ.0) THEN
        WRITE(*,*)
        WRITE(*,*) 'No polars are stored'
        GO TO 500
       ENDIF
C
 52    CONTINUE
       IF(NINPUT.EQ.0) THEN
        IF(NPOL.EQ.1) THEN
         IP = 1
        ELSE
         CALL PLRSUM(1,NPOL,IPACT)
         CALL ASKI('Enter index of polar to modify^',IP)
        ENDIF
       ELSE
        IP = IINPUT(1)
       ENDIF
C
       IF(IP.EQ.0) THEN
        GO TO 500
       ELSEIF(IP.LT.1 .OR. IP.GT.NPOL) THEN
        WRITE(*,*)
        WRITE(*,*) 'Specified polar airfoil does not exist'
        NINPUT = 0
        GO TO 52
       ENDIF
C
       IF(NINPUT.GE.2) THEN
        NREM = NINPUT - 1
       ELSE
        NIPOL = NIPOL0
        IF(LCMINP) THEN
         NIPOL = NIPOL + 1
         IPOL(IMC) = NIPOL
        ENDIF
        IF(LHMOMP) THEN
         NIPOL = NIPOL + 1
         IPOL(ICH) = NIPOL
        ENDIF
C
        WRITE(*,3100) IP
        IA1 = 1
        IA2 = NAPOL(IP)
        CALL POLWRIT(6,' ',ERROR, .TRUE.,
     &                 NAX, IA1,IA2, CPOL(1,1,IP), IPOL,NIPOL,
     &                 REYNP1(IP),MACHP1(IP),ACRITP(IP),XSTRIPP(1,IP),
     &                 PTRATP(IP),ETAPP(IP),
     &                 NAMEPOL(IP), IRETYP(IP),IMATYP(IP),
     &                 ISX,1,CPOLSD(1,1,1,IP), JPOL,NJPOL,
     &                 'XFOIL',VERSION, .FALSE. )
 53     WRITE(*,3220)
 3220   FORMAT(/' Enter alpha(s) of points to be removed:  ', $)
        READ(*,1000) LINE
        NREM = 19
        CALL GETFLT(LINE,RINPUT(2),NREM,ERROR)
        IF(ERROR) GO TO 53
       ENDIF
C
C----- go over specified alphas to be removed
       DO 55 IREM = 1, NREM
C------- check all alpha points in polar IP
         DO IA = 1, NAPOL(IP)
           ADIF = CPOL(IA,IAL,IP) - RINPUT(IREM+1)
           IF(ABS(ADIF) .LT. 0.0005) THEN
C---------- alphas match within 3-digit print tolerance... 
C-             remove point by pulling down all points above it
            DO JA = IA, NAPOL(IP)-1
              DO K = 1, IPTOT
                CPOL(JA,K,IP) = CPOL(JA+1,K,IP)
              ENDDO
              DO K = 1, JPTOT
                CPOLSD(JA,1,K,IP) = CPOLSD(JA+1,1,K,IP)
                CPOLSD(JA,2,K,IP) = CPOLSD(JA+1,2,K,IP)
              ENDDO
            ENDDO
C---------- shrink polar by 1
            NAPOL(IP) = NAPOL(IP) - 1
C
            IF(NAPOL(IP).LE.0) THEN
C----------- last point has been removed... eliminate this polar IP
             DO JP = IP+1, NPOL
               CALL PLRCOP(JP,JP-1)
               IF(IPACT.EQ.JP) IPACT = JP-1
               WRITE(*,1310) JP, JP-1
             ENDDO
             NPOL = NPOL-1
C
             IF(IPACT.GT.0) THEN
              WRITE(*,1320) IPACT
             ENDIF
C
             GO TO 500
            ENDIF
C
C---------- go to next specified alpha to be removed
            GO TO 55
           ENDIF
         ENDDO
 55    CONTINUE
C         
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'PNAM') THEN
       IF(NPOL.EQ.0) THEN
        WRITE(*,*)
        WRITE(*,*) 'No polars are stored'
        GO TO 500
       ENDIF
C
 58    CONTINUE
       IF(NINPUT.EQ.0) THEN
        IF(NPOL.EQ.1) THEN
         IP = 1
        ELSE
         CALL PLRSUM(1,NPOL,IPACT)
         CALL ASKI('Enter index of polar to modify^',IP)
        ENDIF
       ELSE
        IP = IINPUT(1)
       ENDIF
C
       IF(IP.EQ.0) THEN
        GO TO 500
       ELSEIF(IP.LT.1 .OR. IP.GT.NPOL) THEN
        WRITE(*,*)
        WRITE(*,*) 'Specified polar airfoil does not exist'
        NINPUT = 0
        GO TO 58
       ENDIF
C
       NIPOL = NIPOL0
       IF(LCMINP) THEN
        NIPOL = NIPOL + 1
        IPOL(IMC) = NIPOL
       ENDIF
       IF(LHMOMP) THEN
        NIPOL = NIPOL + 1
        IPOL(ICH) = NIPOL
       ENDIF
C
       WRITE(*,3100) IP
       IA1 = 0
       IA2 = -1
       CALL POLWRIT(6,' ',ERROR, .TRUE.,
     &                NAX, IA1,IA2, CPOL(1,1,IP), IPOL,NIPOL,
     &                REYNP1(IP),MACHP1(IP),ACRITP(IP),XSTRIPP(1,IP),
     &                PTRATP(IP),ETAPP(IP),
     &                NAMEPOL(IP), IRETYP(IP),IMATYP(IP),
     &                ISX,1,CPOLSD(1,1,1,IP), JPOL,NJPOL,
     &                'XFOIL',VERSION, .FALSE. )
       NIPOL = NIPOL0
       WRITE(*,3320)
 3320  FORMAT(/' Enter new airfoil name of polar:  ', $)
       READ(*,1000) NAMEPOL(IP)
       CALL STRIP(NAMEPOL(IP),NNP)
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'PPAX') THEN
        CALL POLAXI(CPOLPLF,XCDWID,XALWI,XOCWID)
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'CREF') THEN
       LCPREF = .NOT. LCPREF
       IF(LCPREF) THEN
        WRITE(*,*) 'Reference Cp plotting enabled'
       ELSE
        WRITE(*,*) 'Reference Cp plotting disabled'
       ENDIF
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'FREF') THEN
       LFOREF = .NOT. LFOREF
       IF(LFOREF) THEN
        WRITE(*,*) 'Reference force plotting enabled'
       ELSE
        WRITE(*,*) 'Reference force plotting disabled'
       ENDIF
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'CPX ' .OR.
     &       COMAND.EQ.'CP  '      ) THEN
       CALL CPX
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'UEX ' .OR.
     &       COMAND.EQ.'UE  '      ) THEN
       CALL UEX
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'GRID') THEN
       LCPGRD = .NOT.LCPGRD
       IF(LCPGRD) THEN
        WRITE(*,*) 'Cp grid overlay enabled'
       ELSE
        WRITE(*,*) 'Cp grid overlay disabled'
       ENDIF
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'CPV ') THEN
       CALL CPVEC
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'BL  ') THEN
       IF(.NOT.LVCONV) THEN
        WRITE(*,*) 'Compute valid viscous solution first'
        GO TO 500
       ENDIF
C
       IF(NINPUT.GE.1) THEN
        NPR = MIN( IINPUT(1) , NPRX )
       ELSE
        NPR = 21
        WRITE(*,*) 'Using default number of profiles:', NPR
       ENDIF
C
       IF(NPR.GT.1) THEN
C------ set NPR points along surface, offset slightly for the locating logic
        DOFF = 0.00001*(S(N)-S(1))
        DO IPR = 1, NPR
          FRAC = FLOAT(IPR-1)/FLOAT(NPR-1)
          SPR = S(1) + (S(N)-S(1))*FRAC
          XPR(IPR) = SEVAL(SPR,X,XP,S,N) + DOFF*DEVAL(SPR,Y,YP,S,N)
          YPR(IPR) = SEVAL(SPR,Y,YP,S,N) - DOFF*DEVAL(SPR,X,XP,S,N)
        ENDDO
       ENDIF
C
       CALL CPX
       CALL DPLOT(NPR,XPR,YPR)
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'BLC ') THEN
       IF(.NOT.LVCONV) THEN
        WRITE(*,*) 'Compute valid viscous solution first'
        GO TO 500
       ENDIF
C
       NPR = 0
       CALL CPX
       CALL DPLOT(NPR,XPR,YPR)
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'BLF ') THEN
C
       NPR = 2
       DO IPR = 1, NPR
C
c       WRITE(*,'(1X,A,$)') 'Enter  x/c_BL, delta : '
c       READ(*,1000) LINE
c       NINP = 2
c       CALL GETFLT(LINE,RINPUT,NINP,ERROR)
c       IF(ERROR .OR. NINP.EQ.0) THEN
c        GO TO 500
c       ELSE
c        SGN = SIGN(1.0,RINPUT(1))
c        XOC = ABS(RINPUT(1))
c        DPR = RINPUT(2)
c       ENDIF

        if    (ipr.eq.1) then
         xoc = 0.4
         sgn = 1.0
        elseif(ipr.eq.2) then
         xoc = 0.4
         sgn = -1.0
        endif

        if(ninput .gt. 0) then
         dpr = rinput(1)
        else
         dpr = 0.01
        endif

C       
       IF(SGN .GT. 0.0) THEN
        SPR = SLE + (S(1)-SLE)*XOC
       ELSE
        SPR = SLE + (S(N)-SLE)*XOC
       ENDIF
C
       XPRI = XLE + (XTE-XLE)*XOC
       CALL SINVRT(SPR,XPRI,X,XP,S,N)
C
       DOFF = 0.00001*(S(N)-S(1))
       XPR(IPR) = SEVAL(SPR,X,XP,S,N) + DOFF*DEVAL(SPR,Y,YP,S,N)
       YPR(IPR) = SEVAL(SPR,Y,YP,S,N) - DOFF*DEVAL(SPR,X,XP,S,N)
C
       CALL FBLGET(XPR(IPR),YPR(IPR), DPR,FPR(IPR) )
C
       enddo

       WRITE(*,*)
       DO IPR = 1, NPR
         WRITE(*,7720) 'xBL, Fint =', XPR(IPR), FPR(IPR)*1.0E4
       ENDDO
 7720  FORMAT(1X,A,F7.3,F12.6)

ccc       GO TO 770
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'BLWT') THEN
       IF(NINPUT.GE.1) THEN
        UPRWT = RINPUT(1)
       ELSE
        WRITE(*,*) 'Current u/Qinf profile plot weight =', UPRWT
        CALL ASKR('Enter new plot weight^',UPRWT)
       ENDIF
C
       CALL CPX
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'FMOM') THEN
       CALL MHINGE
       WRITE(*,1500) XOF,YOF,HMOM,HFX,HFY
 1500  FORMAT(/' Flap hinge x,y :', 2F8.4/
     &         '                                           2  2'/
     &         ' Hinge moment/span = ',F8.6,'  x  1/2 rho V  c '/
     &         '                                           2   '/
     &         ' x-Force     /span = ',F8.6,'  x  1/2 rho V  c '/
     &         '                                           2   '/
     &         ' y-Force     /span = ',F8.6,'  x  1/2 rho V  c '/)
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'FNEW') THEN
       IF    (NINPUT.GE.2) THEN
        XOF = RINPUT(1)
        YOF = RINPUT(2)
       ELSEIF(NINPUT.GE.1) THEN
        XOF = RINPUT(1)
        YOF = -999.0
       ELSE
        XOF = -999.0
        YOF = -999.0
       ENDIF
       LFLAP = .FALSE.
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'VELS') THEN
       IF    (NINPUT.GE.2) THEN
        XXX = RINPUT(1)
        YYY = RINPUT(2)
       ELSEIF(NINPUT.GE.1) THEN
        XXX = RINPUT(1)
        CALL ASKR('Enter y^',YYY)
       ELSE
        CALL ASKR('Enter x^',XXX)
        CALL ASKR('Enter y^',YYY)
       ENDIF
       CALL PSILIN(0,XXX,YYY,-1.0,0.0,PSI,VVV,.FALSE.,.TRUE.)
       CALL PSILIN(0,XXX,YYY, 0.0,1.0,PSI,UUU,.FALSE.,.TRUE.)
       QQQ = SQRT(UUU**2 + VVV**2)
       CPP = 1.0 - (UUU**2 + VVV**2)
       WRITE(*,1800) UUU,VVV,QQQ,CPP
 1800  FORMAT(/' u/Uinf = ', F8.4, '   v/Uinf = ', F8.4
     &        /' q/Uinf = ', F8.4, '   Cp     = ', F8.4 /  )
C
       COMOLD = COMAND
       ARGOLD = COMARG
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'DUMP') THEN
       CALL BLDUMP(COMARG)
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'CPWR') THEN
       CALL CPDUMP(COMARG)
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'CPMN') THEN
       IF(LVISC)THEN
        WRITE(*,1769) CPMNI, XCPMNI, CPMNV, XCPMNV
 1769   FORMAT('  Minimum Inviscid Cp =',F8.4,'   at x =',F8.4
     &       / '  Minimum Viscous  Cp =',F8.4,'   at x =',F8.4 )
       ELSE
        WRITE(*,1779) CPMNI, XCPMNI
 1779   FORMAT('  Minimum Inviscid Cp =',F8.4,'   at x =',F8.4)
       ENDIF
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'CINC') THEN
       LCMINP = .NOT.LCMINP
       IF(LCMINP) THEN
        WRITE(*,*) 'Min Cp will be written to polar save file'
       ELSE 
        WRITE(*,*) 'Min Cp won''t be written to polar save file'
       ENDIF
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'HINC') THEN
       LHMOMP = .NOT.LHMOMP
       IF(LHMOMP) THEN
        WRITE(*,*) 'Hinge moment will be written to polar save file'
        IF(.NOT.LFLAP) THEN
         WRITE(*,*)
         WRITE(*,*) 'Note: Flap hinge location not defined'
         WRITE(*,*) '      Set it with FNEW command'
        ENDIF
       ELSE 
        WRITE(*,*) 'Hinge moment won''t be written to polar save file'
       ENDIF
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'ANNO') THEN
       IF(LPLOT) THEN
        CALL ANNOT(CH)
       ELSE
        WRITE(*,*) 'No active plot to annotate'
       ENDIF
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'VPLO' .OR.
     &       COMAND.EQ.'VP  '      ) THEN
       CALL BLPLOT
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'NAME') THEN
       IF(COMARG.EQ.' ') THEN
        CALL NAMMOD(NAME,0,-1)
       ELSE
        NAME = COMARG
       ENDIF
       CALL STRIP(NAME,NNAME)
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'NINC') THEN
       CALL NAMMOD(NAME,1,1)
       CALL STRIP(NAME,NNAME)
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'NDEC') THEN
       CALL NAMMOD(NAME,-1,1)
       CALL STRIP(NAME,NNAME)
C
C--------------------------------------------------------
      ELSEIF(COMAND.EQ.'DAMP') THEN
       IF(IDAMP.EQ.0) THEN
        IDAMP = 1
        WRITE(*,*) 'Modified amplification used'
       ELSE 
        IDAMP = 0
        WRITE(*,*) 'Original amplification used'
       ENDIF
C--------------------------------------------------------
      ELSE
       WRITE(*,8000) COMAND
 8000  FORMAT(1X,A4,' command not recognized.  Type a "?" for list')

      ENDIF
C
C---- go back to top of menu loop
      GO TO 500
C
C--------------------------------------------
 2100 FORMAT(/' * Polar is being accumulated.'
     &       /' * Cannot change its parameters in midstream.')
      END ! OPER


      SUBROUTINE FCPMIN
C------------------------------------------------
C     Finds minimum Cp on dist for cavitation work
C------------------------------------------------
      INCLUDE 'XFOIL.INC'
C
      XCPMNI = X(1)
      XCPMNV = X(1)
      CPMNI = CPI(1)
      CPMNV = CPV(1)
C
      DO I = 2, N + NW
        IF(CPI(I) .LT. CPMNI) THEN
         XCPMNI = X(I)
         CPMNI = CPI(I)
        ENDIF
        IF(CPV(I) .LT. CPMNV) THEN
         XCPMNV = X(I)
         CPMNV = CPV(I)
        ENDIF
      ENDDO
C

      IF(LVISC)THEN
        CPMN = CPMNV
      ELSE
        CPMN = CPMNI
C
        CPMNV = CPMNI
        XCPMNV = XCPMNI
      ENDIF
C
      RETURN
      END ! FCPMIN



      SUBROUTINE MRSHOW(LM,LR)
      INCLUDE 'XFOIL.INC'
      LOGICAL LM, LR
C
      IF(LM .OR. LR) WRITE(*,*)
C
      IF(LM) THEN
       IF(MATYP.EQ.1) WRITE(*,1100) MINF1
       IF(MATYP.EQ.2) WRITE(*,1100) MINF1, ' / sqrt(CL)'
       IF(MATYP.EQ.3) WRITE(*,1100) MINF1, ' / CL'
      ENDIF
C
      IF(LR) THEN
       IF(RETYP.EQ.1) WRITE(*,1200) REINF1
       IF(RETYP.EQ.2) WRITE(*,1200) REINF1, ' / sqrt(CL)'
       IF(RETYP.EQ.3) WRITE(*,1200) REINF1, ' / CL'
      ENDIF
C
      RETURN
C
 1100 FORMAT(1X,'M  =' , F10.4, A)
 1200 FORMAT(1X,'Re =' , G12.4, A)
      END ! MRSHOW



      SUBROUTINE NAMMOD(NAME,KDEL,KMOD0)
      CHARACTER*(*) NAME
C-------------------------------------------
C     Requests new modified NAME with 
C     version number in brackets, e.g.
C            NACA 0012  [5]
C
C     If bracketed index exists in NAME,
C        it is incremented by KDEL.
C     If no bracketed index exists, it 
C        is added with initial value KMOD0,
C        unless KMOD0 is negative in which 
C        case nothing is added.
C-------------------------------------------
      CHARACTER*48 NAMDEF
C
      CALL STRIP(NAME,NNAME)
      KBRACK1 = INDEX(NAME,'[')
      KBRACK2 = INDEX(NAME,']')
C
      NAMDEF = NAME(1:NNAME)
C
      IF(KBRACK1.NE.0 .AND. 
     &   KBRACK2.NE.0 .AND. KBRACK2-KBRACK1.GT.1) THEN
C----- brackets exist... get number, (go get user's input on READ error)
       READ(NAME(KBRACK1+1:KBRACK2-1),*,ERR=40) KMOD
       KMOD = IABS(KMOD)
       KMODP = MOD( KMOD+KDEL , 100 )
       IF(KBRACK1.GE.2) THEN
        NAME = NAME(1:KBRACK1-1)
       ELSE
        NAME = ' '
       ENDIF
       CALL STRIP(NAME,NNAME)
      ELSEIF(KMOD0.GT.0) THEN
       KMODP = MOD( KMOD0 , 100 )
      ELSE
       KMODP = 0
      ENDIF
C
      IF    (KMODP.GE.10) THEN
       NAMDEF = NAME(1:NNAME) // ' [  ]'
       WRITE(NAMDEF(NNAME+3:NNAME+4),1020) KMODP
 1020  FORMAT(I2)
      ELSEIF(KMODP.GE. 1) THEN
       NAMDEF = NAME(1:NNAME) // ' [ ]'
       WRITE(NAMDEF(NNAME+3:NNAME+3),1025) KMODP
 1025  FORMAT(I1)
      ENDIF
C
 40   WRITE(*,1040) NAMDEF
 1040 FORMAT(/' Enter airfoil name or <return> for default:  ',A)
      READ(*,1000) NAME
 1000 FORMAT(A)
      IF(NAME .EQ. ' ') NAME = NAMDEF
C
      RETURN
      END ! NAMMOD



      SUBROUTINE BLDUMP(FNAME1)
      INCLUDE 'XFOIL.INC'
      CHARACTER*(*) FNAME1
C
      CHARACTER*80 FILDEF
C
      CHARACTER*1 DELIM
      CHARACTER*256 LINE
C
      IF    (KDELIM.EQ.0) THEN
       DELIM = ' '
      ELSEIF(KDELIM.EQ.1) THEN
       DELIM = ','
      ELSEIF(KDELIM.EQ.2) THEN
       DELIM = CHAR(9)
      ELSE
       WRITE(*,*) '? Illegal delimiter.  Using blank.'
       DELIM = ' '
      ENDIF
C
 1000 FORMAT(50A)
C
      IF(FNAME1(1:1).NE.' ') THEN
       FNAME = FNAME1
      ELSE
C----- no argument... get it somehow
       IF(NPREFIX.GT.0) THEN
C------ offer default using existing prefix
        FILDEF = PREFIX(1:NPREFIX) // '.bl'
        WRITE(*,1100) FILDEF
 1100   FORMAT(/' Enter filename:  ', A)
        READ(*,1000) FNAME
        CALL STRIP(FNAME,NFN)
        IF(NFN.EQ.0) FNAME = FILDEF
       ELSE
C------ nothing available... just ask for filename
        CALL ASKS('Enter filename^',FNAME)
       ENDIF
      ENDIF
C
      LU = 19
      OPEN(LU,FILE=FNAME,STATUS='UNKNOWN')
      REWIND(LU)
C
      IF(KDELIM.EQ.0) THEN
       WRITE(LU,1000)
     & '#    s        x        y     Ue/Vinf    Dstar     Theta ',
     & '     Cf       H'
C         1.23456  0.23451  0.23451  0.23451  0.012345  0.001234  0.004123  10.512
      ELSE
       WRITE(LU,1000)
     & '#s'     ,DELIM,
     & 'x'      ,DELIM,
     & 'y'      ,DELIM,
     & 'Ue/Vinf',DELIM,
     & 'Dstar'  ,DELIM,
     & 'Theta'  ,DELIM,
     & 'Cf'     ,DELIM,
     & 'H'
      ENDIF
C
      CALL COMSET
      HSTINV = GAMM1*(MINF/QINF)**2 / (1.0 + 0.5*GAMM1*MINF**2)
C
      DO 10 I=1, N
        IS = 1
        IF(GAM(I) .LT. 0.0) IS = 2
C
        IF(LIPAN .AND. LVISC) THEN
          IF(IS.EQ.1) THEN
            IBL = IBLTE(IS) - I + 1
          ELSE
            IBL = IBLTE(IS) + I - N
          ENDIF
          DS = DSTR(IBL,IS)
          TH = THET(IBL,IS)
          CF =  TAU(IBL,IS)/(0.5*QINF**2)
          IF(TH.EQ.0.0) THEN
           H = 1.0
          ELSE
           H = DS/TH
          ENDIF
        ELSE
          DS = 0.
          TH = 0.
          CF = 0.
          H = 1.0
        ENDIF
        UE = (GAM(I)/QINF)*(1.0-TKLAM) / (1.0 - TKLAM*(GAM(I)/QINF)**2)
        AMSQ = UE*UE*HSTINV / (GAMM1*(1.0 - 0.5*UE*UE*HSTINV))
        CALL HKIN( H, AMSQ, HK, DUMMY, DUMMY)
C
        IF(KDELIM.EQ.0) THEN
         WRITE(LU,8500) S(I), X(I), Y(I), UE, DS, TH, CF, HK
 8500    FORMAT(1X, 4F9.5, 3F10.6, F10.3)
C
        ELSE
         WRITE(LINE,8510) 
     &     S(I),DELIM,
     &     X(I),DELIM,
     &     Y(I),DELIM,
     &     UE  ,DELIM,
     &     DS  ,DELIM,
     &     TH  ,DELIM,
     &     CF  ,DELIM,
     &     HK
 8510    FORMAT(1X, 4(F9.5,A), 3(F10.6,A), F10.3)
         CALL BSTRIP(LINE,NLINE)
         WRITE(LU,1000) LINE(1:NLINE)
        ENDIF
C
  10  CONTINUE
C
      IF(LWAKE) THEN
        IS = 2
        DO 20 I=N+1, N+NW
          IBL = IBLTE(IS) + I - N
          DS = DSTR(IBL,IS)
          TH = THET(IBL,IS)
          H = DS/TH
          CF = 0.
          UI = UEDG(IBL,IS)
          UE = (UI/QINF)*(1.0-TKLAM) / (1.0 - TKLAM*(UI/QINF)**2)
          AMSQ = UE*UE*HSTINV / (GAMM1*(1.0 - 0.5*UE*UE*HSTINV))
          CALL HKIN( H, AMSQ, HK, DUMMY, DUMMY)
C
          IF(KDELIM.EQ.0) THEN
           WRITE(LU,8500) S(I), X(I), Y(I), UE, DS, TH, CF, HK
C
          ELSE
           WRITE(LINE,8510) 
     &      S(I),DELIM,
     &      X(I),DELIM,
     &      Y(I),DELIM,
     &      UE  ,DELIM,
     &      DS  ,DELIM,
     &      TH  ,DELIM,
     &      CF  ,DELIM,
     &      HK
           CALL BSTRIP(LINE,NLINE)
           WRITE(LU,1000) LINE(1:NLINE)
          ENDIF
 20     CONTINUE
      ENDIF
C
      CLOSE(LU)
      RETURN
      END ! BLDUMP



      SUBROUTINE CPDUMP(FNAME1)
      INCLUDE 'XFOIL.INC'
      CHARACTER*(*) FNAME1
C
      CHARACTER*80 FILDEF
C
      CHARACTER*1 DELIM
      CHARACTER*128 LINE
C
      IF    (KDELIM.EQ.0) THEN
       DELIM = ' '
      ELSEIF(KDELIM.EQ.1) THEN
       DELIM = ','
      ELSEIF(KDELIM.EQ.2) THEN
       DELIM = CHAR(9)
      ELSE
       WRITE(*,*) '? Illegal delimiter.  Using blank.'
       DELIM = ' '
      ENDIF
C
 1000 FORMAT(8A)
C
      IF(FNAME1(1:1).NE.' ') THEN
       FNAME = FNAME1
      ELSE
C----- no argument... get it somehow
       IF(NPREFIX.GT.0) THEN
C------ offer default using existing prefix
        FILDEF = PREFIX(1:NPREFIX) // '.cp'
        WRITE(*,1100) FILDEF
 1100   FORMAT(/' Enter filename:  ', A)
        READ(*,1000) FNAME
        CALL STRIP(FNAME,NFN)
        IF(NFN.EQ.0) FNAME = FILDEF
       ELSE
C------ nothing available... just ask for filename
        CALL ASKS('Enter filename^',FNAME)
       ENDIF
      ENDIF
C
C
      LU = 19
      OPEN(LU,FILE=FNAME,STATUS='UNKNOWN')
      REWIND(LU)
C
      IF(KDELIM.EQ.0) THEN
       WRITE(LU,1000)
     &  '#      x          Cp  '
C            0.23451    0.23451
      ELSE
       WRITE(LU,1000)
     &  '#x', DELIM,
     &  'Cp'
C
      ENDIF
C
      CALL COMSET
C
      BETA = SQRT(1.0 - MINF**2)
      BFAC = 0.5*MINF**2 / (1.0 + BETA)
C
      DO 10 I=1, N
        CPINC = 1.0 - (GAM(I)/QINF)**2
        DEN = BETA + BFAC*CPINC
        CPCOM = CPINC / DEN
C
        IF(KDELIM.EQ.0) THEN
         WRITE(LU,8500) X(I), CPCOM
 8500    FORMAT(1X,2F11.5)
        ELSE
         WRITE(LINE,8510) 
     &    X(I) , DELIM,
     &    CPCOM
 8510    FORMAT(1X,2(F11.5,A))
         CALL BSTRIP(LINE,NLINE)
         WRITE(LU,1000) LINE(1:NLINE)
        ENDIF
  10  CONTINUE
C
      CLOSE(LU)
      RETURN
      END ! CPDUMP



      SUBROUTINE MHINGE
C----------------------------------------------------
C     Calculates the hinge moment of the flap about
C     (XOF,YOF) by integrating surface pressures.
C----------------------------------------------------
      INCLUDE 'XFOIL.INC'
C
      IF(.NOT.LFLAP) THEN
C
        CALL GETXYF(X,XP,Y,YP,S,N, TOPS,BOTS,XOF,YOF)
        LFLAP = .TRUE.
C
      ELSE
C
C------ find top and bottom y at hinge x location
        TOPS = XOF
        BOTS = S(N) - XOF
        CALL SINVRT(TOPS,XOF,X,XP,S,N)      
        CALL SINVRT(BOTS,XOF,X,XP,S,N)      
C
      ENDIF
C
      TOPX = SEVAL(TOPS,X,XP,S,N)
      TOPY = SEVAL(TOPS,Y,YP,S,N)
      BOTX = SEVAL(BOTS,X,XP,S,N)
      BOTY = SEVAL(BOTS,Y,YP,S,N)
C
C
      HMOM = 0.
      HFX  = 0.
      HFY  = 0.
C
C---- integrate pressures on top and bottom sides of flap
      DO 20 I=2, N
        IF(S(I-1).GE.TOPS .AND. S(I).LE.BOTS) GO TO 20
C
         DX = X(I) - X(I-1)
         DY = Y(I) - Y(I-1)
         XMID = 0.5*(X(I)+X(I-1)) - XOF
         YMID = 0.5*(Y(I)+Y(I-1)) - YOF
         IF(LVISC) THEN
          PMID = 0.5*(CPV(I) + CPV(I-1))
         ELSE
          PMID = 0.5*(CPI(I) + CPI(I-1))
         ENDIF
         HMOM = HMOM + PMID*(XMID*DX + YMID*DY)
         HFX  = HFX  - PMID* DY
         HFY  = HFY  + PMID* DX
   20 CONTINUE
C
C---- find S(I)..S(I-1) interval containing s=TOPS
      DO I=2, N
        IF(S(I).GT.TOPS) GO TO 31
      ENDDO
C
   31 CONTINUE
C---- add on top surface chunk TOPS..S(I-1),  missed in the DO 20 loop.
      DX = TOPX - X(I-1)
      DY = TOPY - Y(I-1)
      XMID = 0.5*(TOPX+X(I-1)) - XOF
      YMID = 0.5*(TOPY+Y(I-1)) - YOF
      IF(S(I) .NE. S(I-1)) THEN
       FRAC = (TOPS-S(I-1))/(S(I)-S(I-1))
      ELSE
       FRAC = 0.
      ENDIF
      IF(LVISC) THEN
       TOPP = CPV(I)*FRAC + CPV(I-1)*(1.0-FRAC)
       PMID = 0.5*(TOPP+CPV(I-1))
      ELSE
       TOPP = CPI(I)*FRAC + CPI(I-1)*(1.0-FRAC)
       PMID = 0.5*(TOPP+CPI(I-1))
      ENDIF
      HMOM = HMOM + PMID*(XMID*DX + YMID*DY)
      HFX  = HFX  - PMID* DY
      HFY  = HFY  + PMID* DX
C
C---- add on inside flap surface contribution from hinge to top surface
      DX = XOF - TOPX
      DY = YOF - TOPY
      XMID = 0.5*(TOPX+XOF) - XOF
      YMID = 0.5*(TOPY+YOF) - YOF
      HMOM = HMOM + PMID*(XMID*DX + YMID*DY)
      HFX  = HFX  - PMID* DY
      HFY  = HFY  + PMID* DX
C
C---- find S(I)..S(I-1) interval containing s=BOTS
      DO I=N, 2, -1
        IF(S(I-1).LT.BOTS) GO TO 41
      ENDDO
C
   41 CONTINUE
C---- add on bottom surface chunk BOTS..S(I),  missed in the DO 20 loop.
      DX = X(I) - BOTX
      DY = Y(I) - BOTY
      XMID = 0.5*(BOTX+X(I)) - XOF
      YMID = 0.5*(BOTY+Y(I)) - YOF
      IF(S(I) .NE. S(I-1)) THEN
       FRAC = (BOTS-S(I-1))/(S(I)-S(I-1))
      ELSE
       FRAC = 0.
      ENDIF
      IF(LVISC) THEN
       BOTP = CPV(I)*FRAC + CPV(I-1)*(1.0-FRAC)
       PMID = 0.5*(BOTP+CPV(I))
      ELSE
       BOTP = CPI(I)*FRAC + CPI(I-1)*(1.0-FRAC)
       PMID = 0.5*(BOTP+CPI(I))
      ENDIF
      HMOM = HMOM + PMID*(XMID*DX + YMID*DY)
      HFX  = HFX  - PMID* DY
      HFY  = HFY  + PMID* DX
C
C---- add on inside flap surface contribution from hinge to bottom surface
      DX = BOTX - XOF
      DY = BOTY - YOF
      XMID = 0.5*(BOTX+XOF) - XOF
      YMID = 0.5*(BOTY+YOF) - YOF
      HMOM = HMOM + PMID*(XMID*DX + YMID*DY)
      HFX  = HFX  - PMID* DY
      HFY  = HFY  + PMID* DX
C
C---- add on TE base thickness contribution
      DX = X(1) - X(N)
      DY = Y(1) - Y(N)
      XMID = 0.5*(X(1)+X(N)) - XOF
      YMID = 0.5*(Y(1)+Y(N)) - YOF
      IF(LVISC) THEN
       PMID = 0.5*(CPV(1)+CPV(N))
      ELSE
       PMID = 0.5*(CPI(1)+CPI(N))
      ENDIF
      HMOM = HMOM + PMID*(XMID*DX + YMID*DY)
      HFX  = HFX  - PMID* DY
      HFY  = HFY  + PMID* DX
C
      RETURN
      END ! MHINGE


      SUBROUTINE VPAR
C---------------------------------------------
C     Viscous parameter change menu routine.
C---------------------------------------------
      INCLUDE 'XFOIL.INC'
      INCLUDE 'BLPAR.INC'
      CHARACTER*4 COMAND
      CHARACTER*128 COMARG
C
      DIMENSION IINPUT(20)
      DIMENSION RINPUT(20)
      LOGICAL ERROR
C
C
 10   TURB = 100.0 * EXP( -(ACRIT + 8.43)/2.4 )
      WRITE(*,1200) XSTRIP(1), XSTRIP(2), ACRIT, TURB, VACCEL,
     &              SCCON, DUXCON, GACON, GBCON, CTCON, CTRCON, CTRCEX
 1200 FORMAT(/' Xtr/c     =', F8.4, '    top    side'
     &       /' Xtr/c     =', F8.4, '    bottom side'
     &       /' Ncrit     =', F8.2, '   (', F6.3, ' % turb. level )'
     &       /' Vacc      =', F8.4,
     &      //' Klag  =', F8.4,'     Uxwt  =', F8.2
     &       /' A     =', F8.4,'     B     =', F8.4,'       KCt =', F8.5
     &       /' CtiniK=', F8.4,'     CtiniX=', F8.4 )
C
C======================================================================
C---- start of user interaction loop
  500 CONTINUE
      CALL ASKC('..VPAR^',COMAND,COMARG)
C
      DO I=1, 20
        IINPUT(I) = 0
        RINPUT(I) = 0.0
      ENDDO
      NINPUT = 20
      CALL GETINT(COMARG,IINPUT,NINPUT,ERROR)
      NINPUT = 20
      CALL GETFLT(COMARG,RINPUT,NINPUT,ERROR)
C
C--------------------------------------------------------------
      IF(COMAND.EQ.'    ') THEN
       RETURN
C
C--------------------------------------------------------------
      ELSEIF(COMAND.EQ.'?   ') THEN
       WRITE(*,1050)
 1050  FORMAT(
     &  /'   <cr>    Return to OPER menu'
     &  /'   SHOW    Display viscous parameters'
     &  /'   XTR  rr Change trip positions Xtr/c'
     &  /'   N    r  Change critical amplification exponent Ncrit'
     &  /'   VACC r  Change Newton solution acceleration parameter'
     &  /'   INIT    BL initialization flag toggle'
     & //'   LAG     change lag equation constants'
     &  /'   GB      change G-beta constants'
     &  /'   CTR     change initial transition-Ctau constants'
     &  /'   REST    restore BL calibration to baseline')
C
C--------------------------------------------------------------
      ELSEIF(COMAND.EQ.'SHOW') THEN
       GO TO 10
C
C--------------------------------------------------------------
      ELSEIF(COMAND.EQ.'XTR ') THEN
       IF(LPACC .AND. LVISC) THEN
        WRITE(*,2100)
        GO TO 500
       ENDIF
       IF(NINPUT.GE.2) THEN
        XSTRIP(1) = RINPUT(1)
        XSTRIP(2) = RINPUT(2)
       ELSE
        CALL ASKR('Enter top    side Xtrip/c^',XSTRIP(1))
        CALL ASKR('Enter bottom side Xtrip/c^',XSTRIP(2))
       ENDIF
       LVCONV = .FALSE.
C
C--------------------------------------------------------------
      ELSEIF(COMAND.EQ.'N   ') THEN
       IF(LPACC .AND. LVISC) THEN
        WRITE(*,2100)
        GO TO 500
       ENDIF
       IF(NINPUT.GE.1) THEN
        ACRIT = RINPUT(1)
       ELSE
        CALL ASKR('Enter critical amplification ratio^',ACRIT)
       ENDIF
       LVCONV = .FALSE.
C
C--------------------------------------------------------------
      ELSEIF(COMAND.EQ.'VACC') THEN
       IF(NINPUT.GE.1) THEN
        VACCEL = RINPUT(1)
       ELSE
        CALL ASKR('Enter viscous acceleration parameter^',VACCEL)
       ENDIF
C
C--------------------------------------------------------------
      ELSEIF(COMAND.EQ.'INIT') THEN
       LBLINI = .NOT.LBLINI
       IF(.NOT.LBLINI) WRITE(*,*)'BLs will be initialized on next point'
       IF(     LBLINI) WRITE(*,*)'BLs are assumed to be initialized'
       IF(.NOT.LBLINI) LIPAN = .FALSE.
C
C--------------------------------------------------------------
      ELSEIF(COMAND.EQ.'LAG ') THEN
       IF(LPACC .AND. LVISC) THEN
        WRITE(*,2100)
        GO TO 500
       ENDIF
       IF(NINPUT.GE.2) THEN
        SCCON  = RINPUT(1)
        DUXCON = RINPUT(2)
       ELSE
        CALL ASKR('Enter shear lag constant^',SCCON)
        CALL ASKR('Enter shear lag UxEQ weight^',DUXCON)
       ENDIF
C
C--------------------------------------------------------------
      ELSEIF(COMAND.EQ.'GB  ') THEN
       IF(LPACC .AND. LVISC) THEN
        WRITE(*,2100)
        GO TO 500
       ENDIF
       IF(NINPUT.GE.2) THEN
        GACON = RINPUT(1)
        GBCON = RINPUT(2)
       ELSE
        CALL ASKR('Enter G-beta constant A^',GACON)
        CALL ASKR('Enter G-beta constant B^',GBCON)
       ENDIF
       CTCON = 0.5/(GACON**2 * GBCON)
C
C--------------------------------------------------------------
      ELSEIF(COMAND.EQ.'CTR ') THEN
       IF(LPACC .AND. LVISC) THEN
        WRITE(*,2100)
        GO TO 500
       ENDIF
       IF(NINPUT.GE.2) THEN
        CTRCON = RINPUT(1)
        CTRCEX = RINPUT(2)
       ELSE
        CALL ASKR('Enter initial-Ctau constant^',CTRCON)
        CALL ASKR('Enter initial-Ctau exponent^',CTRCEX)
       ENDIF
C
C--------------------------------------------------------------
      ELSEIF(COMAND.EQ.'CFAC') THEN
       IF(NINPUT.GE.1) THEN
        CFFAC = RINPUT(1)
       ELSE
        CALL ASKR('Enter Cf scaling factor^',CFFAC)
       ENDIF
C
C--------------------------------------------------------------
      ELSEIF(COMAND.EQ.'REST') THEN
       IF(LPACC .AND. LVISC) THEN
        WRITE(*,2100)
        GO TO 500
       ENDIF
       CALL BLPINI
C
C--------------------------------------------------------------
      ELSE
       WRITE(*,1000) COMAND
 1000  FORMAT(1X,A4,' command not recognized.  Type a "?" for list')
C
      ENDIF
C
      GO TO 500
C--------------------------------------------
 2100 FORMAT(/' * Polar is being accumulated.'
     &       /' * Cannot change its parameters in midstream.')
      END ! VPAR




      SUBROUTINE SPECAL
C-----------------------------------
C     Converges to specified alpha.
C-----------------------------------
      INCLUDE 'XFOIL.INC'
      REAL MINF_CLM, MSQ_CLM
C
C---- calculate surface vorticity distributions for alpha = 0, 90 degrees
      IF(.NOT.LGAMU .OR. .NOT.LQAIJ) CALL GGCALC
C
      COSA = COS(ALFA)
      SINA = SIN(ALFA)
C
C---- superimpose suitably weighted  alpha = 0, 90  distributions
      DO 50 I=1, N
        GAM(I)   =  COSA*GAMU(I,1) + SINA*GAMU(I,2)
        GAM_A(I) = -SINA*GAMU(I,1) + COSA*GAMU(I,2)
   50 CONTINUE
      PSIO = COSA*GAMU(N+1,1) + SINA*GAMU(N+1,2)
C
      CALL TECALC
      CALL QISET
C
C---- set initial guess for the Newton variable CLM
      CLM = 1.0
C
C---- set corresponding  M(CLM), Re(CLM)
      CALL MRCL(CLM,MINF_CLM,REINF_CLM)
      CALL COMSET
C
C---- set corresponding CL(M)
      CALL CLCALC(N,X,Y,GAM,GAM_A,ALFA,MINF,QINF, XCMREF,YCMREF,
     &            CL,CM,CDP, CL_ALF,CL_MSQ)
C
C---- iterate on CLM
      DO 100 ITCL=1, 20
C
        MSQ_CLM = 2.0*MINF*MINF_CLM
        DCLM = (CL - CLM)/(1.0 - CL_MSQ*MSQ_CLM)
C
        CLM1 = CLM
        RLX = 1.0
C
C------ under-relaxation loop to avoid driving M(CL) above 1
        DO 90 IRLX=1, 12
C
          CLM = CLM1 + RLX*DCLM
C
C-------- set new freestream Mach M(CLM)
          CALL MRCL(CLM,MINF_CLM,REINF_CLM)
C
C-------- if Mach is OK, go do next Newton iteration
          IF(MATYP.EQ.1 .OR. MINF.EQ.0.0 .OR. MINF_CLM.NE.0.0) GO TO 91
C
          RLX = 0.5*RLX
   90   CONTINUE
   91   CONTINUE
C
C------ set new CL(M)
        CALL COMSET
        CALL CLCALC(N,X,Y,GAM,GAM_A,ALFA,MINF,QINF, XCMREF,YCMREF,
     &              CL,CM,CDP,CL_ALF,CL_MSQ)
C
        IF(ABS(DCLM).LE.1.0E-6) GO TO 110
C
  100 CONTINUE
      WRITE(*,*) 'SPECAL:  Minf convergence failed'
  110 CONTINUE
C
C---- set final Mach, CL, Cp distributions, and hinge moment
      CALL MRCL(CL,MINF_CL,REINF_CL)
      CALL COMSET
      CALL CLCALC(N,X,Y,GAM,GAM_A,ALFA,MINF,QINF, XCMREF,YCMREF,
     &            CL,CM,CDP, CL_ALF,CL_MSQ)
      CALL CPCALC(N,QINV,QINF,MINF,CPI)
      IF(LVISC) THEN
       CALL CPCALC(N+NW,QVIS,QINF,MINF,CPV)
       CALL CPCALC(N+NW,QINV,QINF,MINF,CPI)
      ELSE
       CALL CPCALC(N,QINV,QINF,MINF,CPI)
      ENDIF
      IF(LFLAP) CALL MHINGE
C
      RETURN
      END ! SPECAL
 
 
      SUBROUTINE SPECCL
C-----------------------------------------
C     Converges to specified inviscid CL.
C-----------------------------------------
      INCLUDE 'XFOIL.INC'
C
C---- calculate surface vorticity distributions for alpha = 0, 90 degrees
      IF(.NOT.LGAMU .OR. .NOT.LQAIJ) CALL GGCALC
C
C---- set freestream Mach from specified CL -- Mach will be held fixed
      CALL MRCL(CLSPEC,MINF_CL,REINF_CL)
      CALL COMSET
C
C---- current alpha is the initial guess for Newton variable ALFA
      COSA = COS(ALFA)
      SINA = SIN(ALFA)
      DO 10 I=1, N
        GAM(I)   =  COSA*GAMU(I,1) + SINA*GAMU(I,2)
        GAM_A(I) = -SINA*GAMU(I,1) + COSA*GAMU(I,2)
   10 CONTINUE
      PSIO = COSA*GAMU(N+1,1) + SINA*GAMU(N+1,2)
C
C---- get corresponding CL, CL_alpha, CL_Mach
      CALL CLCALC(N,X,Y,GAM,GAM_A,ALFA,MINF,QINF, XCMREF,YCMREF,
     &            CL,CM,CDP, CL_ALF,CL_MSQ)
C
C---- Newton loop for alpha to get specified inviscid CL
      DO 100 ITAL=1, 20
C
        DALFA = (CLSPEC - CL) / CL_ALF
        RLX = 1.0
C
        ALFA = ALFA + RLX*DALFA
C
C------ set new surface speed distribution
        COSA = COS(ALFA)
        SINA = SIN(ALFA)
        DO 40 I=1, N
          GAM(I)   =  COSA*GAMU(I,1) + SINA*GAMU(I,2)
          GAM_A(I) = -SINA*GAMU(I,1) + COSA*GAMU(I,2)
   40   CONTINUE
        PSIO = COSA*GAMU(N+1,1) + SINA*GAMU(N+1,2)
C
C------ set new CL(alpha)
        CALL CLCALC(N,X,Y,GAM,GAM_A,ALFA,MINF,QINF, XCMREF,YCMREF,
     &              CL,CM,CDP,CL_ALF,CL_MSQ)
C
        IF(ABS(DALFA).LE.1.0E-6) GO TO 110
  100 CONTINUE
      WRITE(*,*) 'SPECCL:  CL convergence failed'
  110 CONTINUE
C
C---- set final surface speed and Cp distributions
      CALL TECALC
      CALL QISET
      IF(LVISC) THEN
       CALL CPCALC(N+NW,QVIS,QINF,MINF,CPV)
       CALL CPCALC(N+NW,QINV,QINF,MINF,CPI)
      ELSE
       CALL CPCALC(N,QINV,QINF,MINF,CPI)
      ENDIF
      IF(LFLAP) CALL MHINGE
C
      RETURN
      END ! SPECCL


      SUBROUTINE VISCAL(NITER1)
C----------------------------------------
C     Converges viscous operating point
C----------------------------------------
      INCLUDE 'XFOIL.INC'
C
C---- convergence tolerance
      DATA EPS1 / 1.0E-4 /
C
      NITER = NITER1
C
C---- calculate wake trajectory from current inviscid solution if necessary
      IF(.NOT.LWAKE) THEN
       CALL XYWAKE
      ENDIF
C
C---- set velocities on wake from airfoil vorticity for alpha=0, 90
      CALL QWCALC
C
C---- set velocities on airfoil and wake for initial alpha
      CALL QISET
C
      IF(.NOT.LIPAN) THEN
C
       IF(LBLINI) CALL GAMQV
C
C----- locate stagnation point arc length position and panel index
       CALL STFIND
C
C----- set  BL position -> panel position  pointers
       CALL IBLPAN
C
C----- calculate surface arc length array for current stagnation point location
       CALL XICALC
C
C----- set  BL position -> system line  pointers
       CALL IBLSYS
C
      ENDIF
C
C---- set inviscid BL edge velocity UINV from QINV
      CALL UICALC
C
      IF(.NOT.LBLINI) THEN
C
C----- set initial Ue from inviscid Ue
       DO IBL=1, NBL(1)
         UEDG(IBL,1) = UINV(IBL,1)
       ENDDO
C
       DO IBL=1, NBL(2)
         UEDG(IBL,2) = UINV(IBL,2)
       ENDDO
C
      ENDIF
C
      IF(LVCONV) THEN
C----- set correct CL if converged point exists
       CALL QVFUE
       IF(LVISC) THEN
        CALL CPCALC(N+NW,QVIS,QINF,MINF,CPV)
        CALL CPCALC(N+NW,QINV,QINF,MINF,CPI)
       ELSE
        CALL CPCALC(N,QINV,QINF,MINF,CPI)
       ENDIF
       CALL GAMQV
       CALL CLCALC(N,X,Y,GAM,GAM_A,ALFA,MINF,QINF, XCMREF,YCMREF,
     &             CL,CM,CDP, CL_ALF,CL_MSQ)
       CALL CDCALC
      ENDIF
C
C---- set up source influence matrix if it doesn't exist
      IF(.NOT.LWDIJ .OR. .NOT.LADIJ) CALL QDCALC
C
C---- Newton iteration for entire BL solution
      IF(NITER.EQ.0) CALL ASKI('Enter number of iterations^',NITER)
      WRITE(*,*)
      WRITE(*,*) 'Solving BL system ...'
      DO 1000 ITER=1, NITER
C
C------ fill Newton system for BL variables
        CALL SETBL
C
C------ solve Newton system with custom solver
        CALL BLSOLV
C
C------ update BL variables
        CALL UPDATE
C
        IF(LALFA) THEN
C------- set new freestream Mach, Re from new CL
         CALL MRCL(CL,MINF_CL,REINF_CL)
         CALL COMSET
        ELSE
C------- set new inviscid speeds QINV and UINV for new alpha
         CALL QISET
         CALL UICALC
        ENDIF
C
C------ calculate edge velocities QVIS(.) from UEDG(..)
        CALL QVFUE
C
C------ set GAM distribution from QVIS
        CALL GAMQV
C
C------ relocate stagnation point
        CALL STMOVE
C
C------ set updated CL,CD
        CALL CLCALC(N,X,Y,GAM,GAM_A,ALFA,MINF,QINF, XCMREF,YCMREF,
     &              CL,CM,CDP,CL_ALF,CL_MSQ)
        CALL CDCALC
C
C------ display changes and test for convergence
        IF(RLX.LT.1.0) 
     &   WRITE(*,2000) ITER, RMSBL, RMXBL, VMXBL,IMXBL,ISMXBL,RLX
        IF(RLX.EQ.1.0) 
     &   WRITE(*,2010) ITER, RMSBL, RMXBL, VMXBL,IMXBL,ISMXBL
         CDP = CD - CDF
         WRITE(*,2020) ALFA/DTOR, CL, CM, CD, CDF, CDP
C
        IF(RMSBL .LT. EPS1) THEN
         LVCONV = .TRUE.
         AVISC = ALFA
         MVISC = MINF
         GO TO 90
        ENDIF
C
 1000 CONTINUE
      WRITE(*,*) 'VISCAL:  Convergence failed'
C
   90 CONTINUE
      CALL CPCALC(N+NW,QINV,QINF,MINF,CPI)
      CALL CPCALC(N+NW,QVIS,QINF,MINF,CPV)
      IF(LFLAP) CALL MHINGE
      RETURN
C....................................................................
 2000   FORMAT
     &   (/1X,I3,'   rms: ',E10.4,'   max: ',E10.4,3X,A1,' at ',I4,I3,
     &     '   RLX:',F6.3)
 2010   FORMAT
     &   (/1X,I3,'   rms: ',E10.4,'   max: ',E10.4,3X,A1,' at ',I4,I3)
 2020   FORMAT
     &   ( 1X,3X,'   a =', F7.3,'      CL =',F8.4  /
     &     1X,3X,'  Cm =', F8.4, '     CD =',F9.5,
     &           '   =>   CDf =',F9.5,'    CDp =',F9.5)
      END ! VISCAL


      subroutine dcpout
      include 'XFOIL.INC'
c
c     Computes and writes upper and lower-surface 
c     Cp values at two specified x locations
c
c
      x1 = 0.05
      x2 = 0.15
c
      lu = 60
      open(lu,file='dcp.out',status='old',access='append',err=10)
      go to 20
c
 10   continue
      open(lu,file='dcp.out',status='new')
      write(lu,*) '#  ', name
      write(lu,*) '# alpha   CL       ',
     &            ' Cpl05     Cpu05     dCp05    ',
     &            ' Cpl15     Cpu15     dCp15    '
 20   continue
c
      call spline(cpv,w1,s,n)
c
      su1 = sle + x1*(s(1)-sle)
      sl1 = sle + x1*(s(n)-sle)
      su2 = sle + x2*(s(1)-sle)
      sl2 = sle + x2*(s(n)-sle)
c
      call sinvrt(sl1,x1,x,xp,s,n)
      call sinvrt(su1,x1,x,xp,s,n)
      call sinvrt(sl2,x2,x,xp,s,n)
      call sinvrt(su2,x2,x,xp,s,n)
c
      cpl1 = seval(sl1,cpv,w1,s,n)
      cpu1 = seval(su1,cpv,w1,s,n)
      cpl2 = seval(sl2,cpv,w1,s,n)
      cpu2 = seval(su2,cpv,w1,s,n)
c
      write(lu,1200) alfa/dtor, cl, 
     &               cpl1, cpu1, cpl1-cpu1,
     &               cpl2, cpu2, cpl2-cpu2

 1200 format(1x, f7.3, f9.4, 8f10.5)
c
      close(lu)
c
      return
      end