KTH framework for Nek5000 toolboxes; testing version  0.0.1
cnht_tools.f
Go to the documentation of this file.
1 
7 !=======================================================================
11  subroutine cnht_register()
12  implicit none
13 
14  include 'SIZE'
15  include 'INPUT'
16  include 'FRAMELP'
17  include 'CNHTD'
18 
19  ! local variables
20  integer lpmid, il
21  character*2 str
22 !-----------------------------------------------------------------------
23  ! check if the current module was already registered
24  call mntr_mod_is_name_reg(lpmid,cnht_name)
25  if (lpmid.gt.0) then
26  call mntr_warn(lpmid,
27  $ 'module ['//trim(cnht_name)//'] already registered')
28  return
29  endif
30 
31  ! find parent module
32  call mntr_mod_is_name_reg(lpmid,'FRAME')
33  if (lpmid.le.0) then
34  lpmid = 1
35  call mntr_abort(lpmid,
36  $ 'parent module ['//'FRAME'//'] not registered')
37  endif
38 
39  ! register module
40  call mntr_mod_reg(cnht_id,lpmid,cnht_name,
41  $ 'Conjugated heat transfer tools')
42 
43  ! register and set active section
44  call rprm_sec_reg(cnht_sec_id,cnht_id,'_'//adjustl(cnht_name),
45  $ 'Runtime paramere section for conj. heat trans. tool module')
46  call rprm_sec_set_act(.true.,cnht_sec_id)
47 
48  ! register parameters
49  call rprm_rp_reg(cnht_sc_id,cnht_sec_id,'SCLN',
50  $ 'Norm scaling factor',rpar_real,0,3.36558,.false.,' ')
51 
52  call rprm_rp_reg(cnht_sv_id,cnht_sec_id,'SCLV',
53  $ 'Velocity scaling factor (Pareto curve)',
54  $ rpar_real,0,0.5,.false.,' ')
55 
56  call rprm_rp_reg(cnht_st_id,cnht_sec_id,'SCLT',
57  $ 'Temperature scaling factor (Pareto curve)',
58  $ rpar_real,0,0.5,.false.,' ')
59 
60  call rprm_rp_reg(cnht_gx_id,cnht_sec_id,'GRX',
61  $ 'X component of gravitational field',
62  $ rpar_real,0,0.0,.false.,' ')
63 
64  call rprm_rp_reg(cnht_gy_id,cnht_sec_id,'GRY',
65  $ 'Y component of gravitational field',
66  $ rpar_real,0,1.0,.false.,' ')
67 
68  if (if3d) call rprm_rp_reg(cnht_gz_id,cnht_sec_id,'GRZ',
69  $ 'Z component of gravitational field',
70  $ rpar_real,0,0.0,.false.,' ')
71 
72  ! set initialisation flag
73  cnht_ifinit=.false.
74 
75  return
76  end subroutine
77 !=======================================================================
81  subroutine cnht_init()
82  implicit none
83 
84  include 'SIZE'
85  include 'FRAMELP'
86  include 'INPUT'
87  include 'SOLN'
88  include 'ADJOINT'
89  include 'CNHTD'
90 
91  ! local variables
92  integer itmp, il
93  real rtmp
94  logical ltmp
95  character*20 ctmp
96 !-----------------------------------------------------------------------
97  ! check if the module was already initialised
98  if (cnht_ifinit) then
99  call mntr_warn(cnht_id,
100  $ 'module ['//trim(cnht_name)//'] already initiaised.')
101  return
102  endif
103 
104  ! get runtime parameters
105  call rprm_rp_get(itmp,rtmp,ltmp,ctmp,cnht_sc_id,rpar_real)
106  cnht_sc = rtmp
107  call rprm_rp_get(itmp,rtmp,ltmp,ctmp,cnht_sv_id,rpar_real)
108  cnht_sv = rtmp
109  call rprm_rp_get(itmp,rtmp,ltmp,ctmp,cnht_st_id,rpar_real)
110  cnht_st = rtmp
111 
112  call rprm_rp_get(itmp,rtmp,ltmp,ctmp,cnht_gx_id,rpar_real)
113  cnht_gx = rtmp
114  call rprm_rp_get(itmp,rtmp,ltmp,ctmp,cnht_gy_id,rpar_real)
115  cnht_gz = rtmp
116  if (if3d) then
117  call rprm_rp_get(itmp,rtmp,ltmp,ctmp,cnht_gz_id,rpar_real)
118  cnht_gz = rtmp
119  endif
120 
121  ! Rayleight and Prandtl numbers
122  cnht_ra = abs(param(2))
123  cnht_ra = abs(param(1))
124  beta_b = cnht_ra
125 
126  ! gravity
127  g_adj(1) = cnht_gx
128  g_adj(2) = cnht_gy
129  g_adj(3) = cnht_gz
130 
131  ! temperature gradient
132  call gradm1(dtdx,dtdy,dtdz,t)
133 
134  ! everything is initialised
135  cnht_ifinit=.true.
136 
137  return
138  end subroutine
139 !=======================================================================
143  logical function cnht_is_initialised()
144  implicit none
145 
146  include 'SIZE'
147  include 'CNHTD'
148 !-----------------------------------------------------------------------
149  cnht_is_initialised = cnht_ifinit
150 
151  return
152  end function
153 !=======================================================================
159  subroutine cnht_forcing(ffx,ffy,ffz,ix,iy,iz,ieg)
160  implicit none
161 
162  include 'SIZE' !
163  include 'INPUT' ! IF3D, IFHEAT, CPFLD
164  include 'PARALLEL' ! GLLEL
165  include 'TSTEP' ! IFIELD
166  include 'SOLN' ! JP, T, TP
167  include 'ADJOINT' ! IFADJ, G_ADJ, DTD[XYZ]
168  include 'CNHTD' ! CHGR[XYZ]
169 
170 ! argument list
171  real ffx, ffy, ffz
172  integer ix,iy,iz,ieg
173 
174 ! local variables
175  integer iel, ip
176  real rtmp
177 !-----------------------------------------------------------------------
178  if (ifheat) then
179  iel=gllel(ieg)
180  if (jp.eq.0) then
181  rtmp = t(ix,iy,iz,iel,ifield)/cpfld(1,2)
182  ffx = ffx + cnht_gx*rtmp
183  ffy = ffy + cnht_gy*rtmp
184  if (if3d) ffz = ffz + cnht_gz*rtmp
185  else
186  ip=ix+nx1*(iy-1+ny1*(iz-1+nz1*(iel-1)))
187  if (.not.ifadj) then
188  rtmp = tp(ip,ifield,jp)/cpfld(1,2)
189  ffx = ffx + g_adj(1)*rtmp
190  ffy = ffy + g_adj(2)*rtmp
191  if (if3d) ffz = ffz + g_adj(3)*rtmp
192  else
193  ffx = ffx - dtdx(ip)*tp(ip,ifield,jp)
194  ffy = ffy - dtdy(ip)*tp(ip,ifield,jp)
195  if (if3d) ffz = ffz - dtdz(ip)*tp(ip,ifield,jp)
196  end if
197  end if
198  endif
199 
200  return
201  end subroutine
202 !=======================================================================
205  subroutine cnht_cpfld_set()
206  implicit none
207 
208  include 'SIZE' !
209  include 'INPUT' ! CPFLD, PARAM
210  include 'ADJOINT' ! IFADJ
211  include 'CNHTD' ! cnht_Ra, cnht_Ra
212 !-----------------------------------------------------------------------
213  if (ifheat) then
214  if (ifadj) then
215  cpfld(1,1)=cnht_ra/sqrt(cnht_ra)
216  cpfld(1,2)=1.0
217 
218  cpfld(2,1)=1.0/sqrt(cnht_ra)
219  cpfld(2,2)=1.0
220  else
221  cpfld(1,1)=1.0/sqrt(cnht_ra)
222  cpfld(1,2)=1.0/cnht_ra
223 
224  cpfld(2,1)=1.0/sqrt(cnht_ra)
225  cpfld(2,2)=1.0
226  endif
227  else
228  if (param(2).lt.0.0) then
229  cpfld(1,1) = -1.0/param(2)
230  else
231  cpfld(1,1) = param(2)
232  endif
233 
234  if (param(1).lt.0.0) then
235  cpfld(1,2) = -1.0/param(1)
236  else
237  cpfld(1,2) = param(1)
238  endif
239  endif
240 
241  return
242  end subroutine
243 !=======================================================================
248  subroutine cnht_oprzero (a1,a2,a3,a4)
249  implicit none
250 
251  include 'SIZE' ! N[XYZ]1, NEL[VT]
252  include 'INPUT' ! IFFLOW, IFHEAT, IF3D
253 
254 ! argument list
255  real a1(1),a2(1),a3(1),a4(1)
256 
257 ! local variables
258  integer ntotv, ntott
259 !-----------------------------------------------------------------------
260  ntotv = nx1*ny1*nz1*nelv
261  ntott = nx1*ny1*nz1*nelt
262 
263  if (ifflow) then
264  call rzero(a1,ntotv)
265  call rzero(a2,ntotv)
266  if(if3d) call rzero(a3,ntotv)
267  endif
268  if (ifheat) call rzero(a4,ntott)
269  return
270  end subroutine
271 !=======================================================================
278  subroutine cnht_opcopy (a1,a2,a3,a4,b1,b2,b3,b4)
279  implicit none
280 
281  include 'SIZE' ! N[XYZ]1, NEL[VT]
282  include 'INPUT' ! IFFLOW, IFHEAT, IF3D
283 
284 ! argument list
285  real a1(1),a2(1),a3(1),a4(1),b1(1),b2(1),b3(1),b4(1)
286 
287 ! local variables
288  integer ntotv, ntott
289 !-----------------------------------------------------------------------
290  ntotv = nx1*ny1*nz1*nelv
291  ntott = nx1*ny1*nz1*nelt
292 
293  if (ifflow) then
294  call copy(a1,b1,ntotv)
295  call copy(a2,b2,ntotv)
296  if(if3d) call copy(a3,b3,ntotv)
297  endif
298  if (ifheat) call copy(a4,b4,ntott)
299 
300  return
301  end subroutine
302 !=======================================================================
309  subroutine cnht_opadd2 (a1,a2,a3,a4,b1,b2,b3,b4)
310  implicit none
311 
312  include 'SIZE' ! N[XYZ]1, NEL[VT]
313  include 'INPUT' ! IFFLOW, IFHEAT, IF3D
314 
315 ! argument list
316  real a1(1),a2(1),a3(1),a4(1),b1(1),b2(1),b3(1),b4(1)
317 
318 ! local variables
319  integer ntotv, ntott
320 !-----------------------------------------------------------------------
321  ntotv = nx1*ny1*nz1*nelv
322  ntott = nx1*ny1*nz1*nelt
323 
324  if (ifflow) then
325  call add2(a1,b1,ntotv)
326  call add2(a2,b2,ntotv)
327  if(if3d) call add2(a3,b3,ntotv)
328  endif
329  if (ifheat) call add2(a4,b4,ntott)
330 
331  return
332  end subroutine
333 !=======================================================================
340  subroutine cnht_opsub2 (a1,a2,a3,a4,b1,b2,b3,b4)
341  implicit none
342 
343  include 'SIZE' ! N[XYZ]1, NEL[VT]
344  include 'INPUT' ! IFFLOW, IFHEAT, IF3D
345 
346 ! argument list
347  real a1(1),a2(1),a3(1),a4(1),b1(1),b2(1),b3(1),b4(1)
348 
349 ! local variables
350  integer ntotv, ntott
351 !-----------------------------------------------------------------------
352  ntotv = nx1*ny1*nz1*nelv
353  ntott = nx1*ny1*nz1*nelt
354 
355  if (ifflow) then
356  call sub2(a1,b1,ntotv)
357  call sub2(a2,b2,ntotv)
358  if(if3d) call sub2(a3,b3,ntotv)
359  endif
360  if (ifheat) call sub2(a4,b4,ntott)
361 
362  return
363  end subroutine
364 !=======================================================================
373  subroutine cnht_opsub3 (a1,a2,a3,a4,b1,b2,b3,b4,c1,c2,c3,c4)
374  implicit none
375 
376  include 'SIZE' ! N[XYZ]1, NEL[VT]
377  include 'INPUT' ! IFFLOW, IFHEAT, IF3D
378 
379 ! argument list
380  real a1(1),a2(1),a3(1),a4(1),b1(1),b2(1),b3(1),b4(1)
381  real c1(1),c2(1),c3(1),c4(1)
382 
383 ! local variables
384  integer ntotv, ntott
385 !-----------------------------------------------------------------------
386  ntotv = nx1*ny1*nz1*nelv
387  ntott = nx1*ny1*nz1*nelt
388 
389  if (ifflow) then
390  call sub3(a1,b1,c1,ntotv)
391  call sub3(a2,b2,c2,ntotv)
392  if(if3d) call sub3(a3,b3,c3,ntotv)
393  endif
394  if (ifheat) call sub3(a4,b4,c4,ntott)
395  return
396  end subroutine
397 !=======================================================================
404  subroutine cnht_opcmult (a1,a2,a3,a4,const)
405  implicit none
406 
407  include 'SIZE' ! N[XYZ]1, NEL[VT]
408  include 'INPUT' ! IFFLOW, IFHEAT, IF3D
409 
410 ! argument list
411  real a1(1),a2(1),a3(1),a4(1)
412  real const
413 
414 ! local variables
415  integer ntotv, ntott
416 !-----------------------------------------------------------------------
417  ntotv = nx1*ny1*nz1*nelv
418  ntott = nx1*ny1*nz1*nelt
419 
420  if (ifflow) then
421  call cmult(a1,const,ntotv)
422  call cmult(a2,const,ntotv)
423  if(if3d) call cmult(a3,const,ntotv)
424  endif
425  if (ifheat) call cmult(a4,const,ntott)
426  return
427  end subroutine
428 !=======================================================================
436  subroutine cnht_opcmult2c (a1,a2,a3,a4,const1, const2)
437  implicit none
438 
439  include 'SIZE' ! N[XYZ]1, NEL[VT]
440  include 'INPUT' ! IFFLOW, IFHEAT, IF3D
441 
442 ! argument list
443  real a1(1),a2(1),a3(1),a4(1)
444  real const1, const2
445 
446 ! local variables
447  integer ntotv, ntott
448 !-----------------------------------------------------------------------
449  ntotv = nx1*ny1*nz1*nelv
450  ntott = nx1*ny1*nz1*nelt
451 
452  if (ifflow) then
453  call cmult(a1,const1,ntotv)
454  call cmult(a2,const1,ntotv)
455  if(if3d) call cmult(a3,const1,ntotv)
456  endif
457  if (ifheat) call cmult(a4,const2,ntott)
458  return
459  end subroutine
460 !=======================================================================
468  subroutine cnht_opadd2cm (a1,a2,a3,a4,b1,b2,b3,b4,coeff)
469  implicit none
470 
471  include 'SIZE' ! N[XYZ]1, NEL[VT]
472  include 'INPUT' ! IFFLOW, IFHEAT, IF3D
473 
474 ! argument list
475  real a1(1),a2(1),a3(1),a4(1),b1(1),b2(1),b3(1),b4(1)
476  real coeff
477 
478 ! local variables
479  integer ntotv, ntott
480  integer il
481 !-----------------------------------------------------------------------
482  ntotv = nx1*ny1*nz1*nelv
483  ntott = nx1*ny1*nz1*nelt
484 
485  if (ifflow) then
486  if (if3d) then
487  do il=1,ntotv
488  a1(il) = a1(il) + b1(il)*coeff
489  a2(il) = a2(il) + b2(il)*coeff
490  a3(il) = a3(il) + b3(il)*coeff
491  enddo
492  else
493  do il=1,ntotv
494  a1(il) = a1(il) + b1(il)*coeff
495  a2(il) = a2(il) + b2(il)*coeff
496  enddo
497  endif
498  endif
499  if (ifheat) then
500  do il=1,ntott
501  a4(il) = a4(il) + b4(il)*coeff
502  enddo
503  endif
504  return
505  end subroutine
506 !=======================================================================
514  subroutine cnht_opsub2cm (a1,a2,a3,a4,b1,b2,b3,b4,coeff)
515  implicit none
516 
517  include 'SIZE' ! N[XYZ]1, NEL[VT]
518  include 'INPUT' ! IFFLOW, IFHEAT, IF3D
519 
520 ! argument list
521  real a1(1),a2(1),a3(1),a4(1),b1(1),b2(1),b3(1),b4(1)
522  real coeff
523 
524 ! local variables
525  integer ntotv, ntott
526  integer il
527 !-----------------------------------------------------------------------
528  ntotv = nx1*ny1*nz1*nelv
529  ntott = nx1*ny1*nz1*nelt
530 
531  if (ifflow) then
532  if (if3d) then
533  do il=1,ntotv
534  a1(il) = a1(il) - b1(il)*coeff
535  a2(il) = a2(il) - b2(il)*coeff
536  a3(il) = a3(il) - b3(il)*coeff
537  enddo
538  else
539  do il=1,ntotv
540  a1(il) = a1(il) - b1(il)*coeff
541  a2(il) = a2(il) - b2(il)*coeff
542  enddo
543  endif
544  endif
545  if (ifheat) then
546  do il=1,ntott
547  a4(il) = a4(il) - b4(il)*coeff
548  enddo
549  endif
550  return
551  end subroutine
552 !=======================================================================
558  subroutine cnht_weight_fun (lvx,lvy,lvz,lt,coeff)
559  implicit none
560 
561  include 'SIZE' ! N[XYZ]1, NEL[VT]
562  include 'MASS' ! VOLVM1, VOLTM1
563  include 'CNHTD' ! cnht_sc,cnht_sv,cnht_st
564 
565  ! argument list
566  real lvx(1),lvy(1),lvz(1),lt(1)
567  real coeff
568 
569  ! local variables
570  real f1, f2
571 !-----------------------------------------------------------------------
572  f1=cnht_sv/volvm1/coeff
573  f2=cnht_st*cnht_sc/voltm1/coeff
574 
575  !rescale
576  call cnht_opcmult2c (lvx,lvy,lvz,lt,f1,f2)
577 
578  return
579  end subroutine
580 !=======================================================================
589  real function cnht_glsc2_wt (b1,b2,b3,b4,x1,x2,x3,x4,wt)
590  implicit none
591 
592  include 'SIZE' ! N[XYZ]1, NEL[VT]
593  include 'INPUT' ! IFFLOW, IFHEAT, IF3D
594  include 'MASS' ! VOLVM1, VOLTM1
595  include 'CNHTD' ! cnht_sc,cnht_sv,cnht_st
596 
597 ! argument list
598  real b1(1),b2(1),b3(1),b4(1),x1(1),x2(1),x3(1),x4(1),wt(1)
599 
600 ! local variables
601  integer ntotv, ntott
602  real sum, f1, f2
603  integer il
604 ! functions
605  real glsum
606 !-----------------------------------------------------------------------
607  ntotv = nx1*ny1*nz1*nelv
608  ntott = nx1*ny1*nz1*nelt
609 
610  ! scaling factor velocity vs temperature
611  ! veorsion for newton
612  ! f1 = coeff_v
613  ! f2 = coeff_T
614  ! version for oic
615  f1=cnht_sv/volvm1
616  f2=cnht_st*cnht_sc/voltm1
617 
618  sum = 0.
619  if (ifflow) then !if vel
620  if (ifheat) then !if temp & vel
621  if (if3d) then
622  do il=1,ntotv
623  sum = sum + wt(il)*(f1*(b1(il)*x1(il)+b2(il)*x2(il)
624  & +b3(il)*x3(il))+f2*b4(il)*x4(il))
625  end do
626  else
627  do il=1,ntotv
628  sum =sum + wt(il)*(f1*(b1(il)*x1(il)+b2(il)*x2(il))
629  & +f2*b4(il)*x4(il))
630  end do
631  end if
632 
633  ! for conjugate heat transfer
634  if (ntott.gt.ntotv) then
635  do il=ntotv+1,ntott
636  sum = sum + wt(il)*f2*b4(il)*x4(il)
637  end do
638  end if
639  else !just vel
640  if (if3d) then
641  do il=1,ntotv
642  sum = sum + wt(il)*f1*(b1(il)*x1(il)+
643  $ b2(il)*x2(il)+b3(il)*x3(il))
644  end do
645  else
646  do il=1,ntotv
647  sum = sum + wt(il)*f1*(b1(il)*x1(il)+b2(il)*x2(il))
648  end do
649  end if
650  end if
651  else !just temp
652  if (ifheat) then
653  do il=1,ntott
654  sum = sum + wt(il)*(f2*b4(il)*x4(il))
655  end do
656  end if
657  end if
658 
659  cnht_glsc2_wt = glsum(sum,1)
660 
661  return
662  end function
663 !=======================================================================
integer function gllel(ieg)
Definition: dprocmap.f:183
subroutine cnht_opsub2cm(a1, a2, a3, a4, b1, b2, b3, b4, coeff)
Vector subtraction with scaling A = A-c*B (velocity and temperature)
Definition: cnht_tools.f:515
subroutine cnht_opcopy(a1, a2, a3, a4, b1, b2, b3, b4)
Copy vectors A=B (velocity and temperature)
Definition: cnht_tools.f:279
subroutine cnht_opsub3(a1, a2, a3, a4, b1, b2, b3, b4, c1, c2, c3, c4)
Subtract vectors A = B-C (velocity and temperature)
Definition: cnht_tools.f:374
subroutine cnht_register()
Register conjugated heat transfer tools module.
Definition: cnht_tools.f:12
subroutine cnht_opcmult2c(a1, a2, a3, a4, const1, const2)
Multiply vector by constant A = c*A with separate const. for velocity and temperature.
Definition: cnht_tools.f:437
subroutine cnht_opcmult(a1, a2, a3, a4, const)
Multiply vector by constant A = c*A (single coeff. for velocity and temperature)
Definition: cnht_tools.f:405
subroutine cnht_opsub2(a1, a2, a3, a4, b1, b2, b3, b4)
Subtract vectors A = A-B (velocity and temperature)
Definition: cnht_tools.f:341
subroutine cnht_opadd2(a1, a2, a3, a4, b1, b2, b3, b4)
Add velocity and temperature vectors A = A+B.
Definition: cnht_tools.f:310
subroutine cnht_opadd2cm(a1, a2, a3, a4, b1, b2, b3, b4, coeff)
Vector summation with scaling A = A+c*B (velocity and temperature)
Definition: cnht_tools.f:469
subroutine cnht_cpfld_set()
Set cpfld coefficient for given type of simulation.
Definition: cnht_tools.f:206
subroutine cnht_init()
Initilise conjugated heat transfer tools module.
Definition: cnht_tools.f:82
subroutine cnht_oprzero(a1, a2, a3, a4)
Zero velocity and temperature vectors.
Definition: cnht_tools.f:249
real function cnht_glsc2_wt(b1, b2, b3, b4, x1, x2, x3, x4, wt)
Global inner product of velocity and temperature fields.
Definition: cnht_tools.f:590
subroutine cnht_weight_fun(lvx, lvy, lvz, lt, coeff)
Weigth velocity and temperature fields.
Definition: cnht_tools.f:559
subroutine cnht_forcing(ffx, ffy, ffz, ix, iy, iz, ieg)
Calcualte forcing ralted to conjugated heat transfer.
Definition: cnht_tools.f:160
logical function cnht_is_initialised()
Check if module was initialised.
Definition: cnht_tools.f:144
subroutine mntr_warn(mid, logs)
Write warning message.
Definition: mntrlog.f:803
subroutine mntr_mod_is_name_reg(mid, mname)
Check if module name is registered and return its id.
Definition: mntrlog.f:459
subroutine mntr_abort(mid, logs)
Abort simulation.
Definition: mntrlog.f:837
subroutine mntr_mod_reg(mid, pmid, mname, mdscr)
Register new module.
Definition: mntrlog.f:346
subroutine rprm_rp_get(ipval, rpval, lpval, cpval, rpid, ptype)
Get runtime parameter form active section. This operation is performed locally.
Definition: rprm.f:883
subroutine rprm_rp_reg(rpid, mid, pname, pdscr, ptype, ipval, rpval, lpval, cpval)
Register new runtime parameter.
Definition: rprm.f:484
subroutine rprm_sec_set_act(ifact, rpid)
Set section's activation flag. Master value is broadcasted.
Definition: rprm.f:422
subroutine rprm_sec_reg(rpid, mid, pname, pdscr)
Register new parameter section.
Definition: rprm.f:165
subroutine sub2(a, b, n)
Definition: math.f:164
subroutine add2(a, b, n)
Definition: math.f:622
subroutine copy(a, b, n)
Definition: math.f:260
function glsum(x, n)
Definition: math.f:861
subroutine sub3(a, b, c, n)
Definition: math.f:175
subroutine cmult(a, const, n)
Definition: math.f:315
subroutine rzero(a, n)
Definition: math.f:208
subroutine gradm1(ux, uy, uz, u)
Definition: navier5.f:463