KTH framework for Nek5000 toolboxes; testing version  0.0.1
planx.f
Go to the documentation of this file.
1  SUBROUTINE plan3 (IGEOM)
2 C-----------------------------------------------------------------------
3 C
4 C Compute pressure and velocity using consistent approximation spaces.
5 C Operator splitting technique.
6 C
7 C-----------------------------------------------------------------------
8  include 'SIZE'
9  include 'INPUT'
10  include 'EIGEN'
11  include 'SOLN'
12  include 'TSTEP'
13 C
14  COMMON /scrns/ resv1(lx1,ly1,lz1,lelv)
15  $ , resv2(lx1,ly1,lz1,lelv)
16  $ , resv3(lx1,ly1,lz1,lelv)
17  $ , dv1(lx1,ly1,lz1,lelv)
18  $ , dv2(lx1,ly1,lz1,lelv)
19  $ , dv3(lx1,ly1,lz1,lelv)
20  COMMON /scrvh/ h1(lx1,ly1,lz1,lelv)
21  $ , h2(lx1,ly1,lz1,lelv)
22 C
23  IF (igeom.EQ.1) THEN
24 C
25 C Old geometry
26 C
27  CALL makef
28 C
29  ELSE
30 C
31 C New geometry, new b.c.
32 C
33  intype = -1
34  call sethlm (h1,h2,intype)
35  call cresvif (resv1,resv2,resv3,h1,h2)
36 
37  call ophinv (dv1,dv2,dv3,resv1,resv2,resv3,h1,h2,tolhv,nmxv)
38  call opadd2 (vx,vy,vz,dv1,dv2,dv3)
39 c
40  call incomprn(vx,vy,vz,pr)
41 C
42  ENDIF
43 C
44  RETURN
45  END
46 C
47  SUBROUTINE lagpres
48 C--------------------------------------------------------------------
49 C
50 C Keep old pressure values
51 C
52 C--------------------------------------------------------------------
53  include 'SIZE'
54  include 'SOLN'
55  include 'TSTEP'
56 
57  common /cgeom/ igeom
58 
59  IF (nbdinp.EQ.3.and.igeom.le.2) THEN
60  ntot2 = lx2*ly2*lz2*nelv
61  CALL copy (prlag,pr,ntot2)
62  ENDIF
63  RETURN
64  END
65 C
66  subroutine cresvif (resv1,resv2,resv3,h1,h2)
67 C---------------------------------------------------------------------
68 C
69 C Compute startresidual/right-hand-side in the velocity solver
70 C
71 C---------------------------------------------------------------------
72  include 'SIZE'
73  include 'TOTAL'
74  REAL RESV1 (LX1,LY1,LZ1,1)
75  REAL RESV2 (LX1,LY1,LZ1,1)
76  REAL RESV3 (LX1,LY1,LZ1,1)
77  REAL H1 (LX1,LY1,LZ1,1)
78  REAL H2 (LX1,LY1,LZ1,1)
79  COMMON /scruz/ w1(lx1,ly1,lz1,lelv)
80  $ , w2(lx1,ly1,lz1,lelv)
81  $ , w3(lx1,ly1,lz1,lelv)
82 
83  common /cgeom/ igeom
84 
85  ntot1 = lx1*ly1*lz1*nelv
86  ntot2 = lx2*ly2*lz2*nelv
87  if (igeom.eq.2) CALL lagvel
88  CALL bcdirvc (vx,vy,vz,v1mask,v2mask,v3mask)
89  CALL bcneutr
90 C
91  call extrapp (pr,prlag)
92  call opgradt (resv1,resv2,resv3,pr)
93  CALL opadd2 (resv1,resv2,resv3,bfx,bfy,bfz)
94  CALL ophx (w1,w2,w3,vx,vy,vz,h1,h2)
95  CALL opsub2 (resv1,resv2,resv3,w1,w2,w3)
96 C
97  RETURN
98  END
99 c-----------------------------------------------------------------------
subroutine bcdirvc(V1, V2, V3, mask1, mask2, mask3)
Definition: bdry.f:574
subroutine bcneutr
Definition: bdry.f:1200
subroutine ophinv(o1, o2, o3, i1, i2, i3, h1, h2, tolh, nmxhi)
Definition: induct.f:1023
subroutine extrapp(p, plag)
Definition: induct.f:370
subroutine incomprn(ux, uy, uz, up)
Definition: induct.f:282
subroutine copy(a, b, n)
Definition: math.f:260
subroutine opgradt(outx, outy, outz, inpfld)
Definition: navier1.f:4096
subroutine makef
Definition: navier1.f:1426
subroutine lagvel
Definition: navier1.f:1930
subroutine ophx(out1, out2, out3, inp1, inp2, inp3, h1, h2)
Definition: navier1.f:743
subroutine opadd2(a1, a2, a3, b1, b2, b3)
Definition: navier1.f:2350
subroutine opsub2(a1, a2, a3, b1, b2, b3)
Definition: navier1.f:2360
subroutine cresvif(resv1, resv2, resv3, h1, h2)
Definition: planx.f:67
subroutine plan3(IGEOM)
Definition: planx.f:2
subroutine lagpres
Definition: planx.f:48
subroutine sethlm(h1, h2, intloc)
Definition: subs1.f:1014