Nek5000
SEM for Incompressible NS
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
makeq.F90
Go to the documentation of this file.
1 
2 
3 !-----------------------------------------------------------------------
6 subroutine makeq()
7  use kinds, only : dp
8  use input, only : ifmhd, ifaxis, ifadvc, ifchar, iftran, ifcvode
9  use input, only : ifmvbd
10  use tstep, only : ifield
11  use ctimer, only : tmakeq, nmakeq, dnekclock
12  implicit none
13 
14  logical :: if_conv_std
15  real(DP) :: etime
16 
17  nmakeq = nmakeq + 1
18  etime = dnekclock()
19 
20  if_conv_std = .true.
21  if (ifmhd .AND. ifaxis) if_conv_std = .false. ! conv. treated in induct.f
22 
23  call makeq_aux ! nekuq, etc.
24 
25  etime = etime - dnekclock()
26  if (ifadvc(ifield) .AND. .NOT. ifchar .AND. if_conv_std) call convab
27  etime = etime + dnekclock()
28 
29  if (iftran) then
30  if(ifcvode) then
31  write(*,*) "Oops: ifcvode"
32 #if 0
33  ntot = nx1*ny1*nz1*nelfld(ifield)
34  call wlaplacian(w1,t(1,1,1,1,ifield-1),vdiff(1,1,1,1,ifield), &
35  ifield)
36  call add2(bq(1,1,1,1,ifield-1),w1,ntot)
37  if(iftmsh(ifield)) then
38  call dssum(bq,nx1,ny1,nz1)
39  call col2(bq,bintm1,ntot)
40  call col2(bq,bm1,ntot)
41  endif
42 #endif
43  else
44  if (ifmvbd) then ! ifchar is false
45  write(*,*) "Oops: ifmvbd"
46 #if 0
47  call admesht
48  call makeabq
49  call makebdq
50 #endif
51  elseif (ifchar .AND. ifadvc(ifield)) then
52  write(*,*) "Oops: ifchar and ifadvc"
53 #if 0
54  call makeabq
55  call convch
56 #endif
57  else
58  call makeabq
59  call makebdq
60  endif
61  endif
62  endif
63 
64  tmakeq = tmakeq + (dnekclock() - etime)
65 
66  return
67 end subroutine makeq
subroutine dssum(u)
Direct stiffness sum.
Definition: dssum.F90:54
cleaned
Definition: tstep_mod.F90:2
Input parameters from preprocessors.
Definition: input_mod.F90:11
subroutine makeq_aux
Definition: makeq_aux.F90:3
subroutine makeq()
Generate forcing function for the solution of a passive scalar. !! NOTE: Do not change the content of...
Definition: makeq.F90:6
real(dp) function dnekclock()
Definition: ctimer_mod.F90:103
subroutine makeabq
Sum up contributions to 3rd order Adams-Bashforth scheme.
Definition: conduct.F90:240
subroutine makebdq()
Add contributions to F from lagged BD terms.
Definition: conduct.F90:270
subroutine convab()
Eulerian scheme, add convection term to forcing function at current time step.
Definition: conduct.F90:217