source: CPL/oasis3-mct_5.0/lib/mct/testsystem/testall/convertgauss.F90 @ 6328

Last change on this file since 6328 was 6328, checked in by aclsce, 17 months ago

First import of oasis3-mct_5.0 (from oasis git server, branch OASIS3-MCT_5.0)

File size: 15.4 KB
Line 
1!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2!
3!     This program creates a remapping grid file for Gaussian lat/lon
4!     grids (for spectral transform codes).
5!
6!-----------------------------------------------------------------------
7!
8!     CVS:$Id: convertgauss.F90,v 1.3 2002-11-14 17:11:07 eong Exp $
9!     CVS $Name:  $
10!
11!     Copyright (c) 1997, 1998 the Regents of the University of
12!       California.
13!
14!     Unless otherwise indicated, this software has been authored
15!     by an employee or employees of the University of California,
16!     operator of the Los Alamos National Laboratory under Contract
17!     No. W-7405-ENG-36 with the U.S. Department of Energy.  The U.S.
18!     Government has rights to use, reproduce, and distribute this
19!     software.  The public may copy and use this software without
20!     charge, provided that this Notice and any statement of authorship
21!     are reproduced on all copies.  Neither the Government nor the
22!     University makes any warranty, express or implied, or assumes
23!     any liability or responsibility for the use of this software.
24!
25!***********************************************************************
26
27      subroutine convertgauss(GGrid, nx, ny)
28
29!-----------------------------------------------------------------------
30!
31!     This file creates a remapping grid file for a Gaussian grid
32!
33!-----------------------------------------------------------------------
34
35      use m_AttrVect,only    : AttrVect
36!      use m_GeneralGrid,only : MCT_GGrid_init => init
37      use m_GeneralGrid,only : MCT_GGrid_initUnstructured => initUnstructured
38      use m_GeneralGrid,only : MCT_GGrid_indexIA => indexIA
39      use m_GeneralGrid,only : MCT_GGrid_indexRA => indexRA
40      use m_GeneralGrid,only : GeneralGrid
41      use m_die
42      use m_stdio
43
44      implicit none
45
46!-----------------------------------------------------------------------
47!
48!     variables that describe the grid
49!
50!     T42: nx=128 ny=64
51!     T62: nx=192 ny=94
52!
53!-----------------------------------------------------------------------
54
55      type(GeneralGrid),       intent(out) :: GGrid
56      integer,                 intent(in) :: nx
57      integer,                 intent(in) :: ny
58
59      integer :: grid_size
60
61      integer, parameter :: &
62                   grid_rank = 2, &
63                   grid_corners = 4
64
65      integer, dimension(grid_rank) :: &
66                   grid_dims
67
68!-----------------------------------------------------------------------
69!
70!     grid coordinates and masks
71!
72!-----------------------------------------------------------------------
73
74      integer, dimension(:), allocatable :: &
75                   grid_imask
76
77      real, dimension(:), allocatable :: &
78                   grid_area      ,  & ! area weights
79                   grid_center_lat,  & ! lat/lon coordinates for
80                   grid_center_lon     ! each grid center in degrees
81
82      real, dimension(:,:), allocatable :: &
83                   grid_corner_lat,  & ! lat/lon coordinates for
84                   grid_corner_lon     ! each grid corner in degrees
85
86
87!-----------------------------------------------------------------------
88!
89!     defined constants
90!
91!-----------------------------------------------------------------------
92
93      real, parameter ::            &
94            zero   = 0.0,           &
95            one    = 1.0,           &
96            two    = 2.0,           &
97            three  = 3.0,           &
98            four   = 4.0,           &
99            five   = 5.0,           &
100            half   = 0.5,           &
101            quart  = 0.25,          &
102            bignum = 1.e+20,        &
103            tiny   = 1.e-14,        &
104            pi     = 3.14159265359, &
105            pi2    = two*pi,        &
106            pih    = half*pi
107
108!-----------------------------------------------------------------------
109!
110!     other local variables
111!
112!-----------------------------------------------------------------------
113
114      character(len=*),parameter :: myname_= 'convertgauss'
115
116      integer :: i, j, k, p, q, r, ier, atm_add
117
118      integer :: center_lat, center_lon, &
119                 corner_lat, corner_lon, &
120                 imask, area
121
122      real :: dlon, minlon, maxlon, centerlon, &
123              minlat, maxlat, centerlat
124
125      real, dimension(ny) :: gauss_root, gauss_wgt, gauss_lat
126
127      real, dimension(:), pointer :: PointData
128      integer :: offset
129
130!-----------------------------------------------------------------------
131!
132!     compute longitudes of cell centers and corners.  set up alon
133!     array for search routine.
134!
135!-----------------------------------------------------------------------
136
137      grid_size = nx*ny
138
139      allocate(grid_imask(grid_size), &
140               grid_area(grid_size),  &
141               grid_center_lat(grid_size), &
142               grid_center_lon(grid_size), &
143               grid_corner_lat(grid_corners,grid_size), &
144               grid_corner_lon(grid_corners,grid_size), stat=ier)
145
146      if(ier/=0) call die(myname_,"allocate(grid_imask... ", ier)
147
148      grid_dims(1) = nx
149      grid_dims(2) = ny
150
151      dlon = 360./nx
152
153      do i=1,nx
154
155        centerlon = (i-1)*dlon
156        minlon = centerlon - half*dlon
157        maxlon = centerlon + half*dlon
158
159        do j=1,ny
160          atm_add = (j-1)*nx + i
161
162          grid_center_lon(atm_add  ) = centerlon
163          grid_corner_lon(1,atm_add) = minlon
164          grid_corner_lon(2,atm_add) = maxlon
165          grid_corner_lon(3,atm_add) = maxlon
166          grid_corner_lon(4,atm_add) = minlon
167        end do
168
169      end do
170
171!-----------------------------------------------------------------------
172!
173!     compute Gaussian latitudes and store in gauss_wgt.
174!
175!-----------------------------------------------------------------------
176
177      call gquad(ny, gauss_root, gauss_wgt)
178      do j=1,ny
179        gauss_lat(j) = pih - gauss_root(ny+1-j)
180      end do
181
182!-----------------------------------------------------------------------
183!
184!     compute latitudes at cell centers and corners.  set up alat
185!     array for search routine.
186!
187!-----------------------------------------------------------------------
188
189      do j=1,ny
190        centerlat = gauss_lat(j)
191
192        if (j .eq. 1) then
193          minlat = -pih
194        else
195          minlat = ATAN((COS(gauss_lat(j-1)) - &
196                         COS(gauss_lat(j  )))/ &
197                        (SIN(gauss_lat(j  )) - &
198                         SIN(gauss_lat(j-1))))
199        endif
200
201        if (j .eq. ny) then
202          maxlat = pih
203        else
204          maxlat = ATAN((COS(gauss_lat(j  )) - &
205                         COS(gauss_lat(j+1)))/ &
206                        (SIN(gauss_lat(j+1)) - &
207                         SIN(gauss_lat(j  ))))
208        endif
209
210        do i=1,nx
211          atm_add = (j-1)*nx + i
212          grid_center_lat(atm_add  ) = centerlat*360./pi2
213          grid_corner_lat(1,atm_add) = minlat*360./pi2
214          grid_corner_lat(2,atm_add) = minlat*360./pi2
215          grid_corner_lat(3,atm_add) = maxlat*360./pi2
216          grid_corner_lat(4,atm_add) = maxlat*360./pi2
217          grid_area(atm_add) = gauss_wgt(j)*pi2/nx
218        end do
219
220      end do
221
222!-----------------------------------------------------------------------
223!
224!     define mask
225!
226!-----------------------------------------------------------------------
227
228      grid_imask = 1
229
230!-----------------------------------------------------------------------
231!
232!     intialize GeneralGrid
233!
234!-----------------------------------------------------------------------
235
236!      call MCT_GGrid_init(GGrid=GGrid, &
237!                          CoordChars="grid_center_lat:&
238!                                    &grid_center_lon", &
239!                         WeightChars="grid_area", &
240!                         OtherChars="grid_corner_lat_1:&
241!                                    &grid_corner_lat_2:&
242!                                    &grid_corner_lat_3:&
243!                                    &grid_corner_lat_4:&
244!                                    &grid_corner_lon_1:&
245!                                    &grid_corner_lon_2:&
246!                                    &grid_corner_lon_3:&
247!                                    &grid_corner_lon_4", &
248!                         IndexChars="grid_imask", &
249!                          lsize=grid_size)
250
251! Create and fill PointData(:) array for unstructured-style GeneralGrid_init
252
253      allocate(PointData(2*grid_size), stat=ier)
254      if(ier /= 0) then
255         write(stderr,'(2a,i8)') myname_, &
256              ':: allocate(PointData(...) failed with ier=',ier
257         call die(myname_)
258      endif
259
260      do i=1,grid_size
261         offset = 2 * (i-1)
262         PointData(offset+1) = grid_center_lat(i)
263         PointData(offset+2) = grid_center_lon(i)
264      end do
265
266      call MCT_GGrid_initUnstructured(GGrid=GGrid, &
267                                     CoordChars="grid_center_lat:&
268                                     &grid_center_lon", &
269                                     CoordSortOrder="grid_center_lat:&
270                                     &grid_center_lon", &
271                                     WeightChars="grid_area", &
272                                     OtherChars="grid_corner_lat_1:&
273                                     &grid_corner_lat_2:&
274                                     &grid_corner_lat_3:&
275                                     &grid_corner_lat_4:&
276                                     &grid_corner_lon_1:&
277                                     &grid_corner_lon_2:&
278                                     &grid_corner_lon_3:&
279                                     &grid_corner_lon_4", &
280                                     IndexChars="grid_imask", &
281                                     nDims=2, nPoints=grid_size, &
282                                     PointData=PointData)
283
284      deallocate(PointData, stat=ier)
285      if(ier /= 0) then
286         write(stderr,'(2a,i8)') myname_, &
287              ':: deallocate(PointData...) failed with ier=',ier
288         call die(myname_)
289      endif
290
291!      center_lat = MCT_GGrid_indexRA(GGrid,'grid_center_lat')
292!      center_lon = MCT_GGrid_indexRA(GGrid,'grid_center_lon')
293      corner_lat = MCT_GGrid_indexRA(GGrid,'grid_corner_lat_1')
294      corner_lon = MCT_GGrid_indexRA(GGrid,'grid_corner_lon_1')
295      area       = MCT_GGrid_indexRA(GGrid,'grid_area')
296      imask      = MCT_GGrid_indexIA(GGrid,'grid_imask')
297
298!      GGrid%data%rattr(center_lat,1:grid_size) = &
299!          grid_center_lat(1:grid_size)
300!      GGrid%data%rattr(center_lon,1:grid_size) = &
301!          grid_center_lon(1:grid_size)
302      GGrid%data%rattr(area,1:grid_size) = &
303           grid_area(1:grid_size)
304      GGrid%data%iattr(imask,1:grid_size) = &
305           grid_imask(1:grid_size)
306
307      do p = 1,grid_corners
308         GGrid%data%rattr(corner_lat+p-1,1:grid_size) = &
309              grid_corner_lat(p,1:grid_size)
310         GGrid%data%rattr(corner_lon+p-1,1:grid_size) = &
311           grid_corner_lon(p,1:grid_size)
312      enddo
313
314      deallocate(grid_imask, grid_area,            &
315                 grid_center_lat, grid_center_lon, &
316                 grid_corner_lat, grid_corner_lon, &
317                 stat=ier)
318
319      if(ier/=0) call die(myname_,"deallocate(grid_imask... ", ier)
320
321
322!-----------------------------------------------------------------------
323
324    end subroutine convertgauss
325
326!***********************************************************************
327
328      subroutine gquad(l,root,w)
329
330!-----------------------------------------------------------------------
331!
332!     This subroutine finds the l roots (in theta) and gaussian weights
333!     associated with the legendre polynomial of degree l > 1.
334!
335!-----------------------------------------------------------------------
336
337      use m_die
338
339      implicit none
340
341!-----------------------------------------------------------------------
342!
343!     intent(in)
344!
345!-----------------------------------------------------------------------
346
347      integer, intent(in) :: l
348
349!-----------------------------------------------------------------------
350!
351!     intent(out)
352!
353!-----------------------------------------------------------------------
354
355      real, dimension(l), intent(out) :: root, w
356
357!-----------------------------------------------------------------------
358!
359!     defined constants
360!
361!-----------------------------------------------------------------------
362
363      real, parameter ::            &
364            zero   = 0.0,           &
365            one    = 1.0,           &
366            two    = 2.0,           &
367            three  = 3.0,           &
368            four   = 4.0,           &
369            five   = 5.0,           &
370            half   = 0.5,           &
371            quart  = 0.25,          &
372            bignum = 1.e+20,        &
373            tiny   = 1.e-14,        &
374            pi     = 3.14159265359, &
375            pi2    = two*pi,        &
376            pih    = half*pi
377
378!-----------------------------------------------------------------------
379!
380!     local
381!
382!-----------------------------------------------------------------------
383
384      integer :: l1, l2, l22, l3, k, i, j, loop_counter
385
386      real :: del,co,p1,p2,p3,t1,t2,slope,s,c,pp1,pp2,p00
387
388!-----MUST adjust tolerance for newton convergence-----!
389
390      ! Modify tolerance level to the precision of the real numbers:
391      ! Increase for lower precision, decrease for higher precision.
392
393      real, parameter :: RTOL = 1.0e4*epsilon(0.)
394
395!------------------------------------------------------!
396
397!-----------------------------------------------------------------------
398!
399!     Define useful constants.
400!
401!-----------------------------------------------------------------------
402
403      del= pi/float(4*l)
404      l1 = l+1
405      co = float(2*l+3)/float(l1**2)
406      p2 = 1.0
407      t2 = -del
408      l2 = l/2
409      k = 1
410      p00 = one/sqrt(two)
411
412!-----------------------------------------------------------------------
413!
414!     Start search for each root by looking for crossing point.
415!
416!-----------------------------------------------------------------------
417
418      do i=1,l2
419   10    t1 = t2
420         t2 = t1+del
421         p1 = p2
422         s = sin(t2)
423         c = cos(t2)
424         pp1 = 1.0
425         p3 = p00
426         do j=1,l1
427            pp2 = pp1
428            pp1 = p3
429            p3 = 2.0*sqrt((float(j**2)-0.250)/float(j**2))*c*pp1- &
430                 sqrt(float((2*j+1)*(j-1)*(j-1))/ &
431                 float((2*j-3)*j*j))*pp2
432         end do
433         p2 = pp1
434         if ((k*p2).gt.0) goto 10
435
436!-----------------------------------------------------------------------
437!
438!        Now converge using Newton-Raphson.
439!
440!-----------------------------------------------------------------------
441
442         k = -k
443         loop_counter=0
444   20    continue
445            loop_counter=loop_counter+1
446            slope = (t2-t1)/(p2-p1)
447            t1 = t2
448            t2 = t2-slope*p2
449            p1 = p2
450            s = sin(t2)
451            c = cos(t2)
452            pp1 = 1.0
453            p3 = p00
454            do j=1,l1
455               pp2 = pp1
456               pp1 = p3
457               p3 = 2.0*sqrt((float(j**2)-0.250)/float(j**2))*c*pp1- &
458                    sqrt(float((2*j+1)*(j-1)*(j-1))/ &
459                    float((2*j-3)*j*j))*pp2
460            end do
461            p2 = pp1
462
463            if(loop_counter > 1e4) then
464               call die("subroutine gquad",&
465                    "ERROR:: Precision of reals is too low. &
466                    & Increase the magnitude of RTOL.",0)
467            endif
468
469         if (abs(p2).gt.RTOL) goto 20
470         root(i) = t2
471         w(i) = co*(sin(t2)/p3)**2
472      end do
473
474!-----------------------------------------------------------------------
475!
476!     If l is odd, take care of odd point.
477!
478!-----------------------------------------------------------------------
479
480      l22 = 2*l2
481      if (l22 .ne. l) then
482         l2 = l2+1
483         t2 = pi/2.0
484         root(l2) = t2
485         s = sin(t2)
486         c = cos(t2)
487         pp1 = 1.0
488         p3 = p00
489         do j=1,l1
490            pp2 = pp1
491            pp1 = p3
492            p3 = 2.0*sqrt((float(j**2)-0.250)/float(j**2))*c*pp1- &
493                 sqrt(float((2*j+1)*(j-1)*(j-1))/ &
494                 float((2*j-3)*j*j))*pp2
495         end do
496         p2 = pp1
497         w(l2) = co/p3**2
498      endif
499
500!-----------------------------------------------------------------------
501!
502!     Use symmetry to compute remaining roots and weights.
503!
504!-----------------------------------------------------------------------
505
506      l3 = l2+1
507      do i=l3,l
508         root(i) = pi-root(l-i+1)
509         w(i) = w(l-i+1)
510      end do
511
512!-----------------------------------------------------------------------
513
514      end subroutine gquad
515
516!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
Note: See TracBrowser for help on using the repository browser.