Commit 698bfed7978cca17de80435bb75d458f9b2c5638

Authored by daniau
1 parent 36bf795453

git-svn-id: https://lxsd.femto-st.fr/svn/fvn@34 b657c933-2333-4658-acf2-d3c7c2708721

Showing 1 changed file with 26 additions and 12 deletions Side-by-side Diff

... ... @@ -2510,21 +2510,25 @@
2510 2510 double precision, external :: f
2511 2511 double precision, intent(in) :: a,b,epsabs,epsrel
2512 2512 integer, intent(in) :: key
2513   -integer, intent(in) :: limit
  2513 +integer, intent(in),optional :: limit
2514 2514 double precision, intent(out) :: res,abserr
2515 2515 integer, intent(out) :: ier
2516 2516  
2517 2517 double precision, allocatable :: work(:)
2518 2518 integer, allocatable :: iwork(:)
2519 2519 integer :: lenw,neval,last
  2520 +integer :: limitw
2520 2521  
2521 2522 ! imsl value for limit is 500
2522   -lenw=limit*4
  2523 +limitw=500
  2524 +if (present(limit)) limitw=limit
2523 2525  
2524   -allocate(iwork(limit))
  2526 +lenw=limitw*4
  2527 +
  2528 +allocate(iwork(limitw))
2525 2529 allocate(work(lenw))
2526 2530  
2527   -call dqag(f,a,b,epsabs,epsrel,key,res,abserr,neval,ier,limit,lenw,last,iwork,work)
  2531 +call dqag(f,a,b,epsabs,epsrel,key,res,abserr,neval,ier,limitw,lenw,last,iwork,work)
2528 2532  
2529 2533 deallocate(work)
2530 2534 deallocate(iwork)
2531 2535  
2532 2536  
2533 2537  
2534 2538  
... ... @@ -2597,21 +2601,26 @@
2597 2601 implicit none
2598 2602 double precision, external:: f,g,h
2599 2603 double precision, intent(in) :: a,b,epsabs,epsrel
2600   -integer, intent(in) :: key,limit
  2604 +integer, intent(in) :: key
  2605 +integer, intent(in), optional :: limit
2601 2606 integer, intent(out) :: ier
2602 2607 double precision, intent(out) :: res,abserr
2603 2608  
2604 2609  
2605 2610 double precision, allocatable :: work(:)
  2611 +integer :: limitw
2606 2612 integer, allocatable :: iwork(:)
2607 2613 integer :: lenw,neval,last
2608 2614  
2609 2615 ! imsl value for limit is 500
2610   -lenw=limit*4
  2616 +limitw=500
  2617 +if (present(limit)) limitw=limit
  2618 +
  2619 +lenw=limitw*4
2611 2620 allocate(work(lenw))
2612   -allocate(iwork(limit))
  2621 +allocate(iwork(limitw))
2613 2622  
2614   -call dqag_2d_outer(f,a,b,g,h,epsabs,epsrel,key,res,abserr,neval,ier,limit,lenw,last,iwork,work)
  2623 +call dqag_2d_outer(f,a,b,g,h,epsabs,epsrel,key,res,abserr,neval,ier,limitw,lenw,last,iwork,work)
2615 2624  
2616 2625 deallocate(iwork)
2617 2626 deallocate(work)
2618 2627  
2619 2628  
2620 2629  
2621 2630  
... ... @@ -2685,21 +2694,26 @@
2685 2694 implicit none
2686 2695 double precision, external:: f
2687 2696 double precision, intent(in) :: x,a,b,epsabs,epsrel
2688   -integer, intent(in) :: key,limit
  2697 +integer, intent(in) :: key
  2698 +integer, intent(in),optional :: limit
2689 2699 integer, intent(out) :: ier
2690 2700 double precision, intent(out) :: res,abserr
2691 2701  
2692 2702  
2693 2703 double precision, allocatable :: work(:)
  2704 +integer :: limitw
2694 2705 integer, allocatable :: iwork(:)
2695 2706 integer :: lenw,neval,last
2696 2707  
2697 2708 ! imsl value for limit is 500
2698   -lenw=limit*4
  2709 +limitw=500
  2710 +if (present(limit)) limitw=limit
  2711 +
  2712 +lenw=limitw*4
2699 2713 allocate(work(lenw))
2700   -allocate(iwork(limit))
  2714 +allocate(iwork(limitw))
2701 2715  
2702   -call dqag_2d_inner(f,x,a,b,epsabs,epsrel,key,res,abserr,neval,ier,limit,lenw,last,iwork,work)
  2716 +call dqag_2d_inner(f,x,a,b,epsabs,epsrel,key,res,abserr,neval,ier,limitw,lenw,last,iwork,work)
2703 2717  
2704 2718 deallocate(iwork)
2705 2719 deallocate(work)