Nek5000
SEM for Incompressible NS
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
setprop.F90
Go to the documentation of this file.
1 
2 
3 !------------------------------------------------------------------------
5 subroutine setprop
6  use kinds, only : dp
7  use size_m, only :nx1, ny1, nz1, nfield
8  use ctimer, only : icalld, tspro, nspro, etime1, dnekclock
9  use input, only : ifflow, ifmhd
10  use geom, only : bm1
11  use soln, only : vdiff, vtrans
12  use tstep, only : ifield, nelfld, volfld, avdiff, avtran
13  implicit none
14 
15  integer :: nxyz1, mfield, nfldt, ifld, nel, ntot1
16  real(DP) :: vol
17  real(DP), external :: glsc2
18 
19 #ifndef NOTIMER
20  if (icalld == 0) tspro=0.0
21  icalld=icalld+1
22  nspro=icalld
23  etime1=dnekclock()
24 #endif
25 
26  nxyz1 = nx1*ny1*nz1
27  mfield=2
28  IF (ifflow) mfield=1
29  nfldt = nfield
30  if (ifmhd) nfldt = nfield+1
31 
32  ifld = ifield
33 
34  DO ifield=mfield,nfldt
35 !max IF (IFSTRS .AND. IFIELD == 1) CALL STNRINV ! expensive !
36 
37  CALL vprops
38 
39  nel = nelfld(ifield)
40  vol = volfld(ifield)
41  ntot1 = nxyz1*nel
42 
43  avdiff(ifield) = glsc2(bm1,vdiff(1,1,1,1,ifield),ntot1)/vol
44  avtran(ifield) = glsc2(bm1,vtrans(1,1,1,1,ifield),ntot1)/vol
45 
46  ENDDO
47 
48  ifield = ifld
49 
50 #ifndef NOTIMER
51  tspro=tspro+(dnekclock()-etime1)
52 #endif
53 
54 
55  RETURN
56 end subroutine setprop
cleaned
Definition: tstep_mod.F90:2
Input parameters from preprocessors.
Definition: input_mod.F90:11
Definition: soln_mod.F90:1
subroutine vprops
Set material properties Material type: 0 for default (PARAM and PCOND/PRHOCP) 1 for constant props; 2...
Definition: subs1.F90:750
real(dp) function dnekclock()
Definition: ctimer_mod.F90:103
subroutine setprop
Set variable property arrays.
Definition: setprop.F90:5
Geometry arrays.
Definition: geom_mod.F90:2