Nek5000
SEM for Incompressible NS
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
ssolv.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------
3 !-------------------------------------------------------------------
4 SUBROUTINE settolv
5  use kinds, only : dp
6  use size_m, only : nx1, ny1, nz1
7  use eigen, only : eigaa, eigas
8  use input, only : iftran, ifnav
9  use soln, only : vdiff, vtrans
10  use tstep, only : ifield, nelfld, istep, tolabs, dt, vnrml8, vnrmsm
11  use tstep, only : tolps, tolrel, tolhv, tolhr, tolhs, tolpdf
12  implicit none
13 
14  REAL(DP) :: avvisc, avdens, vnorm, factor
15  integer :: ntot
16  real(DP), external :: glmin, glmax
17 
18  ntot = nx1*ny1*nz1*nelfld(ifield)
19  avvisc = glmin(vdiff(1,1,1,1,ifield),ntot)
20  avdens = glmax(vtrans(1,1,1,1,ifield),ntot)
21 
22  vnorm = 0.
23  IF (iftran) THEN
24  IF (istep == 1) vnorm = vnrml8
25  IF (istep > 1) vnorm = vnrmsm
26  IF (vnorm == 0.) vnorm = tolabs
27  factor = 1.+(avdens/(eigaa*avvisc*dt))
28  ELSE
29  vnorm = vnrml8
30  IF (vnorm == 0.) vnorm = tolabs
31  factor = 1.
32  ENDIF
33 
34  tolps = tolrel*vnorm * sqrt(eigas)/(4.*factor)
35  tolhv = tolrel*vnorm * sqrt(eigaa)*avvisc/2.
36  tolhv = tolhv/3.
37  IF ( .NOT. iftran .AND. .NOT. ifnav) tolhv = tolhv/10.
38  tolhr = tolhv
39  tolhs = tolhv
40 
41 ! Non-zero default pressure tolerance
42 ! NOTE: This tolerance may change due to precision problems.
43 ! See subroutine CHKSSV
44 
45  IF (tolpdf /= 0.) tolps = tolpdf
46 
47  RETURN
48 END SUBROUTINE settolv
49 
50 !-------------------------------------------------------------------
52 !-------------------------------------------------------------------
53 SUBROUTINE settolt
54  use kinds, only : dp
55  use size_m, only : nx1, ny1, nz1
56  use eigen, only : eigaa
57  use input, only : iftran
58  use soln, only : vdiff
59  use tstep, only : ifield, nelfld, istep, tnrml8, tnrmsm, tolabs, tolht, tolrel
60  implicit none
61 
62  REAL(DP) :: avcond, tnorm
63  integer :: ntot
64  real(DP), external :: glmin
65 
66  ntot = nx1*ny1*nz1*nelfld(ifield)
67  avcond = glmin(vdiff(1,1,1,1,ifield),ntot)
68 
69  tnorm = 0.
70  IF (iftran) THEN
71  IF (istep == 1) tnorm = tnrml8(ifield-1)
72  IF (istep > 1) tnorm = tnrmsm(ifield-1)
73  IF (tnorm == 0.) tnorm = tolabs
74  ELSE
75  tnorm = tnrml8(ifield-1)
76  IF (tnorm == 0.) tnorm = tolabs
77  ENDIF
78 
79  tolht(ifield) = tolrel*tnorm * sqrt(eigaa)*avcond
80 
81  RETURN
82 END SUBROUTINE settolt
83 
84 !-----------------------------------------------------------------------
89 !----------------------------------------------------------------------
90 SUBROUTINE setchar
91  use kinds, only : dp
92  use size_m, only : nid
93  use input, only : ifchar
94  use tstep, only : ctarg, ntaubd
95  implicit none
96 
97  integer :: ict
98  real(DP) :: rict, dct
99 
100  IF (ifchar) THEN
101  ict = int(ctarg)
102  rict = (ict)
103  dct = ctarg-rict
104  IF (dct == 0.) ntaubd = ict
105  IF (dct > 0.) ntaubd = ict+1
106  ! if (param(78).ne.0) then
107  ! ntaupf=int(param(78))
108  ! if (nid.eq.0) write(6,*) ' new ntaubd:',ntaubd,ntaupf
109  ! ntaubd=max(ntaubd,ntaupf)
110  ! endif
111  ELSE
112  ntaubd = 0
113  IF (ctarg > 0.5) THEN
114  IF (nid == 0) &
115  WRITE (6,*) 'Reset the target Courant number to .5'
116  ctarg = 0.5
117  ENDIF
118  ENDIF
119 
120  RETURN
121 END SUBROUTINE setchar
cleaned
Definition: tstep_mod.F90:2
Input parameters from preprocessors.
Definition: input_mod.F90:11
Definition: soln_mod.F90:1
subroutine setchar
If characteristics, need number of sub-timesteps (DT/DS). Current sub-timeintegration scheme: RK4...
Definition: ssolv.F90:90
subroutine settolt
Set tolerances for temerature/passive scalar solver.
Definition: ssolv.F90:53
subroutine settolv
Set tolerances for velocity solver.
Definition: ssolv.F90:4