~fluidity-core/fluidity/embedded_models

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
!!$BSD Software License
!!$
!!$Pertains to ARPACK and P_ARPACK
!!$
!!$Copyright (c) 1996-2008 Rice University.  
!!$Developed by D.C. Sorensen, R.B. Lehoucq, C. Yang, and K. Maschhoff.
!!$All rights reserved.
!!$
!!$Redistribution and use in source and binary forms, with or without
!!$modification, are permitted provided that the following conditions are
!!$met:
!!$
!!$- Redistributions of source code must retain the above copyright
!!$  notice, this list of conditions and the following disclaimer. 
!!$  
!!$- Redistributions in binary form must reproduce the above copyright
!!$  notice, this list of conditions and the following disclaimer listed
!!$  in this license in the documentation and/or other materials
!!$  provided with the distribution.
!!$  
!!$- Neither the name of the copyright holders nor the names of its
!!$  contributors may be used to endorse or promote products derived from
!!$  this software without specific prior written permission.
!!$  
!!$THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
!!$"AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT  
!!$LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
!!$A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT 
!!$OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
!!$SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
!!$LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
!!$DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
!!$THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT  
!!$(INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
!!$OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 
!!$

#include "fdebug.h"

module snapsvd_module
  use vector_tools 
  implicit none

contains 
  subroutine snapsvd(m,n,snapmatrix,nsvd_total,nsvd,usv,svdval)

    !   
    !     This code shows how to use ARPACK to find a few of the
    !     largest singular values(sigma) and corresponding right singular 
    !     vectors (v) for the the matrix A by solving the symmetric problem:
    !          
    !                        (A'*A)*v = sigma*v
    ! 
    !     where A is an m by n real matrix.
    !
    !     This code may be easily modified to estimate the 2-norm
    !     condition number  largest(sigma)/smallest(sigma) by setting
    !     which = 'BE' below.  This will ask for a few of the smallest
    !     and a few of the largest singular values simultaneously.
    !     The condition number could then be estimated by taking
    !     the ratio of the largest and smallest singular values.
    !
    !     This formulation is appropriate when  m  .ge.  n.
    !     Reverse the roles of A and A' in the case that  m .le. n.
    !
    !     The main points illustrated here are 
    !
    !        1) How to declare sufficient memory to find NEV 
    !           largest singular values of A .  
    !
    !        2) Illustration of the reverse communication interface 
    !           needed to utilize the top level ARPACK routine DSAUPD 
    !           that computes the quantities needed to construct
    !           the desired singular values and vectors(if requested).
    !
    !        3) How to extract the desired singular values and vectors
    !           using the ARPACK routine DSEUPD.
    !
    !        4) How to construct the left singular vectors U from the 
    !           right singular vectors V to obtain the decomposition
    !
    !                        A*V = U*S
    !
    !           where S = diag(sigma_1, sigma_2, ..., sigma_k).
    !
    !     The only thing that must be supplied in order to use this
    !     routine on your problem is to change the array dimensions 
    !     appropriately, to specify WHICH singular values you want to 
    !     compute and to supply a the matrix-vector products 
    !
    !                         w <-  Ax
    !                         y <-  A'w
    !
    !     in place of the calls  to AV( ) and ATV( ) respectively below.  
    !
    !     Further documentation is available in the header of DSAUPD
    !     which may be found in the SRC directory.
    !
    !     This codes implements
    !
    !\Example-1
    !     ... Suppose we want to solve A'A*v = sigma*v in regular mode,
    !         where A is derived from the simplest finite difference 
    !         discretization of the 2-dimensional kernel  K(s,t)dt  where
    !
    !                 K(s,t) =  s(t-1)   if 0 .le. s .le. t .le. 1,
    !                           t(s-1)   if 0 .le. t .lt. s .le. 1. 
    !
    !         See subroutines AV  and ATV for details.
    !     ... OP = A'*A  and  B = I.
    !     ... Assume "call av (n,x,y)" computes y = A*x
    !     ... Assume "call atv (n,y,w)" computes w = A'*y
    !     ... Assume exact shifts are used
    !     ...
    !
    !\BeginLib
    !
    !\Routines called:
    !     dsaupd  ARPACK reverse communication interface routine.
    !     dseupd  ARPACK routine that returns Ritz values and (optionally)
    !             Ritz vectors.
    !     dnrm2   Level 1 BLAS that computes the norm of a vector.
    !     daxpy   Level 1 BLAS that computes y <- alpha*x+y.
    !     dscal   Level 1 BLAS thst computes x <- x*alpha.
    !     dcopy   Level 1 BLAS thst computes y <- x.
    !
    !-----------------------------------------------------------------------
    !
    !     %------------------------------------------------------%
    !     | Storage Declarations:                                |
    !     |                                                      |
    !     | It is assumed that A is M by N with M .ge. N.        |
    !     |                                                      |
    !     | The maximum dimensions for all arrays are            |
    !     | set here to accommodate a problem size of            |
    !     | M .le. MAXM  and  N .le. MAXN                        |
    !     |                                                      |
    !     | The NEV right singular vectors will be computed in   |
    !     | the N by NCV array V.                                |
    !     |                                                      |
    !     | The NEV left singular vectors will be computed in    |
    !     | the M by NEV array U.                                |
    !     |                                                      |
    !     | NEV is the number of singular values requested.      |
    !     |     See specifications for ARPACK usage below.       |
    !     |                                                      |
    !     | NCV is the largest number of basis vectors that will |
    !     |     be used in the Implicitly Restarted Arnoldi      |
    !     |     Process.  Work per major iteration is            |
    !     |     proportional to N*NCV*NCV.                       |
    !     |                                                      |
    !     | You must set:                                        |
    !     |                                                      |
    !     | MAXM:   Maximum number of rows of the A allowed.     |
    !     | MAXN:   Maximum number of columns of the A allowed.  |
    !     | MAXNEV: Maximum NEV allowed                          |
    !     | MAXNCV: Maximum NCV allowed                          |
    !     %------------------------------------------------------%
    !

    use FLDebug

    !include 'paramnew.h'

    integer   maxm, maxn, maxnev, maxncv, ldv, ldu, m, n, i
    integer   ndigit,logfil, msgets, msaitr, msapps, msaupd, msaup2, mseigt, mseupd

    !        parameter (maxm = istate,maxn=nrsnapshots,maxnev=nrsnapshots-1,    &
    !      maxncv=nrsnapshots, ldu = maxm, ldv=maxn )
    !     %--------------%
    !     | Local Arrays |
    !     %--------------%
    !
    !real v(ldv,maxncv), u(ldu, maxnev),            & 
    !                 workl(maxncv*(maxncv+8)), workd(3*maxn),  &
    !                 s(maxncv,2), resid(maxn), ax(maxm)
    !logical          select(maxncv)
    integer          iparam(11), ipntr(11)

    REAL, ALLOCATABLE, DIMENSION(:,:):: v,u,s
    REAL, ALLOCATABLE, DIMENSION(:)::   workl,workd,resid,ax
    logical, ALLOCATABLE, DIMENSION(:)::   select

    !
    !     %---------------%
    !     | Local Scalars |
    !     %---------------%
    !
    character        bmat*1, which*2
    integer          ido, nev, ncv, lworkl, info, ierr,   &
         j, ishfts, maxitr, mode1, nconv
    logical          rvec
    real toll, sigma, temp
    !
    !     %------------%
    !     | Parameters |
    !     %------------%
    !
    real one, zero
    parameter        (one = 1.0D+0, zero = 0.0D+0)


    !ccccc declare input
    real snapmatrix(m,n)
    integer nsvd,nsvd_total
    !cccccccc  declare output
    real usv(m,nsvd),svdval(nsvd) ! left svd, svdval
    !
    !  
    !     %-----------------------------%
    !     | BLAS & LAPACK routines used |
    !     %-----------------------------%
    !
    real dnrm2

    !
    !  %-------------%
    !  | local       |
    !  %-------------%

    REAL TOTAL_E,PARTIAL_E,PERCENT_E

!    external         dnrm2, daxpy, dcopy, dscal
    !
    !     %-----------------------%
    !     | Executable Statements |
    !     %-----------------------%
    !
    !     %-------------------------------------------------%
    !     | The following include statement and assignments |
    !     | initiate trace output from the internal         |
    !     | actions of ARPACK.  See debug.doc in the        |
    !     | DOCUMENTS directory for usage.  Initially, the  |
    !     | most useful information will be a breakdown of  |
    !     | time spent in the various stages of computation |
    !     | given by setting msaupd = 1.                    |
    !     %-------------------------------------------------%
    !
    !      include 'debug1.h'
#ifdef HAVE_LIBARPACK

    maxm = m
    maxn=n
    maxnev=n-1
    maxncv=n
    ldu = maxm
    ldv=maxn

    ALLOCATE(v(ldv,maxncv))
    ALLOCATE(u(ldu,maxnev))
    ALLOCATE(s(maxncv,2))
    ALLOCATE(workl(maxncv*(maxncv+8)))
    ALLOCATE(workd(3*maxn))
    ALLOCATE(resid(maxn))
    ALLOCATE(ax(maxm))
    ALLOCATE(select(maxncv))

    ewrite(3,*) 'in snapsvd.F'
    ewrite(3,*) M,maxm
    !      maxm=m
    !      ldu = maxm
    ndigit = -3
    logfil = 6
    msgets = 0
    msaitr = 0 
    msapps = 0
    msaupd = 1
    msaup2 = 0
    mseigt = 0
    mseupd = 0
    !
    !     %-------------------------------------------------%
    !     | The following sets dimensions for this problem. |
    !     %-------------------------------------------------%
    !
    !      if (n .ne. nrsnapshots) then
    !       write(*,*) 'dimension error in svd routine',n,nrsnapshots
    !       stop
    !      end if
    !
    !     %------------------------------------------------%
    !     | Specifications for ARPACK usage are set        | 
    !     | below:                                         |
    !     |                                                |
    !     |    1) NEV = 4 asks for 4 singular values to be |  
    !     |       computed.                                | 
    !     |                                                |
    !     |    2) NCV = 20 sets the length of the Arnoldi  |
    !     |       factorization                            |
    !     |                                                |
    !     |    3) This is a standard problem               |
    !     |         (indicated by bmat  = 'I')             |
    !     |                                                |
    !     |    4) Ask for the NEV singular values of       |
    !     |       largest magnitude                        |
    !     |         (indicated by which = 'LM')            |
    !     |       See documentation in DSAUPD for the      |
    !     |       other options SM, BE.                    | 
    !     |                                                |
    !     | Note: NEV and NCV must satisfy the following   |
    !     |       conditions:                              |
    !     |                 NEV <= MAXNEV,                 |
    !     |             NEV + 1 <= NCV <= MAXNCV           |
    !     %------------------------------------------------%
    !
    nev   = nsvd_total
    ncv   = n
    bmat  = 'I'
    which = 'LM'
    !
    if ( n .gt. maxn ) then
       ewrite(3,*)  ' ERROR with _SVD: N is greater than MAXN '
       go to 9000
    else if ( m .gt. maxm ) then
       ewrite(3,*)  ' ERROR with _SVD: M is greater than MAXM ',M,MAXM
       go to 9000
    else if ( nev .gt. maxnev ) then
       ewrite(3,*)  ' ERROR with _SVD: NEV is greater than MAXNEV '
       go to 9000
    else if ( ncv .gt. maxncv ) then
       ewrite(3,*)  ' ERROR with _SVD: NCV is greater than MAXNCV '
       go to 9000
    end if
    !
    !     %-----------------------------------------------------%
    !     | Specification of stopping rules and initial         |
    !     | conditions before calling DSAUPD                    |
    !     |                                                     |
    !     |           abs(sigmaC - sigmaT) < TOL*abs(sigmaC)    |
    !     |               computed   true                       |
    !     |                                                     |
    !     |      If TOL .le. 0,  then TOL <- macheps            |
    !     |              (machine precision) is used.           |
    !     |                                                     |
    !     | IDO  is the REVERSE COMMUNICATION parameter         |
    !     |      used to specify actions to be taken on return  |
    !     |      from DSAUPD. (See usage below.)                |
    !     |                                                     |
    !     |      It MUST initially be set to 0 before the first |
    !     |      call to DSAUPD.                                | 
    !     |                                                     |
    !     | INFO on entry specifies starting vector information |
    !     |      and on return indicates error codes            |
    !     |                                                     |
    !     |      Initially, setting INFO=0 indicates that a     | 
    !     |      random starting vector is requested to         |
    !     |      start the ARNOLDI iteration.  Setting INFO to  |
    !     |      a nonzero value on the initial call is used    |
    !     |      if you want to specify your own starting       |
    !     |      vector (This vector must be placed in RESID.)  | 
    !     |                                                     |
    !     | The work array WORKL is used in DSAUPD as           | 
    !     | workspace.  Its dimension LWORKL is set as          |
    !     | illustrated below.                                  |
    !     %-----------------------------------------------------%
    !
    lworkl = ncv*(ncv+8)
    toll = 0.0d-3 ! 1.d-3 
    info = 0
    ido = 0
    !
    !     %---------------------------------------------------%
    !     | Specification of Algorithm Mode:                  |
    !     |                                                   |
    !     | This program uses the exact shift strategy        |
    !     | (indicated by setting IPARAM(1) = 1.)             |
    !     | IPARAM(3) specifies the maximum number of Arnoldi |
    !     | iterations allowed.  Mode 1 of DSAUPD is used     |
    !     | (IPARAM(7) = 1). All these options can be changed |
    !     | by the user. For details see the documentation in |
    !     | DSAUPD.                                           |
    !     %---------------------------------------------------%
    !
    ishfts = 1
    maxitr = n
    mode1 = 1
    !
    iparam(1) = ishfts
    !                
    iparam(3) = maxitr
    !                  
    iparam(7) = mode1
    !
    !     %------------------------------------------------%
    !     | M A I N   L O O P (Reverse communication loop) |
    !     %------------------------------------------------%
    !
10  continue
    !
    !        %---------------------------------------------%
    !        | Repeatedly call the routine DSAUPD and take | 
    !        | actions indicated by parameter IDO until    |
    !        | either convergence is indicated or maxitr   |
    !        | has been exceeded.                          |
    !        %---------------------------------------------%
    !
#ifdef DOUBLEP
    call dsaupd ( ido, bmat, n, which, nev, toll, resid,      &
         ncv, v, ldv, iparam, ipntr, workd, workl,   &
         lworkl, info )
#else
    call ssaupd ( ido, bmat, n, which, nev, toll, resid,      &
         ncv, v, ldv, iparam, ipntr, workd, workl,   &
         lworkl, info )
#endif
    !
    if (ido .eq. -1 .or. ido .eq. 1) then
       !
       !           %---------------------------------------%
       !           | Perform matrix vector multiplications |
       !           |              w <--- A*x       (av())  |
       !           |              y <--- A'*w      (atv()) |
       !           | The user should supply his/her own    |
       !           | matrix vector multiplication routines |
       !           | here that takes workd(ipntr(1)) as    |
       !           | the input, and returns the result in  |
       !           | workd(ipntr(2)).                      |
       !           %---------------------------------------%
       !
       call av (m,n,snapmatrix,workd(ipntr(1)),ax) 
       call atv (m,n,snapmatrix,ax,workd(ipntr(2)))
       !
       !           %-----------------------------------------%
       !           | L O O P   B A C K to call DSAUPD again. |
       !           %-----------------------------------------%
       !
       go to 10
       !
    end if
    !
    !     %----------------------------------------%
    !     | Either we have convergence or there is |
    !     | an error.                              |
    !     %----------------------------------------%
    !
    if ( info .lt. 0 ) then
       !
       !        %--------------------------%
       !        | Error message. Check the |
       !        | documentation in DSAUPD. |
       !        %--------------------------%
       !
       ewrite(3,*)  ' '
       ewrite(3,*)  ' Error with _saupd, info = ', info
       ewrite(3,*)  ' Check documentation in _saupd '
       ewrite(3,*)  ' '
       !
    else 
       !
       !        %--------------------------------------------%
       !        | No fatal errors occurred.                  |
       !        | Post-Process using DSEUPD.                 |
       !        |                                            |
       !        | Computed singular values may be extracted. |  
       !        |                                            |
       !        | Singular vectors may also be computed now  |
       !        | if desired.  (indicated by rvec = .true.)  | 
       !        |                                            |
       !        | The routine DSEUPD now called to do this   |
       !        | post processing                            | 
       !        %--------------------------------------------%
       !           
       rvec = .true.
       !
#ifdef DOUBLEP
       call dseupd ( rvec, 'All', select, s, v, ldv, sigma,   &
            bmat, n, which, nev, toll, resid, ncv, v, ldv,    &
            iparam, ipntr, workd, workl, lworkl, ierr )
#else
       call sseupd ( rvec, 'All', select, s, v, ldv, sigma,   &
            bmat, n, which, nev, toll, resid, ncv, v, ldv,    &
            iparam, ipntr, workd, workl, lworkl, ierr )
#endif
       !
       !        %-----------------------------------------------%
       !        | Singular values are returned in the first     |
       !        | column of the two dimensional array S         |
       !        | and the corresponding right singular vectors  | 
       !        | are returned in the first NEV columns of the  |
       !        | two dimensional array V as requested here.    |
       !        %-----------------------------------------------%
       !

       ewrite(3,*) ipntr(7)
       ewrite(3,*) ipntr(7)+ncv
       ewrite(3,*) (workl(ipntr(7)+ncv+j-1),j=1,ncv)

       TOTAL_E=0.0
       DO J=1,NCV
          IF(workl(ipntr(7)+ncv+j-1).ge.0.0) THEN
             TOTAL_E=TOTAL_E+sqrt(workl(ipntr(7)+ncv+j-1))
          ENDIF
       ENDDO

       if ( ierr .ne. 0) then
          !
          !           %------------------------------------%
          !           | Error condition:                   |
          !           | Check the documentation of DSEUPD. |
          !           %------------------------------------%
          !
          ewrite(3,*)  ' '
          ewrite(3,*)  ' Error with _seupd, info = ', ierr
          ewrite(3,*)  ' Check the documentation of _seupd. '
          ewrite(3,*)  ' '
          !
       else
          !
          nconv =  iparam(5)
          ewrite(3,*) 'nconv',nconv
          PARTIAL_E=0.0
          do 20 j=1, nconv
             s(j,1) = sqrt(s(j,1))
             if(j.le.nsvd) then
                PARTIAL_E=PARTIAL_E+s(j,1)
             endif
             !
             !              %-----------------------------%
             !              | Compute the left singular   |
             !              | vectors from the formula    |
             !              |                             |
             !              |     u = Av/sigma            |
             !              |                             |
             !              | u should have norm 1 so     |
             !              | divide by norm(Av) instead. |
             !              %-----------------------------%
             !
             call av(m, n, snapmatrix, v(1,j), ax)
             u(:,j)=ax
             !call dcopy(m, ax, 1, u(1,j), 1)
             temp = one/norm2(u(:,j))
             !temp = one/dnrm2(m, u(1,j), 1)
             u(:,j)=u(:,j)*temp
             !call dscal(m, temp, u(1,j), 1)
             !
             !              %---------------------------%
             !              |                           |
             !              | Compute the residual norm |
             !              |                           |
             !              |   ||  A*v - sigma*u ||    |
             !              |                           |
             !              | for the NCONV accurately  |
             !              | computed singular values  |
             !              | and vectors.  (iparam(5)  |
             !              | indicates how many are    |
             !              | accurate to the requested |
             !              | tolerance).               |
             !              | Store the result in 2nd   |
             !              | column of array S.        |
             !              %---------------------------%
             !
             !call daxpy(m, -s(j,1), u(1,j), 1, ax, 1)
             ax=ax-s(j,1)*u(:,j)
             s(j,2) = norm2(ax)
             !s(j,2) = dnrm2(m, ax, 1)
             !
20        end do

          PERCENT_E = 100.0*PARTIAL_E/TOTAL_E
          ewrite(3,*)  PARTIAL_E,TOTAL_E
          EWRITE(3,*) 'The energry captured(%)',PERCENT_E
          !
          !           %-------------------------------%
          !           | Display computed residuals    |
          !           %-------------------------------%
#ifdef DOUBLEP
          call dmout(6, nconv, 2, s, maxncv, -6,           &
               'Singular values and direct residuals')
#else
          call smout(6, nconv, 2, s, maxncv, -6,           &
               'Singular values and direct residuals')
#endif
       end if
       !
       !        %------------------------------------------%
       !        | Print additional convergence information |
       !        %------------------------------------------%
       !
       if ( info .eq. 1) then
          ewrite(3,*)  ' '
          ewrite(3,*)  ' Maximum number of iterations reached.'
          ewrite(3,*)  ' '
       else if ( info .eq. 3) then
          ewrite(3,*)  ' ' 
          ewrite(3,*)  ' No shifts could be applied during implicit',   &
               ' Arnoldi update, try increasing NCV.'
          ewrite(3,*)  ' '
       end if
       !
       ewrite(3,*)  ' '
       ewrite(3,*)  ' _SVD '
       ewrite(3,*)  ' ==== '
       ewrite(3,*)  ' '
       ewrite(3,*)  ' Size of the matrix is ', n
       ewrite(3,*)  ' The number of Ritz values requested is ', nev
       ewrite(3,*)  ' The number of Arnoldi vectors generated',     &
            ' (NCV) is ', ncv
       ewrite(3,*)  ' What portion of the spectrum: ', which
       ewrite(3,*)  ' The number of converged Ritz values is ',     &
            nconv 
       ewrite(3,*)  ' The number of Implicit Arnoldi update',       &
            ' iterations taken is ', iparam(3)
       ewrite(3,*)  ' The number of OP*x is ', iparam(9)
       ewrite(3,*)  ' The convergence criterion is ', toll
       ewrite(3,*)  ' '
       !
    end if
    !
    !     %-------------------------%
    !     | Done with program dsvd. |
    !     %-------------------------%
    !
9000 continue

    if (nconv .ne. nsvd_total) then
       write(*,*) 'convergence error in svd computation'
       stop
    end if

    do j = nsvd,1,-1
       svdval(nconv+1-j)= s(j,1)  ! decreasing order
    enddo


    do j = nsvd,1,-1        
       do i=1,m       
          usv(i,nconv+1-j) = u(i,j)
       enddo
    enddo
    ewrite(3,*) 'the end of snapsvd'
1000 format (1X, E24.16, 1X, E24.16)
2000 format (100(1X,E24.16)) 

    DEALLOCATE(v)
    DEALLOCATE(u)
    DEALLOCATE(s)
    DEALLOCATE(workl)
    DEALLOCATE(workd)
    DEALLOCATE(resid)
    DEALLOCATE(ax)
    DEALLOCATE(select)

#else
    FLExit("No ARPACK support, recompile with ARPACK - cannot run reduced model.")
#endif
    return
  end subroutine snapsvd


  !*******************H1***********************
  subroutine snapsvd_H1(m,n,epsilon,snapmatrix,snapmatrix_dx,snapmatrix_dy,snapmatrix_dz,nsvd_total,nsvd,usv,svdval,D3)

    !   
    !     This code shows how to use ARPACK to find a few of the
    !     largest singular values(sigma) and corresponding right singular 
    !     vectors (v) for the the matrix A by solving the symmetric problem:
    !          
    !                        (A'*A)*v = sigma*v
    ! 
    !     where A is an m by n real matrix.
    !
    !     This code may be easily modified to estimate the 2-norm
    !     condition number  largest(sigma)/smallest(sigma) by setting
    !     which = 'BE' below.  This will ask for a few of the smallest
    !     and a few of the largest singular values simultaneously.
    !     The condition number could then be estimated by taking
    !     the ratio of the largest and smallest singular values.
    !
    !     This formulation is appropriate when  m  .ge.  n.
    !     Reverse the roles of A and A' in the case that  m .le. n.
    !
    !     The main points illustrated here are 
    !
    !        1) How to declare sufficient memory to find NEV 
    !           largest singular values of A .  
    !
    !        2) Illustration of the reverse communication interface 
    !           needed to utilize the top level ARPACK routine DSAUPD 
    !           that computes the quantities needed to construct
    !           the desired singular values and vectors(if requested).
    !
    !        3) How to extract the desired singular values and vectors
    !           using the ARPACK routine DSEUPD.
    !
    !        4) How to construct the left singular vectors U from the 
    !           right singular vectors V to obtain the decomposition
    !
    !                        A*V = U*S
    !
    !           where S = diag(sigma_1, sigma_2, ..., sigma_k).
    !
    !     The only thing that must be supplied in order to use this
    !     routine on your problem is to change the array dimensions 
    !     appropriately, to specify WHICH singular values you want to 
    !     compute and to supply a the matrix-vector products 
    !
    !                         w <-  Ax
    !                         y <-  A'w
    !
    !     in place of the calls  to AV( ) and ATV( ) respectively below.  
    !
    !     Further documentation is available in the header of DSAUPD
    !     which may be found in the SRC directory.
    !
    !     This codes implements
    !
    !\Example-1
    !     ... Suppose we want to solve A'A*v = sigma*v in regular mode,
    !         where A is derived from the simplest finite difference 
    !         discretization of the 2-dimensional kernel  K(s,t)dt  where
    !
    !                 K(s,t) =  s(t-1)   if 0 .le. s .le. t .le. 1,
    !                           t(s-1)   if 0 .le. t .lt. s .le. 1. 
    !
    !         See subroutines AV  and ATV for details.
    !     ... OP = A'*A  and  B = I.
    !     ... Assume "call av (n,x,y)" computes y = A*x
    !     ... Assume "call atv (n,y,w)" computes w = A'*y
    !     ... Assume exact shifts are used
    !     ...
    !
    !\BeginLib
    !
    !\Routines called:
    !     dsaupd  ARPACK reverse communication interface routine.
    !     dseupd  ARPACK routine that returns Ritz values and (optionally)
    !             Ritz vectors.
    !     dnrm2   Level 1 BLAS that computes the norm of a vector.
    !     daxpy   Level 1 BLAS that computes y <- alpha*x+y.
    !     dscal   Level 1 BLAS thst computes x <- x*alpha.
    !     dcopy   Level 1 BLAS thst computes y <- x.
    !
    !-----------------------------------------------------------------------
    !
    !     %------------------------------------------------------%
    !     | Storage Declarations:                                |
    !     |                                                      |
    !     | It is assumed that A is M by N with M .ge. N.        |
    !     |                                                      |
    !     | The maximum dimensions for all arrays are            |
    !     | set here to accommodate a problem size of            |
    !     | M .le. MAXM  and  N .le. MAXN                        |
    !     |                                                      |
    !     | The NEV right singular vectors will be computed in   |
    !     | the N by NCV array V.                                |
    !     |                                                      |
    !     | The NEV left singular vectors will be computed in    |
    !     | the M by NEV array U.                                |
    !     |                                                      |
    !     | NEV is the number of singular values requested.      |
    !     |     See specifications for ARPACK usage below.       |
    !     |                                                      |
    !     | NCV is the largest number of basis vectors that will |
    !     |     be used in the Implicitly Restarted Arnoldi      |
    !     |     Process.  Work per major iteration is            |
    !     |     proportional to N*NCV*NCV.                       |
    !     |                                                      |
    !     | You must set:                                        |
    !     |                                                      |
    !     | MAXM:   Maximum number of rows of the A allowed.     |
    !     | MAXN:   Maximum number of columns of the A allowed.  |
    !     | MAXNEV: Maximum NEV allowed                          |
    !     | MAXNCV: Maximum NCV allowed                          |
    !     %------------------------------------------------------%
    !

    use FLDebug

    !IMPLICIT NONE
    !include 'paramnew.h'

    integer   maxm, maxn, maxnev, maxncv, ldv, ldu,m,n,i
    integer   ndigit,logfil, msgets, msaitr, msapps, msaupd, msaup2, mseigt, mseupd

    !        parameter (maxm = istate,maxn=nrsnapshots,maxnev=nrsnapshots-1,    &
    !      maxncv=nrsnapshots, ldu = maxm, ldv=maxn )
    !     %--------------%
    !     | Local Arrays |
    !     %--------------%
    !
    !real v(ldv,maxncv), u(ldu, maxnev),            & 
    !                 workl(maxncv*(maxncv+8)), workd(3*maxn),  &
    !                 s(maxncv,2), resid(maxn), ax(maxm)
    !logical          select(maxncv)
    integer          iparam(11), ipntr(11)

    REAL, ALLOCATABLE, DIMENSION(:,:):: v,u,s
    REAL, ALLOCATABLE, DIMENSION(:)::   workl,workd,resid,ax
    REAL, ALLOCATABLE, DIMENSION(:)::   ay,ay_dx,ay_dy,ay_dz
    logical, ALLOCATABLE, DIMENSION(:)::   select
    logical ::D3

    !
    !     %---------------%
    !     | Local Scalars |
    !     %---------------%
    !
    character        bmat*1, which*2
    integer          ido, nev, ncv, lworkl, info, ierr,   &
         j, ishfts, maxitr, mode1, nconv
    logical          rvec
    real toll, sigma, temp
    !
    !     %------------%
    !     | Parameters |
    !     %------------%
    !
    real one, zero
    parameter        (one = 1.0D+0, zero = 0.0D+0)


    !ccccc declare input
    real snapmatrix(m,n)
    real snapmatrix_dx(m,n),snapmatrix_dy(m,n),snapmatrix_dz(m,n)
    real epsilon
    integer nsvd,nsvd_total
    !cccccccc  declare output
    real usv(m,nsvd),svdval(nsvd) ! left svd, svdval
    !
    !  
    !     %-----------------------------%
    !     | BLAS & LAPACK routines used |
    !     %-----------------------------%
    !
    real dnrm2

    !
    !  %-------------%
    !  | local       |
    !  %-------------%

    REAL TOTAL_E,PARTIAL_E,PERCENT_E

!    external         dnrm2, daxpy, dcopy, dscal
    !
    !     %-----------------------%
    !     | Executable Statements |
    !     %-----------------------%
    !
    !     %-------------------------------------------------%
    !     | The following include statement and assignments |
    !     | initiate trace output from the internal         |
    !     | actions of ARPACK.  See debug.doc in the        |
    !     | DOCUMENTS directory for usage.  Initially, the  |
    !     | most useful information will be a breakdown of  |
    !     | time spent in the various stages of computation |
    !     | given by setting msaupd = 1.                    |
    !     %-------------------------------------------------%
    !
    !      include 'debug1.h'
#ifdef HAVE_LIBARPACK

    maxm = m
    maxn=n
    maxnev=n
    maxncv=n
    ldu = maxm
    ldv=maxn

    ALLOCATE(v(ldv,maxncv))
    ALLOCATE(u(ldu, maxnev))
    ALLOCATE(s(maxncv,2))
    ALLOCATE(workl(maxncv*(maxncv+8)))
    ALLOCATE(workd(3*maxn))
    ALLOCATE(resid(maxn))
    ALLOCATE(ax(maxm))
    ALLOCATE(select(maxncv))
    ALLOCATE(ay(maxn))
    ALLOCATE(ay_dx(maxn))
    ALLOCATE(ay_dy(maxn))
    ALLOCATE(ay_dz(maxn))

    ewrite(3,*) 'in snapsvd.F'
    ewrite(3,*) M,maxm
    !      maxm=m
    !      ldu = maxm
    ndigit = -3
    logfil = 6
    msgets = 0
    msaitr = 0 
    msapps = 0
    msaupd = 1
    msaup2 = 0
    mseigt = 0
    mseupd = 0
    ay=0.0
    ay_dx=0.0
    ay_dy=0.0
    ay_dz=0.0
    !
    !     %-------------------------------------------------%
    !     | The following sets dimensions for this problem. |
    !     %-------------------------------------------------%
    !
    !      if (n .ne. nrsnapshots) then
    !       write(*,*) 'dimension error in svd routine',n,nrsnapshots
    !       stop
    !      end if
    !
    !     %------------------------------------------------%
    !     | Specifications for ARPACK usage are set        | 
    !     | below:                                         |
    !     |                                                |
    !     |    1) NEV = 4 asks for 4 singular values to be |  
    !     |       computed.                                | 
    !     |                                                |
    !     |    2) NCV = 20 sets the length of the Arnoldi  |
    !     |       factorization                            |
    !     |                                                |
    !     |    3) This is a standard problem               |
    !     |         (indicated by bmat  = 'I')             |
    !     |                                                |
    !     |    4) Ask for the NEV singular values of       |
    !     |       largest magnitude                        |
    !     |         (indicated by which = 'LM')            |
    !     |       See documentation in DSAUPD for the      |
    !     |       other options SM, BE.                    | 
    !     |                                                |
    !     | Note: NEV and NCV must satisfy the following   |
    !     |       conditions:                              |
    !     |                 NEV <= MAXNEV,                 |
    !     |             NEV + 1 <= NCV <= MAXNCV           |
    !     %------------------------------------------------%
    !
    nev   = nsvd_total
    ncv   = n
    bmat  = 'I'
    which = 'LM'
    !
    if ( n .gt. maxn ) then
       ewrite(3,*)  ' ERROR with _SVD: N is greater than MAXN '
       go to 9000
    else if ( m .gt. maxm ) then
       ewrite(3,*)  ' ERROR with _SVD: M is greater than MAXM ',M,MAXM
       go to 9000
    else if ( nev .gt. maxnev ) then
       ewrite(3,*)  ' ERROR with _SVD: NEV is greater than MAXNEV '
       go to 9000
    else if ( ncv .gt. maxncv ) then
       ewrite(3,*)  ' ERROR with _SVD: NCV is greater than MAXNCV '
       go to 9000
    end if
    !
    !     %-----------------------------------------------------%
    !     | Specification of stopping rules and initial         |
    !     | conditions before calling DSAUPD                    |
    !     |                                                     |
    !     |           abs(sigmaC - sigmaT) < TOL*abs(sigmaC)    |
    !     |               computed   true                       |
    !     |                                                     |
    !     |      If TOL .le. 0,  then TOL <- macheps            |
    !     |              (machine precision) is used.           |
    !     |                                                     |
    !     | IDO  is the REVERSE COMMUNICATION parameter         |
    !     |      used to specify actions to be taken on return  |
    !     |      from DSAUPD. (See usage below.)                |
    !     |                                                     |
    !     |      It MUST initially be set to 0 before the first |
    !     |      call to DSAUPD.                                | 
    !     |                                                     |
    !     | INFO on entry specifies starting vector information |
    !     |      and on return indicates error codes            |
    !     |                                                     |
    !     |      Initially, setting INFO=0 indicates that a     | 
    !     |      random starting vector is requested to         |
    !     |      start the ARNOLDI iteration.  Setting INFO to  |
    !     |      a nonzero value on the initial call is used    |
    !     |      if you want to specify your own starting       |
    !     |      vector (This vector must be placed in RESID.)  | 
    !     |                                                     |
    !     | The work array WORKL is used in DSAUPD as           | 
    !     | workspace.  Its dimension LWORKL is set as          |
    !     | illustrated below.                                  |
    !     %-----------------------------------------------------%
    !
    lworkl = ncv*(ncv+8)
    toll = 0.0d-3 ! 1.d-3 
    info = 0
    ido = 0
    !
    !     %---------------------------------------------------%
    !     | Specification of Algorithm Mode:                  |
    !     |                                                   |
    !     | This program uses the exact shift strategy        |
    !     | (indicated by setting IPARAM(1) = 1.)             |
    !     | IPARAM(3) specifies the maximum number of Arnoldi |
    !     | iterations allowed.  Mode 1 of DSAUPD is used     |
    !     | (IPARAM(7) = 1). All these options can be changed |
    !     | by the user. For details see the documentation in |
    !     | DSAUPD.                                           |
    !     %---------------------------------------------------%
    !
    ishfts = 1
    maxitr = n
    mode1 = 1
    !
    iparam(1) = ishfts
    !                
    iparam(3) = maxitr
    !                  
    iparam(7) = mode1
    !
    !     %------------------------------------------------%
    !     | M A I N   L O O P (Reverse communication loop) |
    !     %------------------------------------------------%
    !
10  continue
    !
    !        %---------------------------------------------%
    !        | Repeatedly call the routine DSAUPD and take | 
    !        | actions indicated by parameter IDO until    |
    !        | either convergence is indicated or maxitr   |
    !        | has been exceeded.                          |
    !        %---------------------------------------------%
    !
#ifdef DOUBLEP
    call dsaupd ( ido, bmat, n, which, nev, toll, resid,      &
         ncv, v, ldv, iparam, ipntr, workd, workl,   &
         lworkl, info )
#else
    call ssaupd ( ido, bmat, n, which, nev, toll, resid,      &
         ncv, v, ldv, iparam, ipntr, workd, workl,   &
         lworkl, info )
#endif
    !
    if (ido .eq. -1 .or. ido .eq. 1) then
       !
       !           %---------------------------------------%
       !           | Perform matrix vector multiplications |
       !           |              w <--- A*x       (av())  |
       !           |              y <--- A'*w      (atv()) |
       !           | The user should supply his/her own    |
       !           | matrix vector multiplication routines |
       !           | here that takes workd(ipntr(1)) as    |
       !           | the input, and returns the result in  |
       !           | workd(ipntr(2)).                      |
       !           %---------------------------------------%
       !
       call av (m,n,snapmatrix,workd(ipntr(1)),ax(1:m)) 
       call atv (m,n,snapmatrix,ax,ay)

       call av (m,n,snapmatrix_dx,workd(ipntr(1)),ax(1:m)) 
       call atv (m,n,snapmatrix_dx,ax,ay_dx)

       call av (m,n,snapmatrix_dy,workd(ipntr(1)),ax(1:m)) 
       call atv (m,n,snapmatrix_dy,ax,ay_dy)

       if(.false.) then
          call av (m,n,snapmatrix_dz,workd(ipntr(1)),ax(1:m)) 
          call atv (m,n,snapmatrix_dz,ax,ay_dz)
       endif
       open(1,file='left.dat')
       write(1,*) ay(1:n)
       close(1)
       open(1,file='left_dx.dat')
       write(1,*) ay_dx(1:n)
       close(1)
       open(1,file='left_dy.dat')
       write(1,*) ay_dy(1:n)
       close(1)
       open(1,file='left_dz.dat')
       write(1,*) ay_dz(1:n)
       close(1)
       !     stop 1
       workd(ipntr(2):ipntr(2)+n-1)=ay(1:n)+epsilon*(ay_dx(1:n)+ay_dy(1:n)+ay_dz(1:n))
       !
       !           %-----------------------------------------%
       !           | L O O P   B A C K to call DSAUPD again. |
       !           %-----------------------------------------%
       !
       go to 10
       !
    end if
    !
    !     %----------------------------------------%
    !     | Either we have convergence or there is |
    !     | an error.                              |
    !     %----------------------------------------%
    !
    if ( info .lt. 0 ) then
       !
       !        %--------------------------%
       !        | Error message. Check the |
       !        | documentation in DSAUPD. |
       !        %--------------------------%
       !
       ewrite(3,*)  ' '
       ewrite(3,*)  ' Error with _saupd, info = ', info
       ewrite(3,*)  ' Check documentation in _saupd '
       ewrite(3,*)  ' '
       !
    else 
       !
       !        %--------------------------------------------%
       !        | No fatal errors occurred.                  |
       !        | Post-Process using DSEUPD.                 |
       !        |                                            |
       !        | Computed singular values may be extracted. |  
       !        |                                            |
       !        | Singular vectors may also be computed now  |
       !        | if desired.  (indicated by rvec = .true.)  | 
       !        |                                            |
       !        | The routine DSEUPD now called to do this   |
       !        | post processing                            | 
       !        %--------------------------------------------%
       !           
       rvec = .true.
       !
#ifdef DOUBLEP
       call dseupd ( rvec, 'All', select, s, v, ldv, sigma,   &
            bmat, n, which, nev, toll, resid, ncv, v, ldv,    &
            iparam, ipntr, workd, workl, lworkl, ierr )
#else
       call sseupd ( rvec, 'All', select, s, v, ldv, sigma,   &
            bmat, n, which, nev, toll, resid, ncv, v, ldv,    &
            iparam, ipntr, workd, workl, lworkl, ierr )
#endif
       !
       !        %-----------------------------------------------%
       !        | Singular values are returned in the first     |
       !        | column of the two dimensional array S         |
       !        | and the corresponding right singular vectors  | 
       !        | are returned in the first NEV columns of the  |
       !        | two dimensional array V as requested here.    |
       !        %-----------------------------------------------%
       !

       ewrite(3,*) ipntr(7)
       ewrite(3,*) ipntr(7)+ncv
       ewrite(3,*) (workl(ipntr(7)+ncv+j-1),j=1,ncv)

       TOTAL_E=0.0
       DO J=1,NCV
          IF(workl(ipntr(7)+ncv+j-1).ge.0.0) THEN
             TOTAL_E=TOTAL_E+sqrt(workl(ipntr(7)+ncv+j-1))
          ENDIF
       ENDDO
       if ( ierr .ne. 0) then
          !
          !           %------------------------------------%
          !           | Error condition:                   |
          !           | Check the documentation of DSEUPD. |
          !           %------------------------------------%
          !
          ewrite(3,*)  ' '
          ewrite(3,*)  ' Error with _seupd, info = ', ierr
          ewrite(3,*)  ' Check the documentation of _seupd. '
          ewrite(3,*)  ' '
          !
       else
          !
          nconv =  iparam(5)
          ewrite(3,*) 'nconv',nconv
          PARTIAL_E=0.0
          do 20 j=1, nconv
             s(j,1) = sqrt(s(j,1))
             if(j.le.nsvd) then
                PARTIAL_E=PARTIAL_E+s(j,1)
             endif
             !
             !              %-----------------------------%
             !              | Compute the left singular   |
             !              | vectors from the formula    |
             !              |                             |
             !              |     u = Av/sigma            |
             !              |                             |
             !              | u should have norm 1 so     |
             !              | divide by norm(Av) instead. |
             !              %-----------------------------%
             !
             call av(m, n, snapmatrix, v(1,j), ax)
             u(:,j)=ax
             !call dcopy(m, ax, 1, u(1,j), 1)
             temp = one/norm2(u(:,j))
             !temp = one/dnrm2(m, u(1,j), 1)
             u(:,j)=temp*u(:,j)
             !call dscal(m, temp, u(1,j), 1)
             !
             !              %---------------------------%
             !              |                           |
             !              | Compute the residual norm |
             !              |                           |
             !              |   ||  A*v - sigma*u ||    |
             !              |                           |
             !              | for the NCONV accurately  |
             !              | computed singular values  |
             !              | and vectors.  (iparam(5)  |
             !              | indicates how many are    |
             !              | accurate to the requested |
             !              | tolerance).               |
             !              | Store the result in 2nd   |
             !              | column of array S.        |
             !              %---------------------------%
             !
             ax=ax-s(j,1)*u(:,j)
             !call daxpy(m, -s(j,1), u(1,j), 1, ax, 1)
             s(j,2) = norm2(ax)
             !s(j,2) = dnrm2(m, ax, 1)
             !
20        end do

          PERCENT_E = 100.0*PARTIAL_E/TOTAL_E
          ewrite(3,*)  PARTIAL_E,TOTAL_E
          EWRITE(3,*) 'The energry captured(%)',PERCENT_E
          !
          !           %-------------------------------%
          !           | Display computed residuals    |
          !           %-------------------------------%
#ifdef DOUBLEP
          call dmout(6, nconv, 2, s, maxncv, -6,           &
               'Singular values and direct residuals')
#else
          call smout(6, nconv, 2, s, maxncv, -6,           &
               'Singular values and direct residuals')
#endif
       end if
       !
       !        %------------------------------------------%
       !        | Print additional convergence information |
       !        %------------------------------------------%
       !
       if ( info .eq. 1) then
          ewrite(3,*)  ' '
          ewrite(3,*)  ' Maximum number of iterations reached.'
          ewrite(3,*)  ' '
       else if ( info .eq. 3) then
          ewrite(3,*)  ' ' 
          ewrite(3,*)  ' No shifts could be applied during implicit',   &
               ' Arnoldi update, try increasing NCV.'
          ewrite(3,*)  ' '
       end if
       !
       ewrite(3,*)  ' '
       ewrite(3,*)  ' _SVD '
       ewrite(3,*)  ' ==== '
       ewrite(3,*)  ' '
       ewrite(3,*)  ' Size of the matrix is ', n
       ewrite(3,*)  ' The number of Ritz values requested is ', nev
       ewrite(3,*)  ' The number of Arnoldi vectors generated',     &
            ' (NCV) is ', ncv
       ewrite(3,*)  ' What portion of the spectrum: ', which
       ewrite(3,*)  ' The number of converged Ritz values is ',     &
            nconv 
       ewrite(3,*)  ' The number of Implicit Arnoldi update',       &
            ' iterations taken is ', iparam(3)
       ewrite(3,*)  ' The number of OP*x is ', iparam(9)
       ewrite(3,*)  ' The convergence criterion is ', toll
       ewrite(3,*)  ' '
       !
    end if
    !
    !     %-------------------------%
    !     | Done with program dsvd. |
    !     %-------------------------%
    !
9000 continue

    if (nconv .ne. nsvd_total) then
       write(*,*) 'convergence error in svd computation'
       stop
    end if

    do j = nsvd,1,-1
       svdval(nconv+1-j)= s(j,1)  ! decreasing order
    enddo


    do j = nsvd,1,-1        
       do i=1,m       
          usv(i,nconv+1-j) = u(i,j)
       enddo
    enddo
    ewrite(3,*) 'the end of snapsvd'
1000 format (1X, E24.16, 1X, E24.16)
2000 format (100(1X,E24.16)) 

    DEALLOCATE(v)
    DEALLOCATE(u)
    DEALLOCATE(s)
    DEALLOCATE(workl)
    DEALLOCATE(workd)
    DEALLOCATE(resid)
    DEALLOCATE(ax)
    DEALLOCATE(select)
    DEALLOCATE(ay)
    DEALLOCATE(ay_dx)
    DEALLOCATE(ay_dy)
    DEALLOCATE(ay_dz)

#else
    FLAbort("No ARPACK support - cannot run reduced model.")
#endif
    return
  end subroutine snapsvd_H1
  !********************


  ! 
  ! ------------------------------------------------------------------
  !     matrix vector subroutines
  !
  !-------------------------------------------------------------------
  !
  subroutine av (m, n, snapmatrix, x, w)
    !
    !     computes  w <- A*x

    !      include 'paramnew.h'

    real x(n), w(m), psum
    integer i,j,m,n
    real snapmatrix(m,n) 


    !      if (n .ne. nrsnapshots) then
    !         write(*,*)'matrix dimension error: no of snapshots'
    !         stop                     
    !      end if

    !      if(m .ne. istate) then 
    !         write(*,*)'matrix dimension error: state dimension'
    !         stop                     
    !      end if

    do i=1,m
       psum = 0.0d0
       do j=1,n
          psum = psum + snapmatrix(i,j)*x(j)
       enddo
       w(i) = psum
    enddo

    return
  end subroutine av
  !
  !-------------------------------------------------------------------
  !
  subroutine atv (m, n, snapmatrix, w, y)
    !
    !     computes  y <- A'*w

    !      include 'paramnew.h'
    integer         m, n
    real w(m), y(n),psum

    integer i,j
    real snapmatrix(m,n) 

    !      if (n .ne. nrsnapshots .or. m .ne. istate) then
    !         write(*,*)'transpose matrix dimension error'
    !         stop                     
    !      end if

    do j=1,n
       psum = 0.0d0
       do i=1,m
          psum = psum + snapmatrix(i,j)*w(i)
       enddo
       y(j) = psum
    enddo

    return
  end subroutine atv

end module snapsvd_module