Nek5000
SEM for Incompressible NS
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
drive2.F90
Go to the documentation of this file.
1 
3 
4 !-------------------------------------------------------------------
6 !-------------------------------------------------------------------
7 subroutine initdim
8  use size_m
9  use input
10  implicit none
11 
12  nx1=lx1
13  ny1=ly1
14  nz1=lz1
15 
16  nx2=lx2
17  ny2=ly2
18  nz2=lz2
19 
20  nx3=lx3
21  ny3=ly3
22  nz3=lz3
23 
24  nxd=lxd
25  nyd=lyd
26  nzd=lzd
27 
28 
29  nelt=lelt
30  nelv=lelv
31  ndim=ldim
32 
33  RETURN
34 end subroutine initdim
35 
36 !--------------------------------------------------------------------
38 !--------------------------------------------------------------------
39 subroutine initdat
40  use kinds, only : dp
41  use size_m, only : lx1, lx2, lelt, nx1, ny1, nz1, nx2, ny2, nz2
42  use input, only : ifcvode, ifexplvis, ifsplit, param, ccurve, xc, yc, zc
43  use parallel, only : ifgprnt
44  use soln, only : abx1, abx2, aby1, aby2, abz1, abz2, vgradt1, vgradt2
45  use tstep, only : if_full_pres
46  implicit none
47 
48  integer :: nel8, ntot
49 
50 ! Set default logicals
51  ifcvode = .false.
52  ifexplvis = .false.
53 
54  ifsplit = .false.
55  if (lx1 == lx2) ifsplit= .true.
56 
57  if_full_pres = .false.
58 
59 ! Turn off (on) diagnostics for communication
60  ifgprnt= .false.
61 
62  param = 0._dp
63 
64 ! The initialization of CBC is done in READAT
65 
66 ! LCBC = 3*6*LELT*(LDIMT1+1)
67 ! CALL BLANK(CBC,LCBC)
68 
69  CALL blank(ccurve ,12*lelt)
70  nel8 = 8*lelt
71  xc = 0._dp
72  yc = 0._dp
73  zc = 0._dp
74 
75  ntot=nx1*ny1*nz1*lelt
76  abx1 = 0._dp
77  abx2 = 0._dp
78  aby1 = 0._dp
79  aby2 = 0._dp
80  abz1 = 0._dp
81  abz2 = 0._dp
82  vgradt1 = 0._dp
83  vgradt2 = 0._dp
84 
85  ntot=nx2*ny2*nz2*lelt
86 !max CALL RZERO(USRDIV,NTOT)
87 
88  RETURN
89 end subroutine initdat
90 
91 !---------------------------------------------------------------------
93 !---------------------------------------------------------------------
94 subroutine comment
95  use kinds, only : dp
96  use ctimer, only : dnekclock, ttime
97  use geom, only : ifwcno
98  use input, only : ifadvc, iftran
99  use tstep, only : nid, istep, ifield, nfield, lastep, time, dt, courno
100  use tstep, only : re_cell
101  implicit none
102 
103  LOGICAL :: IFCOUR
104  SAVE ifcour
105  REAL(DP) :: EETIME0,EETIME1,EETIME2, TTIME_STP
106  SAVE eetime0,eetime1,eetime2
107  DATA eetime0,eetime1,eetime2 /0.0, 0.0, 0.0/
108 
109 
110 ! Only node zero makes comments.
111  IF (nid /= 0) RETURN
112 
113 
114  IF (eetime0 == 0.0 .AND. istep == 1) eetime0=dnekclock()
115  eetime1=eetime2
116  eetime2=dnekclock()
117 
118  IF (istep == 0) THEN
119  ifcour = .false.
120  DO 10 ifield=1,nfield
121  IF (ifadvc(ifield)) ifcour = .true.
122  10 END DO
123  IF (ifwcno) ifcour = .true.
124  ELSEIF (istep > 0 .AND. lastep == 0 .AND. iftran) THEN
125  ttime_stp = eetime2-eetime1 ! time per timestep
126  ttime = eetime2-eetime0 ! sum of all timesteps
127  IF(istep == 1) THEN
128  ttime_stp = 0
129  ttime = 0
130  ENDIF
131  IF ( ifcour) &
132  WRITE (6,100) istep,time,dt,courno,re_cell,ttime,ttime_stp
133  IF ( .NOT. ifcour) WRITE (6,101) istep,time,dt
134  ELSEIF (lastep == 1) THEN
135  ttime_stp = eetime2-eetime1 ! time per timestep
136  ttime = eetime2-eetime0 ! sum of all timesteps
137  ENDIF
138  100 FORMAT('Step',i7,', t=',1pe14.7,', DT=',1pe10.3 &
139  ,', C=',0pf7.3,', ReC=',0pf8.4,1pe11.4,1pe8.1)
140  101 FORMAT('Step',i7,', time=',1pe12.5,', DT=',1pe11.3)
141 
142  RETURN
143 end subroutine comment
144 
145 !------------------------------------------------------------------------
147 !------------------------------------------------------------------------
148 subroutine setvar
149  use kinds, only : dp
150  use size_m, only : nfield, nid, lorder, nelt, ldimt, lzd, lyd, lxd
151  use size_m, only : nxd, nyd, nzd, nelv
152  use geom, only : ifgmsh3
153  use input, only : param, ifcvode, ifnav, ifnatc, iflomach
154  use input, only : ifadvc, iftmsh, ifmvbd, ifmodel, ifmhd, iftran
155  use input, only : npscal, ifstrs, ifflow, ifsplit, ngeom, ifheat
156  use tstep, only : tolpdf, lastep, iostep, timeio, iocomm, nsteps, fintim
157  use tstep, only : dtinit, dt, gtheta, betag, nmxnl, dtlag, nbd, ifield, tolnl
158  use tstep, only : prelax, nbdinp, tolrel, tolhdf, pi, ctarg, tolabs, nmxe
159  use tstep, only : nelfld, nmxh, nmxp
160 
161  implicit none
162 
163  integer :: NFLDTM, nfldt, mfield, NABMSH, IADV, IFLD1
164  real(DP) :: TLFAC, one
165 
166 ! Enforce splitting/Uzawa according to the way the code was compiled
167  nxd = lxd
168  nyd = lyd
169  nzd = lzd
170 
171 ! Geometry on Mesh 3 or 1?
172  ifgmsh3 = .true.
173  IF ( ifstrs ) ifgmsh3 = .false.
174  IF ( .NOT. ifflow) ifgmsh3 = .false.
175  IF ( ifsplit ) ifgmsh3 = .false.
176 
177  ngeom = 2
178 
179  nfield = 1
180  IF (ifheat) THEN
181  nfield = 2 + npscal
182  nfldtm = 1 + npscal
183  ENDIF
184 
185  nfldt = nfield
186  if (ifmhd) then
187  nfldt = nfield + 1
188  nfldtm = nfldtm + 1
189  endif
190 
191 
192  IF (ifmodel) write(*,*) "Oops: turb"
193 #if 0
194  IF (ifmodel) CALL settmc
195  IF (ifmodel .AND. ifkeps) THEN
196  write(*,*) "Oops: turb"
197  npscal = 1
198  nfldtm = npscal + 1
199  IF (ldimt < nfldtm) THEN
200  WRITE (6,*) 'k-e turbulence model activated'
201  WRITE (6,*) 'Insufficient number of field arrays'
202  WRITE (6,*) 'Rerun through PRE or change SIZE file'
203  call exitt
204  ENDIF
205  nfield = nfield + 2
206  CALL setturb
207  ENDIF
208 #endif
209  mfield = 1
210  IF (ifmvbd) mfield = 0
211 
212  DO 100 ifield=mfield,nfldt+(ldimt-1 - npscal)
213  IF (iftmsh(ifield)) THEN
214  nelfld(ifield) = nelt
215  ELSE
216  nelfld(ifield) = nelv
217  ENDIF
218  100 END DO
219 
220  nmxh = 1000
221  if (iftran) nmxh = 100
222  nmxp = 1000 ! (for testing) 100 ! 2000
223  nmxe = 100 ! 1000
224  nmxnl = 10 ! 100
225 
226  param(86) = 0 ! No skew-symm. convection for now
227 
228  betag = 0 ! PARAM(3)
229  gtheta = 0 ! PARAM(4)
230  dt = abs(param(12))
231  dtinit = dt
232  fintim = param(10)
233  nsteps = int(param(11))
234  iocomm = int(param(13))
235  timeio = param(14)
236  iostep = int(param(15))
237  lastep = 0
238  tolpdf = abs(param(21))
239  tolhdf = abs(param(22))
240  tolrel = abs(param(24))
241  tolabs = abs(param(25))
242  ctarg = param(26)
243  nbdinp = int(param(27))
244  nabmsh = int(param(28))
245 
246  if (nbdinp > lorder) then
247  if (nid == 0) then
248  write(6,*) 'ERROR: torder > lorder.',nbdinp,lorder
249  write(6,*) 'Change SIZEu and recompile entire code.'
250  endif
251  call exitt
252  endif
253 
254  if(abs(param(16)) >= 2) ifcvode = .true.
255 
256 
257 ! Check accuracy requested.
258 
259  IF (tolrel <= 0.) tolrel = 0.01
260 
261 ! Relaxed pressure iteration; maximum decrease in the residual.
262 
263  prelax = 0.1*tolrel
264  IF ( .NOT. iftran .AND. .NOT. ifnav) prelax = 1.e-5
265 
266 ! Tolerance for nonlinear iteration
267 
268  tolnl = 1.e-4
269 
270 ! Fintim overrides nsteps
271  IF (fintim /= 0.) nsteps = 1000000000
272  IF ( .NOT. iftran ) nsteps = 1
273 
274 ! Print interval defaults to 1
275  IF (iocomm == 0) iocomm = nsteps+1
276 
277 
278 ! Set logical for Boussinesq approx (natural convection)
279  ifnatc = .false.
280  IF (betag > 0.) ifnatc= .true.
281  IF(iflomach) ifnatc = .false.
282 
283 ! Set default for mesh integration scheme
284  IF (nabmsh <= 0 .OR. nabmsh > 3) THEN
285  nabmsh = nbdinp
286  param(28) = (nabmsh)
287  ENDIF
288 
289 ! Set default for mixing length factor
290  tlfac = 0.14
291 ! IF (PARAM(49) .LE. 0.0) PARAM(49) = TLFAC
292 
293 ! Courant number only applicable if convection in ANY field.
294  iadv = 0
295  ifld1 = 1
296  IF ( .NOT. ifflow) ifld1 = 2
297  DO 200 ifield=ifld1,nfldt
298  IF (ifadvc(ifield)) iadv = 1
299  200 END DO
300 
301 ! If characteristics, need number of sub-timesteps (DT/DS).
302 ! Current sub-timeintegration scheme: RK4.
303 ! If not characteristics, i.e. standard semi-implicit scheme,
304 ! check user-defined Courant number.
305  IF (iadv == 1) CALL setchar
306 
307 ! Initialize order of time-stepping scheme (BD)
308 ! Initialize time step array.
309  nbd = 0
310  dtlag = 0._dp
311 
312 ! Useful constants
313  one = 1.
314  pi = 4.*atan(one)
315 
316  RETURN
317 end subroutine setvar
318 
320 subroutine echopar
321  use kinds, only : dp
322  use size_m, only : nid, ndim, ldim
323  use input, only : reafle, vnekton, nktonv, param
324  use string, only : ltrunc
325  implicit none
326 
327  CHARACTER(132) :: tmp_string
328  CHARACTER(1) :: tmp_string1(132)
329  equivalence(tmp_string,tmp_string1)
330 
331  real(DP) :: vnekmin
332  integer :: ls, nparam, j, i
333 
334  IF (nid /= 0) RETURN
335 
336  OPEN (unit=9,file=reafle,status='OLD')
337  rewind(unit=9)
338 
339 
340  READ(9,*,err=400)
341  READ(9,*,err=400) vnekton
342  nktonv=vnekton
343  vnekmin=2.5
344  IF(vnekton < vnekmin)THEN
345  print*,' Error: This NEKTON Solver Requires a .rea file'
346  print*,' from prenek version ',vnekmin,' or higher'
347  print*,' Please run the session through the preprocessor'
348  print*,' to bring the .rea file up to date.'
349  call exitt
350  ENDIF
351  READ(9,*,err=400) ndim
352 ! error check
353  IF(ndim /= ldim)THEN
354  WRITE(6,10) ldim,ndim
355  10 FORMAT(//,2x,'Error: This NEKTON Solver has been compiled' &
356  /,2x,' for spatial dimension equal to',i2,'.' &
357  /,2x,' The data file has dimension',i2,'.')
358  CALL exitt
359  ENDIF
360 
361  CALL blank(tmp_string,132)
362  CALL chcopy(tmp_string,reafle,132)
363  ls=ltrunc(tmp_string,132)
364  READ(9,*,err=400) nparam
365  WRITE(6,82) nparam,(tmp_string1(j),j=1,ls)
366 
367  DO 20 i=1,nparam
368  CALL blank(tmp_string,132)
369  READ(9,80,err=400) tmp_string
370  ls=ltrunc(tmp_string,132)
371  IF (param(i) /= 0.0) WRITE(6,81) i,(tmp_string1(j),j=1,ls)
372  20 END DO
373  80 FORMAT(a132)
374  81 FORMAT(i4,3x,132a1)
375  82 FORMAT(i4,3x,'Parameters from file:',132a1)
376  CLOSE (unit=9)
377  write(6,*) ' '
378 
379 ! if(param(2).ne.param(8).and.nid.eq.0) then
380 ! write(6,*) 'Note VISCOS not equal to CONDUCT!'
381 ! write(6,*) 'Note VISCOS =',PARAM(2)
382 ! write(6,*) 'Note CONDUCT =',PARAM(8)
383 ! endif
384 
385  if (param(62) > 0) then
386  if(nid == 0) write(6,*) &
387  'enable byte swap for output'
388  call set_bytesw_write(1)
389  endif
390 
391  return
392 
393 ! Error handling:
394 
395  400 CONTINUE
396  WRITE(6,401)
397  401 FORMAT(2x,'ERROR READING PARAMETER DATA' &
398  ,/,2x,'ABORTING IN ROUTINE ECHOPAR.')
399  CALL exitt
400 
401  WRITE(6,501)
402  501 FORMAT(2x,'ERROR READING LOGICAL DATA' &
403  ,/,2x,'ABORTING IN ROUTINE ECHOPAR.')
404  CALL exitt
405 
406  RETURN
407 end subroutine echopar
408 
409 !----------------------------------------------------------------------
411 !----------------------------------------------------------------------
412 subroutine gengeom (igeom)
413  use size_m, only : nid
414  use geom, only : ifgeom
415  use tstep, only : istep
416  implicit none
417 
418  integer, intent(in) :: igeom
419 
420  if (nid == 0 .AND. istep <= 1) write(6,*) 'generate geometry data'
421 
422  IF (igeom == 1) THEN
423  RETURN
424  ELSEIF (igeom == 2) THEN
425  if (ifgeom) CALL lagmass
426  IF (istep == 0) CALL gencoor()
427 !max IF (ISTEP >= 1) CALL UPDCOOR
428  CALL geom1()!XM3,YM3,ZM3)
429  CALL geom2
430 !max CALL UPDMSYS (1)
431  CALL volume
432  CALL setinvm
433  CALL setdef
434  CALL sfastax
435 !max IF (ISTEP >= 1) CALL EINIT
436  ELSEIF (igeom == 3) THEN
437  ! Take direct stiffness avg of mesh
438  write(*,*) "Oops: igeom"
439 #if 0
440  ifieldo = ifield
441  CALL gencoor(xm3,ym3,zm3)
442  if (ifheat) then
443  ifield = 2
444  CALL dssum(xm3,nx3,ny3,nz3)
445  call col2(xm3,tmult,ntot3)
446  CALL dssum(ym3,nx3,ny3,nz3)
447  call col2(ym3,tmult,ntot3)
448  if (if3d) then
449  CALL dssum(xm3,nx3,ny3,nz3)
450  call col2(xm3,tmult,ntot3)
451  endif
452  else
453  ifield = 1
454  CALL dssum(xm3,nx3,ny3,nz3)
455  call col2(xm3,vmult,ntot3)
456  CALL dssum(ym3,nx3,ny3,nz3)
457  call col2(ym3,vmult,ntot3)
458  if (if3d) then
459  CALL dssum(xm3,nx3,ny3,nz3)
460  call col2(xm3,vmult,ntot3)
461  endif
462  endif
463  CALL geom1(xm3,ym3,zm3)
464  CALL geom2
465 !max CALL UPDMSYS (1)
466  CALL volume
467  CALL setinvm
468  CALL setdef
469  CALL sfastax
470  ifield = ifieldo
471 #endif
472  ENDIF
473 
474  if (nid == 0 .AND. istep <= 1) then
475  write(6,*) 'done :: generate geometry data'
476  write(6,*) ' '
477  endif
478 
479  return
480 end subroutine gengeom
481 
482 !-----------------------------------------------------------------------
484 subroutine files
485  use size_m, only : nid
486  use input, only : session, path
487  use input, only : reafle, re2fle, fldfle, hisfle, schfle
488  use input, only : dmpfle, orefle, nrefle
489  use string, only : ltrunc, indx1
490  implicit none
491 
492  integer :: ls, lpp, lsp, l1, ln, len
493 
494  CHARACTER(132) :: NAME, slash
495  CHARACTER(1) :: DMP(4),FLD(4),REA(4),HIS(4),SCH(4) ,ORE(4), NRE(4)
496  CHARACTER(1) :: RE2(4)
497  CHARACTER(4) :: DMP4 ,FLD4 ,REA4 ,HIS4 ,SCH4 ,ORE4 , NRE4
498  CHARACTER(4) :: RE24
499  equivalence(dmp,dmp4), (fld,fld4), (rea,rea4), (his,his4) &
500  , (sch,sch4), (ore,ore4), (nre,nre4) &
501  , (re2,re24)
502  DATA dmp4,fld4,rea4 /'.dmp','.fld','.rea'/
503  DATA his4,sch4 /'.his','.sch'/
504  DATA ore4,nre4 /'.ore','.nre'/
505  DATA re24 /'.re2' /
506  CHARACTER(78) :: tmp_string
507 
508 ! Find out the session name:
509 
510 ! CALL BLANK(SESSION,132)
511 ! CALL BLANK(PATH ,132)
512 
513 ! ierr = 0
514 ! IF(NID.EQ.0) THEN
515 ! OPEN (UNIT=8,FILE='SESSION.NAME',STATUS='OLD',ERR=24)
516 ! READ(8,10) SESSION
517 ! READ(8,10) PATH
518 ! 10 FORMAT(A132)
519 ! CLOSE(UNIT=8)
520 ! GOTO 23
521 ! 24 ierr = 1
522 ! 23 ENDIF
523 ! call err_chk(ierr,' Cannot open SESSION.NAME!$')
524 
525  slash = '/'
526 
527  len = ltrunc(path,132)
528  if(indx1(path(len:len),slash,1) < 1) then
529  call chcopy(path(len+1:len+1),slash,1)
530  endif
531 
532 ! call bcast(SESSION,132*CSIZE)
533 ! call bcast(PATH,132*CSIZE)
534 
535  CALL blank(reafle,132)
536  CALL blank(re2fle,132)
537  CALL blank(fldfle,132)
538  CALL blank(hisfle,132)
539  CALL blank(schfle,132)
540  CALL blank(dmpfle,132)
541  CALL blank(orefle,132)
542  CALL blank(nrefle,132)
543  CALL blank(name ,132)
544 
545 ! Construct file names containing full path to host:
546 
547  ls=ltrunc(session,132)
548  lpp=ltrunc(path,132)
549  lsp=ls+lpp
550 
551  call chcopy(name(1:1),path,lpp)
552  call chcopy(name(lpp+1:lpp+1),session,ls )
553  l1 = lpp+ls+1
554  ln = lpp+ls+4
555 
556 
557 !.rea file
558  call chcopy(name(l1:l1),rea , 4)
559  call chcopy(reafle ,name,ln)
560 ! write(6,*) 'reafile:',reafle
561 
562 !.re2 file
563  call chcopy(name(l1:l1),re2 , 4)
564  call chcopy(re2fle ,name,ln)
565 
566 !.fld file
567  call chcopy(name(l1:l1),fld , 4)
568  call chcopy(fldfle ,name,ln)
569 
570 !.his file
571  call chcopy(name(l1:l1),his , 4)
572  call chcopy(hisfle ,name,ln)
573 
574 !.sch file
575  call chcopy(name(l1:l1),sch , 4)
576  call chcopy(schfle ,name,ln)
577 
578 
579 !.dmp file
580  call chcopy(name(l1:l1),dmp , 4)
581  call chcopy(dmpfle ,name,ln)
582 
583 !.ore file
584  call chcopy(name(l1:l1),ore , 4)
585  call chcopy(orefle ,name,ln)
586 
587 !.nre file
588  call chcopy(name(l1:l1),nre , 4)
589  call chcopy(nrefle ,name,ln)
590 
591 ! Write the name of the .rea file to the logfile.
592 
593  IF (nid == 0) THEN
594  CALL chcopy(tmp_string,reafle,78)
595  WRITE(6,1000) tmp_string
596  WRITE(6,1001)
597  1000 FORMAT(//,2x,'Beginning session:',/,2x,a78)
598  1001 FORMAT(/,' ')
599  ENDIF
600 
601 
602  RETURN
603 
604 end subroutine files
605 
606 !----------------------------------------------------------------------
611 !----------------------------------------------------------------------
612 subroutine settime
613  use kinds, only : dp
614  use geom, only : ifsurt
615  use input, only : ifmvbd, param, ifprint
616  use tstep, only : dtlag, ab, abmsh, bd, dt, iocomm
617  use tstep, only : istep, nab, nbd, nbdinp, time, timef
618  implicit none
619 
620  integer :: ilag, irst, nabmsh, nbdmsh
621 
622  irst = int(param(46))
623 
624 ! Set time step.
625 
626  DO 10 ilag=10,2,-1
627  dtlag(ilag) = dtlag(ilag-1)
628  10 END DO
629  CALL setdt
630  dtlag(1) = dt
631  IF (istep == 1 .AND. irst <= 0) dtlag(2) = dt
632 
633 ! Set time.
634 
635  timef = time
636  time = time+dt
637 
638 ! Set coefficients in AB/BD-schemes.
639 
640  CALL setordbd
641  if (irst > 0) nbd = nbdinp
642  bd = 0._dp
643  CALL setbd(bd,dtlag,nbd)
644  nab = 3
645  IF (istep <= 2 .AND. irst <= 0) nab = istep
646  ab = 0._dp
647  CALL setabbd(ab,dtlag,nab,nbd)
648  IF (ifmvbd) THEN
649  nbdmsh = 1
650  nabmsh = int(param(28))
651  IF (nabmsh > istep .AND. irst <= 0) nabmsh = istep
652  IF (ifsurt) nabmsh = nbd
653  abmsh = 0._dp
654  CALL setabbd(abmsh,dtlag,nabmsh,nbdmsh)
655  ENDIF
656 
657 
658 ! Set logical for printout to screen/log-file
659 
660  ifprint = .false.
661  IF (iocomm > 0 .AND. mod(istep,iocomm) == 0) ifprint= .true.
662  IF (istep == 1 .OR. istep == 0 ) ifprint= .true.
663 
664  RETURN
665  end subroutine settime
666 
667 !-----------------------------------------------------------------------
673 !-----------------------------------------------------------------------
674 subroutine geneig (igeom)
675  use eigen, only : ifaa, ifae, ifas, ifast, ifga, ifge, ifgs, ifgst
676  use input, only : ifflow, ifheat
677  use tstep, only : ifield, imesh, tolev, tolhdf, tolhe, tolhr, tolhs
678  use tstep, only : tolpdf, tolps
679 
680  implicit none
681 
682  integer, intent(in) :: igeom
683 
684  IF (igeom == 1) RETURN
685 
686 ! Decide which eigenvalues to be computed.
687 
688  IF (ifflow) THEN
689 
690  ifaa = .false.
691  ifae = .false.
692  ifas = .false.
693  ifast = .false.
694  ifga = .true.
695  ifge = .false.
696  ifgs = .false.
697  ifgst = .false.
698 
699  ! For now, only compute eigenvalues during initialization.
700  ! For deforming geometries the eigenvalues should be
701  ! computed every time step (based on old eigenvectors => more memory)
702 
703  imesh = 1
704  ifield = 1
705  tolev = 1.e-3
706  tolhe = tolhdf
707  tolhr = tolhdf
708  tolhs = tolhdf
709  tolps = tolpdf
710  CALL eigenv
711  CALL esteig
712 
713  ELSEIF (ifheat .AND. .NOT. ifflow) THEN
714 
715  CALL esteig
716 
717  ENDIF
718 
719  RETURN
720 end subroutine geneig
721 
722 !-----------------------------------------------------------------------
731 !-----------------------------------------------------------------------
732 subroutine fluid (igeom)
733  use kinds, only : dp
734  use size_m, only : nid
735  use ctimer, only : dnekclock
736  use input, only : ifnav, ifsplit, iftran
737  use tstep, only : ifield, imesh, istep, time
738 
739  implicit none
740 
741  integer, intent(inout) :: igeom
742 
743  real(DP) :: ts
744 
745  ifield = 1
746  imesh = 1
747  call unorm
748  call settolv
749 
750  ts = dnekclock()
751 
752  if(nid == 0 .AND. igeom == 1) &
753  write(6,*) 'Solving for fluid',ifsplit,iftran,ifnav
754 
755  if (ifsplit) then
756 
757  ! PLAN 4: TOMBO SPLITTING
758  ! - Time-dependent Navier-Stokes calculation (Re>>1).
759  ! - Same approximation spaces for pressure and velocity.
760  ! - Incompressibe or Weakly compressible (div u .ne. 0).
761 
762  call plan4
763  igeom = 2
764 #if 0
765  call twalluz(igeom) ! Turbulence model
766 #endif
767 !max call chkptol ! check pressure tolerance
768 !max if (param(55) /= 0) call vol_flow ! check for fixed flow rate
769 
770  elseif (iftran) then
771 #if 0
772 
773  ! call plan1 (igeom) ! Orig. NEKTON time stepper
774 
775  call plan3(igeom) ! Same as PLAN 1 w/o nested iteration
776  ! Std. NEKTON time stepper !
777  if (ifmodel) call twalluz(igeom) ! Turbulence model
778  if (igeom >= 2) call chkptol ! check pressure tolerance
779 !max if (igeom >= 2 and param(55) /= 0) call vol_flow ! check for fixed flow rate
780 
781 #endif
782  else ! steady Stokes, non-split
783 
784  ! - Steady/Unsteady Stokes/Navier-Stokes calculation.
785  ! - Consistent approximation spaces for velocity and pressure.
786  ! - Explicit treatment of the convection term.
787  ! - Velocity/stress formulation.
788 #if 0
789  call plan1(igeom) ! The NEKTON "Classic".
790 #endif
791  endif
792 
793  if(nid == 0 .AND. igeom >= 2) &
794  write(*,'(4x,i7,1x,1p2e12.4,a)') &
795  istep,time,dnekclock()-ts,' Fluid done'
796 
797  return
798 end subroutine fluid
799 
800 !-----------------------------------------------------------------------
811 !-----------------------------------------------------------------------
812 subroutine heat (igeom)
813  use kinds, only : dp
814  use size_m, only : nid, nfield
815  use ctimer, only : dnekclock
816  use input, only : ifcvode, ifsplit, iftmsh
817  use tstep, only : ifield, imesh, istep, time
818 
819  implicit none
820 
821  integer, intent(inout) :: igeom
822  real(DP) :: ts
823  integer :: intype, igeo
824 
825  ts = dnekclock()
826 
827  if (nid == 0 .AND. igeom == 1) &
828  write(*,'(13x,a)') 'Solving for heat'
829 
830  if (ifcvode) then
831 #if 0
832 
833  call cdscal_cvode(igeom)
834  igeom = 2
835 #endif
836  elseif (ifsplit) then
837 
838  do igeo=1,2
839  do ifield=2,nfield
840  intype = -1
841  if ( .NOT. iftmsh(ifield)) imesh = 1
842  if ( iftmsh(ifield)) imesh = 2
843  call unorm
844  call settolt
845  call cdscal(igeo)
846  enddo
847  enddo
848  igeom = 2
849 
850  else ! PN-PN-2
851 
852  do ifield=2,nfield
853  intype = -1
854  if ( .NOT. iftmsh(ifield)) imesh = 1
855  if ( iftmsh(ifield)) imesh = 2
856  call unorm
857  call settolt
858  call cdscal(igeom)
859  enddo
860 
861  endif
862 
863  if (nid == 0 .AND. igeom >= 2) &
864  write(*,'(4x,i7,1x,1p2e12.4,a)') &
865  istep,time,dnekclock()-ts,' Heat done'
866 
867  return
868 end subroutine heat
869 
870 !-----------------------------------------------------------------------
872 subroutine time00
873  use ctimer
874  implicit none
875 
876  nmxmf=0
877  nmxms=0
878  ndsum=0
879  nvdss=0
880  nsett=0
881  ncdtp=0
882  npres=0
883  nmltd=0
884  ngsum=0
885  nprep=0
886  ndsnd=0
887  ndadd=0
888  nhmhz=0
889  ngop =0
890  nusbc=0
891  ncopy=0
892  ninvc=0
893  ninv3=0
894  nsolv=0
895  nslvb=0
896  nddsl=0
897  ncrsl=0
898  ndott=0
899  nbsol=0
900  nadvc=0
901  nspro=0
902 
903  tmxmf=0.0
904  tmxms=0.0
905  tdsum=0.0
906  tvdss=0.0
907  tvdss=0.0
908  tdsmn=9.9e9
909  tdsmx=0.0
910  tsett=0.0
911  tcdtp=0.0
912  tpres=0.0
913  teslv=0.0
914  tmltd=0.0
915  tgsum=0.0
916  tgsmn=9.9e9
917  tgsmx=0.0
918  tprep=0.0
919  tdsnd=0.0
920  tdadd=0.0
921  thmhz=0.0
922  tgop =0.0
923  tusbc=0.0
924  tcopy=0.0
925  tinvc=0.0
926  tinv3=0.0
927  tsolv=0.0
928  tslvb=0.0
929  tddsl=0.0
930  tcrsl=0.0
931  tdott=0.0
932  tbsol=0.0
933  tbso2=0.0
934  tspro=0.0
935  tadvc=0.0
936  ttime=0.0
937 
938  return
939 end subroutine time00
940 
941 !-----------------------------------------------------------------------
943 subroutine runstat
944 #ifndef NOTIMER
945  use kinds, only : dp
946  use ctimer, only : ifsync, nadvc, ncdtp, ncrsl, ndadd, nddsl
947  use ctimer, only : pinvc, pinv3, phmhz, peslv, pdsum, pddsl, pdadd
948  use ctimer, only : pvdss, pusbc, pspro, psolv, ppres, pmltd, pcrsl
949  use ctimer, only : pcdtp
950  use ctimer, only : tgop_sync, tgop, teslv, tdsum, tddsl, tdadd, tcrsl, tcdtp
951  use ctimer, only : tspro, tsolv, tpres, tprep, tmltd, tinvc, tinv3, thmhz
952  use ctimer, only : tadvc, twal, tvdss, tusbc, tttstp, ttime, tsyc
953  use ctimer, only : nwal, nvdss, nusbc, nsyc, nspro, nsolv, npres
954  use ctimer, only : ninvc, ninv3, nhmhz, ngop, neslv, nmltd, ndsum
955  use ctimer, only : dnekclock
956  use ctimer, only : nproj, tproj, proj_flop, proj_mop
957  use ctimer, only : nhconj, thconj, hconj_flop, hconj_mop
958  use ctimer, only : ncggo, tcggo, cggo_flop, cggo_mop
959  use ctimer, only : naxhm, taxhm, axhelm_flop, axhelm_mop
960  use ctimer, only : nintp, tintp, intp_flop, intp_mop
961  use ctimer, only : ngrst, tgrst, grst_flop, grst_mop
962  use ctimer, only : nscn, tscn
963  use ctimer, only : nmakef, tmakef
964  use ctimer, only : nmakeq, tmakeq
965  use ctimer, only : nsetfast, tsetfast
966  use ctimer, only : ndpc, tdpc
967  use ctimer, only : nmg_mask, tmg_mask
968  use ctimer, only : np4misc, tp4misc
969  use ctimer, only : total_flop, total_mop, time_flop, sum_flops
970  use ctimer, only : ngmres, tgmres, gmres_flop, gmres_mop
971  use ctimer, only : nh1mg, th1mg, h1mg_flop, h1mg_mop
972  use ctimer, only : nschw, tschw, schw_flop, schw_mop
973  use ctimer, only : nscps, tscps
974  use ctimer, only : nnmsc, tnmsc
975  use ctimer, only : nnmvc, tnmvc
976  use ctimer, only : ncrespsp, tcrespsp
977  use ctimer, only : ncresvsp, tcresvsp
978  use ctimer, only : nheat2, theat2
979  use ctimer, only : nprep, tprep
980  use ctimer, only : nfoo, tfoo
981  use size_m, only : nid
982  use parallel, only : np
983  implicit none
984 
985  real(DP) :: min_dsum, max_dsum, avg_dsum
986  real(DP) :: min_vdss, max_vdss, avg_vdss
987  real(DP) :: min_gop, max_gop, avg_gop
988  real(DP) :: min_gop_sync, max_gop_sync, avg_gop_sync
989  real(DP) :: min_crsl, max_crsl, avg_crsl
990  real(DP) :: min_usbc, max_usbc, avg_usbc
991  real(DP) :: min_syc, max_syc, avg_syc
992  real(DP) :: min_wal, max_wal, avg_wal
993  real(DP) :: min_irc, max_irc, avg_irc
994  real(DP) :: min_isd, max_isd, avg_isd
995  real(DP) :: min_comm, max_comm, avg_comm
996  real(DP) :: tstop, tirc, tisd, trc, tsd, tcomm, wwork, padvc
997  real(DP) :: total_share
998  integer :: nirc, nisd
999 
1000  real(DP) :: comm_timers(8)
1001  integer :: comm_counters(8)
1002  character(132) :: s132
1003 
1004  tstop=dnekclock()
1005  tttstp=ttime ! sum over all timesteps
1006 
1007 ! call opcount(3) ! print op-counters
1008 
1009  call nek_comm_getstat(comm_timers,comm_counters)
1010  tgop = comm_timers(1)
1011  tgop_sync = comm_timers(2)
1012  twal = comm_timers(3)
1013  tsyc = comm_timers(4)
1014  tirc = comm_timers(5)
1015  tisd = comm_timers(6)
1016  trc = comm_timers(7)
1017  tsd = comm_timers(8)
1018  ngop = comm_counters(1)
1019  nwal = comm_counters(3)
1020  nsyc = comm_counters(4)
1021  nirc = comm_counters(5)
1022  nisd = comm_counters(6)
1023 
1024  tcomm = tisd + tirc + tsyc + tgop + twal + trc + tsd
1025  min_comm = tcomm
1026  call gop(min_comm,wwork,'m ',1)
1027  max_comm = tcomm
1028  call gop(max_comm,wwork,'M ',1)
1029  avg_comm = tcomm
1030  call gop(avg_comm,wwork,'+ ',1)
1031  avg_comm = avg_comm/np
1032 
1033  min_isd = tisd
1034  call gop(min_isd,wwork,'m ',1)
1035  max_isd = tisd
1036  call gop(max_isd,wwork,'M ',1)
1037  avg_isd = tisd
1038  call gop(avg_isd,wwork,'+ ',1)
1039  avg_isd = avg_isd/np
1040 
1041  min_irc = tirc
1042  call gop(min_irc,wwork,'m ',1)
1043  max_irc = tirc
1044  call gop(max_irc,wwork,'M ',1)
1045  avg_irc = tirc
1046  call gop(avg_irc,wwork,'+ ',1)
1047  avg_irc = avg_irc/np
1048 
1049  min_syc = tsyc
1050  call gop(min_syc,wwork,'m ',1)
1051  max_syc = tsyc
1052  call gop(max_syc,wwork,'M ',1)
1053  avg_syc = tsyc
1054  call gop(avg_syc,wwork,'+ ',1)
1055  avg_syc = avg_syc/np
1056 
1057  min_wal = twal
1058  call gop(min_wal,wwork,'m ',1)
1059  max_wal = twal
1060  call gop(max_wal,wwork,'M ',1)
1061  avg_wal = twal
1062  call gop(avg_wal,wwork,'+ ',1)
1063  avg_wal = avg_wal/np
1064 
1065  min_gop = tgop
1066  call gop(min_gop,wwork,'m ',1)
1067  max_gop = tgop
1068  call gop(max_gop,wwork,'M ',1)
1069  avg_gop = tgop
1070  call gop(avg_gop,wwork,'+ ',1)
1071  avg_gop = avg_gop/np
1072 
1073  min_gop_sync = tgop_sync
1074  call gop(min_gop_sync,wwork,'m ',1)
1075  max_gop_sync = tgop_sync
1076  call gop(max_gop_sync,wwork,'M ',1)
1077  avg_gop_sync = tgop_sync
1078  call gop(avg_gop_sync,wwork,'+ ',1)
1079  avg_gop_sync = avg_gop_sync/np
1080 
1081  min_vdss = tvdss
1082  call gop(min_vdss,wwork,'m ',1)
1083  max_vdss = tvdss
1084  call gop(max_vdss,wwork,'M ',1)
1085  avg_vdss = tvdss
1086  call gop(avg_vdss,wwork,'+ ',1)
1087  avg_vdss = avg_vdss/np
1088 
1089  min_dsum = tdsum
1090  call gop(min_dsum,wwork,'m ',1)
1091  max_dsum = tdsum
1092  call gop(max_dsum,wwork,'M ',1)
1093  avg_dsum = tdsum
1094  call gop(avg_dsum,wwork,'+ ',1)
1095  avg_dsum = avg_dsum/np
1096 
1097 
1098  min_crsl = tcrsl
1099  call gop(min_crsl,wwork,'m ',1)
1100  max_crsl = tcrsl
1101  call gop(max_crsl,wwork,'M ',1)
1102  avg_crsl = tcrsl
1103  call gop(avg_crsl,wwork,'+ ',1)
1104  avg_crsl = avg_crsl/np
1105 
1106  min_usbc = tusbc
1107  call gop(min_usbc,wwork,'m ',1)
1108  max_usbc = tusbc
1109  call gop(max_usbc,wwork,'M ',1)
1110  avg_usbc = tusbc
1111  call gop(avg_usbc,wwork,'+ ',1)
1112  avg_usbc = avg_usbc/np
1113 
1114  tttstp = tttstp + 1e-7
1115  if (nid == 0) then
1116  write(6,'(A)') 'runtime statistics:'
1117  write(6,*) 'total time',tttstp
1118 
1119  ! Projection timings
1120  total_share = 0._dp
1121 
1122  call print_times('cggo time', ncggo , tcggo , tttstp, total_share)
1123  call print_times('axhm time', naxhm , taxhm , tttstp, total_share)
1124  call print_times('h1mg time', nh1mg , th1mg , tttstp, total_share)
1125  call print_times('schw time', nschw , tschw , tttstp, total_share)
1126  call print_times('gmrs time', ngmres , tgmres , tttstp, total_share)
1127  call print_times('prhs time', ncrespsp, tcrespsp, tttstp, total_share)
1128  call print_times('vrhs time', ncresvsp, tcresvsp, tttstp, total_share)
1129  !call print_times('grst time', ngrst , tgrst , tttstp, total_share)
1130  !call print_times('intp time', nintp , tintp , tttstp, total_share)
1131  call print_times('proj time', nproj , tproj , tttstp, total_share)
1132  call print_times('hcoj time', nhconj , thconj , tttstp, total_share)
1133  call print_times('dpc time', ndpc , tdpc , tttstp, total_share)
1134  call print_times('scps time', nscps , tscps , tttstp, total_share)
1135  call print_times('stft time', nsetfast, tsetfast, tttstp, total_share)
1136  call print_times('scn time', nscn , tscn , tttstp, total_share)
1137  call print_times('het2 time', nheat2 , theat2 , tttstp, total_share)
1138  call print_times('p4mc time', np4misc , tp4misc , tttstp, total_share)
1139  call print_times('makf time', nmakef , tmakef , tttstp, total_share)
1140  call print_times('makq time', nmakeq , tmakeq , tttstp, total_share)
1141  call print_times('prep time', nprep , tprep , tttstp, total_share)
1142  call print_times('mmsk time', nmg_mask, tmg_mask, tttstp, total_share)
1143  call print_times('nmsc time', nnmsc , tnmsc , tttstp, total_share)
1144  call print_times('nmvc time', nnmvc , tnmvc , tttstp, total_share)
1145  call print_times('advc time', nadvc , tadvc , tttstp, total_share)
1146  call print_times('foo time', nfoo , tfoo , tttstp, total_share)
1147 
1148  write(6,'(A,F8.4)') "Still missing ", (tttstp - total_share) / tttstp
1149 
1150  call print_flops('cggo flop/s', cggo_flop, cggo_mop, tcggo)
1151  call print_flops('gmres flop/s', gmres_flop, gmres_mop, tgmres)
1152  call print_flops('h1mg flop/s', h1mg_flop, h1mg_mop, th1mg)
1153  call print_flops('schw flop/s', schw_flop, schw_mop, tschw)
1154  call print_flops('axhm flop/s', axhelm_flop, axhelm_mop, taxhm)
1155  !call print_flops('grst flop/s', grst_flop, grst_mop, tgrst)
1156  !call print_flops('intp flop/s', intp_flop, intp_mop, tintp)
1157  call print_flops('proj flop/s', proj_flop, proj_mop, tproj)
1158  call print_flops('hcoj flop/s', hconj_flop, hconj_mop, thconj)
1159  write(6,*) ""
1160  call sum_flops()
1161  call print_flops('Subset FLOPS/s', total_flop, total_mop, time_flop)
1162  call print_flops('Total FLOPS/s', total_flop, total_mop, tttstp)
1163  write(6,*) ""
1164 
1165 
1166  ! pcopy=tcopy/tttstp
1167  ! write(6,*) 'copy time',ncopy,tcopy,pcopy
1168  ! pmxmf=tmxmf/tttstp
1169  ! write(6,*) 'mxmf time',nmxmf,tmxmf,pmxmf
1170 
1171  pinv3=tinv3/tttstp
1172  write(6,*) 'inv3 time',ninv3,tinv3,pinv3
1173  pinvc=tinvc/tttstp
1174  write(6,*) 'invc time',ninvc,tinvc,pinvc
1175  pmltd=tmltd/tttstp
1176  write(6,*) 'mltd time',nmltd,tmltd,pmltd
1177  pcdtp=tcdtp/tttstp
1178  write(6,*) 'cdtp time',ncdtp,tcdtp,pcdtp
1179  peslv=teslv/tttstp
1180  write(6,*) 'eslv time',neslv,teslv,peslv
1181 
1182  ! Pressure solver timings
1183  ppres=tpres/tttstp
1184  write(6,*) 'pres time',npres,tpres,ppres
1185 
1186  ! Coarse grid solver timings
1187  pcrsl=tcrsl/tttstp
1188  write(6,*) 'crsl time',ncrsl,tcrsl,pcrsl
1189  write(6,*) 'crsl min ',min_crsl
1190  write(6,*) 'crsl max ',max_crsl
1191  write(6,*) 'crsl avg ',avg_crsl
1192 
1193  ! Helmholz solver timings
1194  phmhz=thmhz/tttstp
1195  write(6,*) 'hmhz time',nhmhz,thmhz,phmhz
1196 
1197 
1198  pspro=tspro/tttstp
1199  write(6,*) 'spro time',nspro,tspro,pspro
1200 
1201  ! USERBC timings
1202  pusbc=tusbc/tttstp
1203  write(6,*) 'usbc time',nusbc,tusbc,pusbc
1204  write(6,*) 'usbc min ',min_usbc
1205  write(6,*) 'usbc max ',max_usbc
1206  write(6,*) 'usb avg ',avg_usbc
1207 
1208 
1209  ! Convection timings
1210  padvc=tadvc/tttstp
1211  write(6,*) 'advc time',nadvc,tadvc,padvc
1212 
1213  ! Vector direct stiffness summuation timings
1214  pvdss=tvdss/tttstp
1215  write(6,*) 'vdss time',nvdss,tvdss,pvdss
1216  write(6,*) 'vdss min ',min_vdss
1217  write(6,*) 'vdss max ',max_vdss
1218  write(6,*) 'vdss avg ',avg_vdss
1219 
1220  ! Direct stiffness summuation timings
1221  pdsum=tdsum/tttstp
1222  write(6,*) 'dsum time',ndsum,tdsum,pdsum
1223  write(6,*) 'dsum min ',min_dsum
1224  write(6,*) 'dsum max ',max_dsum
1225  write(6,*) 'dsum avg ',avg_dsum
1226 
1227  ! pgsum=tgsum/tttstp
1228  ! write(6,*) 'gsum time',ngsum,tgsum,pgsum
1229 
1230  ! pdsnd=tdsnd/tttstp
1231  ! write(6,*) 'dsnd time',ndsnd,tdsnd,pdsnd
1232 
1233  pdadd=tdadd/tttstp
1234  write(6,*) 'dadd time',ndadd,tdadd,pdadd
1235 
1236  ! pdsmx=tdsmx/tttstp
1237  ! write(6,*) 'dsmx time',ndsmx,tdsmx,pdsmx
1238  ! pdsmn=tdsmn/tttstp
1239  ! write(6,*) 'dsmn time',ndsmn,tdsmn,pdsmn
1240  ! pslvb=tslvb/tttstp
1241  ! write(6,*) 'slvb time',nslvb,tslvb,pslvb
1242  pddsl=tddsl/tttstp
1243  write(6,*) 'ddsl time',nddsl,tddsl,pddsl
1244 
1245  psolv=tsolv/tttstp
1246  write(6,*) 'solv time',nsolv,tsolv,psolv
1247 
1248  ! psett=tsett/tttstp
1249  ! write(6,*) 'sett time',nsett,tsett,psett
1250 
1251  ! pbsol=tbsol/tttstp
1252  ! write(6,*) 'bsol time',nbsol,tbsol,pbsol
1253  ! pbso2=tbso2/tttstp
1254  ! write(6,*) 'bso2 time',nbso2,tbso2,pbso2
1255 
1256 #ifdef MPITIMER
1257  write(6,'(/,A)') 'MPI timings'
1258  ! MPI timings
1259  write(6,*) 'total comm time',tcomm, max_comm/ttime
1260  write(6,*) 'comm min ',min_comm
1261  write(6,*) 'comm max ',max_comm
1262  write(6,*) 'comm avg ',avg_comm
1263 
1264  ! MPI_Barrier timings
1265  psyc=tsyc/tcomm
1266  write(6,*) 'barrier time',nsyc,tsyc,psyc
1267  write(6,*) 'barrier min ',min_syc
1268  write(6,*) 'barrier max ',max_syc
1269  write(6,*) 'barrier avg ',avg_syc
1270 
1271  ! MPI_Waitall timings
1272  pwal=twal/tcomm
1273  write(6,*) 'waitall time',nwal,twal,pwal
1274  write(6,*) 'waitall min ',min_wal
1275  write(6,*) 'waitall max ',max_wal
1276  write(6,*) 'waitall avg ',avg_wal
1277 
1278  ! MPI_Allreduce timings
1279  pgop=tgop/tcomm
1280  write(6,*) 'allreduce time',ngop,tgop,pgop
1281  write(6,*) 'allreduce min ',min_gop
1282  write(6,*) 'allreduce max ',max_gop
1283  write(6,*) 'allreduce avg ',avg_gop
1284 
1285  ! MPI_Allreduce(sync) timings
1286  pgop_sync=tgop_sync/tcomm
1287  write(6,*) 'allreduce_sync time',tgop_sync,pgop_sync
1288  write(6,*) 'allreduce_sync min ',min_gop_sync
1289  write(6,*) 'allreduce_sync max ',max_gop_sync
1290  write(6,*) 'allreduce_sync avg ',avg_gop_sync
1291 #endif
1292  endif
1293 
1294  if (np > 1024) return
1295  if (nid == 0) & ! header for timing
1296  write(6,1) 'tusbc','tdadd','tcrsl','tvdss','tdsum',' tgop',ifsync
1297  1 format(/,'#',2x,'nid',6(7x,a5),4x,'qqq',1x,l4)
1298 
1299  call blank(s132,132)
1300  write(s132,132) nid,tusbc,tdadd,tcrsl,tvdss,tdsum,tgop
1301  132 format(i12,1p6e12.4,' qqq')
1302  call pprint_all(s132,132,6)
1303 
1304 #endif
1305 
1306  return
1307 end subroutine runstat
1308 
1309 !-----------------------------------------------------------------------
1311 subroutine pprint_all(s,n_in,io)
1312  use size_m, only : nid
1313  use parallel, only : np
1314  use string, only : ltrunc
1315  implicit none
1316 
1317  integer, intent(in) :: n_in, io
1318  character(n_in), intent(in) :: s
1319  character(132) :: w
1320  integer :: n, mtag, m, i, k, l
1321 
1322  n = min(132,n_in)
1323 
1324  mtag = 999
1325  m = 1
1326  call nekgsync()
1327 
1328  if (nid == 0) then
1329  l = ltrunc(s,n)
1330  write(io,1) (s(k:k),k=1,l)
1331  1 format(132a1)
1332 
1333  do i=1,np-1
1334  call csend(mtag,s,1,i,0) ! send handshake
1335  m = 132
1336  call blank(w,m)
1337  call crecv(i,w,m)
1338  if (m <= 132) then
1339  l = ltrunc(w,m)
1340  write(io,1) (w(k:k),k=1,l)
1341  else
1342  write(io,*) 'pprint long message: ',i,m
1343  l = ltrunc(w,132)
1344  write(io,1) (w(k:k),k=1,l)
1345  endif
1346  enddo
1347  else
1348  call crecv(mtag,w,m) ! wait for handshake
1349  l = ltrunc(s,n)
1350  call csend(nid,s,l,0,0) ! send data to node 0
1351  endif
1352  return
1353 end subroutine pprint_all
1354 
1355 !-----------------------------------------------------------------------
1357 subroutine opcount(ICALL)
1358  use opctr, only : nrout, dcount, ncall, dct, maxrts
1359  implicit none
1360 
1361  integer, intent(in) :: icall
1362  integer :: i
1363 
1364  if (icall == 1) then
1365  nrout=0
1366  endif
1367  if (icall == 1 .OR. icall == 2) then
1368  dcount = 0.0
1369  do 100 i=1,maxrts
1370  ncall(i) = 0
1371  dct(i) = 0.0
1372  100 END DO
1373  endif
1374  if (icall == 3) then
1375  write(*,*) "Oops: icall = 3"
1376 #if 0
1377  ! Sort and print out diagnostics
1378 
1379  if (nid == 0) then
1380  write(6,*) nid,' opcount',dcount
1381  do i = 1,np-1
1382  call csend(i,idum,4,i,0)
1383  call crecv(i,ddcount,wdsize)
1384  write(6,*) i,' opcount',ddcount
1385  enddo
1386  else
1387  call crecv(nid,idum,4)
1388  call csend(nid,dcount,wdsize,0,0)
1389  endif
1390 
1391  dhc = dcount
1392  call gop(dhc,dwork,'+ ',1)
1393  if (nid == 0) then
1394  write(6,*) ' TOTAL OPCOUNT',dhc
1395  endif
1396 
1397  CALL drcopy(rct,dct,nrout)
1398  CALL sort(rct,ind,nrout)
1399  CALL chswapr(rname,6,ind,nrout,sname)
1400  call iswap(ncall,ind,nrout,idum)
1401 
1402  if (nid == 0) then
1403  do 200 i=1,nrout
1404  write(6,201) nid,rname(i),rct(i),ncall(i)
1405  200 END DO
1406  201 format(2x,' opnode',i4,2x,a6,g18.7,i12)
1407  endif
1408 #endif
1409  endif
1410  return
1411 end subroutine opcount
1412 
1413 !-----------------------------------------------------------------------
1415 subroutine dofcnt
1416  use kinds, only : dp, i8
1417  use size_m, only : nx1, ny1, nz1, nelv, nid
1418  use size_m, only : nx2, ny2, nz2
1419  use input, only : ifflow, ifsplit
1420  use parallel, only : nvtot
1421  use soln, only : vmult, tmult, tmask, v1mask, pmask
1422  implicit none
1423 
1424  real(DP) :: WORK(1)!LCTMP1)
1425 
1426  integer :: nxyz, nel, ntot1, ntot2
1427  real(DP) :: vpts, ppts
1428  real(DP), external :: glsum
1429  integer, external :: glsc2, i8glsum
1430  integer(i8) :: ntot,ntotp,ntotv
1431 
1432  nxyz = nx1*ny1*nz1
1433  nel = nelv
1434 
1435 ! unique points on v-mesh
1436  vpts = glsum(vmult,nel*nxyz) + .1
1437  nvtot=int(vpts)
1438 
1439 ! unique points on pressure mesh
1440  work(1)=nel*nxyz
1441  ppts = glsum(work,1) + .1
1442  ntot=int(ppts)
1443 
1444  if (nid == 0) write(6,'(A,2i13)') &
1445  'gridpoints unique/tot: ',nvtot,ntot
1446 
1447  ntot1=nx1*ny1*nz1*nelv
1448  ntot2=nx2*ny2*nz2*nelv
1449 
1450  ntotv = glsc2(tmult,tmask,ntot1)
1451  ntotp = i8glsum(ntot2,1)
1452 
1453  if (ifflow) ntotv = glsc2(vmult,v1mask,ntot1)
1454  if (ifsplit) ntotp = glsc2(vmult,pmask ,ntot1)
1455  if (nid == 0) write(6,*) ' dofs:',ntotv,ntotp
1456 
1457  return
1458 end subroutine dofcnt
1459 !-----------------------------------------------------------------------
1460 subroutine print_times(label, num_calls, time, total_time, cum_time)
1461  use kinds, only : dp
1462  implicit none
1463  character(*), intent(in) :: label
1464  integer, intent(in) :: num_calls
1465  real(DP), intent(in) :: time, total_time
1466  real(DP), intent(inout) :: cum_time
1467 
1468  cum_time = cum_time + time
1469  write(6,42) label,num_calls, time, time/total_time
1470 
1471  42 format(a,i12,f10.2,f8.4)
1472 end subroutine print_times
1473 !-----------------------------------------------------------------------
1474 subroutine print_flops(label, flops, mops, time)
1475  use kinds, only : i8, dp
1476  implicit none
1477  character(*) :: label
1478  integer(i8) :: flops, mops
1479  real(DP) :: time
1480  real(DP), parameter :: bandwidth = 30.*1024 / 64
1481  !real(DP), parameter :: bandwidth = 59.7*1024 / 4
1482  real(DP), parameter :: compute = 204.*1024 / 64
1483  !real(DP), parameter :: compute = 3400*8
1484  real(DP) :: peak
1485 
1486 
1487  if (time > 1.e-6 .and. flops > 1 .and. mops > 1) then
1488  peak = min(compute, flops * bandwidth / (mops*8))
1489  write(6,'(A,2F14.1,F8.4)') label, &
1490  real(flops, kind=DP) / (1.e6 * time), &
1491  peak, real(flops, kind=DP) / (1.e6 * time * peak)
1492  endif
1493 
1494 end subroutine print_flops
subroutine setdt
Set the new time step. All cases covered.
Definition: subs1.F90:3
subroutine plan4()
Splitting scheme .
Definition: plan4.F90:14
subroutine dssum(u)
Direct stiffness sum.
Definition: dssum.F90:54
cleaned
Definition: tstep_mod.F90:2
subroutine sum_flops()
Definition: ctimer_mod.F90:124
Input parameters from preprocessors.
Definition: input_mod.F90:11
subroutine sort(a, ind, n)
Use Heap Sort (p 231 Num. Rec., 1st Ed.)
Definition: math.F90:476
Definition: soln_mod.F90:1
subroutine setvar
Initialize variables.
Definition: drive2.F90:148
subroutine files
Defines machine specific input and output file names.
Definition: drive2.F90:484
subroutine print_flops(label, flops, mops, time)
Definition: drive2.F90:1474
subroutine esteig
Estimate eigenvalues.
Definition: eigsolv.F90:9
subroutine initdim
Transfer array dimensions to common.
Definition: drive2.F90:7
subroutine iswap(b, ind, n, temp)
SORT ASSOCIATED ELEMENTS BY PUTTING ITEM(JJ) into item(i) where JJ = ind(i)
Definition: math.F90:152
subroutine runstat
print run statistics from ctimer
Definition: drive2.F90:943
subroutine crecv(mtype, buf, lenm)
Definition: comm_mpi.F90:223
void exitt()
Definition: comm_mpi.F90:411
integer function indx1(S1, S2, L2)
Definition: string_mod.F90:43
subroutine opcount(ICALL)
init opcounter
Definition: drive2.F90:1357
subroutine dofcnt
count degrees of freedom
Definition: drive2.F90:1415
subroutine time00
zero the ctimer
Definition: drive2.F90:872
subroutine pprint_all(s, n_in, io)
pretty-print from runstat
Definition: drive2.F90:1311
real(dp) function dnekclock()
Definition: ctimer_mod.F90:103
subroutine setchar
If characteristics, need number of sub-timesteps (DT/DS). Current sub-timeintegration scheme: RK4...
Definition: ssolv.F90:90
subroutine sfastax()
For undeformed elements, set up appropriate elemental matrices and geometric factors for fast evaluat...
Definition: hmholtz.F90:387
subroutine unorm
Norm calculation.
Definition: subs1.F90:163
subroutine fluid(igeom)
Driver for solving the incompressible Navier-Stokes equations.
Definition: drive2.F90:732
subroutine gengeom(igeom)
Generate geometry data.
Definition: drive2.F90:412
#define set_bytesw_write
Definition: byte.c:40
subroutine geneig(igeom)
Compute eigenvalues.
Definition: drive2.F90:674
cleaned
Definition: parallel_mod.F90:2
subroutine eigenv
Compute the following eigenvalues:. EIGAA = minimum eigenvalue of the matrix A (=Laplacian) EIGAE = m...
Definition: eigsolv.F90:107
subroutine blank(A, N)
blank a string
Definition: math.F90:38
subroutine nekgsync()
Definition: comm_mpi.F90:318
Geometry arrays.
Definition: geom_mod.F90:2
subroutine lagmass()
Lag the mass matrix (matrices)
Definition: coef.F90:1144
subroutine setinvm()
Invert the mass matrix. 1) Copy BM1 to BINVM1 2) Perform direct stiffness summation on BINVM1 3) Comp...
Definition: coef.F90:1168
subroutine gop(x, w, op, n)
Global vector commutative operation.
Definition: comm_mpi.F90:104
subroutine volume
Compute the volume based on mesh M1 and mesh M2.
Definition: coef.F90:914
subroutine geom2
Routine to generate all elemental geometric data for mesh 2 (Gauss-Legendre mesh). RXM2, RYM2, RZM2 - dr/dx, dr/dy, dr/dz SXM2, SYM2, SZM2 - ds/dx, ds/dy, ds/dz TXM2, TYM2, TZM2 - dt/dx, dt/dy, dt/dz JACM2 - Jacobian BM2 - Mass matrix.
Definition: coef.F90:725
subroutine chcopy(a, b, n)
Definition: math.F90:63
subroutine gencoor()
Generate xyz coordinates for all elements. Velocity formulation : mesh 3 is used Stress formulation :...
Definition: genxyz.F90:196
subroutine echopar
Echo the nonzero parameters from the readfile to the logfile.
Definition: drive2.F90:320
#define nek_comm_getstat
Definition: nek_comm.c:35
subroutine print_times(label, num_calls, time, total_time, cum_time)
Definition: drive2.F90:1460
static uint np
Definition: findpts_test.c:63
subroutine settolt
Set tolerances for temerature/passive scalar solver.
Definition: ssolv.F90:53
integer function ltrunc(string, l)
Definition: string_mod.F90:260
subroutine cdscal(igeom)
Solve the convection-diffusion equation for passive scalar IPSCAL.
Definition: conduct.F90:3
subroutine geom1()
Routine to generate all elemental geometric data for mesh 1. Velocity formulation : global-to-local m...
Definition: coef.F90:383
subroutine heat(igeom)
Driver for temperature or passive scalar.
Definition: drive2.F90:812
subroutine setdef()
Set up deformed element logical switches.
Definition: genxyz.F90:4
subroutine settolv
Set tolerances for velocity solver.
Definition: ssolv.F90:4
subroutine settime
setup time-stepping
Definition: drive2.F90:612
subroutine comment
No need to comment !!
Definition: drive2.F90:94
subroutine csend(mtype, buf, len, jnid, jpid)
Definition: comm_mpi.F90:209
subroutine initdat
Initialize and set default values.
Definition: drive2.F90:39
Definition: io_mod.F90:9