source: LMDZ5/trunk/tools/Max_diff_nc_with_lib/Jumble/compare.f90 @ 5223

Last change on this file since 5223 was 1907, checked in by lguez, 11 years ago

Added a copyright property to every file of the distribution, except
for the fcm files (which have their own copyright). Use svn propget on
a file to see the copyright. For instance:

$ svn propget copyright libf/phylmd/physiq.F90
Name of program: LMDZ
Creation date: 1984
Version: LMDZ5
License: CeCILL version 2
Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
See the license file in the root directory

Also added the files defining the CeCILL version 2 license, in French
and English, at the top of the LMDZ tree.

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 24.6 KB
Line 
1module compare_m
2
3  use prt_cmp_m, only: prt_cmp
4  use avg_mag_m, only: avg_mag
5
6  implicit none
7
8  interface compare
9     ! Makes a numerical comparison between two arrays of rank 1 to 4,
10     ! of types real or double precision.
11     module procedure compare1, compare1_dble, compare2, compare2_dble, &
12          compare3, compare3_dble, compare4, compare4_dble
13  end interface
14
15  character(len=*), parameter:: dashes &
16       = "---------------------------------------------"
17
18  private
19  public compare
20
21contains
22
23  subroutine compare1(data_old, data_new, tag, comp_mag, report_id, quiet)
24
25    ! Rank 1, real
26
27    use nr_util, only: assert_eq
28
29    real, intent(in):: data_old(:), data_new(:)
30    character(len=*), intent(in):: tag
31    logical, intent(in):: comp_mag
32    logical, intent(in):: report_id ! report identical variables
33    logical, intent(in):: quiet
34
35    ! Variables local to the subprogram:
36
37    logical zero(size(data_old))
38
39    real rel_diff(size(data_old))
40    ! (absolute value of relative difference)
41
42    real abs_diff(size(data_old))
43    ! (absolute value of absolute difference)
44
45    integer data_size
46    integer location(1)
47    character(len=len(dashes)+len(tag)+20) tag_fmt
48
49    !------------------------------------------------------
50
51    data_size = assert_eq(size(data_old), size(data_new), &
52         "compare1: sizes differ -- " // tag)
53    if (quiet) then
54       tag_fmt = '(a, ":")'
55    else
56       tag_fmt = '(/, "' // dashes // '", /, a, ":")'
57    end if
58
59    if (data_size == 0) then
60       print tag_fmt, tag
61       print *, "Zero-sized array"
62    else
63       if (all(data_old == data_new)) then
64          if (report_id) then
65             write(unit=*, fmt=tag_fmt, advance="no") tag
66             print *, "arrays are identical."
67          end if
68       else
69          if (quiet) then
70             write(unit=*, fmt=tag_fmt, advance="no") tag
71             print *, "arrays are different."
72          else
73             zero = data_old == 0. .or. data_new == 0.
74             where (.not. zero) rel_diff = abs(data_new / data_old - 1.)
75             abs_diff = abs(data_new - data_old)
76
77             print tag_fmt, tag
78
79             if (comp_mag) then
80                print *
81                print '("Average value of log10(abs(data_old)): ", 1pg8.1)', &
82                     avg_mag(data_old)
83             end if
84
85             if (any(.not. zero)) then
86                print *
87                print *, 'Relative difference for non-zero values:'
88                location = maxloc(rel_diff, mask=.not. zero)
89                call prt_cmp(data_old(location(1)), data_new(location(1)), &
90                     rel_diff(location(1)), location)
91             end if
92
93             if (any(zero)) then
94                print *
95                print *, 'Absolute difference when there is a zero:'
96                location = maxloc(abs_diff, mask=zero)
97                call prt_cmp(data_old(location(1)), data_new(location(1)), &
98                     abs_diff(location(1)), location)
99             end if
100
101             print *
102             print *, 'Absolute difference:'
103             location = maxloc(abs_diff)
104             call prt_cmp(data_old(location(1)), data_new(location(1)), &
105                  abs_diff(location(1)), location)
106          end if
107       end if
108    end if
109
110  end subroutine compare1
111
112  !***********************************************************
113
114  subroutine compare1_dble(data_old, data_new, tag, comp_mag, report_id, quiet)
115
116    ! Rank 1, double precision
117
118    use nr_util, only: assert_eq
119
120    double precision, intent(in):: data_old(:), data_new(:)
121    character(len=*), intent(in):: tag
122    logical, intent(in):: comp_mag
123    logical, intent(in):: report_id ! report identical variables
124    logical, intent(in):: quiet
125
126    ! Variables local to the subprogram:
127
128    logical zero(size(data_old))
129
130    double precision rel_diff(size(data_old))
131    ! (absolute value of relative difference)
132
133    double precision abs_diff(size(data_old))
134    ! (absolute value of absolute difference)
135
136    integer data_size
137    integer location(1)
138    character(len=len(dashes)+len(tag)+20) tag_fmt
139
140    !------------------------------------------------------
141
142    data_size = assert_eq(size(data_old), size(data_new), &
143         "compare1_dble: sizes differ -- " // tag)
144    if (quiet) then
145       tag_fmt = '(a, ":")'
146    else
147       tag_fmt = '(/, "' // dashes // '", /, a, ":")'
148    end if
149
150    if (data_size == 0) then
151       print tag_fmt, tag
152       print *, "Zero-sized array"
153    else
154       if (all(data_old == data_new)) then
155          if (report_id) then
156             write(unit=*, fmt=tag_fmt, advance="no") tag
157             print *, "arrays are identical."
158          end if
159       else
160          if (quiet) then
161             write(unit=*, fmt=tag_fmt, advance="no") tag
162             print *, "arrays are different."
163          else
164             zero = data_old == 0d0 .or. data_new == 0d0
165             where (.not. zero) rel_diff = abs(data_new / data_old - 1d0)
166             abs_diff = abs(data_new - data_old)
167
168             print tag_fmt, tag
169
170             if (comp_mag) then
171                print *
172                print '("Average value of log10(abs(data_old)): ", 1pg8.1)', &
173                     avg_mag(data_old)
174             end if
175
176             if (any(.not. zero)) then
177                print *
178                print *, 'Relative difference for non-zero values:'
179                location = maxloc(rel_diff, mask=.not. zero)
180                call prt_cmp(data_old(location(1)), data_new(location(1)), &
181                     rel_diff(location(1)), location)
182             end if
183
184             if (any(zero)) then
185                print *
186                print *, 'Absolute difference when there is a zero:'
187                location = maxloc(abs_diff, mask=zero)
188                call prt_cmp(data_old(location(1)), data_new(location(1)), &
189                     abs_diff(location(1)), location)
190             end if
191
192             print *
193             print *, 'Absolute difference:'
194             location = maxloc(abs_diff)
195             call prt_cmp(data_old(location(1)), data_new(location(1)), &
196                  abs_diff(location(1)), location)
197          end if
198       end if
199    end if
200
201  end subroutine compare1_dble
202
203  !***********************************************************
204
205  subroutine compare2(data_old, data_new, tag, comp_mag, report_id, quiet)
206
207    ! Rank 2, real
208
209    use nr_util, only: assert
210    use point_m, only: point
211
212    real, intent(in):: data_old(:,:), data_new(:,:)
213    character(len=*), intent(in):: tag
214    logical, intent(in):: comp_mag
215    logical, intent(in):: report_id ! report identical variables
216    logical, intent(in):: quiet
217
218    ! Variables local to the subprogram:
219
220    logical zero(size(data_old,1), size(data_old,2))
221
222    real rel_diff(size(data_old,1), size(data_old,2))
223    ! (absolute value of relative difference)
224
225    real abs_diff(size(data_old,1), size(data_old,2))
226    ! (absolute value of absolute difference)
227
228    integer location(2)
229    character(len=len(dashes)+len(tag)+20) tag_fmt
230
231    !------------------------------------------------------
232
233    call assert(shape(data_old) == shape(data_new), &
234         "compare2: shapes differ -- " // tag)
235    if (quiet) then
236       tag_fmt = '(a, ":")'
237    else
238       tag_fmt = '(/, "' // dashes // '", /, a, ":")'
239    end if
240
241    if (size(data_old) == 0) then
242       print tag_fmt, tag
243       print *, "Zero-sized array"
244    else
245       if (all(data_old == data_new)) then
246          if (report_id) then
247             write(unit=*, fmt=tag_fmt, advance="no") tag
248             print *, "arrays are identical."
249          end if
250       else
251          if (quiet) then
252             write(unit=*, fmt=tag_fmt, advance="no") tag
253             print *, "arrays are different."
254          else
255             zero = data_old == 0. .or. data_new == 0.
256             where (.not. zero) rel_diff = abs(data_new / data_old - 1.)
257             abs_diff = abs(data_new - data_old)
258
259             print tag_fmt, tag
260
261             if (comp_mag) then
262                print *
263                print '("Average value of log10(abs(data_old)): ", 1pg8.1)', &
264                     avg_mag(data_old)
265             end if
266
267             if (any(.not. zero)) then
268                print *
269                print *, 'Relative difference for non-zero values:'
270                location = maxloc(rel_diff, mask=.not. zero)
271                call prt_cmp(point(data_old, location), &
272                     point(data_new, location), point(rel_diff, location), &
273                     location)
274             end if
275
276             if (any(zero)) then
277                print *
278                print *, 'Absolute difference when there is a zero:'
279                location = maxloc(abs_diff, mask=zero)
280                call prt_cmp(point(data_old, location), &
281                     point(data_new, location), point(abs_diff, location), &
282                     location)
283             end if
284
285             print *
286             print *, 'Absolute difference:'
287             location = maxloc(abs_diff)
288             call prt_cmp(point(data_old, location), &
289                  point(data_new, location), point(abs_diff, location), &
290                  location)
291          end if
292       end if
293    end if
294
295  end subroutine compare2
296
297  !***********************************************************
298
299  subroutine compare2_dble(data_old, data_new, tag, comp_mag, report_id, quiet)
300
301    ! Rank 2, double precision
302
303    use nr_util, only: assert
304    use point_m, only: point
305
306    double precision, intent(in):: data_old(:,:), data_new(:,:)
307    character(len=*), intent(in):: tag
308    logical, intent(in):: comp_mag
309    logical, intent(in):: report_id ! report identical variables
310    logical, intent(in):: quiet
311
312    ! Variables local to the subprogram:
313
314    logical zero(size(data_old,1), size(data_old,2))
315
316    double precision rel_diff(size(data_old,1), size(data_old,2))
317    ! (absolute value of relative difference)
318
319    double precision abs_diff(size(data_old,1), size(data_old,2))
320    ! (absolute value of absolute difference)
321
322    integer location(2)
323    character(len=len(dashes)+len(tag)+20) tag_fmt
324
325    !------------------------------------------------------
326
327    call assert(shape(data_old) == shape(data_new), &
328         "compare2_dble: shapes differ -- "  // tag)
329    if (quiet) then
330       tag_fmt = '(a, ":")'
331    else
332       tag_fmt = '(/, "' // dashes // '", /, a, ":")'
333    end if
334
335    if (size(data_old) == 0) then
336       print tag_fmt, tag
337       print *, "Zero-sized array"
338    else
339       if (all(data_old == data_new)) then
340          if (report_id) then
341             write(unit=*, fmt=tag_fmt, advance="no") tag
342             print *, "arrays are identical."
343          end if
344       else
345          if (quiet) then
346             write(unit=*, fmt=tag_fmt, advance="no") tag
347             print *, "arrays are different."
348          else
349             zero = data_old == 0d0 .or. data_new == 0d0
350             where (.not. zero) rel_diff = abs(data_new / data_old - 1d0)
351             abs_diff = abs(data_new - data_old)
352
353             print tag_fmt, tag
354
355             if (comp_mag) then
356                print *
357                print '("Average value of log10(abs(data_old)): ", 1pg8.1)', &
358                     avg_mag(data_old)
359             end if
360
361             if (any(.not. zero)) then
362                print *
363                print *, 'Relative difference for non-zero values:'
364                location = maxloc(rel_diff, mask=.not. zero)
365                call prt_cmp(point(data_old, location), &
366                     point(data_new, location), point(rel_diff, location), &
367                     location)
368             end if
369
370             if (any(zero)) then
371                print *
372                print *, 'Absolute difference when there is a zero:'
373                location = maxloc(abs_diff, mask=zero)
374                call prt_cmp(point(data_old, location), &
375                     point(data_new, location), point(abs_diff, location), &
376                     location)
377             end if
378
379             print *
380             print *, 'Absolute difference:'
381             location = maxloc(abs_diff)
382             call prt_cmp(point(data_old, location), &
383                  point(data_new, location), point(abs_diff, location), &
384                  location)
385          end if
386       end if
387    end if
388
389  end subroutine compare2_dble
390
391  !***********************************************************
392
393  subroutine compare3(data_old, data_new, tag, comp_mag, report_id, quiet)
394
395    ! Rank 3, real
396
397    use nr_util, only: assert
398    use point_m, only: point
399
400    real, intent(in):: data_old(:, :, :), data_new(:, : ,:)
401    character(len=*), intent(in):: tag
402    logical, intent(in):: comp_mag
403    logical, intent(in):: report_id ! report identical variables
404    logical, intent(in):: quiet
405
406    ! Variables local to the subprogram:
407
408    logical zero(size(data_old,1), size(data_old,2), size(data_old,3))
409
410    real rel_diff(size(data_old,1), size(data_old,2), size(data_old,3))
411    ! (absolute value of relative difference)
412
413    real abs_diff(size(data_old,1), size(data_old,2), size(data_old,3))
414    ! (absolute value of absolute difference)
415
416    integer location(3)
417    character(len=len(dashes)+len(tag)+20) tag_fmt
418
419    !------------------------------------------------------
420
421    call assert(shape(data_old) == shape(data_new), &
422         "compare3: shapes differ -- " // tag)
423    if (quiet) then
424       tag_fmt = '(a, ":")'
425    else
426       tag_fmt = '(/, "' // dashes // '", /, a, ":")'
427    end if
428
429    if (size(data_old) == 0) then
430       print tag_fmt, tag
431       print *, "Zero-sized array"
432    else
433       if (all(data_old == data_new)) then
434          if (report_id) then
435             write(unit=*, fmt=tag_fmt, advance="no") tag
436             print *, "arrays are identical."
437          end if
438       else
439          if (quiet) then
440             write(unit=*, fmt=tag_fmt, advance="no") tag
441             print *, "arrays are different."
442          else
443             zero = data_old == 0. .or. data_new == 0.
444             where (.not. zero) rel_diff = abs(data_new / data_old - 1.)
445             abs_diff = abs(data_new - data_old)
446
447             print tag_fmt, tag
448
449             if (comp_mag) then
450                print *
451                print '("Average value of log10(abs(data_old)): ", 1pg8.1)', &
452                     avg_mag(data_old)
453             end if
454
455             if (any(.not. zero)) then
456                print *
457                print *, 'Relative difference for non-zero values:'
458                location = maxloc(rel_diff, mask=.not. zero)
459                call prt_cmp(point(data_old, location), &
460                     point(data_new, location), point(rel_diff, location), &
461                     location)
462             end if
463
464             if (any(zero)) then
465                print *
466                print *, 'Absolute difference when there is a zero:'
467                location = maxloc(abs_diff, mask=zero)
468                call prt_cmp(point(data_old, location), &
469                     point(data_new, location), point(abs_diff, location), &
470                     location)
471             end if
472
473             print *
474             print *, 'Absolute difference:'
475             location = maxloc(abs_diff)
476             call prt_cmp(point(data_old, location), &
477                  point(data_new, location), point(abs_diff, location), &
478                  location)
479          end if
480       end if
481    end if
482
483  end subroutine compare3
484
485  !***********************************************************
486
487  subroutine compare3_dble(data_old, data_new, tag, comp_mag, report_id, quiet)
488
489    ! Rank 3, double precision
490
491    use nr_util, only: assert
492    use point_m, only: point
493
494    double precision, intent(in):: data_old(:, :, :), data_new(:, : ,:)
495    character(len=*), intent(in):: tag
496    logical, intent(in):: comp_mag
497    logical, intent(in):: report_id ! report identical variables
498    logical, intent(in):: quiet
499
500    ! Variables local to the subprogram:
501
502    logical zero(size(data_old,1), size(data_old,2), size(data_old,3))
503
504    double precision rel_diff(size(data_old,1), size(data_old,2), &
505         size(data_old,3))
506    ! (absolute value of relative difference)
507
508    double precision abs_diff(size(data_old,1), size(data_old,2), &
509         size(data_old,3))
510    ! (absolute value of absolute difference)
511
512    integer location(3)
513    character(len=len(dashes)+len(tag)+20) tag_fmt
514
515    !------------------------------------------------------
516
517    call assert(shape(data_old) == shape(data_new), &
518         "compare3_dble: shapes differ -- " // tag)
519    if (quiet) then
520       tag_fmt = '(a, ":")'
521    else
522       tag_fmt = '(/, "' // dashes // '", /, a, ":")'
523    end if
524
525    if (size(data_old) == 0) then
526       print tag_fmt, tag
527       print *, "Zero-sized array"
528    else
529       if (all(data_old == data_new)) then
530          if (report_id) then
531             write(unit=*, fmt=tag_fmt, advance="no") tag
532             print *, "arrays are identical."
533          end if
534       else
535          if (quiet) then
536             write(unit=*, fmt=tag_fmt, advance="no") tag
537             print *, "arrays are different."
538          else
539             zero = data_old == 0. .or. data_new == 0.
540             where (.not. zero) rel_diff = abs(data_new / data_old - 1.)
541             abs_diff = abs(data_new - data_old)
542
543             print tag_fmt, tag
544
545             if (comp_mag) then
546                print *
547                print '("Average value of log10(abs(data_old)): ", 1pg8.1)', &
548                     avg_mag(data_old)
549             end if
550
551             if (any(.not. zero)) then
552                print *
553                print *, 'Relative difference for non-zero values:'
554                location = maxloc(rel_diff, mask=.not. zero)
555                call prt_cmp(point(data_old, location), &
556                     point(data_new, location), point(rel_diff, location), &
557                     location)
558             end if
559
560             if (any(zero)) then
561                print *
562                print *, 'Absolute difference when there is a zero:'
563                location = maxloc(abs_diff, mask=zero)
564                call prt_cmp(point(data_old, location), &
565                     point(data_new, location), point(abs_diff, location), &
566                     location)
567             end if
568
569             print *
570             print *, 'Absolute difference:'
571             location = maxloc(abs_diff)
572             call prt_cmp(point(data_old, location), &
573                  point(data_new, location), point(abs_diff, location), &
574                  location)
575          end if
576       end if
577    end if
578
579  end subroutine compare3_dble
580
581  !***********************************************************
582
583  subroutine compare4(data_old, data_new, tag, comp_mag, report_id, quiet)
584
585    ! Rank 4, real
586
587    use nr_util, only: assert
588    use point_m, only: point
589
590    real, intent(in):: data_old(:, :, :, :), data_new(:, : ,:, :)
591    character(len=*), intent(in):: tag
592    logical, intent(in):: comp_mag
593    logical, intent(in):: report_id ! report identical variables
594    logical, intent(in):: quiet
595
596    ! Variables local to the subprogram:
597
598    logical zero(size(data_old,1), size(data_old,2), size(data_old,3), &
599         size(data_old,4))
600
601    real rel_diff(size(data_old,1), size(data_old,2), size(data_old,3), &
602         size(data_old,4))
603    ! (absolute value of relative difference)
604
605    real abs_diff(size(data_old,1), size(data_old,2), size(data_old,3), &
606         size(data_old,4))
607    ! (absolute value of absolute difference)
608
609    integer location(4)
610    character(len=len(dashes)+len(tag)+20) tag_fmt
611
612    !------------------------------------------------------
613
614    call assert(shape(data_old) == shape(data_new), &
615         "compare4: shapes differ -- " // tag)
616    if (quiet) then
617       tag_fmt = '(a, ":")'
618    else
619       tag_fmt = '(/, "' // dashes // '", /, a, ":")'
620    end if
621
622    if (size(data_old) == 0) then
623       print tag_fmt, tag
624       print *, "Zero-sized array"
625    else
626       if (all(data_old == data_new)) then
627          if (report_id) then
628             write(unit=*, fmt=tag_fmt, advance="no") tag
629             print *, "arrays are identical."
630          end if
631       else
632          if (quiet) then
633             write(unit=*, fmt=tag_fmt, advance="no") tag
634             print *, "arrays are different."
635          else
636             zero = data_old == 0. .or. data_new == 0.
637             where (.not. zero) rel_diff = abs(data_new / data_old - 1.)
638             abs_diff = abs(data_new - data_old)
639
640             print tag_fmt, tag
641
642             if (comp_mag) then
643                print *
644                print '("Average value of log10(abs(data_old)): ", 1pg8.1)', &
645                     avg_mag(data_old)
646             end if
647
648             if (any(.not. zero)) then
649                print *
650                print *, 'Relative difference for non-zero values:'
651                location = maxloc(rel_diff, mask=.not. zero)
652                call prt_cmp(point(data_old, location), &
653                     point(data_new, location), point(rel_diff, location), &
654                     location)
655             end if
656
657             if (any(zero)) then
658                print *
659                print *, 'Absolute difference when there is a zero:'
660                location = maxloc(abs_diff, mask=zero)
661                call prt_cmp(point(data_old, location), &
662                     point(data_new, location), point(abs_diff, location), &
663                     location)
664             end if
665
666             print *
667             print *, 'Absolute difference:'
668             location = maxloc(abs_diff)
669             call prt_cmp(point(data_old, location), &
670                  point(data_new, location), point(abs_diff, location), &
671                  location)
672          end if
673       end if
674    end if
675
676  end subroutine compare4
677
678  !***********************************************************
679
680  subroutine compare4_dble(data_old, data_new, tag, comp_mag, report_id, quiet)
681
682    ! Rank 4, double precision
683
684    use nr_util, only: assert
685    use point_m, only: point
686
687    double precision, intent(in):: data_old(:, :, :, :), data_new(:, : ,:, :)
688    character(len=*), intent(in):: tag
689    logical, intent(in):: comp_mag
690    logical, intent(in):: report_id ! report identical variables
691    logical, intent(in):: quiet
692
693    ! Variables local to the subprogram:
694
695    logical zero(size(data_old,1), size(data_old,2), size(data_old,3), &
696         size(data_old,4))
697
698    double precision rel_diff(size(data_old,1), size(data_old,2), &
699         size(data_old,3), size(data_old,4))
700    ! (absolute value of relative difference)
701
702    double precision abs_diff(size(data_old,1), size(data_old,2), &
703         size(data_old,3), size(data_old,4))
704    ! (absolute value of absolute difference)
705
706    integer location(4)
707    character(len=len(dashes)+len(tag)+20) tag_fmt
708
709    !------------------------------------------------------
710
711    call assert(shape(data_old) == shape(data_new), &
712         "compare4_dble: shapes differ -- " // tag)
713    if (quiet) then
714       tag_fmt = '(a, ":")'
715    else
716       tag_fmt = '(/, "' // dashes // '", /, a, ":")'
717    end if
718
719    if (size(data_old) == 0) then
720       print tag_fmt, tag
721       print *, "Zero-sized array"
722    else
723       if (all(data_old == data_new)) then
724          if (report_id) then
725             write(unit=*, fmt=tag_fmt, advance="no") tag
726             print *, "arrays are identical."
727          end if
728       else
729          if (quiet) then
730             write(unit=*, fmt=tag_fmt, advance="no") tag
731             print *, "arrays are different."
732          else
733             zero = data_old == 0. .or. data_new == 0.
734             where (.not. zero) rel_diff = abs(data_new / data_old - 1.)
735             abs_diff = abs(data_new - data_old)
736
737             print tag_fmt, tag
738
739             if (comp_mag) then
740                print *
741                print '("Average value of log10(abs(data_old)): ", 1pg8.1)', &
742                     avg_mag(data_old)
743             end if
744
745             if (any(.not. zero)) then
746                print *
747                print *, 'Relative difference for non-zero values:'
748                location = maxloc(rel_diff, mask=.not. zero)
749                call prt_cmp(point(data_old, location), &
750                     point(data_new, location), point(rel_diff, location), &
751                     location)
752             end if
753
754             if (any(zero)) then
755                print *
756                print *, 'Absolute difference when there is a zero:'
757                location = maxloc(abs_diff, mask=zero)
758                call prt_cmp(point(data_old, location), &
759                     point(data_new, location), point(abs_diff, location), &
760                     location)
761             end if
762
763             print *
764             print *, 'Absolute difference:'
765             location = maxloc(abs_diff)
766             call prt_cmp(point(data_old, location), &
767                  point(data_new, location), point(abs_diff, location), &
768                  location)
769          end if
770       end if
771    end if
772
773  end subroutine compare4_dble
774
775end module compare_m
Note: See TracBrowser for help on using the repository browser.