Neko  0.9.99
A portable framework for high-order spectral element flow simulations
speclib.f90
Go to the documentation of this file.
1 ! Copyright (c) 2008-2020, UCHICAGO ARGONNE, LLC.
2 !
3 ! The UChicago Argonne, LLC as Operator of Argonne National
4 ! Laboratory holds copyright in the Software. The copyright holder
5 ! reserves all rights except those expressly granted to licensees,
6 ! and U.S. Government license rights.
7 !
8 ! Redistribution and use in source and binary forms, with or without
9 ! modification, are permitted provided that the following conditions
10 ! are met:
11 !
12 ! 1. Redistributions of source code must retain the above copyright
13 ! notice, this list of conditions and the disclaimer below.
14 !
15 ! 2. Redistributions in binary form must reproduce the above copyright
16 ! notice, this list of conditions and the disclaimer (as noted below)
17 ! in the documentation and/or other materials provided with the
18 ! distribution.
19 !
20 ! 3. Neither the name of ANL nor the names of its contributors
21 ! may be used to endorse or promote products derived from this software
22 ! without specific prior written permission.
23 !
24 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
25 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
26 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
27 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL
28 ! UCHICAGO ARGONNE, LLC, THE U.S. DEPARTMENT OF
29 ! ENERGY OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
30 ! SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
31 ! TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
32 ! DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
33 ! THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
34 ! (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
35 ! OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 !
37 ! Additional BSD Notice
38 ! ---------------------
39 ! 1. This notice is required to be provided under our contract with
40 ! the U.S. Department of Energy (DOE). This work was produced at
41 ! Argonne National Laboratory under Contract
42 ! No. DE-AC02-06CH11357 with the DOE.
43 !
44 ! 2. Neither the United States Government nor UCHICAGO ARGONNE,
45 ! LLC nor any of their employees, makes any warranty,
46 ! express or implied, or assumes any liability or responsibility for the
47 ! accuracy, completeness, or usefulness of any information, apparatus,
48 ! product, or process disclosed, or represents that its use would not
49 ! infringe privately-owned rights.
50 !
51 ! 3. Also, reference herein to any specific commercial products, process,
52 ! or services by trade name, trademark, manufacturer or otherwise does
53 ! not necessarily constitute or imply its endorsement, recommendation,
54 ! or favoring by the United States Government or UCHICAGO ARGONNE LLC.
55 ! The views and opinions of authors expressed
56 ! herein do not necessarily state or reflect those of the United States
57 ! Government or UCHICAGO ARGONNE, LLC, and shall
58 ! not be used for advertising or product endorsement purposes.
59 !
60 !==============================================================================
61 !
75 !
76 !------------------------------------------------------------------------------
77 !
78 ! ABBRIVIATIONS:
79 !
80 ! M - Set of mesh points
81 ! Z - Set of collocation/quadrature points
82 ! W - Set of quadrature weights
83 ! H - Lagrangian interpolant
84 ! D - Derivative operator
85 ! I - Interpolation operator
86 ! GL - Gauss Legendre
87 ! GLL - Gauss-Lobatto Legendre
88 ! GJ - Gauss Jacobi
89 ! GLJ - Gauss-Lobatto Jacobi
90 !
91 !
92 ! MAIN ROUTINES:
93 !
94 ! Points and weights:
95 !
96 ! ZWGL Compute Gauss Legendre points and weights
97 ! ZWGLL Compute Gauss-Lobatto Legendre points and weights
98 ! ZWGJ Compute Gauss Jacobi points and weights (general)
99 ! ZWGLJ Compute Gauss-Lobatto Jacobi points and weights (general)
100 !
101 ! Lagrangian interpolants:
102 !
103 ! HGL Compute Gauss Legendre Lagrangian interpolant
104 ! HGLL Compute Gauss-Lobatto Legendre Lagrangian interpolant
105 ! HGJ Compute Gauss Jacobi Lagrangian interpolant (general)
106 ! HGLJ Compute Gauss-Lobatto Jacobi Lagrangian interpolant (general)
107 !
108 ! Derivative operators:
109 !
110 ! DGLL Compute Gauss-Lobatto Legendre derivative matrix
111 ! DGLLGL Compute derivative matrix for a staggered mesh (GLL->GL)
112 ! DGJ Compute Gauss Jacobi derivative matrix (general)
113 ! DGLJ Compute Gauss-Lobatto Jacobi derivative matrix (general)
114 ! DGLJGJ Compute derivative matrix for a staggered mesh (GLJ->GJ) (general)
115 !
116 ! Interpolation operators:
117 !
118 ! IGLM Compute interpolation operator GL -> M
119 ! IGLLM Compute interpolation operator GLL -> M
120 ! IGJM Compute interpolation operator GJ -> M (general)
121 ! IGLJM Compute interpolation operator GLJ -> M (general)
122 !
123 ! Other:
124 !
125 ! PNLEG Compute Legendre polynomial of degree N
126 ! legendre_poly Compute Legendre polynomial of degree 0-N
127 ! PNDLEG Compute derivative of Legendre polynomial of degree N
128 !
129 ! Comments:
130 !
131 ! Note that many of the above routines exist in both single and
132 ! double precision. If the name of the single precision routine is
133 ! SUB, the double precision version is called SUBD. In most cases
134 ! all the "low-level" arithmetic is done in double precision, even
135 ! for the single precsion versions.
136 !
137 ! Useful references:
138 !
139 ! [1] Gabor Szego: Orthogonal Polynomials, American Mathematical Society,
140 ! Providence, Rhode Island, 1939.
141 ! [2] Abramowitz & Stegun: Handbook of Mathematical Functions,
142 ! Dover, New York, 1972.
143 ! [3] Canuto, Hussaini, Quarteroni & Zang: Spectral Methods in Fluid
144 ! Dynamics, Springer-Verlag, 1988.
145 !
146 !
147 !==============================================================================
148 module speclib
149  use num_types, only : rp, xp
150  use utils, only: neko_error
151 
152 
153 contains
160  SUBROUTINE zwgl (Z,W,NP)
161  REAL(KIND=rp) z(1),w(1), alpha, beta
162  alpha = 0.
163  beta = 0.
164  CALL zwgj (z,w,np,alpha,beta)
165  RETURN
166  end subroutine zwgl
167 
168  SUBROUTINE zwgll (Z,W,NP)
169 !--------------------------------------------------------------------
170 !
171 ! Generate NP Gauss-Lobatto Legendre points (Z) and weights (W)
172 ! associated with Jacobi polynomial P(N)(alpha=0,beta=0).
173 ! The polynomial degree N=NP-1.
174 ! Z and W are in single precision, but all the arithmetic
175 ! operations are done in double precision.
176 !
177 !--------------------------------------------------------------------
178  REAL(KIND=rp) z(1),w(1), alpha, beta
179  alpha = 0.
180  beta = 0.
181  CALL zwglj (z,w,np,alpha,beta)
182  RETURN
183  end subroutine zwgll
184 
185  SUBROUTINE zwgj (Z,W,NP,ALPHA,BETA)
186 !--------------------------------------------------------------------
187 !
188 ! Generate NP GAUSS JACOBI points (Z) and weights (W)
189 ! associated with Jacobi polynomial P(N)(alpha>-1,beta>-1).
190 ! The polynomial degree N=NP-1.
191 ! Single precision version.
192 !
193 !--------------------------------------------------------------------
194  parameter(nmax=84)
195  parameter(nzd = nmax)
196  REAL(KIND=xp) zd(nzd),wd(nzd),alphad,betad
197  REAL(KIND=rp) z(1),w(1),alpha,beta
198 
199  npmax = nzd
200  IF (np.GT.npmax) THEN
201  WRITE (6,*) 'Too large polynomial degree in ZWGJ'
202  WRITE (6,*) 'Maximum polynomial degree is',nmax
203  WRITE (6,*) 'Here NP=',np
204  call neko_error
205  ENDIF
206  alphad = alpha
207  betad = beta
208  CALL zwgjd (zd,wd,np,alphad,betad)
209  DO 100 i=1,np
210  z(i) = zd(i)
211  w(i) = wd(i)
212 100 CONTINUE
213  RETURN
214  end subroutine zwgj
215 
216  SUBROUTINE zwgjd (Z,W,NP,ALPHA,BETA)
217 !--------------------------------------------------------------------
218 !
219 ! Generate NP GAUSS JACOBI points (Z) and weights (W)
220 ! associated with Jacobi polynomial P(N)(alpha>-1,beta>-1).
221 ! The polynomial degree N=NP-1.
222 ! Double precision version.
223 !
224 !--------------------------------------------------------------------
225  IMPLICIT REAL(KIND=XP) (a-h,o-z)
226  REAL(KIND=xp) z(1),w(1),alpha,beta
227 
228  n = np-1
229  dn = ((n))
230  one = 1.
231  two = 2.
232  apb = alpha+beta
233 
234  IF (np.LE.0) THEN
235  WRITE (6,*) 'ZWGJD: Minimum number of Gauss points is 1',np
236  call neko_error
237  ENDIF
238  IF ((alpha.LE.-one).OR.(beta.LE.-one)) THEN
239  WRITE (6,*) 'ZWGJD: Alpha and Beta must be greater than -1'
240  call neko_error
241  ENDIF
242 
243  IF (np.EQ.1) THEN
244  z(1) = (beta-alpha)/(apb+two)
245  w(1) = gammaf(alpha+one)*gammaf(beta+one)/gammaf(apb+two) &
246  * two**(apb+one)
247  RETURN
248  ENDIF
249 
250  CALL jacg (z,np,alpha,beta)
251 
252  np1 = n+1
253  np2 = n+2
254  dnp1 = ((np1))
255  dnp2 = ((np2))
256  fac1 = dnp1+alpha+beta+one
257  fac2 = fac1+dnp1
258  fac3 = fac2+one
259  fnorm = pnormj(np1,alpha,beta)
260  rcoef = (fnorm*fac2*fac3)/(two*fac1*dnp2)
261  DO 100 i=1,np
262  CALL jacobf (p,pd,pm1,pdm1,pm2,pdm2,np2,alpha,beta,z(i))
263  w(i) = -rcoef/(p*pdm1)
264 100 CONTINUE
265  RETURN
266  end subroutine zwgjd
267 
268  SUBROUTINE zwglj (Z,W,NP,ALPHA,BETA)
269 !--------------------------------------------------------------------
270 !
271 ! Generate NP GAUSS LOBATTO JACOBI points (Z) and weights (W)
272 ! associated with Jacobi polynomial P(N)(alpha>-1,beta>-1).
273 ! The polynomial degree N=NP-1.
274 ! Single precision version.
275 !
276 !--------------------------------------------------------------------
277  parameter(nmax=84)
278  parameter(nzd = nmax)
279  REAL(KIND=xp) zd(nzd),wd(nzd),alphad,betad
280  REAL(KIND=rp) z(1),w(1),alpha,beta
281 
282  npmax = nzd
283  IF (np.GT.npmax) THEN
284  WRITE (6,*) 'Too large polynomial degree in ZWGLJ'
285  WRITE (6,*) 'Maximum polynomial degree is',nmax
286  WRITE (6,*) 'Here NP=',np
287  call neko_error
288  ENDIF
289  alphad = alpha
290  betad = beta
291  CALL zwgljd (zd,wd,np,alphad,betad)
292  DO 100 i=1,np
293  z(i) = zd(i)
294  w(i) = wd(i)
295 100 CONTINUE
296  RETURN
297  end subroutine zwglj
298 
299  SUBROUTINE zwgljd (Z,W,NP,ALPHA,BETA)
300 !--------------------------------------------------------------------
301 !
302 ! Generate NP GAUSS LOBATTO JACOBI points (Z) and weights (W)
303 ! associated with Jacobi polynomial P(N)(alpha>-1,beta>-1).
304 ! The polynomial degree N=NP-1.
305 ! Double precision version.
306 !
307 !--------------------------------------------------------------------
308  IMPLICIT REAL(KIND=XP) (a-h,o-z)
309  REAL(KIND=xp) z(np),w(np),alpha,beta
310 
311  n = np-1
312  nm1 = n-1
313  one = 1.
314  two = 2.
315 
316  IF (np.LE.1) THEN
317  WRITE (6,*) 'ZWGLJD: Minimum number of Gauss-Lobatto points is 2'
318  WRITE (6,*) 'ZWGLJD: alpha,beta:',alpha,beta,np
319  call neko_error
320  ENDIF
321  IF ((alpha.LE.-one).OR.(beta.LE.-one)) THEN
322  WRITE (6,*) 'ZWGLJD: Alpha and Beta must be greater than -1'
323  call neko_error
324  ENDIF
325 
326  IF (nm1.GT.0) THEN
327  alpg = alpha+one
328  betg = beta+one
329  CALL zwgjd (z(2),w(2),nm1,alpg,betg)
330  ENDIF
331  z(1) = -one
332  z(np) = one
333  DO 100 i=2,np-1
334  w(i) = w(i)/(one-z(i)**2)
335 100 CONTINUE
336  CALL jacobf (p,pd,pm1,pdm1,pm2,pdm2,n,alpha,beta,z(1))
337  w(1) = endw1(n,alpha,beta)/(two*pd)
338  CALL jacobf (p,pd,pm1,pdm1,pm2,pdm2,n,alpha,beta,z(np))
339  w(np) = endw2(n,alpha,beta)/(two*pd)
340 
341 ! RETURN
342  end subroutine zwgljd
343 
344  REAL(kind=xp) FUNCTION endw1 (N,ALPHA,BETA)
345  IMPLICIT REAL(KIND=XP) (a-h,o-z)
346  REAL(kind=xp) alpha,beta
347  zero = 0.
348  one = 1.
349  two = 2.
350  three = 3.
351  four = 4.
352  apb = alpha+beta
353  IF (n.EQ.0) THEN
354  endw1 = zero
355  RETURN
356  ENDIF
357  f1 = gammaf(alpha+two)*gammaf(beta+one)/gammaf(apb+three)
358  f1 = f1*(apb+two)*two**(apb+two)/two
359  IF (n.EQ.1) THEN
360  endw1 = f1
361  RETURN
362  ENDIF
363  fint1 = gammaf(alpha+two)*gammaf(beta+one)/gammaf(apb+three)
364  fint1 = fint1*two**(apb+two)
365  fint2 = gammaf(alpha+two)*gammaf(beta+two)/gammaf(apb+four)
366  fint2 = fint2*two**(apb+three)
367  f2 = (-two*(beta+two)*fint1 + (apb+four)*fint2) &
368  * (apb+three)/four
369  IF (n.EQ.2) THEN
370  endw1 = f2
371  RETURN
372  ENDIF
373  DO 100 i=3,n
374  di = ((i-1))
375  abn = alpha+beta+di
376  abnn = abn+di
377  a1 = -(two*(di+alpha)*(di+beta))/(abn*abnn*(abnn+one))
378  a2 = (two*(alpha-beta))/(abnn*(abnn+two))
379  a3 = (two*(abn+one))/((abnn+two)*(abnn+one))
380  f3 = -(a2*f2+a1*f1)/a3
381  f1 = f2
382  f2 = f3
383 100 CONTINUE
384  endw1 = f3
385  RETURN
386  end function endw1
387 
388  REAL(kind=xp) FUNCTION endw2 (N,ALPHA,BETA)
389  IMPLICIT REAL(KIND=XP) (a-h,o-z)
390  REAL(kind=xp) alpha,beta
391  zero = 0.
392  one = 1.
393  two = 2.
394  three = 3.
395  four = 4.
396  apb = alpha+beta
397  IF (n.EQ.0) THEN
398  endw2 = zero
399  RETURN
400  ENDIF
401  f1 = gammaf(alpha+one)*gammaf(beta+two)/gammaf(apb+three)
402  f1 = f1*(apb+two)*two**(apb+two)/two
403  IF (n.EQ.1) THEN
404  endw2 = f1
405  RETURN
406  ENDIF
407  fint1 = gammaf(alpha+one)*gammaf(beta+two)/gammaf(apb+three)
408  fint1 = fint1*two**(apb+two)
409  fint2 = gammaf(alpha+two)*gammaf(beta+two)/gammaf(apb+four)
410  fint2 = fint2*two**(apb+three)
411  f2 = (two*(alpha+two)*fint1 - (apb+four)*fint2) &
412  * (apb+three)/four
413  IF (n.EQ.2) THEN
414  endw2 = f2
415  RETURN
416  ENDIF
417  DO 100 i=3,n
418  di = ((i-1))
419  abn = alpha+beta+di
420  abnn = abn+di
421  a1 = -(two*(di+alpha)*(di+beta))/(abn*abnn*(abnn+one))
422  a2 = (two*(alpha-beta))/(abnn*(abnn+two))
423  a3 = (two*(abn+one))/((abnn+two)*(abnn+one))
424  f3 = -(a2*f2+a1*f1)/a3
425  f1 = f2
426  f2 = f3
427 100 CONTINUE
428  endw2 = f3
429  RETURN
430  end function endw2
431 
432  REAL(kind=xp) FUNCTION gammaf (X)
433  IMPLICIT REAL(KIND=XP) (a-h,o-z)
434  REAL(kind=xp) x
435  zero = 0.0
436  half = 0.5
437  one = 1.0
438  two = 2.0
439  four = 4.0
440  pi = four*atan(one)
441  gammaf = one
442  IF (x.EQ.-half) gammaf = -two*sqrt(pi)
443  IF (x.EQ. half) gammaf = sqrt(pi)
444  IF (x.EQ. one ) gammaf = one
445  IF (x.EQ. two ) gammaf = one
446  IF (x.EQ. 1.5 ) gammaf = sqrt(pi)/2.
447  IF (x.EQ. 2.5) gammaf = 1.5*sqrt(pi)/2.
448  IF (x.EQ. 3.5) gammaf = 0.5*(2.5*(1.5*sqrt(pi)))
449  IF (x.EQ. 3. ) gammaf = 2.
450  IF (x.EQ. 4. ) gammaf = 6.
451  IF (x.EQ. 5. ) gammaf = 24.
452  IF (x.EQ. 6. ) gammaf = 120.
453  RETURN
454  end function gammaf
455 
456  REAL(kind=xp) FUNCTION pnormj (N,ALPHA,BETA)
457  IMPLICIT REAL(KIND=XP) (a-h,o-z)
458  REAL(kind=xp) alpha,beta
459  one = 1.
460  two = 2.
461  dn = ((n))
462  const = alpha+beta+one
463  IF (n.LE.1) THEN
464  prod = gammaf(dn+alpha)*gammaf(dn+beta)
465  prod = prod/(gammaf(dn)*gammaf(dn+alpha+beta))
466  pnormj = prod * two**const/(two*dn+const)
467  RETURN
468  ENDIF
469  prod = gammaf(alpha+one)*gammaf(beta+one)
470  prod = prod/(two*(one+const)*gammaf(const+one))
471  prod = prod*(one+alpha)*(two+alpha)
472  prod = prod*(one+beta)*(two+beta)
473  DO 100 i=3,n
474  dindx = ((i))
475  frac = (dindx+alpha)*(dindx+beta)/(dindx*(dindx+alpha+beta))
476  prod = prod*frac
477 100 CONTINUE
478  pnormj = prod * two**const/(two*dn+const)
479  RETURN
480  end function pnormj
481 
482  SUBROUTINE jacg (XJAC,NP,ALPHA,BETA)
483 !--------------------------------------------------------------------
484 !
485 ! Compute NP Gauss points XJAC, which are the zeros of the
486 ! Jacobi polynomial J(NP) with parameters ALPHA and BETA.
487 ! ALPHA and BETA determines the specific type of Gauss points.
488 ! Examples:
489 ! ALPHA = BETA = 0.0 -> Legendre points
490 ! ALPHA = BETA = -0.5 -> Chebyshev points
491 !
492 !--------------------------------------------------------------------
493  IMPLICIT REAL(KIND=XP) (a-h,o-z)
494  REAL(KIND=xp) xjac(1)
495  DATA kstop /10/
496  DATA eps/1.0e-12_rp/
497  n = np-1
498  one = 1.
499  dth = 4.*atan(one)/(2.*((n))+2.)
500  DO 40 j=1,np
501  IF (j.EQ.1) THEN
502  x = cos((2.*(((j))-1.)+1.)*dth)
503  ELSE
504  x1 = cos((2.*(((j))-1.)+1.)*dth)
505  x2 = xlast
506  x = (x1+x2)/2.
507  ENDIF
508  DO 30 k=1,kstop
509  CALL jacobf (p,pd,pm1,pdm1,pm2,pdm2,np,alpha,beta,x)
510  recsum = 0.
511  jm = j-1
512  DO 29 i=1,jm
513  recsum = recsum+1./(x-xjac(np-i+1))
514 29 CONTINUE
515  delx = -p/(pd-recsum*p)
516  x = x+delx
517  IF (abs(delx) .LT. eps) GOTO 31
518 30 CONTINUE
519 31 CONTINUE
520  xjac(np-j+1) = x
521  xlast = x
522 40 CONTINUE
523  DO 200 i=1,np
524  xmin = 2.
525  DO 100 j=i,np
526  IF (xjac(j).LT.xmin) THEN
527  xmin = xjac(j)
528  jmin = j
529  ENDIF
530 100 CONTINUE
531  IF (jmin.NE.i) THEN
532  swap = xjac(i)
533  xjac(i) = xjac(jmin)
534  xjac(jmin) = swap
535  ENDIF
536 200 CONTINUE
537  RETURN
538  end subroutine jacg
539 
540  SUBROUTINE jacobf (POLY,PDER,POLYM1,PDERM1,POLYM2,PDERM2,N,ALP,BET,X)
541 !--------------------------------------------------------------------
542 !
543 ! Computes the Jacobi polynomial (POLY) and its derivative (PDER)
544 ! of degree N at X.
545 !
546 !--------------------------------------------------------------------
547  IMPLICIT REAL(KIND=XP) (a-h,o-z)
548  apb = alp+bet
549  poly = 1.
550  pder = 0.
551  IF (n .EQ. 0) RETURN
552  polyl = poly
553  pderl = pder
554  poly = (alp-bet+(apb+2.)*x)/2.
555  pder = (apb+2.)/2.
556  IF (n .EQ. 1) RETURN
557  DO 20 k=2,n
558  dk = ((k))
559  a1 = 2.*dk*(dk+apb)*(2.*dk+apb-2.)
560  a2 = (2.*dk+apb-1.)*(alp**2-bet**2)
561  b3 = (2.*dk+apb-2.)
562  a3 = b3*(b3+1.)*(b3+2.)
563  a4 = 2.*(dk+alp-1.)*(dk+bet-1.)*(2.*dk+apb)
564  polyn = ((a2+a3*x)*poly-a4*polyl)/a1
565  pdern = ((a2+a3*x)*pder-a4*pderl+a3*poly)/a1
566  psave = polyl
567  pdsave = pderl
568  polyl = poly
569  poly = polyn
570  pderl = pder
571  pder = pdern
572 20 CONTINUE
573  polym1 = polyl
574  pderm1 = pderl
575  polym2 = psave
576  pderm2 = pdsave
577  RETURN
578  end subroutine jacobf
579 
580  REAL(kind=xp) FUNCTION hgj (II,Z,ZGJ,NP,ALPHA,BETA)
581 !---------------------------------------------------------------------
582 !
583 ! Compute the value of the Lagrangian interpolant HGJ through
584 ! the NP Gauss Jacobi points ZGJ at the point Z.
585 ! Single precision version.
586 !
587 !---------------------------------------------------------------------
588  parameter(nmax=84)
589  parameter(nzd = nmax)
590  REAL(kind=xp) zd,zgjd(nzd),alphad,betad
591  REAL(kind=xp) z,zgj(1),alpha,beta
592  npmax = nzd
593  IF (np.GT.npmax) THEN
594  WRITE (6,*) 'Too large polynomial degree in HGJ'
595  WRITE (6,*) 'Maximum polynomial degree is',nmax
596  WRITE (6,*) 'Here NP=',np
597  call neko_error
598  ENDIF
599  zd = z
600  DO 100 i=1,np
601  zgjd(i) = zgj(i)
602 100 CONTINUE
603  alphad = alpha
604  betad = beta
605  hgj = hgjd(ii,zd,zgjd,np,alphad,betad)
606  RETURN
607  end function hgj
608 
609  REAL(kind=xp) FUNCTION hgjd (II,Z,ZGJ,NP,ALPHA,BETA)
610 !---------------------------------------------------------------------
611 !
612 ! Compute the value of the Lagrangian interpolant HGJD through
613 ! the NZ Gauss-Lobatto Jacobi points ZGJ at the point Z.
614 ! Double precision version.
615 !
616 !---------------------------------------------------------------------
617  IMPLICIT REAL(KIND=XP) (a-h,o-z)
618  REAL(kind=xp) z,zgj(1),alpha,beta
619  eps = 1.e-5
620  one = 1.
621  zi = zgj(ii)
622  dz = z-zi
623  IF (abs(dz).LT.eps) THEN
624  hgjd = one
625  RETURN
626  ENDIF
627  CALL jacobf (pzi,pdzi,pm1,pdm1,pm2,pdm2,np,alpha,beta,zi)
628  CALL jacobf (pz,pdz,pm1,pdm1,pm2,pdm2,np,alpha,beta,z)
629  hgjd = pz/(pdzi*(z-zi))
630  RETURN
631  end function hgjd
632 
633  REAL(kind=xp) FUNCTION hglj (II,Z,ZGLJ,NP,ALPHA,BETA)
634 !---------------------------------------------------------------------
635 !
636 ! Compute the value of the Lagrangian interpolant HGLJ through
637 ! the NZ Gauss-Lobatto Jacobi points ZGLJ at the point Z.
638 ! Single precision version.
639 !
640 !---------------------------------------------------------------------
641  parameter(nmax=84)
642  parameter(nzd = nmax)
643  REAL(kind=xp) zd,zgljd(nzd),alphad,betad
644  REAL(kind=xp) z,zglj(1),alpha,beta
645  npmax = nzd
646  IF (np.GT.npmax) THEN
647  WRITE (6,*) 'Too large polynomial degree in HGLJ'
648  WRITE (6,*) 'Maximum polynomial degree is',nmax
649  WRITE (6,*) 'Here NP=',np
650  call neko_error
651  ENDIF
652  zd = z
653  DO 100 i=1,np
654  zgljd(i) = zglj(i)
655 100 CONTINUE
656  alphad = alpha
657  betad = beta
658  hglj = hgljd(ii,zd,zgljd,np,alphad,betad)
659  RETURN
660  end function hglj
661 
662  REAL(kind=xp) FUNCTION hgljd (I,Z,ZGLJ,NP,ALPHA,BETA)
663 !---------------------------------------------------------------------
664 !
665 ! Compute the value of the Lagrangian interpolant HGLJD through
666 ! the NZ Gauss-Lobatto Jacobi points ZJACL at the point Z.
667 ! Double precision version.
668 !
669 !---------------------------------------------------------------------
670  IMPLICIT REAL(KIND=XP) (a-h,o-z)
671  REAL(kind=xp) z,zglj(1),alpha,beta
672  eps = 1.e-5
673  one = 1.
674  zi = zglj(i)
675  dz = z-zi
676  IF (abs(dz).LT.eps) THEN
677  hgljd = one
678  RETURN
679  ENDIF
680  n = np-1
681  dn = ((n))
682  eigval = -dn*(dn+alpha+beta+one)
683  CALL jacobf (pi,pdi,pm1,pdm1,pm2,pdm2,n,alpha,beta,zi)
684  const = eigval*pi+alpha*(one+zi)*pdi-beta*(one-zi)*pdi
685  CALL jacobf (p,pd,pm1,pdm1,pm2,pdm2,n,alpha,beta,z)
686  hgljd = (one-z**2)*pd/(const*(z-zi))
687  RETURN
688  end function hgljd
689 
690  SUBROUTINE dgj (D,DT,Z,NZ,NZD,ALPHA,BETA)
691 !-----------------------------------------------------------------
692 !
693 ! Compute the derivative matrix D and its transpose DT
694 ! associated with the Nth order Lagrangian interpolants
695 ! through the NZ Gauss Jacobi points Z.
696 ! Note: D and DT are square matrices.
697 ! Single precision version.
698 !
699 !-----------------------------------------------------------------
700  parameter(nmax=84)
701  parameter(nzdd = nmax)
702  REAL(KIND=xp) dd(nzdd,nzdd),dtd(nzdd,nzdd),zd(nzdd),alphad,betad
703  REAL(KIND=xp) d(nzd,nzd),dt(nzd,nzd),z(1),alpha,beta
704 
705  IF (nz.LE.0) THEN
706  WRITE (6,*) 'DGJ: Minimum number of Gauss points is 1'
707  call neko_error
708  ENDIF
709  IF (nz .GT. nmax) THEN
710  WRITE (6,*) 'Too large polynomial degree in DGJ'
711  WRITE (6,*) 'Maximum polynomial degree is',nmax
712  WRITE (6,*) 'Here Nz=',nz
713  call neko_error
714  ENDIF
715  IF ((alpha.LE.-1.).OR.(beta.LE.-1.)) THEN
716  WRITE (6,*) 'DGJ: Alpha and Beta must be greater than -1'
717  call neko_error
718  ENDIF
719  alphad = alpha
720  betad = beta
721  DO 100 i=1,nz
722  zd(i) = z(i)
723 100 CONTINUE
724  CALL dgjd (dd,dtd,zd,nz,nzdd,alphad,betad)
725  DO i=1,nz
726  DO j=1,nz
727  d(i,j) = dd(i,j)
728  dt(i,j) = dtd(i,j)
729  END DO
730  END DO
731  RETURN
732  end subroutine dgj
733 
734  SUBROUTINE dgjd (D,DT,Z,NZ,NZD,ALPHA,BETA)
735 !-----------------------------------------------------------------
736 !
737 ! Compute the derivative matrix D and its transpose DT
738 ! associated with the Nth order Lagrangian interpolants
739 ! through the NZ Gauss Jacobi points Z.
740 ! Note: D and DT are square matrices.
741 ! Double precision version.
742 !
743 !-----------------------------------------------------------------
744  IMPLICIT REAL(KIND=XP) (a-h,o-z)
745  REAL(KIND=xp) d(nzd,nzd),dt(nzd,nzd),z(1),alpha,beta
746  n = nz-1
747  dn = ((n))
748  one = 1.
749  two = 2.
750 
751  IF (nz.LE.1) THEN
752  WRITE (6,*) 'DGJD: Minimum number of Gauss-Lobatto points is 2'
753  call neko_error
754  ENDIF
755  IF ((alpha.LE.-one).OR.(beta.LE.-one)) THEN
756  WRITE (6,*) 'DGJD: Alpha and Beta must be greater than -1'
757  call neko_error
758  ENDIF
759 
760  DO i=1,nz
761  DO j=1,nz
762  CALL jacobf (pi,pdi,pm1,pdm1,pm2,pdm2,nz,alpha,beta,z(i))
763  CALL jacobf (pj,pdj,pm1,pdm1,pm2,pdm2,nz,alpha,beta,z(j))
764  IF (i.NE.j) d(i,j) = pdi/(pdj*(z(i)-z(j)))
765  IF (i.EQ.j) d(i,j) = ((alpha+beta+two)*z(i)+alpha-beta)/ &
766  (two*(one-z(i)**2))
767  dt(j,i) = d(i,j)
768  END DO
769  END DO
770  RETURN
771  end subroutine dgjd
772 
773  SUBROUTINE dglj (D,DT,Z,NZ,NZD,ALPHA,BETA)
774 !-----------------------------------------------------------------
775 !
776 ! Compute the derivative matrix D and its transpose DT
777 ! associated with the Nth order Lagrangian interpolants
778 ! through the NZ Gauss-Lobatto Jacobi points Z.
779 ! Note: D and DT are square matrices.
780 ! Single precision version.
781 !
782 !-----------------------------------------------------------------
783  parameter(nmax=84)
784  parameter(nzdd = nmax)
785  REAL(KIND=xp) dd(nzdd,nzdd),dtd(nzdd,nzdd),zd(nzdd),alphad,betad
786  REAL(KIND=xp) d(nzd,nzd),dt(nzd,nzd),z(1),alpha,beta
787 
788  IF (nz.LE.1) THEN
789  WRITE (6,*) 'DGLJ: Minimum number of Gauss-Lobatto points is 2'
790  call neko_error
791  ENDIF
792  IF (nz .GT. nmax) THEN
793  WRITE (6,*) 'Too large polynomial degree in DGLJ'
794  WRITE (6,*) 'Maximum polynomial degree is',nmax
795  WRITE (6,*) 'Here NZ=',nz
796  call neko_error
797  ENDIF
798  IF ((alpha.LE.-1.).OR.(beta.LE.-1.)) THEN
799  WRITE (6,*) 'DGLJ: Alpha and Beta must be greater than -1'
800  call neko_error
801  ENDIF
802  alphad = alpha
803  betad = beta
804  DO 100 i=1,nz
805  zd(i) = z(i)
806 100 CONTINUE
807  CALL dgljd (dd,dtd,zd,nz,nzdd,alphad,betad)
808  DO i=1,nz
809  DO j=1,nz
810  d(i,j) = dd(i,j)
811  dt(i,j) = dtd(i,j)
812  END DO
813  END DO
814  RETURN
815  end subroutine dglj
816 
817  SUBROUTINE dgljd (D,DT,Z,NZ,NZD,ALPHA,BETA)
818 !-----------------------------------------------------------------
819 !
820 ! Compute the derivative matrix D and its transpose DT
821 ! associated with the Nth order Lagrangian interpolants
822 ! through the NZ Gauss-Lobatto Jacobi points Z.
823 ! Note: D and DT are square matrices.
824 ! Double precision version.
825 !
826 !-----------------------------------------------------------------
827  IMPLICIT REAL(KIND=XP) (a-h,o-z)
828  REAL(KIND=xp) d(nzd,nzd),dt(nzd,nzd),z(1),alpha,beta
829  n = nz-1
830  dn = ((n))
831  one = 1.
832  two = 2.
833  eigval = -dn*(dn+alpha+beta+one)
834 
835  IF (nz.LE.1) THEN
836  WRITE (6,*) 'DGLJD: Minimum number of Gauss-Lobatto points is 2'
837  call neko_error
838  ENDIF
839  IF ((alpha.LE.-one).OR.(beta.LE.-one)) THEN
840  WRITE (6,*) 'DGLJD: Alpha and Beta must be greater than -1'
841  call neko_error
842  ENDIF
843 
844  DO i=1,nz
845  DO j=1,nz
846  CALL jacobf (pi,pdi,pm1,pdm1,pm2,pdm2,n,alpha,beta,z(i))
847  CALL jacobf (pj,pdj,pm1,pdm1,pm2,pdm2,n,alpha,beta,z(j))
848  ci = eigval*pi-(beta*(one-z(i))-alpha*(one+z(i)))*pdi
849  cj = eigval*pj-(beta*(one-z(j))-alpha*(one+z(j)))*pdj
850  IF (i.NE.j) d(i,j) = ci/(cj*(z(i)-z(j)))
851  IF ((i.EQ.j).AND.(i.NE.1).AND.(i.NE.nz)) &
852  d(i,j) = (alpha*(one+z(i))-beta*(one-z(i)))/ &
853  (two*(one-z(i)**2))
854  IF ((i.EQ.j).AND.(i.EQ.1)) &
855  d(i,j) = (eigval+alpha)/(two*(beta+two))
856  IF ((i.EQ.j).AND.(i.EQ.nz)) &
857  d(i,j) = -(eigval+beta)/(two*(alpha+two))
858  dt(j,i) = d(i,j)
859  END DO
860  END DO
861  RETURN
862  end subroutine dgljd
863 
864  SUBROUTINE dgll (D,DT,Z,NZ,NZD)
865 !-----------------------------------------------------------------
866 !
867 ! Compute the derivative matrix D and its transpose DT
868 ! associated with the Nth order Lagrangian interpolants
869 ! through the NZ Gauss-Lobatto Legendre points Z.
870 ! Note: D and DT are square matrices.
871 !
872 !-----------------------------------------------------------------
873  IMPLICIT REAL(KIND=XP) (a-h,o-z)
874  parameter(nmax=84)
875  REAL(KIND=rp) d(nzd,nzd),dt(nzd,nzd),z(1)
876  n = nz-1
877  IF (nz .GT. nmax) THEN
878  WRITE (6,*) 'Subroutine DGLL'
879  WRITE (6,*) 'Maximum polynomial degree =',nmax
880  WRITE (6,*) 'Polynomial degree =',nz
881  ENDIF
882  IF (nz .EQ. 1) THEN
883  d(1,1) = 0.
884  RETURN
885  ENDIF
886  fn = (n)
887  d0 = fn*(fn+1.)/4.
888  DO i=1,nz
889  DO j=1,nz
890  d(i,j) = 0.
891  IF (i.NE.j) d(i,j) = pnleg(real(z(i),xp),n)/ &
892  (pnleg(real(z(j),xp),n)*(z(i)-z(j)))
893  IF ((i.EQ.j).AND.(i.EQ.1)) d(i,j) = -d0
894  IF ((i.EQ.j).AND.(i.EQ.nz)) d(i,j) = d0
895  dt(j,i) = d(i,j)
896  END DO
897  END DO
898  RETURN
899  end subroutine dgll
900 
901  REAL(kind=xp) FUNCTION hgll (I,Z,ZGLL,NZ)
902 !---------------------------------------------------------------------
903 !
904 ! Compute the value of the Lagrangian interpolant L through
905 ! the NZ Gauss-Lobatto Legendre points ZGLL at the point Z.
906 !
907 !---------------------------------------------------------------------
908  IMPLICIT REAL(KIND=XP) (a-h,o-z)
909  REAL(kind=xp) zgll(1), eps, dz, z
910  eps = 1.e-5
911  dz = z - zgll(i)
912  IF (abs(dz) .LT. eps) THEN
913  hgll = 1.
914  RETURN
915  ENDIF
916  n = nz - 1
917  alfan = (n)*((n)+1.)
918  hgll = - (1.-z*z)*pndleg(z,n)/ &
919  (alfan*pnleg(zgll(i),n)*(z-zgll(i)))
920  RETURN
921  end function hgll
922 
923  REAL(kind=xp) FUNCTION hgl (I,Z,ZGL,NZ)
924 !---------------------------------------------------------------------
925 !
926 ! Compute the value of the Lagrangian interpolant HGL through
927 ! the NZ Gauss Legendre points ZGL at the point Z.
928 !
929 !---------------------------------------------------------------------
930  REAL(kind=xp) zgl(1), z, eps, dz
931  eps = 1.e-5
932  dz = z - zgl(i)
933  IF (abs(dz) .LT. eps) THEN
934  hgl = 1.
935  RETURN
936  ENDIF
937  n = nz-1
938  hgl = pnleg(z,nz)/(pndleg(zgl(i),nz)*(z-zgl(i)))
939  RETURN
940  end function hgl
941 
942  REAL(kind=xp) FUNCTION pnleg (Z,N)
943 !---------------------------------------------------------------------
944 !
945 ! Compute the value of the Nth order Legendre polynomial at Z.
946 ! (Simpler than JACOBF)
947 ! Based on the recursion formula for the Legendre polynomials.
948 !
949 !---------------------------------------------------------------------
950 !
951 ! This next statement is to overcome the underflow bug in the i860.
952 ! It can be removed at a later date. 11 Aug 1990 pff.
953 !
954  IMPLICIT REAL(KIND=XP) (a-h,o-z)
955  REAL(kind=xp) z, p1, p2, p3
956  IF(abs(z) .LT. 1.0e-25) z = 0.0
957 
958 
959  p1 = 1.
960  IF (n.EQ.0) THEN
961  pnleg = p1
962  RETURN
963  ENDIF
964  p2 = z
965  p3 = p2
966  DO 10 k = 1, n-1
967  fk = (k)
968  p3 = ((2.*fk+1.)*z*p2 - fk*p1)/(fk+1.)
969  p1 = p2
970  p2 = p3
971 10 CONTINUE
972  pnleg = p3
973  if (n.eq.0) pnleg = 1.
974  RETURN
975  end function pnleg
976 
979  subroutine legendre_poly(L, x, N)
980  real(kind=rp), intent(inout):: l(0:n)
981  real(kind=rp) :: x
982  integer :: N, j
983 
984  l(0) = 1.0_xp
985  if (n .eq. 0) return
986  l(1) = x
987 
988  do j=1, n-1
989  l(j+1) = ( (2.0_xp * real(j, xp) + 1.0_xp) * x * l(j) &
990  - real(j, xp) * l(j-1) ) / (real(j, xp) + 1.0_xp)
991  end do
992  end subroutine legendre_poly
993 
994  REAL(kind=xp) FUNCTION pndleg (Z,N)
995 !----------------------------------------------------------------------
996 !
997 ! Compute the derivative of the Nth order Legendre polynomial at Z.
998 ! (Simpler than JACOBF)
999 ! Based on the recursion formula for the Legendre polynomials.
1000 !
1001 !----------------------------------------------------------------------
1002  IMPLICIT REAL(KIND=XP) (a-h,o-z)
1003  REAL(kind=xp) p1, p2, p1d, p2d, p3d, z
1004  p1 = 1.
1005  p2 = z
1006  p1d = 0.
1007  p2d = 1.
1008  p3d = 1.
1009  DO 10 k = 1, n-1
1010  fk = (k)
1011  p3 = ((2.*fk+1.)*z*p2 - fk*p1)/(fk+1.)
1012  p3d = ((2.*fk+1.)*p2 + (2.*fk+1.)*z*p2d - fk*p1d)/(fk+1.)
1013  p1 = p2
1014  p2 = p3
1015  p1d = p2d
1016  p2d = p3d
1017 10 CONTINUE
1018  pndleg = p3d
1019  IF (n.eq.0) pndleg = 0.
1020  RETURN
1021  end function pndleg
1022 
1023  SUBROUTINE dgllgl (D,DT,ZM1,ZM2,IM12,NZM1,NZM2,ND1,ND2)
1024 !-----------------------------------------------------------------------
1025 !
1026 ! Compute the (one-dimensional) derivative matrix D and its
1027 ! transpose DT associated with taking the derivative of a variable
1028 ! expanded on a Gauss-Lobatto Legendre mesh (M1), and evaluate its
1029 ! derivative on a Guass Legendre mesh (M2).
1030 ! Need the one-dimensional interpolation operator IM12
1031 ! (see subroutine IGLLGL).
1032 ! Note: D and DT are rectangular matrices.
1033 !
1034 !-----------------------------------------------------------------------
1035  REAL(KIND=xp) d(nd2,nd1), dt(nd1,nd2), zm1(nd1), zm2(nd2), im12(nd2,nd1)
1036  REAL(KIND=xp) eps, zp, zq
1037  IF (nzm1.EQ.1) THEN
1038  d(1,1) = 0.
1039  dt(1,1) = 0.
1040  RETURN
1041  ENDIF
1042  eps = 1.e-6
1043  nm1 = nzm1-1
1044  DO ip = 1, nzm2
1045  DO jq = 1, nzm1
1046  zp = zm2(ip)
1047  zq = zm1(jq)
1048  IF ((abs(zp) .LT. eps).AND.(abs(zq) .LT. eps)) THEN
1049  d(ip,jq) = 0.
1050  ELSE
1051  d(ip,jq) = (pnleg(zp,nm1)/pnleg(zq,nm1) &
1052  -im12(ip,jq))/(zp-zq)
1053  ENDIF
1054  dt(jq,ip) = d(ip,jq)
1055  END DO
1056  END DO
1057  RETURN
1058  end subroutine dgllgl
1059 
1060  SUBROUTINE dgljgj (D,DT,ZGL,ZG,IGLG,NPGL,NPG,ND1,ND2,ALPHA,BETA)
1061 !-----------------------------------------------------------------------
1062 !
1063 ! Compute the (one-dimensional) derivative matrix D and its
1064 ! transpose DT associated with taking the derivative of a variable
1065 ! expanded on a Gauss-Lobatto Jacobi mesh (M1), and evaluate its
1066 ! derivative on a Guass Jacobi mesh (M2).
1067 ! Need the one-dimensional interpolation operator IM12
1068 ! (see subroutine IGLJGJ).
1069 ! Note: D and DT are rectangular matrices.
1070 ! Single precision version.
1071 !
1072 !-----------------------------------------------------------------------
1073  REAL(KIND=xp) d(nd2,nd1), dt(nd1,nd2), zgl(nd1), zg(nd2), iglg(nd2,nd1)
1074  parameter(nmax=84)
1075  parameter(ndd = nmax)
1076  REAL(KIND=xp) dd(ndd,ndd), dtd(ndd,ndd)
1077  REAL(KIND=xp) zgd(ndd), zgld(ndd), iglgd(ndd,ndd)
1078  REAL(KIND=xp) alphad, betad
1079 
1080  IF (npgl.LE.1) THEN
1081  WRITE(6,*) 'DGLJGJ: Minimum number of Gauss-Lobatto points is 2'
1082  call neko_error
1083  ENDIF
1084  IF (npgl.GT.nmax) THEN
1085  WRITE(6,*) 'Polynomial degree too high in DGLJGJ'
1086  WRITE(6,*) 'Maximum polynomial degree is',nmax
1087  WRITE(6,*) 'Here NPGL=',npgl
1088  call neko_error
1089  ENDIF
1090  IF ((alpha.LE.-1.).OR.(beta.LE.-1.)) THEN
1091  WRITE(6,*) 'DGLJGJ: Alpha and Beta must be greater than -1'
1092  call neko_error
1093  ENDIF
1094 
1095  alphad = alpha
1096  betad = beta
1097  DO i=1,npg
1098  zgd(i) = zg(i)
1099  DO j=1,npgl
1100  iglgd(i,j) = iglg(i,j)
1101  END DO
1102  END DO
1103  DO 200 i=1,npgl
1104  zgld(i) = zgl(i)
1105 200 CONTINUE
1106  CALL dgljgjd (dd,dtd,zgld,zgd,iglgd,npgl,npg,ndd,ndd,alphad,betad)
1107  DO i=1,npg
1108  DO j=1,npgl
1109  d(i,j) = dd(i,j)
1110  dt(j,i) = dtd(j,i)
1111  END DO
1112  END DO
1113  RETURN
1114  end subroutine dgljgj
1115 
1116  SUBROUTINE dgljgjd (D,DT,ZGL,ZG,IGLG,NPGL,NPG,ND1,ND2,ALPHA,BETA)
1117 !-----------------------------------------------------------------------
1118 !
1119 ! Compute the (one-dimensional) derivative matrix D and its
1120 ! transpose DT associated with taking the derivative of a variable
1121 ! expanded on a Gauss-Lobatto Jacobi mesh (M1), and evaluate its
1122 ! derivative on a Guass Jacobi mesh (M2).
1123 ! Need the one-dimensional interpolation operator IM12
1124 ! (see subroutine IGLJGJ).
1125 ! Note: D and DT are rectangular matrices.
1126 ! Double precision version.
1127 !
1128 !-----------------------------------------------------------------------
1129  IMPLICIT REAL(KIND=XP) (a-h,o-z)
1130  REAL(KIND=xp) d(nd2,nd1), dt(nd1,nd2), zgl(nd1), zg(nd2)
1131  REAL(KIND=xp) iglg(nd2,nd1), alpha, beta
1132 
1133  IF (npgl.LE.1) THEN
1134  WRITE(6,*) 'DGLJGJD: Minimum number of Gauss-Lobatto points is 2'
1135  call neko_error
1136  ENDIF
1137  IF ((alpha.LE.-1.).OR.(beta.LE.-1.)) THEN
1138  WRITE(6,*) 'DGLJGJD: Alpha and Beta must be greater than -1'
1139  call neko_error
1140  ENDIF
1141 
1142  eps = 1.e-6
1143  one = 1.
1144  two = 2.
1145  ngl = npgl-1
1146  dn = ((ngl))
1147  eigval = -dn*(dn+alpha+beta+one)
1148 
1149  DO i=1,npg
1150  DO j=1,npgl
1151  dz = abs(zg(i)-zgl(j))
1152  IF (dz.LT.eps) THEN
1153  d(i,j) = (alpha*(one+zg(i))-beta*(one-zg(i)))/ &
1154  (two*(one-zg(i)**2))
1155  ELSE
1156  CALL jacobf (pi,pdi,pm1,pdm1,pm2,pdm2,ngl,alpha,beta,zg(i))
1157  CALL jacobf (pj,pdj,pm1,pdm1,pm2,pdm2,ngl,alpha,beta,zgl(j))
1158  faci = alpha*(one+zg(i))-beta*(one-zg(i))
1159  facj = alpha*(one+zgl(j))-beta*(one-zgl(j))
1160  const = eigval*pj+facj*pdj
1161  d(i,j) = ((eigval*pi+faci*pdi)*(zg(i)-zgl(j)) &
1162  -(one-zg(i)**2)*pdi)/(const*(zg(i)-zgl(j))**2)
1163  ENDIF
1164  dt(j,i) = d(i,j)
1165  END DO
1166  END DO
1167  RETURN
1168  end subroutine dgljgjd
1169 
1170  SUBROUTINE iglm (I12,IT12,Z1,Z2,NZ1,NZ2,ND1,ND2)
1171 !----------------------------------------------------------------------
1172 !
1173 ! Compute the one-dimensional interpolation operator (matrix) I12
1174 ! ands its transpose IT12 for interpolating a variable from a
1175 ! Gauss Legendre mesh (1) to a another mesh M (2).
1176 ! Z1 : NZ1 Gauss Legendre points.
1177 ! Z2 : NZ2 points on mesh M.
1178 !
1179 !--------------------------------------------------------------------
1180  REAL(KIND=xp) i12(nd2,nd1),it12(nd1,nd2),z1(nd1),z2(nd2), zi
1181  IF (nz1 .EQ. 1) THEN
1182  i12(1,1) = 1.
1183  it12(1,1) = 1.
1184  RETURN
1185  ENDIF
1186  DO i=1,nz2
1187  zi = z2(i)
1188  DO j=1,nz1
1189  i12(i,j) = hgl(j,zi,z1,nz1)
1190  it12(j,i) = i12(i,j)
1191  END DO
1192  END DO
1193  RETURN
1194  end subroutine iglm
1195 
1196  SUBROUTINE igllm (I12,IT12,Z1,Z2,NZ1,NZ2,ND1,ND2)
1197 !----------------------------------------------------------------------
1198 !
1199 ! Compute the one-dimensional interpolation operator (matrix) I12
1200 ! ands its transpose IT12 for interpolating a variable from a
1201 ! Gauss-Lobatto Legendre mesh (1) to a another mesh M (2).
1202 ! Z1 : NZ1 Gauss-Lobatto Legendre points.
1203 ! Z2 : NZ2 points on mesh M.
1204 !
1205 !--------------------------------------------------------------------
1206  REAL(KIND=xp) i12(nd2,nd1),it12(nd1,nd2),z1(nd1),z2(nd2),zi
1207  IF (nz1 .EQ. 1) THEN
1208  i12(1,1) = 1.
1209  it12(1,1) = 1.
1210  RETURN
1211  ENDIF
1212  DO i=1,nz2
1213  zi = z2(i)
1214  DO j=1,nz1
1215  i12(i,j) = hgll(j,zi,z1,nz1)
1216  it12(j,i) = i12(i,j)
1217  END DO
1218  END DO
1219  RETURN
1220  end subroutine igllm
1221 
1222  SUBROUTINE igjm (I12,IT12,Z1,Z2,NZ1,NZ2,ND1,ND2,ALPHA,BETA)
1223 !----------------------------------------------------------------------
1224 !
1225 ! Compute the one-dimensional interpolation operator (matrix) I12
1226 ! ands its transpose IT12 for interpolating a variable from a
1227 ! Gauss Jacobi mesh (1) to a another mesh M (2).
1228 ! Z1 : NZ1 Gauss Jacobi points.
1229 ! Z2 : NZ2 points on mesh M.
1230 ! Single precision version.
1231 !
1232 !--------------------------------------------------------------------
1233  REAL(KIND=xp) i12(nd2,nd1),it12(nd1,nd2),z1(nd1),z2(nd2),zi,alpha,beta
1234  IF (nz1 .EQ. 1) THEN
1235  i12(1,1) = 1.
1236  it12(1,1) = 1.
1237  RETURN
1238  ENDIF
1239  DO i=1,nz2
1240  zi = z2(i)
1241  DO j=1,nz1
1242  i12(i,j) = hgj(j,zi,z1,nz1,alpha,beta)
1243  it12(j,i) = i12(i,j)
1244  END DO
1245  END DO
1246  RETURN
1247  end subroutine igjm
1248 
1249  SUBROUTINE igljm (I12,IT12,Z1,Z2,NZ1,NZ2,ND1,ND2,ALPHA,BETA)
1250 !----------------------------------------------------------------------
1251 !
1252 ! Compute the one-dimensional interpolation operator (matrix) I12
1253 ! ands its transpose IT12 for interpolating a variable from a
1254 ! Gauss-Lobatto Jacobi mesh (1) to a another mesh M (2).
1255 ! Z1 : NZ1 Gauss-Lobatto Jacobi points.
1256 ! Z2 : NZ2 points on mesh M.
1257 ! Single precision version.
1258 !
1259 !--------------------------------------------------------------------
1260  REAL(KIND=xp) i12(nd2,nd1),it12(nd1,nd2),z1(nd1),z2(nd2),zi,alpha,beta
1261  IF (nz1 .EQ. 1) THEN
1262  i12(1,1) = 1.
1263  it12(1,1) = 1.
1264  RETURN
1265  ENDIF
1266  DO i=1,nz2
1267  zi = z2(i)
1268  DO j=1,nz1
1269  i12(i,j) = hglj(j,zi,z1,nz1,alpha,beta)
1270  it12(j,i) = i12(i,j)
1271  END DO
1272  END DO
1273  RETURN
1274  end subroutine igljm
1275 end module speclib
double real
Definition: device_config.h:12
integer, parameter, public xp
Definition: num_types.f90:14
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
LIBRARY ROUTINES FOR SPECTRAL METHODS.
Definition: speclib.f90:148
real(kind=xp) function hgljd(I, Z, ZGLJ, NP, ALPHA, BETA)
Definition: speclib.f90:663
real(kind=xp) function hgl(I, Z, ZGL, NZ)
Definition: speclib.f90:924
subroutine zwgll(Z, W, NP)
Definition: speclib.f90:169
real(kind=xp) function hgjd(II, Z, ZGJ, NP, ALPHA, BETA)
Definition: speclib.f90:610
real(kind=xp) function endw2(N, ALPHA, BETA)
Definition: speclib.f90:389
subroutine dgljgj(D, DT, ZGL, ZG, IGLG, NPGL, NPG, ND1, ND2, ALPHA, BETA)
Definition: speclib.f90:1061
real(kind=xp) function pndleg(Z, N)
Definition: speclib.f90:995
subroutine zwgj(Z, W, NP, ALPHA, BETA)
Definition: speclib.f90:186
subroutine legendre_poly(L, x, N)
Evaluate Legendre polynomials of degrees 0-N at point x and store in array L.
Definition: speclib.f90:980
subroutine dgj(D, DT, Z, NZ, NZD, ALPHA, BETA)
Definition: speclib.f90:691
subroutine jacg(XJAC, NP, ALPHA, BETA)
Definition: speclib.f90:483
subroutine zwgjd(Z, W, NP, ALPHA, BETA)
Definition: speclib.f90:217
subroutine dgll(D, DT, Z, NZ, NZD)
Definition: speclib.f90:865
subroutine zwglj(Z, W, NP, ALPHA, BETA)
Definition: speclib.f90:269
subroutine dglj(D, DT, Z, NZ, NZD, ALPHA, BETA)
Definition: speclib.f90:774
subroutine igjm(I12, IT12, Z1, Z2, NZ1, NZ2, ND1, ND2, ALPHA, BETA)
Definition: speclib.f90:1223
real(kind=xp) function hglj(II, Z, ZGLJ, NP, ALPHA, BETA)
Definition: speclib.f90:634
subroutine dgljgjd(D, DT, ZGL, ZG, IGLG, NPGL, NPG, ND1, ND2, ALPHA, BETA)
Definition: speclib.f90:1117
subroutine igljm(I12, IT12, Z1, Z2, NZ1, NZ2, ND1, ND2, ALPHA, BETA)
Definition: speclib.f90:1250
subroutine iglm(I12, IT12, Z1, Z2, NZ1, NZ2, ND1, ND2)
Definition: speclib.f90:1171
subroutine zwgl(Z, W, NP)
Generate NP Gauss Legendre points Z and weights W associated with Jacobi polynomial ....
Definition: speclib.f90:161
subroutine dgljd(D, DT, Z, NZ, NZD, ALPHA, BETA)
Definition: speclib.f90:818
real(kind=xp) function endw1(N, ALPHA, BETA)
Definition: speclib.f90:345
real(kind=xp) function hgll(I, Z, ZGLL, NZ)
Definition: speclib.f90:902
real(kind=xp) function hgj(II, Z, ZGJ, NP, ALPHA, BETA)
Definition: speclib.f90:581
subroutine zwgljd(Z, W, NP, ALPHA, BETA)
Definition: speclib.f90:300
subroutine dgjd(D, DT, Z, NZ, NZD, ALPHA, BETA)
Definition: speclib.f90:735
subroutine jacobf(POLY, PDER, POLYM1, PDERM1, POLYM2, PDERM2, N, ALP, BET, X)
Definition: speclib.f90:541
real(kind=xp) function gammaf(X)
Definition: speclib.f90:433
subroutine dgllgl(D, DT, ZM1, ZM2, IM12, NZM1, NZM2, ND1, ND2)
Definition: speclib.f90:1024
subroutine igllm(I12, IT12, Z1, Z2, NZ1, NZ2, ND1, ND2)
Definition: speclib.f90:1197
real(kind=xp) function pnleg(Z, N)
Definition: speclib.f90:943
real(kind=xp) function pnormj(N, ALPHA, BETA)
Definition: speclib.f90:457
Utilities.
Definition: utils.f90:35