source: LMDZ5/trunk/libf/phylmd/cosp/array_lib.F90 @ 2850

Last change on this file since 2850 was 2428, checked in by idelkadi, 9 years ago

Mise a jour du simulateur COSP (passage de la version v3.2 a la version v1.4) :

  • mise a jour des sources pour ISCCP, CALIPSO et PARASOL
  • prise en compte des changements de phases pour les nuages (Calipso)
  • rajout de plusieurs diagnostiques (fraction nuageuse en fonction de la temperature, ...)

http://lmdz.lmd.jussieu.fr/Members/aidelkadi/cosp

  • 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: 4.1 KB
Line 
1! $Revision: 23 $, $Date: 2011-03-31 15:41:37 +0200 (jeu. 31 mars 2011) $
2! $URL: http://cfmip-obs-sim.googlecode.com/svn/stable/v1.4.0/quickbeam/array_lib.f90 $
3! ARRAY_LIB: Array procedures for F90
4! Compiled/Modified:
5!   07/01/06  John Haynes (haynes@atmos.colostate.edu)
6!
7! infind (function)
8! lin_interpolate (function)
9 
10  module array_lib
11  implicit none
12
13  contains
14
15! ----------------------------------------------------------------------------
16! function INFIND
17! ----------------------------------------------------------------------------
18  function infind(list,val,sort,dist)
19  use m_mrgrnk
20  implicit none
21!
22! Purpose:
23!   Finds the index of an array that is closest to a value, plus the
24!   difference between the value found and the value specified
25!
26! Inputs:
27!   [list]   an array of sequential values
28!   [val]    a value to locate
29! Optional input:
30!   [sort]   set to 1 if [list] is in unknown/non-sequential order
31!
32! Returns:
33!   index of [list] that is closest to [val]
34!
35! Optional output:
36!   [dist]   set to variable containing [list([result])] - [val]
37!
38! Requires:
39!   mrgrnk library
40!
41! Created:
42!   10/16/03  John Haynes (haynes@atmos.colostate.edu)
43! Modified:
44!   01/31/06  IDL to Fortran 90
45 
46! ----- INPUTS -----
47  real*8, dimension(:), intent(in) :: list
48  real*8, intent(in) :: val 
49  integer, intent(in), optional :: sort
50 
51! ----- OUTPUTS -----
52  integer*4 :: infind
53  real*8, intent(out), optional :: dist
54
55! ----- INTERNAL -----
56  real*8, dimension(size(list)) :: lists
57  integer*4 :: nlist, result, tmp(1), sort_list
58  integer*4, dimension(size(list)) :: mask, idx
59
60  if (present(sort)) then
61    sort_list = sort
62  else
63    sort_list = 0
64  endif 
65
66  nlist = size(list)
67  if (sort_list == 1) then
68    call mrgrnk(list,idx)
69    lists = list(idx)
70  else
71    lists = list
72  endif
73
74  if (val >= lists(nlist)) then
75    result = nlist
76  else if (val <= lists(1)) then
77    result = 1
78  else
79    mask(:) = 0
80    where (lists < val) mask = 1
81      tmp = minloc(mask,1)
82      if (abs(lists(tmp(1)-1)-val) < abs(lists(tmp(1))-val)) then
83        result = tmp(1) - 1
84      else
85        result = tmp(1)
86      endif
87  endif
88  if (present(dist)) dist = lists(result)-val
89  if (sort_list == 1) then
90    infind = idx(result)
91  else
92    infind = result
93  endif
94
95  end function infind
96
97! ----------------------------------------------------------------------------
98! function LIN_INTERPOLATE
99! ---------------------------------------------------------------------------- 
100  subroutine lin_interpolate(yarr,xarr,yyarr,xxarr,tol)
101  use m_mrgrnk
102  implicit none
103!
104! Purpose:
105!   linearly interpolate a set of y2 values given a set of y1,x1,x2
106!
107! Inputs:
108!   [yarr]    an array of y1 values
109!   [xarr]    an array of x1 values
110!   [xxarr]   an array of x2 values
111!   [tol]     maximum distance for a match
112!
113! Output:
114!   [yyarr]   interpolated array of y2 values
115!
116! Requires:
117!   mrgrnk library
118!
119! Created:
120!   06/07/06  John Haynes (haynes@atmos.colostate.edu)
121
122! ----- INPUTS -----
123  real*8, dimension(:), intent(in) :: yarr, xarr, xxarr
124  real*8, intent(in) :: tol
125
126! ----- OUTPUTS -----
127  real*8, dimension(size(xxarr)), intent(out) :: yyarr
128
129! ----- INTERNAL -----
130  real*8, dimension(size(xarr)) :: ysort, xsort
131  integer*4, dimension(size(xarr)) :: ist
132  integer*4 :: nx, nxx, i, iloc
133  real*8 :: d, m
134
135  nx = size(xarr)
136  nxx = size(xxarr)
137
138! // xsort, ysort are sorted versions of xarr, yarr 
139  call mrgrnk(xarr,ist)
140  ysort = yarr(ist)
141  xsort = xarr(ist)
142 
143  do i=1,nxx
144    iloc = infind(xsort,xxarr(i),dist=d)
145    if (d > tol) then
146      print *, 'interpolation error'
147      stop
148    endif
149    if (iloc == nx) then
150!     :: set to the last value
151      yyarr(i) = ysort(nx)
152    else
153!     :: is there another closeby value?
154      if (abs(xxarr(i)-xsort(iloc+1)) < 2*tol) then
155!       :: yes, do a linear interpolation     
156        m = (ysort(iloc+1)-ysort(iloc))/(xsort(iloc+1)-xsort(iloc))
157        yyarr(i) = ysort(iloc) + m*(xxarr(i)-xsort(iloc))
158      else
159!       :: no, set to the only nearby value
160        yyarr(i) = ysort(iloc)
161      endif
162    endif
163  enddo
164 
165  end subroutine lin_interpolate
166
167  end module array_lib
Note: See TracBrowser for help on using the repository browser.