KTH framework for Nek5000 toolboxes; testing version  0.0.1
math.f
Go to the documentation of this file.
1 c-----------------------------------------------------------------------
2  subroutine icopy48(a,b,n)
3  integer*8 a(1)
4  integer*4 b(1)
5  do 100 i = 1, n
6  100 a(i) = b(i)
7  return
8  end
9 c-----------------------------------------------------------------------
10  subroutine icopy84(a,b,n)
11  integer*4 a(1)
12  integer*8 b(1)
13  do 100 i = 1, n
14  100 a(i) = b(i)
15  return
16  end
17 c-----------------------------------------------------------------------
18  SUBROUTINE blank(A,N)
19  CHARACTER*1 A(1)
20  CHARACTER*1 BLNK
21  SAVE blnk
22  DATA blnk /' '/
23 C
24  DO 10 i=1,n
25  a(i)=blnk
26  10 CONTINUE
27  RETURN
28  END
29 c-----------------------------------------------------------------------
30  SUBROUTINE vsq (A,N)
31  dimension a(1)
32 C
33  include 'OPCTR'
34 C
35  DO 100 i = 1, n
36  100 a(i) = a(i)**2
37  RETURN
38  END
39 c-----------------------------------------------------------------------
40  SUBROUTINE vsqrt(A,N)
41  dimension a(1)
42 C
43  include 'OPCTR'
44 C
45  DO 100 i = 1, n
46  100 a(i) = sqrt(a(i))
47  RETURN
48  END
49 c-----------------------------------------------------------------------
50  subroutine invers2(a,b,n)
51  REAL A(1),B(1)
52 C
53  include 'OPCTR'
54 C
55  DO 100 i=1,n
56  a(i)=1./b(i)
57  100 CONTINUE
58  return
59  END
60 c-----------------------------------------------------------------------
61  subroutine invcol1(a,n)
62  REAL A(1)
63 C
64  include 'OPCTR'
65 C
66  DO 100 i=1,n
67  a(i)=1./a(i)
68  100 CONTINUE
69  return
70  END
71 c-----------------------------------------------------------------------
72  subroutine invcol2(a,b,n)
73 C
74  REAL A(1),B(1)
75  include 'CTIMER'
76  include 'OPCTR'
77 
78 #ifdef TIMER2
79  if (icalld.eq.0) then
80  tinvc=0.0
81  icalld=icalld+1
82  endif
83  etime1=dnekclock()
84 #endif
85 
86  DO 100 i=1,n
87  a(i)=a(i)/b(i)
88  100 CONTINUE
89 
90 #ifdef TIMER2
91  tinvc=tinvc+(dnekclock()-etime1)
92 #endif
93 
94  return
95  END
96 c-----------------------------------------------------------------------
97  subroutine invcol3(a,b,c,n)
98  REAL A(1),B(1),C(1)
99 C
100  include 'OPCTR'
101  include 'CTIMER'
102 
103 #ifdef TIMER2
104  if (icalld.eq.0) tinv3=0.0
105  icalld=icalld+1
106  ninv3=icalld
107  etime1=dnekclock()
108 #endif
109 C
110  DO 100 i=1,n
111  a(i)=b(i)/c(i)
112  100 CONTINUE
113 #ifdef TIMER2
114  tinv3=tinv3+(dnekclock()-etime1)
115 #endif
116  return
117  END
118 c-----------------------------------------------------------------------
119  subroutine col4(a,b,c,d,n)
120  REAL A(1),B(1),C(1),D(1)
121 C
122  include 'OPCTR'
123 C
124  DO 100 i=1,n
125  a(i)=b(i)*c(i)*d(i)
126  100 CONTINUE
127  return
128  END
129 c-----------------------------------------------------------------------
130  subroutine xaddcol3(a,b,c,n)
131  REAL A(1),B(1),C(1)
132 C
133  include 'OPCTR'
134 C
135  DO 100 i=1,n
136  a(i)=a(i)+b(i)*c(i)
137  100 CONTINUE
138  return
139  END
140 c-----------------------------------------------------------------------
141  subroutine addcol4(a,b,c,d,n)
142  REAL A(1),B(1),C(1),D(1)
143 C
144  include 'OPCTR'
145 C
146  DO 100 i=1,n
147  a(i)=a(i)+b(i)*c(i)*d(i)
148  100 CONTINUE
149  return
150  END
151 c-----------------------------------------------------------------------
152  subroutine ascol5 (a,b,c,d,e,n)
153  REAL A(1),B(1),C(1),D(1),E(1)
154 C
155  include 'OPCTR'
156 C
157  DO 100 i=1,n
158  a(i) = b(i)*c(i)-d(i)*e(i)
159  100 CONTINUE
160  return
161  END
162 c-----------------------------------------------------------------------
163  subroutine sub2(a,b,n)
164  REAL A(1),B(1)
165 C
166  include 'OPCTR'
167 C
168  DO 100 i=1,n
169  a(i)=a(i)-b(i)
170  100 CONTINUE
171  return
172  END
173 c-----------------------------------------------------------------------
174  subroutine sub3(a,b,c,n)
175  REAL A(1),B(1),C(1)
176 C
177  include 'OPCTR'
178 C
179  DO 100 i=1,n
180  a(i)=b(i)-c(i)
181  100 CONTINUE
182  return
183  END
184 c-----------------------------------------------------------------------
185  subroutine subcol3(a,b,c,n)
186  REAL A(1),B(1),C(1)
187 C
188  include 'OPCTR'
189 C
190  DO 100 i=1,n
191  a(i)=a(i)-b(i)*c(i)
192  100 CONTINUE
193  return
194  END
195 c-----------------------------------------------------------------------
196  subroutine subcol4(a,b,c,d,n)
197  REAL A(1),B(1),C(1),D(1)
198 C
199  include 'OPCTR'
200 C
201  DO 100 i=1,n
202  a(i)=a(i)-b(i)*c(i)*d(i)
203  100 CONTINUE
204  return
205  END
206 c-----------------------------------------------------------------------
207  subroutine rzero(a,n)
208  dimension a(1)
209  DO 100 i = 1, n
210  100 a(i ) = 0.0
211  return
212  END
213 c-----------------------------------------------------------------------
214  subroutine izero(a,n)
215  INTEGER A(1)
216 C
217  DO 100 i = 1, n
218  100 a(i ) = 0
219  return
220  END
221 c-----------------------------------------------------------------------
222  subroutine ione(a,n)
223  INTEGER A(1)
224  DO 100 i = 1, n
225  100 a(i ) = 1
226  return
227  END
228 c-----------------------------------------------------------------------
229  subroutine rone(a,n)
230  dimension a(1)
231  DO 100 i = 1, n
232  100 a(i ) = 1.0
233  return
234  END
235 c-----------------------------------------------------------------------
236  subroutine ltrue(a,n)
237  LOGICAL A(1)
238  DO 100 i = 1, n
239  100 a(i ) = .true.
240  return
241  END
242 c-----------------------------------------------------------------------
243  subroutine cfill(a,b,n)
244  dimension a(1)
245 C
246  DO 100 i = 1, n
247  100 a(i) = b
248  return
249  END
250 c-----------------------------------------------------------------------
251  subroutine ifill(ia,ib,n)
252  dimension ia(1)
253 C
254  DO 100 i = 1, n
255  100 ia(i) = ib
256  return
257  END
258 c-----------------------------------------------------------------------
259  subroutine copy(a,b,n)
260  real a(1),b(1)
261 
262  do i=1,n
263  a(i)=b(i)
264  enddo
265 
266  return
267  end
268 c-----------------------------------------------------------------------
269  subroutine copyi4(a,b,n)
270  integer a(1)
271  real b(1)
272 
273  do i=1,n
274  a(i)=b(i)
275  enddo
276 
277  return
278  end
279 c-----------------------------------------------------------------------
280  subroutine chcopy(a,b,n)
281  CHARACTER*1 A(1), B(1)
282 C
283  DO 100 i = 1, n
284  100 a(i) = b(i)
285  return
286  END
287 C
288  subroutine icopy(a,b,n)
289  INTEGER A(1), B(1)
290 C
291  DO 100 i = 1, n
292  100 a(i) = b(i)
293  return
294  END
295 c-----------------------------------------------------------------------
296  subroutine i8copy(a,b,n)
297  INTEGER*8 A(1), B(1)
298 C
299  DO 100 i = 1, n
300  100 a(i) = b(i)
301  return
302  END
303 c-----------------------------------------------------------------------
304  subroutine chsign(a,n)
305  REAL A(1)
306 C
307  DO 100 i=1,n
308  a(i) = -a(i)
309  100 CONTINUE
310  return
311  END
312 C
313 c-----------------------------------------------------------------------
314  subroutine cmult(a,const,n)
315  REAL A(1)
316 C
317  include 'OPCTR'
318 C
319  DO 100 i=1,n
320  a(i)=a(i)*const
321  100 CONTINUE
322  return
323  END
324 c-----------------------------------------------------------------------
325  subroutine cadd(a,const,n)
326  REAL A(1)
327 C
328  include 'OPCTR'
329 C
330  DO 100 i=1,n
331  a(i)=a(i)+const
332  100 CONTINUE
333  return
334  END
335 c-----------------------------------------------------------------------
336  subroutine iadd(i1,iscal,n)
337  dimension i1(1)
338 C
339  DO 10 i=1,n
340  i1(i)=i1(i)+iscal
341  10 CONTINUE
342  return
343  END
344 c-----------------------------------------------------------------------
345  subroutine cadd2(a,b,const,n)
346  REAL A(1),B(1)
347 C
348  include 'OPCTR'
349 C
350  DO 100 i=1,n
351  a(i)=b(i)+const
352  100 CONTINUE
353  return
354  END
355 c-----------------------------------------------------------------------
356  real function vlmin(vec,n)
357  REAL vec(1)
358  tmin = 99.0e20
359 C
360  DO 100 i=1,n
361  tmin = min(tmin,vec(i))
362  100 CONTINUE
363  vlmin = tmin
364  return
365  END
366 c-----------------------------------------------------------------------
367  integer function ivlmin(vec,n)
368  integer vec(1),tmin
369  if (n.eq.0) then
370  ivlmin=0
371  return
372  endif
373  tmin = 8888888
374  do i=1,n
375  tmin = min(tmin,vec(i))
376  enddo
377  ivlmin = tmin
378  return
379  end
380 c-----------------------------------------------------------------------
381  integer function ivlmax(vec,n)
382  integer vec(1),tmax
383  if (n.eq.0) then
384  ivlmax=0
385  return
386  endif
387  tmax =-8888888
388  do i=1,n
389  tmax = max(tmax,vec(i))
390  enddo
391  ivlmax = tmax
392  return
393  end
394 c-----------------------------------------------------------------------
395  real function vlmax(vec,n)
396  REAL vec(1)
397  tmax =-99.0e20
398  do i=1,n
399  tmax = max(tmax,vec(i))
400  enddo
401  vlmax = tmax
402  return
403  END
404 c-----------------------------------------------------------------------
405  real function vlamax(vec,n)
406  REAL vec(1)
407  tamax = 0.0
408 C
409  DO 100 i=1,n
410  tamax = max(tamax,abs(vec(i)))
411  100 CONTINUE
412  vlamax = tamax
413  return
414  END
415 c-----------------------------------------------------------------------
416  real function vlsum(vec,n)
417  REAL vec(1)
418  include 'OPCTR'
419 C
420  sum = 0.
421 C
422  DO 100 i=1,n
423  sum=sum+vec(i)
424  100 CONTINUE
425  vlsum = sum
426  return
427  END
428 c-----------------------------------------------------------------------
429  subroutine vcross (u1,u2,u3,v1,v2,v3,w1,w2,w3,n)
430 C
431 C Compute a Cartesian vector cross product.
432 C
433  dimension u1(1),u2(1),u3(1)
434  dimension v1(1),v2(1),v3(1)
435  dimension w1(1),w2(1),w3(1)
436 C
437 C
438  DO 100 i=1,n
439  u1(i) = v2(i)*w3(i) - v3(i)*w2(i)
440  u2(i) = v3(i)*w1(i) - v1(i)*w3(i)
441  u3(i) = v1(i)*w2(i) - v2(i)*w1(i)
442  100 CONTINUE
443  return
444  END
445 c-----------------------------------------------------------------------
446  subroutine vdot2 (dot,u1,u2,v1,v2,n)
447 C
448 C Compute a Cartesian vector dot product. 2-d version
449 C
450  dimension dot(1)
451  dimension u1(1),u2(1)
452  dimension v1(1),v2(1)
453 C
454 C
455  DO 100 i=1,n
456  dot(i) = u1(i)*v1(i) + u2(i)*v2(i)
457  100 CONTINUE
458  return
459  END
460 c-----------------------------------------------------------------------
461  subroutine vdot3 (dot,u1,u2,u3,v1,v2,v3,n)
462 C
463 C Compute a Cartesian vector dot product. 3-d version
464 C
465  dimension dot(1)
466  dimension u1(1),u2(1),u3(1)
467  dimension v1(1),v2(1),v3(1)
468 C
469 C
470  DO 100 i=1,n
471  dot(i) = u1(i)*v1(i) + u2(i)*v2(i) + u3(i)*v3(i)
472  100 CONTINUE
473  return
474  END
475 c-----------------------------------------------------------------------
476  subroutine addtnsr(s,h1,h2,h3,nx,ny,nz)
477 C
478 C Map and add to S a tensor product form of the three functions H1,H2,H3.
479 C This is a single element routine used for deforming geometry.
480 C
481  dimension h1(1),h2(1),h3(1)
482  dimension s(nx,ny,nz)
483 C
484  DO 200 iz=1,nz
485  DO 200 iy=1,ny
486  hh = h2(iy)*h3(iz)
487  DO 100 ix=1,nx
488  s(ix,iy,iz)=s(ix,iy,iz)+hh*h1(ix)
489  100 CONTINUE
490  200 CONTINUE
491  return
492  END
493  function ltrunc(string,l)
494  CHARACTER*1 string(l)
495  CHARACTER*1 blnk
496  DATA blnk/' '/
497 C
498  DO 100 i=l,1,-1
499  l1=i
500  IF (string(i).NE.blnk) GOTO 200
501  100 CONTINUE
502  l1=0
503  200 CONTINUE
504  ltrunc=l1
505  return
506  END
507 c-----------------------------------------------------------------------
508  function mod1(i,n)
509 C
510 C Yields MOD(I,N) with the exception that if I=K*N, result is N.
511 C
512  mod1=0
513  IF (i.EQ.0) THEN
514  return
515  ENDIF
516  IF (n.EQ.0) THEN
517  WRITE(6,*)
518  $ 'WARNING: Attempt to take MOD(I,0) in function mod1.'
519  return
520  ENDIF
521  ii = i+n-1
522  mod1 = mod(ii,n)+1
523  return
524  END
525 c-----------------------------------------------------------------------
526  integer function log2(k)
527  rk=(k)
528  rlog=log10(rk)
529  rlog2=log10(2.0)
530  rlog=rlog/rlog2+0.5
531  log2=int(rlog)
532  return
533  END
534 c-----------------------------------------------------------------------
535  subroutine iflip(i1,n)
536  dimension i1(1)
537  n1=n+1
538  n2=n/2
539  DO 10 i=1,n2
540  ilast=n1-i
541  itmp=i1(ilast)
542  i1(ilast)=i1(i)
543  i1(i)=itmp
544  10 CONTINUE
545  return
546  END
547 c-----------------------------------------------------------------------
548  subroutine iswap(b,ind,n,temp)
549  INTEGER B(1),IND(1),TEMP(1)
550 C***
551 C*** SORT ASSOCIATED ELEMENTS BY PUTTING ITEM(JJ)
552 C*** INTO ITEM(I), WHERE JJ=IND(I).
553 C***
554  DO 20 i=1,n
555  jj=ind(i)
556  temp(i)=b(jj)
557  20 CONTINUE
558  DO 30 i=1,n
559  30 b(i)=temp(i)
560  return
561  END
562 c-----------------------------------------------------------------------
563  subroutine col2(a,b,n)
564  real a(1),b(1)
565  include 'OPCTR'
566  include 'CTIMER'
567 
568 #ifdef TIMER2
569  if (icalld.eq.0) then
570  icalld=1
571  tcol2=0
572  endif
573  etime1=dnekclock()
574 #endif
575 
576  do i=1,n
577  a(i)=a(i)*b(i)
578  enddo
579 
580 #ifdef TIMER2
581  tcol2=tcol2+(dnekclock()-etime1)
582 #endif
583 
584  return
585  end
586 c-----------------------------------------------------------------------
587  subroutine col2c(a,b,c,n)
588  real a(1),b(1),c
589 
590  do i=1,n
591  a(i)=a(i)*b(i)*c
592  enddo
593 
594  return
595  end
596 c-----------------------------------------------------------------------
597  subroutine col3(a,b,c,n)
598  real a(1),b(1),c(1)
599  include 'OPCTR'
600  include 'CTIMER'
601 
602 #ifdef TIMER2
603  if (icalld.eq.0) then
604  icalld=1
605  tcol3=0
606  endif
607  etime1=dnekclock()
608 #endif
609 
610  do i=1,n
611  a(i)=b(i)*c(i)
612  enddo
613 
614 #ifdef TIMER2
615  tcol3=tcol3+(dnekclock()-etime1)
616 #endif
617 
618  return
619  end
620 c-----------------------------------------------------------------------
621  subroutine add2(a,b,n)
622  real a(1),b(1)
623  include 'OPCTR'
624  include 'CTIMER'
625 
626 #ifdef TIMER2
627  if (icalld.eq.0) then
628  icalld=1
629  tadd2=0
630  endif
631  etime1=dnekclock()
632 #endif
633  do i=1,n
634  a(i)=a(i)+b(i)
635  enddo
636 
637 #ifdef TIMER2
638  tadd2=tadd2+(dnekclock()-etime1)
639 #endif
640  return
641  end
642 c-----------------------------------------------------------------------
643  subroutine add3(a,b,c,n)
644  real a(1),b(1),c(1)
645  include 'OPCTR'
646 
647  do i=1,n
648  a(i)=b(i)+c(i)
649  enddo
650  return
651  end
652 c-----------------------------------------------------------------------
653  subroutine addcol3(a,b,c,n)
654  real a(1),b(1),c(1)
655  include 'OPCTR'
656  include 'CTIMER'
657 
658 #ifdef TIMER2
659  if (icalld.eq.0) then
660  icalld=1
661  tadc3=0
662  endif
663  etime1=dnekclock()
664 #endif
665 
666  do i=1,n
667  a(i)=a(i)+b(i)*c(i)
668  enddo
669 
670 #ifdef TIMER2
671  tadc3=tadc3+(dnekclock()-etime1)
672 #endif
673 
674  return
675  end
676 c-----------------------------------------------------------------------
677  subroutine add2s1(a,b,c1,n)
678  real a(1),b(1)
679 C
680  include 'OPCTR'
681 C
682  DO 100 i=1,n
683  a(i)=c1*a(i)+b(i)
684  100 CONTINUE
685  return
686  END
687 C
688 c-----------------------------------------------------------------------
689  subroutine add2s2(a,b,c1,n)
690  real a(1),b(1)
691 C
692  include 'OPCTR'
693  include 'CTIMER'
694 C
695 #ifdef TIMER2
696  if (icalld.eq.0) then
697  icalld=1
698  ta2s2=0
699  endif
700  etime1=dnekclock()
701 #endif
702 C
703  DO 100 i=1,n
704  a(i)=a(i)+c1*b(i)
705  100 CONTINUE
706 
707 #ifdef TIMER2
708  ta2s2=ta2s2+(dnekclock()-etime1)
709 #endif
710 
711  return
712  END
713 C
714 c-----------------------------------------------------------------------
715  subroutine add3s2(a,b,c,c1,c2,n)
716  real a(1),b(1),c(1)
717 C
718  include 'OPCTR'
719 C
720  DO 100 i=1,n
721  a(i)=c1*b(i)+c2*c(i)
722  100 CONTINUE
723  return
724  END
725 C
726 c-----------------------------------------------------------------------
727  subroutine add4(a,b,c,d,n)
728  REAL A(1),B(1),C(1),D(1)
729 C
730  include 'OPCTR'
731 C
732  DO 100 i=1,n
733  a(i)=b(i)+c(i)+d(i)
734  100 CONTINUE
735  return
736  END
737 c-----------------------------------------------------------------------
738  real function vlsc2(x,y,n)
739  REAL x(1),y(1)
740  include 'SIZE'
741  include 'OPCTR'
742  include 'PARALLEL'
743 C
744  s = 0.
745  do i=1,n
746  s = s + x(i)*y(i)
747  enddo
748  vlsc2=s
749  return
750  end
751 c-----------------------------------------------------------------------
752  real function vlsc21(x,y,n)
753  real x(1),y(1)
754  include 'SIZE'
755  include 'OPCTR'
756  include 'PARALLEL'
757 C
758  s = 0.
759  do i=1,n
760  s = s + x(i)*x(i)*y(i)
761  enddo
762  vlsc21=s
763  return
764  end
765 
766 C----------------------------------------------------------------------------
767 C
768 C Vector reduction routines which require communication
769 C on a parallel machine. These routines must be substituted with
770 C appropriate routines which take into account the specific architecture.
771 C
772 C----------------------------------------------------------------------------
773 
774 
775  function glsc3(a,b,mult,n)
776 C
777 C Perform inner-product in double precision
778 C
779  REAL a(1),b(1),mult(1)
780  REAL tmp,work(1)
781 C
782  include 'OPCTR'
783 C
784  tmp = 0.0
785  DO 10 i=1,n
786  tmp = tmp + a(i)*b(i)*mult(i)
787  10 CONTINUE
788  CALL gop(tmp,work,'+ ',1)
789  glsc3 = tmp
790  return
791  END
792 c-----------------------------------------------------------------------
793  function glsc2(x,y,n)
794 C
795 C Perform inner-product in double precision
796 C
797  include 'OPCTR'
798 c
799  real x(1), y(1)
800  real tmp,work(1)
801 C
802  tmp=0.0
803  do 10 i=1,n
804  tmp = tmp+ x(i)*y(i)
805  10 continue
806  CALL gop(tmp,work,'+ ',1)
807  glsc2 = tmp
808  return
809  END
810 c-----------------------------------------------------------------------
811  function glsc23(x,y,z,n)
812 c
813 C Perform inner-product x*x*y*z
814 c
815  real x(1), y(1),z(1)
816  real tmp,work(1)
817 
818  ds = 0.0
819  do 10 i=1,n
820  ds=ds+x(i)*x(i)*y(i)*z(i)
821  10 continue
822  tmp=ds
823  call gop(tmp,work,'+ ',1)
824  glsc23 = tmp
825  return
826  end
827 c-----------------------------------------------------------------------
828  real function gl2norm(a,n)
829 
830  include 'SIZE'
831  include 'MASS'
832 
833  real a(1)
834 
835  common /scrsf/ w1(lx1,ly1,lz1,lelt)
836 
837  call col3 (w1,a,a,n)
838  call col2 (w1,bm1,n)
839  gl2norm = sqrt(glsum(w1,n)/volvm1)
840 
841  return
842  end
843 c-----------------------------------------------------------------------
844  real function gl2norm2(a,n)
845 
846  include 'SIZE'
847  include 'MASS'
848 
849  real a(n)
850 
851  common /scrsf/ w1(lx2*ly2*lz2*lelt)
852 
853  call col3 (w1,a,a,n)
854  call col2 (w1,bm2,n)
855  gl2norm2 = sqrt(glsum(w1,n)/volvm2)
856 
857  return
858  end
859 c-----------------------------------------------------------------------
860  function glsum (x,n)
861  dimension x(1)
862  dimension tmp(1),work(1)
863  tsum = 0.
864  DO 100 i=1,n
865  tsum = tsum+x(i)
866  100 CONTINUE
867  tmp(1)=tsum
868  CALL gop(tmp,work,'+ ',1)
869  glsum = tmp(1)
870  return
871  END
872 c-----------------------------------------------------------------------
873  real function glamax(a,n)
874  REAL a(1)
875  dimension tmp(1),work(1)
876  tmax = 0.0
877  DO 100 i=1,n
878  tmax = max(tmax,abs(a(i)))
879  100 CONTINUE
880  tmp(1)=tmax
881  CALL gop(tmp,work,'M ',1)
882  glamax=abs(tmp(1))
883  return
884  END
885 c-----------------------------------------------------------------------
886  real function glamin(a,n)
887  real a(1)
888  dimension tmp(1),work(1)
889  tmin = 9.e28
890  do 100 i=1,n
891  tmin = min(tmin,abs(a(i)))
892  100 continue
893  tmp(1)=tmin
894  call gop(tmp,work,'m ',1)
895  glamin=abs(tmp(1))
896  return
897  end
898 c-----------------------------------------------------------------------
899  function iglmin(a,n)
900  integer a(1),tmin
901  integer tmp(1),work(1)
902  tmin= 999999999
903  do i=1,n
904  tmin=min(tmin,a(i))
905  enddo
906  tmp(1)=tmin
907  call igop(tmp,work,'m ',1)
908  iglmin=tmp(1)
909  return
910  end
911 c-----------------------------------------------------------------------
912  function iglmax(a,n)
913  integer a(1),tmax
914  integer tmp(1),work(1)
915  tmax= -999999999
916  do i=1,n
917  tmax=max(tmax,a(i))
918  enddo
919  tmp(1)=tmax
920  call igop(tmp,work,'M ',1)
921  iglmax=tmp(1)
922  return
923  end
924 c-----------------------------------------------------------------------
925  function iglsum(a,n)
926  integer a(1),tsum
927  integer tmp(1),work(1)
928  tsum= 0
929  do i=1,n
930  tsum=tsum+a(i)
931  enddo
932  tmp(1)=tsum
933  call igop(tmp,work,'+ ',1)
934  iglsum=tmp(1)
935  return
936  end
937 c-----------------------------------------------------------------------
938  subroutine i8zero(a,n)
939  integer*8 a(1)
940  do i=1,n
941  a(i)=0
942  enddo
943  return
944  end
945 c-----------------------------------------------------------------------
946  integer*8 function i8glsum(a,n)
947  integer*8 a(1),tsum
948  integer*8 tmp(1),work(1)
949  tsum= 0
950  do i=1,n
951  tsum=tsum+a(i)
952  enddo
953  tmp(1)=tsum
954  call i8gop(tmp,work,'+ ',1)
955  i8glsum=tmp(1)
956  return
957  end
958 C-----------------------------------------------------------------------
959  function glmax(a,n)
960  REAL a(1)
961  dimension tmp(1),work(1)
962  tmax=-99.0e20
963  DO 100 i=1,n
964  tmax=max(tmax,a(i))
965  100 CONTINUE
966  tmp(1)=tmax
967  CALL gop(tmp,work,'M ',1)
968  glmax=tmp(1)
969  return
970  END
971 c-----------------------------------------------------------------------
972  function glmin(a,n)
973  REAL a(1)
974  dimension tmp(1),work(1)
975  tmin=99.0e20
976  DO 100 i=1,n
977  tmin=min(tmin,a(i))
978  100 CONTINUE
979  tmp(1)=tmin
980  CALL gop(tmp,work,'m ',1)
981  glmin = tmp(1)
982  return
983  END
984 c-----------------------------------------------------------------------
985  subroutine gllog(la,lb)
986 C
987 C If ANY LA=LB, then ALL LA=LB.
988 C
989  LOGICAL LA,LB
990  dimension tmp(1),work(1)
991 C
992  tmp(1)=1.0
993  IF (lb) THEN
994  IF (la) tmp(1)=0.0
995  ELSE
996  IF (.NOT.la) tmp(1)=0.0
997  ENDIF
998  CALL gop(tmp,work,'* ',1)
999  IF (tmp(1).EQ.0.0) la=lb
1000  return
1001  END
1002 c-----------------------------------------------------------------------
1003  function fmdian(a,n,ifok)
1004 C find the Median of the (global) set A
1005  include 'SIZE'
1006  dimension a(1)
1007  dimension work1(5),work2(5)
1008  dimension gues(100)
1009  LOGICAL ifok
1010 C
1011  amp =1.5
1012  afac =1.5
1013  gmin =glmin(a,n)
1014  gmax =glmax(a,n)
1015  gmin0=glmin(a,n)
1016  gmax0=glmax(a,n)
1017  guess=(gmax+gmin)/2.0
1018  eps =(gmax-gmin)
1019  IF (eps.EQ.0.0) THEN
1020  fmdian=gmax
1021  return
1022  ENDIF
1023  work1(1)=n
1024  CALL gop(work1,work2,'+ ',1)
1025  ntot=work1(1)
1026  n2 = (ntot+1)/2
1027  IF (.NOT.ifok) THEN
1028  WRITE(6,8) nid,n,(a(i),i=1,n)
1029  WRITE(6,9) nid,ntot,n2,n,gmin,gmax
1030  8 FORMAT(i5,'N,A:',i5,10(6f10.5,/))
1031  9 FORMAT(i5,'mnx:',3i6,2f10.5)
1032  ENDIF
1033 C
1034 C This is the trial loop
1035 C
1036  itry=-1
1037  10 CONTINUE
1038  itry=itry+1
1039  ii=itry+1
1040  IF (ii.LE.100) gues(ii)=guess
1041 C error check for infinite loop
1042  IF (itry.GT.2*ntot) GOTO 9000
1043  CALL rzero(work1,5)
1044  nlt=0
1045  ngt=0
1046  clt=gmin0
1047  cgt=gmax0
1048  DO 100 i=1,n
1049  aa=a(i)
1050  IF (aa.NE.guess) THEN
1051  IF (aa.LT.guess) THEN
1052  nlt=nlt+1
1053 C CLT - closest value to GUESS Less Than GUESS
1054  IF (aa.GT.clt) clt=aa
1055  ENDIF
1056  IF (aa.GT.guess) THEN
1057  ngt=ngt+1
1058 C CGT - closest value to GUESS Greater Than GUESS
1059  IF (aa.LT.cgt) cgt=aa
1060  ENDIF
1061  dum=1./(eps+abs(aa-guess))
1062  work1(1)=work1(1)+dum
1063  work1(2)=work1(2)+dum*aa
1064  ELSE
1065 C detected values equaling the guess.
1066  work1(5)=work1(5)+1.0
1067  ENDIF
1068  100 CONTINUE
1069 C Invoke vector reduction across processors:
1070  work2(1)=clt
1071  clt=glmax(work2,1)
1072  work2(1)=cgt
1073  cgt=glmin(work2,1)
1074  work1(3)=nlt
1075  work1(4)=ngt
1076  CALL gop(work1,work2,'+ ',5)
1077  nlt=work1(3)
1078  ngt=work1(4)
1079  IF (.NOT.ifok) THEN
1080  WRITE(6,101) nid,guess,clt,cgt
1081  WRITE(6,102) nid,(work1(i),i=1,5)
1082  101 FORMAT(i5,'Glg:',3f12.5)
1083  102 FORMAT(i5,'WORK1:',5f12.5)
1084  ENDIF
1085 C
1086 C Done?
1087 C
1088  IF (nlt.GT.n2.OR.ngt.GT.n2) THEN
1089 C we're not done.....
1090  IF (ngt.GT.nlt) THEN
1091 C guess is too low
1092  gmin=cgt
1093  g2=cgt+max(0.,work1(2)/work1(1)-guess)*amp
1094  IF (g2.GT.gmax) g2=0.5*(guess+gmax)
1095  eps=afac*abs(g2-guess)
1096 C see that we move at least as far as the next closest value.
1097  guess=max(g2,cgt)
1098  GOTO 10
1099  ELSE IF (nlt.GT.ngt) THEN
1100 C guess is too high
1101  gmax=clt
1102  g2=clt+min(0.,work1(2)/work1(1)-guess)*amp
1103  IF (g2.LT.gmin) g2=0.5*(guess+gmin)
1104  eps=afac*abs(g2-guess)
1105 C see that we move at least as far as the next closest value.
1106  guess=min(g2,clt)
1107  GOTO 10
1108  ENDIF
1109  ELSE
1110 C
1111 C we're done....
1112  IF (work1(5).NE.0) THEN
1113 C the median is (usually) one of the values
1114  fmdian=guess
1115  IF (work1(5).EQ.1.0) THEN
1116  IF (mod(ntot,2).EQ.0) THEN
1117  IF (ngt.GT.nlt) THEN
1118  fmdian=0.5*(guess+cgt)
1119  ELSE
1120  fmdian=0.5*(guess+clt)
1121  ENDIF
1122  ELSE
1123  IF (ngt.EQ.nlt) THEN
1124  fmdian=guess
1125  ELSE IF(ngt.GT.nlt) THEN
1126  fmdian=cgt
1127  ELSE
1128  fmdian=clt
1129  ENDIF
1130  ENDIF
1131  ENDIF
1132  ELSE
1133  IF (mod(ntot,2).EQ.0) THEN
1134  IF (ngt.EQ.nlt) THEN
1135  fmdian=0.5*(clt+cgt)
1136  ELSE IF(ngt.GT.nlt) THEN
1137  fmdian=0.5*(guess+cgt)
1138  ELSE
1139  fmdian=0.5*(guess+clt)
1140  ENDIF
1141  ELSE
1142  IF (ngt.EQ.nlt) THEN
1143  fmdian=guess
1144  ELSE IF(ngt.GT.nlt) THEN
1145  fmdian=cgt
1146  ELSE
1147  fmdian=clt
1148  ENDIF
1149  ENDIF
1150  ENDIF
1151 C
1152  ENDIF
1153  IF (.NOT.ifok) WRITE(6,*) nid,'FMDIAN2',fmdian,(a(i),i=1,n)
1154  return
1155 C
1156 C Error handling
1157 C
1158  9000 CONTINUE
1159  WRITE(6,11) ntot,gmin0,gmax0,guess
1160  11 FORMAT('ABORTING IN FMDIAN: N,AMIN,AMAX:',i6,3g14.6)
1161  DO 13 i1=1,n,5
1162  in=i1+5
1163  in=min(in,n)
1164  WRITE(6,12) nid,(a(i),i=i1,in)
1165  12 FORMAT(i4,' FMA:',5g14.6)
1166  13 CONTINUE
1167  DO 15 i1=1,itry,5
1168  in=i1+5
1169  in=min(in,itry)
1170  WRITE(6,14) nid,(gues(i),i=i1,in)
1171  14 FORMAT(i4,' FMG:',5g14.6)
1172  15 CONTINUE
1173  call exitt
1174  END
1175 
1176 C========================================================================
1177 C Double precision matrix and vector routines
1178 C========================================================================
1179 
1180 c-----------------------------------------------------------------------
1181  subroutine dcadd(a,const,n)
1182  real*8 a(1),const
1183 C
1184  DO 100 i=1,n
1185  a(i)=a(i)+const
1186  100 CONTINUE
1187  return
1188  END
1189 c-----------------------------------------------------------------------
1190  subroutine dsub2(a,b,n)
1191  real*8 a(1), b(1)
1192 C
1193  DO 100 i=1,n
1194  a(i)=a(i)-b(i)
1195  100 CONTINUE
1196  return
1197  END
1198 C
1199 c-----------------------------------------------------------------------
1200  subroutine dadd2(a,b,n)
1201  real*8 a(1), b(1)
1202 C
1203  DO 100 i=1,n
1204  a(i)=a(i)+b(i)
1205  100 CONTINUE
1206  return
1207  END
1208 c-----------------------------------------------------------------------
1209  subroutine chswapr(b,L,ind,n,temp)
1210  INTEGER IND(1)
1211  CHARACTER*6 B(1),TEMP(1)
1212 C***
1213 C*** SORT ASSOCIATED ELEMENTS BY PUTTING ITEM(JJ)
1214 C*** INTO ITEM(I), WHERE JJ=IND(I).
1215 C***
1216  DO 20 i=1,n
1217  jj=ind(i)
1218  temp(i)=b(jj)
1219  20 CONTINUE
1220  DO 30 i=1,n
1221  30 b(i)=temp(i)
1222  return
1223  END
1224 c-----------------------------------------------------------------------
1225  subroutine drcopy(r,d,N)
1226  real*8 d(1)
1227  dimension r(1)
1228  do 10 i=1,n
1229  r(i)=d(i)
1230  10 continue
1231  return
1232  end
1233 c-----------------------------------------------------------------------
1234  subroutine rrcopy(r,d,N)
1235  real*4 d(1)
1236  real*4 r(1)
1237  do 10 i=1,n
1238  r(i)=d(i)
1239  10 continue
1240  return
1241  end
1242 c-----------------------------------------------------------------------
1243  subroutine sorts(xout,xin,work,n)
1244  real xout(1),xin(1),work(1)
1245  call copy(xout,xin,n)
1246  call sort(xout,work,n)
1247  return
1248  end
1249 C
1250 c-----------------------------------------------------------------------
1251  function ivlsum(a,n)
1252  INTEGER a(1)
1253  INTEGER tsum
1254  if (n.eq.0) then
1255  ivlsum = 0
1256  return
1257  endif
1258  tsum=a(1)
1259  DO 100 i=2,n
1260  tsum=tsum+a(i)
1261  100 CONTINUE
1262  ivlsum=tsum
1263  return
1264  END
1265 c-----------------------------------------------------------------------
1266  subroutine icadd(a,c,n)
1267  INTEGER A(1),C
1268  DO 100 i = 1, n
1269  100 a(i) = a(i) + c
1270  return
1271  END
1272  subroutine isort(a,ind,n)
1273 C
1274 C Use Heap Sort (p 231 Num. Rec., 1st Ed.)
1275 C
1276  integer a(1),ind(1)
1277  integer aa
1278 C
1279  dO 10 j=1,n
1280  ind(j)=j
1281  10 continue
1282 C
1283  if (n.le.1) return
1284  l=n/2+1
1285  ir=n
1286  100 continue
1287  if (l.gt.1) then
1288  l=l-1
1289  aa = a(l)
1290  ii = ind(l)
1291  else
1292  aa = a(ir)
1293  ii = ind(ir)
1294  a(ir) = a( 1)
1295  ind(ir) = ind( 1)
1296  ir=ir-1
1297  if (ir.eq.1) then
1298  a(1) = aa
1299  ind(1) = ii
1300  return
1301  endif
1302  endif
1303  i=l
1304  j=l+l
1305  200 continue
1306  if (j.le.ir) then
1307  if (j.lt.ir) then
1308  if ( a(j).lt.a(j+1) ) j=j+1
1309  endif
1310  if (aa.lt.a(j)) then
1311  a(i) = a(j)
1312  ind(i) = ind(j)
1313  i=j
1314  j=j+j
1315  else
1316  j=ir+1
1317  endif
1318  GOTO 200
1319  endif
1320  a(i) = aa
1321  ind(i) = ii
1322  GOTO 100
1323  end
1324  subroutine sort(a,ind,n)
1325 C
1326 C Use Heap Sort (p 231 Num. Rec., 1st Ed.)
1327 C
1328  real a(1),aa
1329  integer ind(1)
1330 C
1331  dO 10 j=1,n
1332  ind(j)=j
1333  10 continue
1334 C
1335  if (n.le.1) return
1336  l=n/2+1
1337  ir=n
1338  100 continue
1339  if (l.gt.1) then
1340  l=l-1
1341  aa = a(l)
1342  ii = ind(l)
1343  else
1344  aa = a(ir)
1345  ii = ind(ir)
1346  a(ir) = a( 1)
1347  ind(ir) = ind( 1)
1348  ir=ir-1
1349  if (ir.eq.1) then
1350  a(1) = aa
1351  ind(1) = ii
1352  return
1353  endif
1354  endif
1355  i=l
1356  j=l+l
1357  200 continue
1358  if (j.le.ir) then
1359  if (j.lt.ir) then
1360  if ( a(j).lt.a(j+1) ) j=j+1
1361  endif
1362  if (aa.lt.a(j)) then
1363  a(i) = a(j)
1364  ind(i) = ind(j)
1365  i=j
1366  j=j+j
1367  else
1368  j=ir+1
1369  endif
1370  GOTO 200
1371  endif
1372  a(i) = aa
1373  ind(i) = ii
1374  GOTO 100
1375  end
1376 c-----------------------------------------------------------------------
1377  subroutine iswap_ip(x,p,n)
1378  integer x(1),xstart
1379  integer p(1)
1380 c
1381 c In-place permutation: x' = x(p)
1382 c
1383  do k=1,n
1384  if (p(k).gt.0) then ! not swapped
1385  xstart = x(k)
1386  loop_start = k
1387  last = k
1388  do j=k,n
1389  next = p(last)
1390  if (next.lt.0) then
1391  write(6,*) 'Hey! iswap_ip problem.',j,k,n,next
1392  call exitt
1393  elseif (next.eq.loop_start) then
1394  x(last) = xstart
1395  p(last) = -p(last)
1396  goto 10
1397  else
1398  x(last) = x(next)
1399  p(last) = -p(last)
1400  last = next
1401  endif
1402  enddo
1403  10 continue
1404  endif
1405  enddo
1406 c
1407  do k=1,n
1408  p(k) = -p(k)
1409  enddo
1410  return
1411  end
1412 c-----------------------------------------------------------------------
1413  subroutine iswapt_ip(x,p,n)
1414  integer x(1),t1,t2
1415  integer p(1)
1416 c
1417 c In-place permutation: x'(p) = x
1418 c
1419 
1420  do k=1,n
1421  if (p(k).gt.0) then ! not swapped
1422  loop_start = k
1423  next = p(loop_start)
1424  t1 = x(loop_start)
1425  do j=1,n
1426  if (next.lt.0) then
1427  write(6,*) 'Hey! iswapt_ip problem.',j,k,n,next
1428  call exitt
1429  elseif (next.eq.loop_start) then
1430  x(next) = t1
1431  p(next) = -p(next)
1432  goto 10
1433  else
1434  t2 = x(next)
1435  x(next) = t1
1436  t1 = t2
1437  nextp = p(next)
1438  p(next) = -p(next)
1439  next = nextp
1440  endif
1441  enddo
1442  10 continue
1443  endif
1444  enddo
1445 c
1446  do k=1,n
1447  p(k) = -p(k)
1448  enddo
1449  return
1450  end
1451 c-----------------------------------------------------------------------
1452  subroutine swap_ip(x,p,n)
1453  real x(1),xstart
1454  integer p(1)
1455 c
1456 c In-place permutation: x' = x(p)
1457 c
1458  do k=1,n
1459  if (p(k).gt.0) then ! not swapped
1460  xstart = x(k)
1461  loop_start = k
1462  last = k
1463  do j=k,n
1464  next = p(last)
1465  if (next.lt.0) then
1466  write(6,*) 'Hey! swap_ip problem.',j,k,n,next
1467  call exitt
1468  elseif (next.eq.loop_start) then
1469  x(last) = xstart
1470  p(last) = -p(last)
1471  goto 10
1472  else
1473  x(last) = x(next)
1474  p(last) = -p(last)
1475  last = next
1476  endif
1477  enddo
1478  10 continue
1479  endif
1480  enddo
1481 c
1482  do k=1,n
1483  p(k) = -p(k)
1484  enddo
1485  return
1486  end
1487 c-----------------------------------------------------------------------
1488  subroutine swapt_ip(x,p,n)
1489  real x(1),t1,t2
1490  integer p(1)
1491 c
1492 c In-place permutation: x'(p) = x
1493 c
1494 
1495  do k=1,n
1496  if (p(k).gt.0) then ! not swapped
1497  loop_start = k
1498  next = p(loop_start)
1499  t1 = x(loop_start)
1500  do j=1,n
1501  if (next.lt.0) then
1502  write(6,*) 'Hey! swapt_ip problem.',j,k,n,next
1503  call exitt
1504  elseif (next.eq.loop_start) then
1505  x(next) = t1
1506  p(next) = -p(next)
1507  goto 10
1508  else
1509  t2 = x(next)
1510  x(next) = t1
1511  t1 = t2
1512  nextp = p(next)
1513  p(next) = -p(next)
1514  next = nextp
1515  endif
1516  enddo
1517  10 continue
1518  endif
1519  enddo
1520 c
1521  do k=1,n
1522  p(k) = -p(k)
1523  enddo
1524  return
1525  end
1526 c-----------------------------------------------------------------------
1527  subroutine glvadd(x,w,n)
1528  real x(1),w(1)
1529  call gop(x,w,'+ ',n)
1530  return
1531  end
1532 c-----------------------------------------------------------------------
1533  subroutine add3s12(x,y,z,c1,c2,n)
1534  real x(1),y(1),z(1),c1,c2
1535  do i=1,n
1536  x(i) = c1*y(i)+c2*z(i)
1537  enddo
1538  return
1539  end
1540 c-----------------------------------------------------------------------
1541  integer*8 function i8glmax(a,n)
1542  integer*8 a(1),tmax
1543  integer*8 tmp(1),work(1)
1544  tmax= -999999
1545  do i=1,n
1546  tmax=max(tmax,a(i))
1547  enddo
1548  tmp(1)=tmax
1549  call i8gop(tmp,work,'M ',1)
1550  i8glmax=tmp(1)
1551  if (i8glmax .eq. -999999) i8glmax=0
1552  return
1553  end
1554 c-----------------------------------------------------------------------
1555  subroutine admcol3(a,b,c,d,n)
1556  REAL A(1),B(1),C(1),D
1557 C
1558  DO 100 i=1,n
1559  a(i)=a(i)+b(i)*c(i)*d
1560  100 CONTINUE
1561  return
1562  END
1563 c-----------------------------------------------------------------------
1564  subroutine add2col2(a,b,c,n)
1565  real a(1),b(1),c(1)
1566 c
1567  do i=1,n
1568  a(i) = a(i) + b(i)*c(i)
1569  enddo
1570  return
1571  end
1572 c-----------------------------------------------------------------------
1573  subroutine add2sxy(x,a,y,b,n)
1574  real x(1),y(1)
1575 c
1576  do i=1,n
1577  x(i) = a*x(i) + b*y(i)
1578  enddo
1579 c
1580  return
1581  end
1582 c-----------------------------------------------------------------------
1583  subroutine col2s2(x,y,s,n)
1584  real x(n),y(n)
1585 c
1586  do i=1,n
1587  x(i)=s*x(i)*y(i)
1588  enddo
1589 c
1590  return
1591  end
1592 cc-----------------------------------------------------------------------
1593  subroutine transpose(a,lda,b,ldb)
1594  real a(lda,1),b(ldb,1)
1595 c
1596  do j=1,ldb
1597  do i=1,lda
1598  a(i,j) = b(j,i)
1599  enddo
1600  enddo
1601  return
1602  end
1603 c-----------------------------------------------------------------------
1604  subroutine transpose1(a,n)
1605  real a(n,n)
1606 c
1607  do j=1,n
1608  do i=j+1,n
1609  ta = a(i,j)
1610  a(i,j) = a(j,i)
1611  a(j,i) = ta
1612  enddo
1613  enddo
1614  return
1615  end
1616 c-----------------------------------------------------------------------
1617  real function difmax(a,b,n)
1618  real a(1),b(1)
1619 
1620  d=0
1621  do i=1,n
1622  diff = abs(a(i)-b(i))
1623  d = max(d,diff)
1624  enddo
1625  difmax = glamax(d,1)
1626 
1627  return
1628  end
1629 c-----------------------------------------------------------------------
1630 ccc Nek-Nek routines
1631 c-----------------------------------------------------------------------
1632  function glmin_ms(a,n)
1633  include 'SIZE'
1634  include 'PARALLEL'
1635  real a(1)
1636 
1637  call setnekcomm(iglobalcomm)
1638  glmin_ms = glmin(a,n)
1639  call setnekcomm(intracomm)
1640 
1641  return
1642  end
1643 c-----------------------------------------------------------------------
1644  function glamin_ms(a,n)
1645  include 'SIZE'
1646  include 'PARALLEL'
1647  real a(1)
1648 
1649  call setnekcomm(iglobalcomm)
1650  glamin_ms = glamin(a,n)
1651  call setnekcomm(intracomm)
1652 
1653  return
1654  end
1655 c-----------------------------------------------------------------------
1656  function glmax_ms(a,n)
1657  include 'SIZE'
1658  include 'PARALLEL'
1659  real a(1)
1660 
1661  call setnekcomm(iglobalcomm)
1662  glmax_ms = glmax(a,n)
1663  call setnekcomm(intracomm)
1664 
1665  return
1666  end
1667 c------------------------------------------------------------------------
1668  function glamax_ms(a,n)
1669  include 'SIZE'
1670  include 'PARALLEL'
1671  real a(1)
1672 
1673  call setnekcomm(iglobalcomm)
1674  glamax_ms = glamax(a,n)
1675  call setnekcomm(intracomm)
1676 
1677  return
1678  end
1679 c------------------------------------------------------------------------
1680  function glsum_ms(a,n)
1681  include 'SIZE'
1682  include 'PARALLEL'
1683  real a(1)
1684 
1685  call setnekcomm(iglobalcomm)
1686  glsum_ms = glsum(a,n)
1687  call setnekcomm(intracomm)
1688 
1689  return
1690  end
1691 c-----------------------------------------------------------------------
1692  function iglsum_ms(a,n)
1693  include 'SIZE'
1694  include 'PARALLEL'
1695  integer a(1),n
1696 
1697  call setnekcomm(iglobalcomm)
1698  iglsum_ms = iglsum(a,n)
1699  call setnekcomm(intracomm)
1700 
1701  return
1702  end
1703 c-----------------------------------------------------------------------
1704  function glsc3_ms(a,b,c,n)
1705  include 'SIZE'
1706  include 'PARALLEL'
1707  real a(1),b(1),c(1)
1708 
1709  call setnekcomm(iglobalcomm)
1710  glsc3_ms = glsc3(a,b,c,n)
1711  call setnekcomm(intracomm)
1712 
1713  return
1714  end
1715 c-----------------------------------------------------------------------
1716  function glsc2_ms(a,b,n)
1717  include 'SIZE'
1718  include 'PARALLEL'
1719  real a(1),b(1)
1720 
1721  call setnekcomm(iglobalcomm)
1722  glsc2_ms = glsc2(a,b,n)
1723  call setnekcomm(intracomm)
1724 
1725  return
1726  end
1727 c-----------------------------------------------------------------------
1728  function iglmax_ms(a,n)
1729  include 'SIZE'
1730  include 'PARALLEL'
1731  integer a(1)
1732 
1733  call setnekcomm(iglobalcomm)
1734  iglmax_ms = iglmax(a,n)
1735  call setnekcomm(intracomm)
1736 
1737  return
1738  end
1739 c-----------------------------------------------------------------------
1740  function iglmin_ms(a,n)
1741  include 'SIZE'
1742  include 'PARALLEL'
1743  integer a(1)
1744 
1745  call setnekcomm(iglobalcomm)
1746  iglmin_ms = iglmin(a,n)
1747  call setnekcomm(intracomm)
1748 
1749  return
1750  end
1751 c-----------------------------------------------------------------------
void exitt()
Definition: comm_mpi.f:604
subroutine setnekcomm(comm_in)
Definition: comm_mpi.f:1411
subroutine igop(x, w, op, n)
Definition: comm_mpi.f:247
subroutine gop(x, w, op, n)
Definition: comm_mpi.f:201
real *8 function dnekclock()
Definition: comm_mpi.f:393
subroutine i8gop(x, w, op, n)
Definition: comm_mpi.f:275
real function dot(V1, V2, N)
Definition: genxyz.f:885
function glsc3_ms(a, b, c, n)
Definition: math.f:1705
subroutine iswap(b, ind, n, temp)
Definition: math.f:549
subroutine dadd2(a, b, n)
Definition: math.f:1201
subroutine cadd(a, const, n)
Definition: math.f:326
subroutine col3(a, b, c, n)
Definition: math.f:598
subroutine iadd(i1, iscal, n)
Definition: math.f:337
subroutine ascol5(a, b, c, d, e, n)
Definition: math.f:153
subroutine sub2(a, b, n)
Definition: math.f:164
subroutine invers2(a, b, n)
Definition: math.f:51
function glmin(a, n)
Definition: math.f:973
subroutine invcol2(a, b, n)
Definition: math.f:73
integer function ivlmin(vec, n)
Definition: math.f:368
subroutine ifill(ia, ib, n)
Definition: math.f:252
subroutine transpose(a, lda, b, ldb)
Definition: math.f:1594
subroutine add3s12(x, y, z, c1, c2, n)
Definition: math.f:1534
function mod1(i, n)
Definition: math.f:509
subroutine i8zero(a, n)
Definition: math.f:939
subroutine add2col2(a, b, c, n)
Definition: math.f:1565
real function difmax(a, b, n)
Definition: math.f:1618
subroutine icopy(a, b, n)
Definition: math.f:289
subroutine addcol3(a, b, c, n)
Definition: math.f:654
function glmin_ms(a, n)
Definition: math.f:1633
subroutine dsub2(a, b, n)
Definition: math.f:1191
subroutine admcol3(a, b, c, d, n)
Definition: math.f:1556
subroutine ione(a, n)
Definition: math.f:223
function glsc2(x, y, n)
Definition: math.f:794
subroutine copyi4(a, b, n)
Definition: math.f:270
function iglsum(a, n)
Definition: math.f:926
subroutine isort(a, ind, n)
Definition: math.f:1273
real function vlsc21(x, y, n)
Definition: math.f:753
function glsc3(a, b, mult, n)
Definition: math.f:776
function glsum_ms(a, n)
Definition: math.f:1681
real function vlamax(vec, n)
Definition: math.f:406
subroutine add2(a, b, n)
Definition: math.f:622
subroutine icopy84(a, b, n)
Definition: math.f:11
real function gl2norm2(a, n)
Definition: math.f:845
subroutine col2(a, b, n)
Definition: math.f:564
subroutine add2s2(a, b, c1, n)
Definition: math.f:690
real function vlmax(vec, n)
Definition: math.f:396
subroutine copy(a, b, n)
Definition: math.f:260
subroutine addcol4(a, b, c, d, n)
Definition: math.f:142
subroutine rrcopy(r, d, N)
Definition: math.f:1235
subroutine addtnsr(s, h1, h2, h3, nx, ny, nz)
Definition: math.f:477
real function glamin(a, n)
Definition: math.f:887
subroutine add2sxy(x, a, y, b, n)
Definition: math.f:1574
subroutine iswapt_ip(x, p, n)
Definition: math.f:1414
subroutine sorts(xout, xin, work, n)
Definition: math.f:1244
subroutine swapt_ip(x, p, n)
Definition: math.f:1489
subroutine xaddcol3(a, b, c, n)
Definition: math.f:131
subroutine vdot2(dot, u1, u2, v1, v2, n)
Definition: math.f:447
subroutine iflip(i1, n)
Definition: math.f:536
subroutine glvadd(x, w, n)
Definition: math.f:1528
subroutine blank(A, N)
Definition: math.f:19
subroutine add3(a, b, c, n)
Definition: math.f:644
integer function log2(k)
Definition: math.f:527
subroutine iswap_ip(x, p, n)
Definition: math.f:1378
subroutine col4(a, b, c, d, n)
Definition: math.f:120
real function vlsum(vec, n)
Definition: math.f:417
subroutine dcadd(a, const, n)
Definition: math.f:1182
subroutine col2s2(x, y, s, n)
Definition: math.f:1584
function iglmax(a, n)
Definition: math.f:913
function glamin_ms(a, n)
Definition: math.f:1645
subroutine vsqrt(A, N)
Definition: math.f:41
subroutine add4(a, b, c, d, n)
Definition: math.f:728
integer *8 function i8glsum(a, n)
Definition: math.f:947
subroutine icadd(a, c, n)
Definition: math.f:1267
subroutine col2c(a, b, c, n)
Definition: math.f:588
integer *8 function i8glmax(a, n)
Definition: math.f:1542
subroutine izero(a, n)
Definition: math.f:215
function ivlsum(a, n)
Definition: math.f:1252
function iglsum_ms(a, n)
Definition: math.f:1693
function iglmin_ms(a, n)
Definition: math.f:1741
function glmax_ms(a, n)
Definition: math.f:1657
subroutine ltrue(a, n)
Definition: math.f:237
subroutine chswapr(b, L, ind, n, temp)
Definition: math.f:1210
real function vlsc2(x, y, n)
Definition: math.f:739
subroutine subcol3(a, b, c, n)
Definition: math.f:186
function ltrunc(string, l)
Definition: math.f:494
subroutine cadd2(a, b, const, n)
Definition: math.f:346
function fmdian(a, n, ifok)
Definition: math.f:1004
function iglmin(a, n)
Definition: math.f:900
integer function ivlmax(vec, n)
Definition: math.f:382
subroutine gllog(la, lb)
Definition: math.f:986
subroutine transpose1(a, n)
Definition: math.f:1605
function glsum(x, n)
Definition: math.f:861
subroutine drcopy(r, d, N)
Definition: math.f:1226
subroutine add3s2(a, b, c, c1, c2, n)
Definition: math.f:716
subroutine subcol4(a, b, c, d, n)
Definition: math.f:197
subroutine sub3(a, b, c, n)
Definition: math.f:175
subroutine swap_ip(x, p, n)
Definition: math.f:1453
subroutine cmult(a, const, n)
Definition: math.f:315
real function glamax(a, n)
Definition: math.f:874
subroutine chcopy(a, b, n)
Definition: math.f:281
subroutine icopy48(a, b, n)
Definition: math.f:3
function iglmax_ms(a, n)
Definition: math.f:1729
subroutine vsq(A, N)
Definition: math.f:31
function glamax_ms(a, n)
Definition: math.f:1669
subroutine cfill(a, b, n)
Definition: math.f:244
subroutine rone(a, n)
Definition: math.f:230
subroutine i8copy(a, b, n)
Definition: math.f:297
subroutine invcol1(a, n)
Definition: math.f:62
subroutine invcol3(a, b, c, n)
Definition: math.f:98
subroutine add2s1(a, b, c1, n)
Definition: math.f:678
function glsc2_ms(a, b, n)
Definition: math.f:1717
real function vlmin(vec, n)
Definition: math.f:357
subroutine vdot3(dot, u1, u2, u3, v1, v2, v3, n)
Definition: math.f:462
subroutine vcross(u1, u2, u3, v1, v2, v3, w1, w2, w3, n)
Definition: math.f:430
subroutine rzero(a, n)
Definition: math.f:208
function glmax(a, n)
Definition: math.f:960
subroutine chsign(a, n)
Definition: math.f:305
real function gl2norm(a, n)
Definition: math.f:829
function glsc23(x, y, z, n)
Definition: math.f:812
subroutine sort(a, ind, n)
Definition: math.f:1325